Top | ![]() |
![]() |
![]() |
![]() |
#define G_PARAM_STATIC_STRINGS (G_PARAM_STATIC_NAME | G_PARAM_STATIC_NICK | G_PARAM_STATIC_BLURB)
void
gstlal_fftw_lock (void
);
Aquire the lock to protect the global shared state in the FFTW wisdom. This function is a proxy for LAL's FFTW wisdom lock.
See also: gstlal_fftw_unlock()
void
gstlal_fftw_unlock (void
);
Release the lock to protect the global shared state in the FFTW wisdom. This function is a proxy for LAL's FFTW wisdom unlock.
See also: gstlal_fftw_lock()
void
gstlal_load_fftw_wisdom (void
);
Attempt to load double-precision and single-precision FFTW wisdom files. The names for these files are taken from the environment variables GSTLAL_FFTW_WISDOM and GSTLAL_FFTWF_WISDOM, respectively, or if one or the other isn't set then the respective FFTW default is used. This function acquires and releases the GstLAL FFTW locks and is thread-safe.
GValueArray * gstlal_g_value_array_from_ints (const gint *src
,gint n
);
Build a GValueArray from an array of ints.
gint * gstlal_ints_from_g_value_array (GValueArray *va
,gint *dest
,gint *n
);
Convert a GValueArray of ints to an array of ints. If dest
is NULL
then new memory will be allocated otherwise the ints are copied into the
memory pointed to by dest
, which must be large enough to hold them. If
memory is allocated by this function, free with g_free()
. If n
is not
NULL
it is set to the number of elements in the array.
va |
the GValueArray from which to copy the elements |
|
dest |
address of memory large enough to hold elements, or |
|
n |
address of integer that will be set to the number of elements, or
|
gsl_vector_int *
gstlal_gsl_vector_int_from_g_value_array
(GValueArray *va
);
Build a gsl_vector_int from a GValueArray of ints.
GValueArray *
gstlal_g_value_array_from_gsl_vector_int
(const gsl_vector_int *vector
);
Build a GValueArray of ints from a gsl_vector_int.
gsl_matrix_int *
gstlal_gsl_matrix_int_from_g_value_array
(GValueArray *va
);
Build a gsl_matrix_int from a GValueArray of GValueArrays of ints.
GValueArray *
gstlal_g_value_array_from_gsl_matrix_int
(const gsl_matrix_int *matrix
);
Build a GValueArray of GValueArrays of ints from a gsl_matrix_int.
guint64 * gstlal_uint64s_from_g_value_array (GValueArray *va
,guint64 *dest
,gint *n
);
Convert a GValueArray of guint64 to an array of guint64. If dest
is
NULL
then new memory will be allocated otherwise the doubles are copied
into the memory pointed to by dest
, which must be large enough to hold
them. If memory is allocated by this function, free with g_free()
. If
n
is not NULL
it is set to the number of elements in the array.
va |
the GValueArray from which to copy the elements |
|
dest |
address of memory large enough to hold elements, or |
|
n |
address of integer that will be set to the number of elements, or
|
gsl_matrix_ulong *
gstlal_gsl_matrix_ulong_from_g_value_array
(GValueArray *va
);
Build a gsl_matrix_ulong from a GValueArray of GValueArrays of guint64.
GValueArray *
gstlal_g_value_array_from_gsl_matrix_ulong
(const gsl_matrix_ulong *matrix
);
Build a GValueArray of GValueArrays of guin64 from a gsl_matrix_ulong.
the newly-allocated GValueArray of newly-allocated
GValueArrays of guint64s or NULL
on failure.
GValueArray * gstlal_g_value_array_from_uint64s (const guint64 *src
,gint n
);
Build a GValueArray from an array of guint64.
GValueArray * gstlal_g_value_array_from_doubles (const gdouble *src
,gint n
);
Build a GValueArray from an array of doubles.
gdouble * gstlal_doubles_from_g_value_array (GValueArray *va
,gdouble *dest
,gint *n
);
Convert a GValueArray of doubles to an array of doubles. If dest
is
NULL
then new memory will be allocated otherwise the doubles are copied
into the memory pointed to by dest
, which must be large enough to hold
them. If memory is allocated by this function, free with g_free()
. If
n
is not NULL
it is set to the number of elements in the array.
va |
the GValueArray from which to copy the elements |
|
dest |
address of memory large enough to hold elements, or |
|
n |
address of integer that will be set to the number of elements, or
|
gsl_vector *
gstlal_gsl_vector_from_g_value_array (GValueArray *va
);
Build a gsl_vector from a GValueArray of doubles.
GValueArray *
gstlal_g_value_array_from_gsl_vector (const gsl_vector *vector
);
Build a GValueArray of doubles from a gsl_vector.
gsl_matrix *
gstlal_gsl_matrix_from_g_value_array (GValueArray *va
);
Build a gsl_matrix from a GValueArray of GValueArrays of doubles.
GValueArray *
gstlal_g_value_array_from_gsl_matrix (const gsl_matrix *matrix
);
Build a GValueArray of GValueArrays of doubles from a gsl_matrix.
the newly-allocated GValueArray of newly-allocated
GValueArrays of doubles or NULL
on failure.
gsl_vector_complex *
gstlal_gsl_vector_complex_from_g_value_array
(GValueArray *va
);
Build a gsl_vector_complex from a GValueArray of doubles packed as real,imag,real,imag,...
GValueArray *
gstlal_g_value_array_from_gsl_vector_complex
(const gsl_vector_complex *vector
);
Build a GValueArray of doubles from a gsl_vector_complex. The complex numbers are packed into the result as real,imag,real,imag,...
gsl_matrix_complex *
gstlal_gsl_matrix_complex_from_g_value_array
(GValueArray *va
);
Build a gsl_matrix_complex from a GValueArray of GValueArrays of doubles. The complex numbers are unpacked from the doubles as real,imag,real,imag,...
GValueArray *
gstlal_g_value_array_from_gsl_matrix_complex
(const gsl_matrix_complex *matrix
);
Build a GValueArray of GValueArrays of doubles from a gsl_matrix_complex. The complex numbers are packed into the doubles as real,imag,real,imag,...
the newly-allocated GValueArray of newly-allocated
GValueArrays of doubles or NULL
on failure.
char * gstlal_build_full_channel_name (const char *instrument
,const char *channel_name
);
Prefix a channel name with the instrument name. I.e., turn "H1" and
"LSC-STRAIN" into "H1:LSC-STRAIN". If either instrument or channel_name
is NULL
, then the corresponding part of the result is left blank and
the colon is omited. If both are NULL
an empty string is returned.
REAL8TimeSeries * gstlal_REAL8TimeSeries_from_buffer (GstBuffer *buf
,const char *instrument
,const char *channel_name
,const char *units
);
Wrap a GstBuffer in a REAL8TimeSeries. The time series's data->data
pointer points to the GstBuffer's own data, so this pointer must not be
free()
ed --- it must be set to NULL
before passing the time series to
XLALDestroyREAL8TimeSeries()
. Only GstBuffers containing
single-channel time series data are supported. The instrument
and
channel_name
are used to build the REAL8TimeSeries' name using
gstlal_build_full_channel_name()
.
Example 1. Create a REAL8TimeSeries and free it
1 2 3 4 5 6 7 8 9 10 |
REAL8TimeSeries *series; series = gstlal_REAL8TimeSeries_from_buffer(buf, "H1", "LSC-STRAIN", "strain"); if(!series) handle_error(); blah_blah_blah(); series->data->data = NULL; XLALDestroyREAL8TimeSeries(series); |
newly-allocated REAL8TimeSeries. Free with
XLALDestroyREAL8TimeSeries()
after setting the internal data->data
pointer to NULL
.
GstDateTime *
gstlal_datetime_new_from_gps (GstClockTime gps
);
Converts a GstClockTime containing a GPS time (in integer nanoseconds) to a GstDateTime object representing the same time.