Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALInference 4.1.9.1-5e288d3
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
LALInferenceGenerateROQ.h File Reference

Prototypes

REAL8 LALInferenceGenerateREAL8OrthonormalBasis (REAL8Array **RB, const REAL8Vector *delta, REAL8 tolerance, REAL8Array **TS, UINT4Vector **greedypoints)
 Create a orthonormal basis set from a training set of real waveforms. More...
 
REAL8 LALInferenceGenerateCOMPLEX16OrthonormalBasis (COMPLEX16Array **RB, const REAL8Vector *delta, REAL8 tolerance, COMPLEX16Array **TS, UINT4Vector **greedypoints)
 Create a orthonormal basis set from a training set of complex waveforms. More...
 
void LALInferenceValidateREAL8OrthonormalBasis (REAL8Vector **projerr, const REAL8Vector *delta, const REAL8Array *RB, REAL8Array **testmodels)
 Validate the real reduced basis against another set of waveforms. More...
 
void LALInferenceValidateCOMPLEX16OrthonormalBasis (REAL8Vector **projerr, const REAL8Vector *delta, const COMPLEX16Array *RB, COMPLEX16Array **testmodels)
 Validate the complex reduced basis against another set of waveforms. More...
 
INT4 LALInferenceTestREAL8OrthonormalBasis (const REAL8Vector *delta, REAL8 tolerance, const REAL8Array *RB, REAL8Array **testmodels)
 Test the reduced basis against another set of waveforms. More...
 
INT4 LALInferenceTestCOMPLEX16OrthonormalBasis (const REAL8Vector *delta, REAL8 tolerance, const COMPLEX16Array *RB, COMPLEX16Array **testmodels)
 Test the reduced basis against another set of waveforms. More...
 
REAL8 LALInferenceEnrichREAL8Basis (const REAL8Vector *delta, const REAL8 tolerance, REAL8Array **RB, UINT4Vector **greedypoints, const REAL8Array *testmodels, REAL8Array **testmodelsnew)
 Expand the real training waveforms with ones from a set of new training waveforms. More...
 
REAL8 LALInferenceEnrichCOMPLEX16Basis (const REAL8Vector *delta, const REAL8 tolerance, COMPLEX16Array **RB, UINT4Vector **greedypoints, const COMPLEX16Array *testmodels, COMPLEX16Array **testmodelsnew)
 Expand the complex training waveforms with ones from a set of new training waveforms. More...
 
LALInferenceREALROQInterpolantLALInferenceGenerateREALROQInterpolant (REAL8Array *RB)
 Create a real empirical interpolant from a set of orthonormal basis functions. More...
 
LALInferenceCOMPLEXROQInterpolantLALInferenceGenerateCOMPLEXROQInterpolant (COMPLEX16Array *RB)
 Create a complex empirical interpolant from a set of orthonormal basis functions. More...
 
REAL8VectorLALInferenceGenerateREAL8LinearWeights (REAL8Array *B, REAL8Vector *data, REAL8Vector *vars)
 Create the weights for the ROQ interpolant for the linear data and model dot product <d|h> More...
 
REAL8VectorLALInferenceGenerateQuadraticWeights (REAL8Array *B, REAL8Vector *vars)
 Create the weights for the ROQ interpolant for the model quadratic model term real(<h|h>). More...
 
COMPLEX16VectorLALInferenceGenerateCOMPLEX16LinearWeights (COMPLEX16Array *B, COMPLEX16Vector *data, REAL8Vector *vars)
 Create the weights for the ROQ interpolant for the complex data and model dot product <d|h> More...
 
REAL8 LALInferenceROQREAL8DotProduct (REAL8Vector *weights, REAL8Vector *model)
 Calculate the dot product of two vectors using the ROQ iterpolant. More...
 
COMPLEX16 LALInferenceROQCOMPLEX16DotProduct (COMPLEX16Vector *weights, COMPLEX16Vector *model)
 Calculate the dot product of two complex vectors using the ROQ iterpolant. More...
 
void LALInferenceRemoveREALROQInterpolant (LALInferenceREALROQInterpolant *a)
 Free memory for a LALInferenceREALROQInterpolant. More...
 
void LALInferenceRemoveCOMPLEXROQInterpolant (LALInferenceCOMPLEXROQInterpolant *a)
 Free memory for a LALInferenceCOMPLEXROQInterpolant. More...
 

Go to the source code of this file.

Data Structures

struct  LALInferenceREALROQInterpolant
 A structure to hold a real (double precision) interpolant matrix and interpolation node indices. More...
 
struct  LALInferenceCOMPLEXROQInterpolant
 A structure to hold a complex (double precision) interpolant matrix and interpolation node indices. More...
 

Function Documentation

◆ LALInferenceGenerateREAL8OrthonormalBasis()

REAL8 LALInferenceGenerateREAL8OrthonormalBasis ( REAL8Array **  RBin,
const REAL8Vector delta,
REAL8  tolerance,
REAL8Array **  TS,
UINT4Vector **  greedypoints 
)

Create a orthonormal basis set from a training set of real waveforms.

Given a gsl_matrix containing a training set of real waveforms (where the waveforms are created at time or frequency steps seperated by delta) an orthonormal basis will be generated using the greedy binning Algorithm 1 of [1] . The stopping criteria for the algorithm is controlled by the tolerance value, which defined the maximum residual between the current basis set (at a given iteration) and the training set (and example tolerance is \(10^{-12}\)). In this function the training set will be normalised, so the input TS will be modified.

If the RBin value is NULL then a new reduced basis will be formed from the given training set. However, if RBin already contains a previously produced basis, then this basis will be enriched with bases if possible using the new training set. NOTE: when using small tolerances enriching the basis in this way can lead to numerical precision issues, so in general you should use LALInferenceEnrichREAL8Basis for enrichment.

Parameters
[out]RBinA REAL8Array to return the reduced basis.
[in]deltaThe time/frequency step(s) in the training set used to normalise the models. This can be a vector containing just one value.
[in]toleranceThe tolerance used as a stopping criteria for the basis generation.
[in]TSA REAL8Array matrix containing the complex training set, where the number of waveforms in the training set is given by the rows and the waveform points by the columns. This will be modified by this function, as the training set will be normalised.
[out]greedypointsA UINT4Vector to return the indices of the training set rows that have been used to form the reduced basis.
Returns
A REAL8 with the maximum projection error for the final reduced basis.
See also
LALInferenceEnrichREAL8Basis

Definition at line 762 of file LALInferenceGenerateROQ.c.

◆ LALInferenceGenerateCOMPLEX16OrthonormalBasis()

REAL8 LALInferenceGenerateCOMPLEX16OrthonormalBasis ( COMPLEX16Array **  RBin,
const REAL8Vector delta,
REAL8  tolerance,
COMPLEX16Array **  TS,
UINT4Vector **  greedypoints 
)

Create a orthonormal basis set from a training set of complex waveforms.

Given a gsl_matrix containing a training set TS of complex waveforms (where the waveforms are created at time or frequency steps seperated by delta) an orthonormal basis will be generated using the greedy binning Algorithm 1 of [1] . The stopping criteria for the algorithm is controlled by the tolerance value, which is defined the maximum residual between the current basis set (at a given iteration) and the training set (and example tolerance is \(10^{-12}\)). In this function the training set will be normalised, so the input TS will be modified.

If the RBin value is NULL then a new reduced basis will be formed from the given training set. However, if RBin already contains a previously produced basis, then this basis will be enriched with bases if possible using the new training set. NOTE: when using small tolerances enriching the basis in this way can lead to numerical precision issues, so in general you should use LALInferenceEnrichCOMPLEX16Basis for enrichment.

Note that in this function we have to cast the COMPLEX16 array as a double to use gsl_matrix_view_array, which assume that the data is passed as a double array with memory laid out so that adjacent double memory blocks hold the corresponding real and imaginary parts.

Parameters
[out]RBinA COMPLEX16Array to return the reduced basis.
[in]deltaThe time/frequency step(s) in the training set used to normalise the models. This can be a vector containing just one value.
[in]toleranceThe tolerance used as a stopping criteria for the basis generation.
[in]TSA COMPLEX16Array matrix containing the complex training set, where the number of waveforms in the training set is given by the rows and the waveform points by the columns. This will be modified by this function, as the training set will be normalised.
[out]greedypointsA UINT4Vector to return the indices of the training set rows that have been used to form the reduced basis.
Returns
A REAL8 with the maximum projection error for the final reduced basis.
See also
LALInferenceEnrichCOMPLEX16Basis

Definition at line 968 of file LALInferenceGenerateROQ.c.

◆ LALInferenceValidateREAL8OrthonormalBasis()

void LALInferenceValidateREAL8OrthonormalBasis ( REAL8Vector **  projerr,
const REAL8Vector delta,
const REAL8Array RB,
REAL8Array **  testmodels 
)

Validate the real reduced basis against another set of waveforms.

This function projects a set of waveforms onto the reduced basis and checks that the residuals are within a given tolerance. It returns the projection errors.

Note that the projection error returned are the square of the residual errors, as is used as the criterion for adding new bases in the reduced basis generation function LALInferenceGenerateREAL8OrthonormalBasis. This is different to the validation.cpp code in greedycpp, which returns the square root of the value we return.

Parameters
[out]projerrThe projection errors for each test waveform.
[in]deltaThe time/frequency step(s) in the training set used to normalise the models. This can be a vector containing just one value.
[in]RBThe reduced basis set
[in]testmodelsThe set of waveform models to project onto the basis (these will be changed by this function, as the waveforms will get normalised).
See also
LALInferenceTestREAL8OrthonormalBasis

Definition at line 1159 of file LALInferenceGenerateROQ.c.

◆ LALInferenceValidateCOMPLEX16OrthonormalBasis()

void LALInferenceValidateCOMPLEX16OrthonormalBasis ( REAL8Vector **  projerr,
const REAL8Vector delta,
const COMPLEX16Array RB,
COMPLEX16Array **  testmodels 
)

Validate the complex reduced basis against another set of waveforms.

This function projects a set of waveforms onto the reduced basis and checks that the residuals are within a given tolerance. It returns the projection errors.

Note that the projection error returned are the square of the residual errors, as is used as the criterion for adding new bases in the reduced basis generation function LALInferenceGenerateCOMPLEX16OrthonormalBasis. This is different to the validation.cpp code in greedycpp, which returns the square root of the value we return.

Parameters
[out]projerrThe projection errors (square of the residual) for each test waveform.
[in]deltaThe time/frequency step(s) in the training set used to normalise the models. This can be a vector containing just one value.
[in]RBThe reduced basis set
[in]testmodelsThe set of waveform models to project onto the basis (these will be changed by this function, as the waveforms will get normalised).
See also
LALInferenceTestCOMPLEX16OrthonormalBasis

Definition at line 1416 of file LALInferenceGenerateROQ.c.

◆ LALInferenceTestREAL8OrthonormalBasis()

INT4 LALInferenceTestREAL8OrthonormalBasis ( const REAL8Vector delta,
REAL8  tolerance,
const REAL8Array RB,
REAL8Array **  testmodels 
)

Test the reduced basis against another set of waveforms.

This function projects a set of waveforms onto the reduced basis and checks that the residuals are within a given tolerance

Parameters
[in]deltaThe time/frequency step(s) in the training set used to normalise the models. This can be a vector containing just one value.
[in]toleranceThe allowed (squared) residual tolerence for the test waveforms
[in]RBThe reduced basis set
[in]testmodelsThe set of waveform models to project onto the basis
Returns
Returns XLAL_SUCCESS if all test waveforms meet the tolerance
See also
LALInferenceValidateREAL8OrthonormalBasis

Definition at line 1240 of file LALInferenceGenerateROQ.c.

◆ LALInferenceTestCOMPLEX16OrthonormalBasis()

INT4 LALInferenceTestCOMPLEX16OrthonormalBasis ( const REAL8Vector delta,
REAL8  tolerance,
const COMPLEX16Array RB,
COMPLEX16Array **  testmodels 
)

Test the reduced basis against another set of waveforms.

This function projects a set of waveforms onto the reduced basis and checks that the residuals are within a given tolerance

Parameters
[in]deltaThe time/frequency step(s) in the training set used to normalise the models. This can be a vector containing just one value.
[in]toleranceThe allowed residual (squared) tolerence for the test waveforms
[in]RBThe reduced basis set
[in]testmodelsThe set of waveform models to project onto the basis
Returns
Returns XLAL_SUCCESS if all test waveforms meet the tolerance
See also
LALInferenceValidateCOMPLEX16OrthonormalBasis

Definition at line 1508 of file LALInferenceGenerateROQ.c.

◆ LALInferenceEnrichREAL8Basis()

REAL8 LALInferenceEnrichREAL8Basis ( const REAL8Vector delta,
REAL8  tolerance,
REAL8Array **  RB,
UINT4Vector **  greedypoints,
const REAL8Array testmodels,
REAL8Array **  testmodelsnew 
)

Expand the real training waveforms with ones from a set of new training waveforms.

Expands an original set of training waveforms with waveforms from a new set a training waveforms if, when projected onto the reduced basis, those new waveforms are greater than the given tolerance. The reduced basis will then be recomputed using the expanded training set.

Parameters
[in]deltaThe time/frequency step(s) in the training set used to normalise the models. This can be a vector containing just one value.
[in]toleranceThe allowed residual tolerence for the test waveforms
[in]RBThe reduced basis set
[in]greedypointsThe rows in testmodels that formed the reduced basis
[in]testmodelsThe set of waveform models used to produce the reduced basis (already normalised)
[in]testmodelsnewA new set of waveforms to project onto the current reduced basis (these will be changed by this function, as the waveforms will get normalised when pass to LALInferenceValidateREAL8OrthonormalBasis, and editted to contain the new set of training waveforms from which the enriched basis was formed).
Returns
Returns the maximum projection error of the new basis set

Definition at line 1291 of file LALInferenceGenerateROQ.c.

◆ LALInferenceEnrichCOMPLEX16Basis()

REAL8 LALInferenceEnrichCOMPLEX16Basis ( const REAL8Vector delta,
REAL8  tolerance,
COMPLEX16Array **  RB,
UINT4Vector **  greedypoints,
const COMPLEX16Array testmodels,
COMPLEX16Array **  testmodelsnew 
)

Expand the complex training waveforms with ones from a set of new training waveforms.

Expands an original set of training waveforms with waveforms from a new set a training waveforms if, when projected onto the reduced basis, those new waveforms are greater than the given tolerance. The reduced basis will then be recomputed using the expanded training set.

Parameters
[in]deltaThe time/frequency step(s) in the training set used to normalise the models. This can be a vector containing just one value.
[in]toleranceThe allowed residual tolerence for the test waveforms
[in]RBThe reduced basis set
[in]greedypointsThe rows in testmodels that formed the reduced basis
[in]testmodelsThe set of waveform models used to produce the reduced basis
[in]testmodelsnewA new set of waveforms to project onto the current reduced basis (these will be changed by this function, as the waveforms will get normalised when pass to LALInferenceValidateCOMPLEX16OrthonormalBasis, and editted to contain the new set of training waveforms from which the enriched basis was formed)
Returns
Returns the maximum projection error of the new basis set

Definition at line 1560 of file LALInferenceGenerateROQ.c.

◆ LALInferenceGenerateREALROQInterpolant()

LALInferenceREALROQInterpolant * LALInferenceGenerateREALROQInterpolant ( REAL8Array RB)

Create a real empirical interpolant from a set of orthonormal basis functions.

Given a real REAL8Array containing a set of orthonormal basis functions generate an empirical intopolant, and set of interpolation points, using Algorithm 2 of [1] .

Parameters
[in]RBThe set of basis functions
Returns
A LALInferenceREALROQInterpolant structure containing the interpolant and its nodes

Definition at line 1679 of file LALInferenceGenerateROQ.c.

◆ LALInferenceGenerateCOMPLEXROQInterpolant()

LALInferenceCOMPLEXROQInterpolant * LALInferenceGenerateCOMPLEXROQInterpolant ( COMPLEX16Array RB)

Create a complex empirical interpolant from a set of orthonormal basis functions.

Given a complex gsl_matrix_complex containing a set of orthonormal basis functions generate an empirical intopolant, and set of interpolation points, using Algorithm 2 of [1] .

Parameters
[in]RBThe set of basis functions
Returns
A LALInferenceCOMPLEXROQInterpolant structure containing the interpolant and its nodes

Definition at line 1771 of file LALInferenceGenerateROQ.c.

◆ LALInferenceGenerateREAL8LinearWeights()

REAL8Vector * LALInferenceGenerateREAL8LinearWeights ( REAL8Array B,
REAL8Vector data,
REAL8Vector vars 
)

Create the weights for the ROQ interpolant for the linear data and model dot product <d|h>

Parameters
[in]BThe interpolant matrix
[in]dataThe real data vector
[in]varsA vector of data noise variance values (or a single value) to weight the "weights"
Returns
The vector of weights

Definition at line 1860 of file LALInferenceGenerateROQ.c.

◆ LALInferenceGenerateQuadraticWeights()

REAL8Vector * LALInferenceGenerateQuadraticWeights ( REAL8Array B,
REAL8Vector vars 
)

Create the weights for the ROQ interpolant for the model quadratic model term real(<h|h>).

Parameters
[in]BThe interpolant matrix
[in]varsA vector of data noise variance values (or a single value) to weight the "weights"
Returns
The vector of weights

Definition at line 1899 of file LALInferenceGenerateROQ.c.

◆ LALInferenceGenerateCOMPLEX16LinearWeights()

COMPLEX16Vector * LALInferenceGenerateCOMPLEX16LinearWeights ( COMPLEX16Array B,
COMPLEX16Vector data,
REAL8Vector vars 
)

Create the weights for the ROQ interpolant for the complex data and model dot product <d|h>

Parameters
[in]BThe interpolant matrix
[in]dataThe complex data vector
[in]varsA vector of data noise variance values (or a single value) to weight the "weights"
Returns
The vector of weights

Definition at line 1936 of file LALInferenceGenerateROQ.c.

◆ LALInferenceROQREAL8DotProduct()

REAL8 LALInferenceROQREAL8DotProduct ( REAL8Vector weights,
REAL8Vector model 
)

Calculate the dot product of two vectors using the ROQ iterpolant.

This function calculates the dot product of two real vectors using the ROQ interpolant. This required the interpolant weights computed with LALInferenceCreateREAL8DataModelWeights and the waveform model defined at the interpolation node.

Parameters
[in]weightsThe ROQ interpolation weights
[in]modelThe waveform model (or quadratic model) defined at the interpolation points
Returns
The dot product weights and model

Definition at line 1982 of file LALInferenceGenerateROQ.c.

◆ LALInferenceROQCOMPLEX16DotProduct()

COMPLEX16 LALInferenceROQCOMPLEX16DotProduct ( COMPLEX16Vector weights,
COMPLEX16Vector model 
)

Calculate the dot product of two complex vectors using the ROQ iterpolant.

This function calculates the dot product of two complex vectors using the ROQ interpolant. This required the interpolant weights computed with LALInferenceCreateCOMPLEX16DataModelWeights and the waveform model defined at the interolation node.

Parameters
[in]weightsThe ROQ interpolation weights
[in]modelThe waveform model defined at the interpolation points
Returns
The dot product of the data and model

Definition at line 2005 of file LALInferenceGenerateROQ.c.

◆ LALInferenceRemoveREALROQInterpolant()

void LALInferenceRemoveREALROQInterpolant ( LALInferenceREALROQInterpolant a)

Free memory for a LALInferenceREALROQInterpolant.

Parameters
[in]aA pointer to a LALInferenceREALROQInterpolant

Definition at line 2022 of file LALInferenceGenerateROQ.c.

◆ LALInferenceRemoveCOMPLEXROQInterpolant()

void LALInferenceRemoveCOMPLEXROQInterpolant ( LALInferenceCOMPLEXROQInterpolant a)

Free memory for a LALInferenceCOMPLEXROQInterpolant.

Parameters
[in]aA pointer to a LALInferenceCOMPLEXROQInterpolant

Definition at line 2036 of file LALInferenceGenerateROQ.c.