LALInference  4.1.6.1-b72065a
Header LALInference.h

Detailed Description

Main header file for LALInference common routines and structures.

LALInference is a Bayesian analysis toolkit for use with LAL. It contains common requirements for Bayesian codes such as Likelihood functions, data handling routines, MCMC and Nested Sampling algorithms and a template generation interface to the LALInspiral package.

This file contains the basic structures for the algorithm state, interferometer data, manipulation of variables and type declarations for the standard function types.

Prototypes

char ** LALInferenceGetHeaderLine (FILE *inp)
 Returns an array of header strings (terminated by NULL) from a common-format output file. More...
 
const charLALInferenceTranslateInternalToExternalParamName (const char *inName)
 Converts between internally used parameter names and those external (e.g. More...
 
void LALInferenceTranslateExternalToInternalParamName (char *outName, const char *inName)
 Converts between externally used parameter names and those internal. More...
 
int LALInferenceFprintParameterHeaders (FILE *out, LALInferenceVariables *params)
 Print the parameter names to a file as a tab-separated ASCII line. More...
 
INT4 LALInferenceFprintParameterNonFixedHeaders (FILE *out, LALInferenceVariables *params)
 Print the parameters which do not vary to a file as a tab-separated ASCII line. More...
 
INT4 LALInferenceFprintParameterNonFixedHeadersWithSuffix (FILE *out, LALInferenceVariables *params, const char *suffix)
 Print the parameters which do not vary to a file as a tab-separated ASCII line, adding the given suffix. More...
 
UINT4 LALInferencePrintNVariableItem (char *out, UINT4 N, const LALInferenceVariableItem *const ptr)
 Prints a variable item to a string. More...
 
void * LALInferenceGetVariable (const LALInferenceVariables *vars, const char *name)
 Return a pointer to the memory the variable vars is stored in specified by name User must cast this pointer to the expected type before dereferencing it to get the value of the variable. More...
 
INT4 LALInferenceGetVariableDimension (LALInferenceVariables *vars)
 Get number of dimensions in variable vars. More...
 
INT4 LALInferenceGetVariableDimensionNonFixed (LALInferenceVariables *vars)
 Get number of dimensions in vars which are not fixed to a certain value. More...
 
INT4 LALInferenceGetVariableDimensionNonFixedChooseVectors (LALInferenceVariables *vars, INT4 count_vectors)
 Get number of dimensions in vars which are not fixed to a certain value, with a flag for skipping counting vectors. More...
 
INT4 LALInferenceGetVariableTypeByIndex (LALInferenceVariables *vars, int idx)
 Get the LALInferenceVariableType of the idx -th item in the vars Indexing starts at 1. More...
 
LALInferenceVariableType LALInferenceGetVariableType (const LALInferenceVariables *vars, const char *name)
 Get the LALInferenceVariableType of the parameter named name in vars. More...
 
LALInferenceParamVaryType LALInferenceGetVariableVaryType (LALInferenceVariables *vars, const char *name)
 Get the LALInferenceParamVaryType of the parameter named name in vars see the declaration of LALInferenceParamVaryType for possibilities. More...
 
void LALInferenceSetParamVaryType (LALInferenceVariables *vars, const char *name, LALInferenceParamVaryType vary)
 Set the LALInferenceParamVaryType of the parameter named name in vars, see the declaration of LALInferenceParamVaryType for possibilities. More...
 
charLALInferenceGetVariableName (LALInferenceVariables *vars, int idx)
 Get the name of the idx-th variable Indexing starts at 1. More...
 
void LALInferenceSetVariable (LALInferenceVariables *vars, const char *name, const void *value)
 Set a variable named name in vars with a value. More...
 
void LALInferenceAddVariable (LALInferenceVariables *vars, const char *name, const void *value, LALInferenceVariableType type, LALInferenceParamVaryType vary)
 Add a variable named name to vars with initial value referenced by value. More...
 
void LALInferenceRemoveVariable (LALInferenceVariables *vars, const char *name)
 Remove name from vars Frees the memory for the name structure and its contents. More...
 
int LALInferenceCheckVariable (LALInferenceVariables *vars, const char *name)
 Checks for name being present in vars returns 1(==true) or 0. More...
 
int LALInferenceCheckVariableNonFixed (LALInferenceVariables *vars, const char *name)
 Checks for name being present in vars and having type LINEAR or CIRCULAR. More...
 
int LALInferenceCheckVariableToPrint (LALInferenceVariables *vars, const char *name)
 
void LALInferenceClearVariables (LALInferenceVariables *vars)
 Delete the variables in this structure. More...
 
void LALInferenceCopyVariables (LALInferenceVariables *origin, LALInferenceVariables *target)
 Deep copy the variables from one to another LALInferenceVariables structure. More...
 
void LALInferenceCopyUnsetREAL8Variables (LALInferenceVariables *origin, LALInferenceVariables *target, ProcessParamsTable *commandLine)
 
void LALInferencePrintVariables (LALInferenceVariables *var)
 Print variables to stdout. More...
 
int LALInferenceCompareVariables (LALInferenceVariables *var1, LALInferenceVariables *var2)
 Check for equality in two variables. More...
 
int LALInferenceSplineCalibrationFactor (REAL8Vector *freqs, REAL8Vector *deltaAmps, REAL8Vector *deltaPhases, COMPLEX16FrequencySeries *calFactor)
 Computes the factor relating the physical waveform to a measured waveform for a spline-fit calibration model in amplitude and phase. More...
 
int LALInferenceSplineCalibrationFactorROQ (REAL8Vector *logfreqs, REAL8Vector *deltaAmps, REAL8Vector *deltaPhases, REAL8Sequence *freqNodesLin, COMPLEX16Sequence **calFactorROQLin, REAL8Sequence *freqNodesQuad, COMPLEX16Sequence **calFactorROQQuad)
 Modified version of LALInferenceSplineCalibrationFactor to compute the calibration factors for the specific frequency nodes used for Reduced Order Quadrature likelihoods. More...
 
LALInferenceThreadStateLALInferenceInitThread (LALInferenceThreadState *thread)
 Structure to contain data-related Reduced Order Quadrature quantities. More...
 
LALInferenceThreadStateLALInferenceInitThreads (INT4 nthreads)
 
ProcessParamsTableLALInferenceGetProcParamVal (ProcessParamsTable *procparams, const char *name)
 Returns the element of the process params table with "name". More...
 
void LALInferenceParseCharacterOptionString (char *input, char **strings[], UINT4 *n)
 parses a character string (passed as one of the options) and decomposes it into individual parameter character strings. More...
 
ProcessParamsTableLALInferenceParseCommandLine (int argc, char **argv)
 Return a ProcessParamsTable from the command line arguments. More...
 
ProcessParamsTableLALInferenceParseStringVector (LALStringVector *arglist)
 Return a ProcessParamsTrable from a string vector. More...
 
ProcessParamsTableLALInferenceParseCommandLineStringVector (LALStringVector *args)
 Return a ProcessParamsTable from the command line arguments (SWIG version) More...
 
charLALInferencePrintCommandLine (ProcessParamsTable *procparams)
 Output the command line based on the ProcessParamsTable procparams. More...
 
void LALInferenceExecuteFT (LALInferenceModel *model)
 Execute FFT for data in IFOdata. More...
 
void LALInferenceExecuteInvFT (LALInferenceModel *model)
 Execute Inverse FFT for data in IFOdata. More...
 
LALInferenceVariableItemLALInferenceGetItem (const LALInferenceVariables *vars, const char *name)
 Return the list node for "name" - do not rely on this. More...
 
LALInferenceVariableItemLALInferenceGetItemNr (LALInferenceVariables *vars, int idx)
 Return the list node for the idx-th item - do not rely on this Indexing starts at 1. More...
 
LALInferenceVariableItemLALInferencePopVariableItem (LALInferenceVariables *vars, const char *name)
 Pop the list node for "name". More...
 
void LALInferencePrintSample (FILE *fp, LALInferenceVariables *sample)
 Output the sample to file *fp, in ASCII format. More...
 
void LALInferencePrintSampleNonFixed (FILE *fp, LALInferenceVariables *sample)
 Output only non-fixed parameters. More...
 
void LALInferencePrintSplineCalibration (FILE *fp, LALInferenceThreadState *thread)
 Output spline calibration parameters. More...
 
void LALInferenceReadSampleNonFixed (FILE *fp, LALInferenceVariables *sample)
 Read in the non-fixed parameters from the given file (position in the file must be arranged properly before calling this function). More...
 
REAL8LALInferenceParseDelimitedAscii (FILE *input, INT4 nCols, INT4 *wantedCols, INT4 *nLines)
 Utility for readling in delimited ASCII files. More...
 
void parseLine (char *record, const char *delim, char arr[][VARNAME_MAX], INT4 *cnt)
 Parse a single line of delimited ASCII. More...
 
void LALInferenceDiscardPTMCMCHeader (FILE *filestream)
 Discard the standard header of a PTMCMC chain file. More...
 
void LALInferenceBurninPTMCMC (FILE *filestream, INT4 logl_idx, INT4 nPar)
 Determine burnin cycle from delta-logl criteria. More...
 
void LALInferenceBurninStream (FILE *filestream, INT4 burnin)
 Burn-in a generic ASCII stream. More...
 
void LALInferenceReadAsciiHeader (FILE *input, char params[][VARNAME_MAX], INT4 *nCols)
 Read column names from an ASCII file. More...
 
REAL8 ** LALInferenceSelectColsFromArray (REAL8 **inarray, INT4 nRows, INT4 nCols, INT4 nSelCols, INT4 *selCols)
 Utility for selecting columns from an array, in the specified order. More...
 
int LALInferencePrintProposalStatsHeader (FILE *fp, LALInferenceProposalCycle *cycle)
 Output proposal statistics header to file *fp. More...
 
void LALInferencePrintProposalStats (FILE *fp, LALInferenceProposalCycle *cycle)
 Output proposal statistics to file *fp. More...
 
void LALInferenceProcessParamLine (FILE *inp, char **headers, LALInferenceVariables *vars)
 Reads one line from the given file and stores the values there into the variable structure, using the given header array to name the columns. More...
 
void LALInferenceSortVariablesByName (LALInferenceVariables *vars)
 Sorts the variable structure by name. More...
 
INT4 LALInferenceThinnedBufferToArray (LALInferenceThreadState *thread, REAL8 **DEarray, INT4 step)
 LALInferenceVariable buffer to array and vica versa. More...
 
INT4 LALInferenceBufferToArray (LALInferenceThreadState *thread, REAL8 **DEarray)
 
void LALInferenceCopyVariablesToArray (LALInferenceVariables *origin, REAL8 *target)
 LALInference variables to an array, and vica versa. More...
 
void LALInferenceCopyArrayToVariables (REAL8 *origin, LALInferenceVariables *target)
 
void LALInferenceLogSampleToFile (LALInferenceVariables *algorithmParams, LALInferenceVariables *vars)
 Append the sample to a file. More...
 
void LALInferenceLogSampleToArray (LALInferenceVariables *algorithmParams, LALInferenceVariables *vars)
 Append the sample to an array which can be later processed by the user. More...
 
void LALInferenceMcEta2Masses (double mc, double eta, double *m1, double *m2)
 Convert from Mc, eta space to m1, m2 space (note m1 > m2). More...
 
void LALInferenceMcQ2Masses (double mc, double q, double *m1, double *m2)
 Convert from Mc, q space to m1, m2 space (q = m2/m1, with m1 > m2). More...
 
void LALInferenceQ2Eta (double q, double *eta)
 Convert from q to eta (q = m2/m1, with m1 > m2). More...
 
void LALInferencedQuadMonSdQuadMonA (REAL8 dQuadMonS, REAL8 dQuadMonA, REAL8 *dQuadMon1, REAL8 *dQuadMon2)
 Convert from dQuadMonS, dQuadMonA to dQuadMon1, dQuadMon2. More...
 
void LALInferenceLambdaTsEta2Lambdas (REAL8 lambdaT, REAL8 dLambdaT, REAL8 eta, REAL8 *lambda1, REAL8 *lambda2)
 Convert from lambdaT, dLambdaT, and eta to lambda1 and lambda2. More...
 
void LALInferenceLogp1GammasMasses2Lambdas (REAL8 logp1, REAL8 gamma1, REAL8 gamma2, REAL8 gamma3, REAL8 mass1, REAL8 mass2, REAL8 *lambda1, REAL8 *lambda2)
 Calculate lambda1,2(m1,2|eos(logp1,gamma1,gamma2,gamma3)) More...
 
void LALInferenceSDGammasMasses2Lambdas (REAL8 gamma[], REAL8 mass1, REAL8 mass2, REAL8 *lambda1, REAL8 *lambda2, int size)
 Convert from spectral parameters to lambda1, lambda2. More...
 
int LALInferenceEOSPhysicalCheck (LALInferenceVariables *params, ProcessParamsTable *commandLine)
 Check for causality violation and mass conflict given masses and eos. More...
 
double AdiabaticIndex (double gamma[], double x, int size)
 Specral decomposition of eos's adiabatic index. More...
 
int LALInferenceSDGammaCheck (double gamma[], int size)
 Determine if the Adiabatic index is within 'prior' range. More...
 
void LALInferenceKDTreeDelete (LALInferenceKDTree *tree)
 Delete a kD-tree. More...
 
LALInferenceKDTreeLALInferenceKDEmpty (REAL8 *lowerLeft, REAL8 *upperRight, size_t ndim)
 Constructs a fresh, empty kD tree. More...
 
int LALInferenceKDAddPoint (LALInferenceKDTree *tree, REAL8 *pt)
 Adds a point to the kD-tree, returns 0 on successful exit. More...
 
LALInferenceKDTreeLALInferenceKDFindCell (LALInferenceKDTree *tree, REAL8 *pt, size_t Npts)
 Returns the first cell that contains the given point that also contains fewer than Npts points, if possible. More...
 
double LALInferenceKDLogCellVolume (LALInferenceKDTree *cell)
 Returns the log of the volume of the given cell, which is part of the given tree. More...
 
double LALInferenceKDLogCellEigenVolume (LALInferenceKDTree *cell)
 Returns the log of the volume of the box aligned with the principal axes of the points in the given cell that tightly encloses those points. More...
 
void LALInferenceKDVariablesToREAL8 (LALInferenceVariables *params, REAL8 *pt, LALInferenceVariables *templt)
 Fills in the given REAL8 array with the parameter values from params; the ordering of the variables is taken from the order of the non-fixed variables in templt. More...
 
void LALInferenceKDREAL8ToVariables (LALInferenceVariables *params, REAL8 *pt, LALInferenceVariables *templt)
 Fills in the non-fixed variables in params from the given REAL8 array. More...
 
void LALInferenceKDDrawEigenFrame (gsl_rng *rng, LALInferenceKDTree *tree, REAL8 *pt, size_t Npts)
 Draws a pt uniformly from a randomly chosen cell of tree. More...
 
REAL8 LALInferenceKDLogProposalRatio (LALInferenceKDTree *tree, REAL8 *current, REAL8 *proposed, size_t Npts)
 Returns the log of the jump proposal probability ratio for the LALInferenceKDDrawEigenFrame() proposal to propose the point proposed given the current position current , where Npts is the parameter used to select the box to draw from in LALInferenceKDDrawEigenFrame(). More...
 
UINT4 LALInferenceCheckPositiveDefinite (gsl_matrix *matrix, UINT4 dim)
 Check matrix is positive definite. More...
 
void XLALMultiNormalDeviates (REAL4Vector *vector, gsl_matrix *matrix, UINT4 dim, RandomParams *randParam)
 Draw a random multivariate vector from Gaussian distr given covariance matrix. More...
 
void XLALMultiStudentDeviates (REAL4Vector *vector, gsl_matrix *matrix, UINT4 dim, UINT4 n, RandomParams *randParam)
 Draw a random multivariate vector from student-t distr given covariance matrix. More...
 
REAL8 LALInferenceAngularDistance (REAL8 a1, REAL8 a2)
 Calculate shortest angular distance between a1 and a2 (modulo 2PI) More...
 
REAL8 LALInferenceAngularVariance (LALInferenceVariables **list, const char *pname, int N)
 Calculate the variance of a distribution on an angle (modulo 2PI) More...
 
INT4 LALInferenceSanityCheck (LALInferenceRunState *state)
 Sanity check the structures in the given state. More...
 
void LALInferenceDumpWaveforms (LALInferenceModel *model, const char *basefilename)
 Dump all waveforms from the ifodata structure. More...
 
int LALInferenceWriteVariablesBinary (FILE *file, LALInferenceVariables *vars)
 Write a LALInferenceVariables as binary to a given FILE pointer, returns the number of items written (should be the dimension of the variables) or -1 on error. More...
 
LALInferenceVariablesLALInferenceReadVariablesBinary (FILE *stream)
 Read from the given FILE * a LALInferenceVariables, which was previously written with LALInferenceWriteVariablesBinary() Returns a new LALInferenceVariables. More...
 
int LALInferenceWriteVariablesArrayBinary (FILE *file, LALInferenceVariables **vars, UINT4 N)
 Write an array N of LALInferenceVariables to the given FILE * using LALInferenceWriteVariablesBinary(). More...
 
int LALInferenceReadVariablesArrayBinary (FILE *file, LALInferenceVariables **vars, UINT4 N)
 Read N LALInferenceVariables from the binary FILE *file, previously written with LALInferenceWriteVariablesArrayBinary() returns the number read. More...
 
int LALInferenceWriteRunStateBinary (FILE *file, LALInferenceRunState *state)
 Write the LALInferenceVariables contents of runState to a file in binary format. More...
 
int LALInferenceReadRunStateBinary (FILE *file, LALInferenceRunState *state)
 Reads the file and populates LALInferenceVariables in the runState that were saved with LALInferenceReadVariablesArrayBinary() More...
 
void LALInferenceAddINT4Variable (LALInferenceVariables *vars, const char *name, INT4 value, LALInferenceParamVaryType vary)
 
INT4 LALInferenceGetINT4Variable (LALInferenceVariables *vars, const char *name)
 
void LALInferenceSetINT4Variable (LALInferenceVariables *vars, const char *name, INT4 value)
 
void LALInferenceAddINT8Variable (LALInferenceVariables *vars, const char *name, INT8 value, LALInferenceParamVaryType vary)
 
INT8 LALInferenceGetINT8Variable (LALInferenceVariables *vars, const char *name)
 
void LALInferenceSetINT8Variable (LALInferenceVariables *vars, const char *name, INT8 value)
 
void LALInferenceAddUINT4Variable (LALInferenceVariables *vars, const char *name, UINT4 value, LALInferenceParamVaryType vary)
 
UINT4 LALInferenceGetUINT4Variable (LALInferenceVariables *vars, const char *name)
 
void LALInferenceSetUINT4Variable (LALInferenceVariables *vars, const char *name, UINT4 value)
 
void LALInferenceAddREAL4Variable (LALInferenceVariables *vars, const char *name, REAL4 value, LALInferenceParamVaryType vary)
 
REAL4 LALInferenceGetREAL4Variable (LALInferenceVariables *vars, const char *name)
 
void LALInferenceSetREAL4Variable (LALInferenceVariables *vars, const char *name, REAL4 value)
 
void LALInferenceAddREAL8Variable (LALInferenceVariables *vars, const char *name, REAL8 value, LALInferenceParamVaryType vary)
 
REAL8 LALInferenceGetREAL8Variable (LALInferenceVariables *vars, const char *name)
 
void LALInferenceSetREAL8Variable (LALInferenceVariables *vars, const char *name, REAL8 value)
 
void LALInferenceAddCOMPLEX8Variable (LALInferenceVariables *vars, const char *name, COMPLEX8 value, LALInferenceParamVaryType vary)
 
COMPLEX8 LALInferenceGetCOMPLEX8Variable (LALInferenceVariables *vars, const char *name)
 
void LALInferenceSetCOMPLEX8Variable (LALInferenceVariables *vars, const char *name, COMPLEX8 value)
 
void LALInferenceAddCOMPLEX16Variable (LALInferenceVariables *vars, const char *name, COMPLEX16 value, LALInferenceParamVaryType vary)
 
COMPLEX16 LALInferenceGetCOMPLEX16Variable (LALInferenceVariables *vars, const char *name)
 
void LALInferenceSetCOMPLEX16Variable (LALInferenceVariables *vars, const char *name, COMPLEX16 value)
 
void LALInferenceAddgslMatrixVariable (LALInferenceVariables *vars, const char *name, gsl_matrix *value, LALInferenceParamVaryType vary)
 
gsl_matrix * LALInferenceGetgslMatrixVariable (LALInferenceVariables *vars, const char *name)
 
void LALInferenceSetgslMatrixVariable (LALInferenceVariables *vars, const char *name, gsl_matrix *value)
 
void LALInferenceAddREAL8VectorVariable (LALInferenceVariables *vars, const char *name, REAL8Vector *value, LALInferenceParamVaryType vary)
 
REAL8VectorLALInferenceGetREAL8VectorVariable (LALInferenceVariables *vars, const char *name)
 
void LALInferenceSetREAL8VectorVariable (LALInferenceVariables *vars, const char *name, REAL8Vector *value)
 
void LALInferenceAddCOMPLEX16VectorVariable (LALInferenceVariables *vars, const char *name, COMPLEX16Vector *value, LALInferenceParamVaryType vary)
 
COMPLEX16VectorLALInferenceGetCOMPLEX16VectorVariable (LALInferenceVariables *vars, const char *name)
 
void LALInferenceSetCOMPLEX16VectorVariable (LALInferenceVariables *vars, const char *name, COMPLEX16Vector *value)
 
void LALInferenceAddINT4VectorVariable (LALInferenceVariables *vars, const char *name, INT4Vector *value, LALInferenceParamVaryType vary)
 
void LALInferenceAddUINT4VectorVariable (LALInferenceVariables *vars, const char *name, UINT4Vector *value, LALInferenceParamVaryType vary)
 
INT4VectorLALInferenceGetINT4VectorVariable (LALInferenceVariables *vars, const char *name)
 
UINT4VectorLALInferenceGetUINT4VectorVariable (LALInferenceVariables *vars, const char *name)
 
void LALInferenceSetINT4VectorVariable (LALInferenceVariables *vars, const char *name, INT4Vector *value)
 
void LALInferenceSetUINT4VectorVariable (LALInferenceVariables *vars, const char *name, UINT4Vector *value)
 
void LALInferenceAddMCMCrunphase_ptrVariable (LALInferenceVariables *vars, const char *name, LALInferenceMCMCRunPhase *value, LALInferenceParamVaryType vary)
 
LALInferenceMCMCRunPhaseLALInferenceGetMCMCrunphase_ptrVariable (LALInferenceVariables *vars, const char *name)
 
void LALInferenceSetMCMCrunphase_ptrVariable (LALInferenceVariables *vars, const char *name, LALInferenceMCMCRunPhase *value)
 
void LALInferenceAddstringVariable (LALInferenceVariables *vars, const char *name, const CHAR *value, LALInferenceParamVaryType vary)
 
const CHARLALInferenceGetstringVariable (LALInferenceVariables *vars, const char *name)
 
void LALInferenceSetstringVariable (LALInferenceVariables *vars, const char *name, const CHAR *value)
 
void LALInferenceFprintSplineCalibrationHeader (FILE *output, LALInferenceThreadState *thread)
 Print spline calibration parameter names as tab-separated ASCII. More...
 
void LALInferenceBinaryLove (LALInferenceVariables *vars, REAL8 *lambda1, REAL8 *lambda2)
 Compute Tidal deformabilities following BinaryLove Universal relations. More...
 
void LALInferenceDetFrameToEquatorial (LALDetector *det0, LALDetector *det1, REAL8 t0, REAL8 alpha, REAL8 theta, REAL8 *tg, REAL8 *ra, REAL8 *dec)
 Conversion routines between Equatorial (RA,DEC) and detector-based coordinate systems, where new "north pole" points along vector from det0 to det1. More...
 
void LALInferenceEquatorialToDetFrame (LALDetector *det0, LALDetector *det1, REAL8 tg, REAL8 ra, REAL8 dec, REAL8 *t0, REAL8 *alpha, REAL8 *theta)
 

Data Structures

struct  LALInferenceVariableItem
 The LALInferenceVariableItem list node structure This should only be accessed using the accessor functions below Implementation may change to hash table so please use only the accessor functions below. More...
 
struct  LALInferenceVariables
 The LALInferenceVariables structure to contain a set of parameters Implemented as a linked list of LALInferenceVariableItems. More...
 
struct  LALInferenceIFOModel
 Detector-dependent buffers and parameters A linked list meant for parameters and buffers that are separately specified for each detector. More...
 
struct  LALInferenceModel
 Structure to constain a model and its parameters. More...
 
struct  LALInferenceProposal
 Structure for holding a LALInference proposal, along with name and stats. More...
 
struct  LALInferenceProposalCycle
 Structure for holding a proposal cycle. More...
 
struct  LALInferenceThreadState
 Structure containing chain-specific variables. More...
 
struct  LALInferenceRunState
 Structure containing inference run state This includes pointers to the function types required to run the algorithm, and data structures as required. More...
 
struct  LALInferenceIFOData
 Structure to contain IFO data. More...
 
struct  LALInferenceROQData
 
struct  LALInferenceROQSplineWeights
 
struct  LALInferenceROQModel
 
struct  LALInferenceKDTree
 The kD trees in LALInference are composed of cells. More...
 

Typedefs

typedef void(* LALInferenceTemplateFunction) (struct tagLALInferenceModel *model)
 Type declaration for template function, which operates on a LALInferenceIFOData structure *data. More...
 
typedef REAL8(* LALInferenceProposalFunction) (struct tagLALInferenceThreadState *thread, LALInferenceVariables *currentParams, LALInferenceVariables *proposedParams)
 Jump proposal distribution Computes proposedParams based on currentParams and additional variables stored as proposalArgs inside runState, which could include correlation matrix, etc., as well as forward and reverse proposal probability. More...
 
typedef REAL8(* LALInferencePriorFunction) (struct tagLALInferenceRunState *runState, LALInferenceVariables *params, struct tagLALInferenceModel *model)
 Type declaration for prior function which returns p(params) Can depend on runState ->priorArgs. More...
 
typedef UINT4(* LALInferenceCubeToPriorFunction) (struct tagLALInferenceRunState *runState, LALInferenceVariables *params, struct tagLALInferenceModel *model, double *cube, void *context)
 Type declaration for CubeToPrior function which converts parameters in unit hypercube to their corresponding physical values according to the prior. More...
 
typedef LALInferenceModel *(* LALInferenceInitModelFunction) (struct tagLALInferenceRunState *runState)
 Type declaration for variables init function, can be user-declared. More...
 
typedef REAL8(* LALInferenceLikelihoodFunction) (LALInferenceVariables *currentParams, struct tagLALInferenceIFOData *data, LALInferenceModel *model)
 Type declaration for likelihood function Computes p(data | currentParams, templt ) templt is a LALInferenceTemplateFunction defined below. More...
 
typedef INT4(* LALInferenceEvolveOneStepFunction) (struct tagLALInferenceRunState *runState)
 Perform one step of an algorithm, replaces runState ->currentParams. More...
 
typedef void(* LALInferenceSwapRoutine) (struct tagLALInferenceRunState *runState, FILE *)
 Propose a swap between chain locations. More...
 
typedef void(* LALInferenceAlgorithm) (struct tagLALInferenceRunState *runState)
 Type declaration for an algorithm function which is called by the driver code The user must initialise runState before use. More...
 
typedef void(* LALInferenceLogFunction) (LALInferenceVariables *algorithmParams, LALInferenceVariables *vars)
 Type declaration for output logging function, can be user-declared. More...
 

Enumerations

enum  LALInferenceVariableType {
  LALINFERENCE_INT4_t , LALINFERENCE_INT8_t , LALINFERENCE_UINT4_t , LALINFERENCE_REAL4_t ,
  LALINFERENCE_REAL8_t , LALINFERENCE_COMPLEX8_t , LALINFERENCE_COMPLEX16_t , LALINFERENCE_gslMatrix_t ,
  LALINFERENCE_REAL8Vector_t , LALINFERENCE_INT4Vector_t , LALINFERENCE_UINT4Vector_t , LALINFERENCE_COMPLEX16Vector_t ,
  LALINFERENCE_string_t , LALINFERENCE_MCMCrunphase_ptr_t , LALINFERENCE_void_ptr_t
}
 An enumerated type for denoting the type of a variable. More...
 
enum  LALInferenceParamVaryType { LALINFERENCE_PARAM_LINEAR , LALINFERENCE_PARAM_CIRCULAR , LALINFERENCE_PARAM_FIXED , LALINFERENCE_PARAM_OUTPUT }
 An enumerated type for denoting the topology of a parameter. More...
 
enum  LALInferenceMCMCRunPhase {
  LALINFERENCE_ONLY_PT , LALINFERENCE_TEMP_PT , LALINFERENCE_ANNEALING , LALINFERENCE_SINGLE_CHAIN ,
  LALINFERENCE_LADDER_UPDATE
}
 Phase of MCMC run (depending on burn-in status, different actions are performed during the run, and this tag controls the activity). More...
 

Macros

#define VARNAME_MAX   40
 
#define VARVALSTRINGSIZE_MAX   128
 
#define DETNAMELEN   8
 

Variables

size_t LALInferenceTypeSize [15]
 

Function Documentation

◆ LALInferenceGetHeaderLine()

char** LALInferenceGetHeaderLine ( FILE *  inp)

Returns an array of header strings (terminated by NULL) from a common-format output file.

Definition at line 2087 of file LALInference.c.

◆ LALInferenceTranslateInternalToExternalParamName()

const char* LALInferenceTranslateInternalToExternalParamName ( const char inName)

Converts between internally used parameter names and those external (e.g.

in SimInspiralTable?)

Definition at line 1287 of file LALInference.c.

◆ LALInferenceTranslateExternalToInternalParamName()

void LALInferenceTranslateExternalToInternalParamName ( char outName,
const char inName 
)

Converts between externally used parameter names and those internal.

Definition at line 1319 of file LALInference.c.

◆ LALInferenceFprintParameterHeaders()

int LALInferenceFprintParameterHeaders ( FILE *  out,
LALInferenceVariables params 
)

Print the parameter names to a file as a tab-separated ASCII line.

Parameters
out[in] pointer to output file
params[in] LALInferenceVaraibles structure to print

Definition at line 1356 of file LALInference.c.

◆ LALInferenceFprintParameterNonFixedHeaders()

INT4 LALInferenceFprintParameterNonFixedHeaders ( FILE *  out,
LALInferenceVariables params 
)

Print the parameters which do not vary to a file as a tab-separated ASCII line.

Parameters
out[in] pointer to output file
params[in] LALInferenceVaraibles structure to print

Definition at line 1398 of file LALInference.c.

◆ LALInferenceFprintParameterNonFixedHeadersWithSuffix()

INT4 LALInferenceFprintParameterNonFixedHeadersWithSuffix ( FILE *  out,
LALInferenceVariables params,
const char suffix 
)

Print the parameters which do not vary to a file as a tab-separated ASCII line, adding the given suffix.

Parameters
out[in] pointer to output file
params[in] LALInferenceVaraibles structure to print
suffix[in] Suffix string to add to each parameter name

Definition at line 1449 of file LALInference.c.

◆ LALInferencePrintNVariableItem()

UINT4 LALInferencePrintNVariableItem ( char out,
UINT4  strsize,
const LALInferenceVariableItem *const  ptr 
)

Prints a variable item to a string.

Print at most N characters. Returns the number of characters actually required to store the output (can be more or less than N)

Prints a variable item to a string.

Definition at line 687 of file LALInference.c.

◆ LALInferenceGetVariable()

void* LALInferenceGetVariable ( const LALInferenceVariables vars,
const char name 
)

Return a pointer to the memory the variable vars is stored in specified by name User must cast this pointer to the expected type before dereferencing it to get the value of the variable.

Definition at line 238 of file LALInference.c.

◆ LALInferenceGetVariableDimension()

INT4 LALInferenceGetVariableDimension ( LALInferenceVariables vars)

Get number of dimensions in variable vars.

Definition at line 250 of file LALInference.c.

◆ LALInferenceGetVariableDimensionNonFixed()

INT4 LALInferenceGetVariableDimensionNonFixed ( LALInferenceVariables vars)

Get number of dimensions in vars which are not fixed to a certain value.

Definition at line 255 of file LALInference.c.

◆ LALInferenceGetVariableDimensionNonFixedChooseVectors()

INT4 LALInferenceGetVariableDimensionNonFixedChooseVectors ( LALInferenceVariables vars,
INT4  count_vectors 
)

Get number of dimensions in vars which are not fixed to a certain value, with a flag for skipping counting vectors.

Definition at line 260 of file LALInference.c.

◆ LALInferenceGetVariableTypeByIndex()

INT4 LALInferenceGetVariableTypeByIndex ( LALInferenceVariables vars,
int  idx 
)

Get the LALInferenceVariableType of the idx -th item in the vars Indexing starts at 1.

Definition at line 326 of file LALInference.c.

◆ LALInferenceGetVariableType()

LALInferenceVariableType LALInferenceGetVariableType ( const LALInferenceVariables vars,
const char name 
)

Get the LALInferenceVariableType of the parameter named name in vars.

Definition at line 321 of file LALInference.c.

◆ LALInferenceGetVariableVaryType()

LALInferenceParamVaryType LALInferenceGetVariableVaryType ( LALInferenceVariables vars,
const char name 
)

Get the LALInferenceParamVaryType of the parameter named name in vars see the declaration of LALInferenceParamVaryType for possibilities.

Definition at line 226 of file LALInference.c.

◆ LALInferenceSetParamVaryType()

void LALInferenceSetParamVaryType ( LALInferenceVariables vars,
const char name,
LALInferenceParamVaryType  vary 
)

Set the LALInferenceParamVaryType of the parameter named name in vars, see the declaration of LALInferenceParamVaryType for possibilities.

Definition at line 231 of file LALInference.c.

◆ LALInferenceGetVariableName()

char* LALInferenceGetVariableName ( LALInferenceVariables vars,
int  idx 
)

Get the name of the idx-th variable Indexing starts at 1.

Definition at line 339 of file LALInference.c.

◆ LALInferenceSetVariable()

void LALInferenceSetVariable ( LALInferenceVariables vars,
const char name,
const void *  value 
)

Set a variable named name in vars with a value.

Pass a void * in value to the value you wish to set, i.e. LALInferenceSetVariable(vars, "mu", (void *)&mu);

Definition at line 352 of file LALInference.c.

◆ LALInferenceAddVariable()

void LALInferenceAddVariable ( LALInferenceVariables vars,
const char name,
const void *  value,
LALInferenceVariableType  type,
LALInferenceParamVaryType  vary 
)

Add a variable named name to vars with initial value referenced by value.

Parameters
typeis a LALInferenceVariableType (enumerated above)
varyis a LALInferenceParamVaryType (enumerated above)
varsUNDOCUMENTED
nameUNDOCUMENTED
valueUNDOCUMENTED If the variable already exists it will be over-written UNLESS IT HAS A CONFLICTING TYPE

Definition at line 395 of file LALInference.c.

◆ LALInferenceRemoveVariable()

void LALInferenceRemoveVariable ( LALInferenceVariables vars,
const char name 
)

Remove name from vars Frees the memory for the name structure and its contents.

Definition at line 450 of file LALInference.c.

◆ LALInferenceCheckVariable()

int LALInferenceCheckVariable ( LALInferenceVariables vars,
const char name 
)

Checks for name being present in vars returns 1(==true) or 0.

Definition at line 525 of file LALInference.c.

◆ LALInferenceCheckVariableNonFixed()

int LALInferenceCheckVariableNonFixed ( LALInferenceVariables vars,
const char name 
)

Checks for name being present in vars and having type LINEAR or CIRCULAR.

returns 1 or 0

Definition at line 503 of file LALInference.c.

◆ LALInferenceCheckVariableToPrint()

int LALInferenceCheckVariableToPrint ( LALInferenceVariables vars,
const char name 
)

Definition at line 514 of file LALInference.c.

◆ LALInferenceClearVariables()

void LALInferenceClearVariables ( LALInferenceVariables vars)

Delete the variables in this structure.

Does not free the LALInferenceVariables itself

Parameters
varswill have its dimension set to 0

Definition at line 532 of file LALInference.c.

◆ LALInferenceCopyVariables()

void LALInferenceCopyVariables ( LALInferenceVariables origin,
LALInferenceVariables target 
)

Deep copy the variables from one to another LALInferenceVariables structure.

Definition at line 558 of file LALInference.c.

◆ LALInferenceCopyUnsetREAL8Variables()

void LALInferenceCopyUnsetREAL8Variables ( LALInferenceVariables origin,
LALInferenceVariables target,
ProcessParamsTable commandLine 
)

Definition at line 659 of file LALInference.c.

◆ LALInferencePrintVariables()

void LALInferencePrintVariables ( LALInferenceVariables var)

Print variables to stdout.

Print variables to stdout.

Definition at line 798 of file LALInference.c.

◆ LALInferenceCompareVariables()

int LALInferenceCompareVariables ( LALInferenceVariables var1,
LALInferenceVariables var2 
)

Check for equality in two variables.

Definition at line 1476 of file LALInference.c.

◆ LALInferenceSplineCalibrationFactor()

int LALInferenceSplineCalibrationFactor ( REAL8Vector freqs,
REAL8Vector deltaAmps,
REAL8Vector deltaPhases,
COMPLEX16FrequencySeries calFactor 
)

Computes the factor relating the physical waveform to a measured waveform for a spline-fit calibration model in amplitude and phase.

The spline points can be arbitrary frequencies, and the values of the calibration offset at these points can be arbitary, too. For amplitude offset, \(\delta A\), and phase offset, \(\delta \phi\), the measured waveform is related to the physical waveform via

\[ h_\mathrm{meas} = h_\mathrm{phys} \left(1 + \delta A \right) \frac{2 + i \delta \phi}{2 - i \delta \phi} \]

The phase factor takes the form above rather than the more obvious \(\exp(i \delta \phi)\) or \(1 + \delta \phi\) because it is faster to compute than the former (but, for small \(\phi\), equivalent through second order in \(\delta \phi\)) and, unlike the latter, it ensures that the complex amplitude of the correction is always 1. (A similar technique is used when using finite-difference approximations to solve the multi-dimensional Schrodinger equation to ensure that the evolution remains unitary.) Note that this implies that the phase calibration parameter ranges over \(\delta \phi \in [-\infty, \infty]\).

The values of \(\delta A\) and \(\delta \phi\) at arbitrary frequency are obtained by a spline curve that passes through the given values at the given frequencies.

Definition at line 4255 of file LALInference.c.

◆ LALInferenceSplineCalibrationFactorROQ()

int LALInferenceSplineCalibrationFactorROQ ( REAL8Vector logfreqs,
REAL8Vector deltaAmps,
REAL8Vector deltaPhases,
REAL8Sequence freqNodesLin,
COMPLEX16Sequence **  calFactorROQLin,
REAL8Sequence freqNodesQuad,
COMPLEX16Sequence **  calFactorROQQuad 
)

Modified version of LALInferenceSplineCalibrationFactor to compute the calibration factors for the specific frequency nodes used for Reduced Order Quadrature likelihoods.

Definition at line 4333 of file LALInference.c.

◆ LALInferenceInitThread()

LALInferenceThreadState* LALInferenceInitThread ( LALInferenceThreadState thread)

Structure to contain data-related Reduced Order Quadrature quantities.

Definition at line 105 of file LALInference.c.

◆ LALInferenceInitThreads()

LALInferenceThreadState* LALInferenceInitThreads ( INT4  nthreads)

Definition at line 138 of file LALInference.c.

◆ LALInferenceGetProcParamVal()

ProcessParamsTable* LALInferenceGetProcParamVal ( ProcessParamsTable procparams,
const char name 
)

Returns the element of the process params table with "name".

Definition at line 1762 of file LALInference.c.

◆ LALInferenceParseCharacterOptionString()

void LALInferenceParseCharacterOptionString ( char input,
char **  strings[],
UINT4 n 
)

parses a character string (passed as one of the options) and decomposes it into individual parameter character strings.

input is of the form input : "[one,two,three]" and the resulting output strings is strings : {"one", "two", "three"} length of parameter names is for now limited to 512 characters. (should 'theoretically' (untested) be able to digest white space as well. Irrelevant for command line options, though.) n UNDOCUMENTED

Definition at line 1783 of file LALInference.c.

◆ LALInferenceParseCommandLine()

ProcessParamsTable* LALInferenceParseCommandLine ( int  argc,
char **  argv 
)

Return a ProcessParamsTable from the command line arguments.

◆ LALInferenceParseStringVector()

ProcessParamsTable* LALInferenceParseStringVector ( LALStringVector arglist)

Return a ProcessParamsTrable from a string vector.

Definition at line 1844 of file LALInference.c.

◆ LALInferenceParseCommandLineStringVector()

ProcessParamsTable* LALInferenceParseCommandLineStringVector ( LALStringVector args)

Return a ProcessParamsTable from the command line arguments (SWIG version)

Definition at line 1911 of file LALInference.c.

◆ LALInferencePrintCommandLine()

char* LALInferencePrintCommandLine ( ProcessParamsTable procparams)

Output the command line based on the ProcessParamsTable procparams.

Definition at line 1916 of file LALInference.c.

◆ LALInferenceExecuteFT()

void LALInferenceExecuteFT ( LALInferenceModel model)

Execute FFT for data in IFOdata.

Definition at line 1945 of file LALInference.c.

◆ LALInferenceExecuteInvFT()

void LALInferenceExecuteInvFT ( LALInferenceModel model)

Execute Inverse FFT for data in IFOdata.

Definition at line 2038 of file LALInference.c.

◆ LALInferenceGetItem()

LALInferenceVariableItem* LALInferenceGetItem ( const LALInferenceVariables vars,
const char name 
)

Return the list node for "name" - do not rely on this.

< [in] Pointer to hash table

< [in] Hash element to match

< [out] Pointer to matched hash element, or NULL if not found

Definition at line 175 of file LALInference.c.

◆ LALInferenceGetItemNr()

LALInferenceVariableItem* LALInferenceGetItemNr ( LALInferenceVariables vars,
int  idx 
)

Return the list node for the idx-th item - do not rely on this Indexing starts at 1.

Definition at line 207 of file LALInference.c.

◆ LALInferencePopVariableItem()

LALInferenceVariableItem* LALInferencePopVariableItem ( LALInferenceVariables vars,
const char name 
)

Pop the list node for "name".

Returns a pointer to the node, which is removed from vars

Definition at line 2187 of file LALInference.c.

◆ LALInferencePrintSample()

void LALInferencePrintSample ( FILE *  fp,
LALInferenceVariables sample 
)

Output the sample to file *fp, in ASCII format.

Definition at line 923 of file LALInference.c.

◆ LALInferencePrintSampleNonFixed()

void LALInferencePrintSampleNonFixed ( FILE *  fp,
LALInferenceVariables sample 
)

Output only non-fixed parameters.

Definition at line 948 of file LALInference.c.

◆ LALInferencePrintSplineCalibration()

void LALInferencePrintSplineCalibration ( FILE *  fp,
LALInferenceThreadState thread 
)

Output spline calibration parameters.

Definition at line 4456 of file LALInference.c.

◆ LALInferenceReadSampleNonFixed()

void LALInferenceReadSampleNonFixed ( FILE *  fp,
LALInferenceVariables sample 
)

Read in the non-fixed parameters from the given file (position in the file must be arranged properly before calling this function).

Definition at line 974 of file LALInference.c.

◆ LALInferenceParseDelimitedAscii()

REAL8* LALInferenceParseDelimitedAscii ( FILE *  input,
INT4  nCols,
INT4 wantedCols,
INT4 nLines 
)

Utility for readling in delimited ASCII files.

Reads in an ASCII (delimited) file, and returns the results in a REAL8 array.

Parameters
[in]inputInput stream to be parsed.
[in]nColsNumber of total columns in the input stream.
[in]wantedColsArray of 0/1 flags (should be nCols long), indicating desired columns.
nLinesIf 0, read until end of file and set to be number of lines read. If > 0, only read nLines.
Returns
A REAL8 array containing the parsed data.

Definition at line 1017 of file LALInference.c.

◆ parseLine()

void parseLine ( char record,
const char delim,
char  arr[][VARNAME_MAX],
INT4 cnt 
)

Parse a single line of delimited ASCII.

Splits a line using delim into an array of strings.

Parameters
[in]recordLine to be parsed.
[in]delimDelimiter to split on.
[out]arrParsed string.
[out]cntTotal number of fields read.

Definition at line 1075 of file LALInference.c.

◆ LALInferenceDiscardPTMCMCHeader()

void LALInferenceDiscardPTMCMCHeader ( FILE *  filestream)

Discard the standard header of a PTMCMC chain file.

Reads a PTMCMC input stream, moving foward until it finds a line beginning with "cycle", which is typically the line containing the column names of a PTMCMC output file. The final stream points to the line containing the column names.

Parameters
filestreamThe PTMCMC input stream to discard the header of.

Definition at line 1101 of file LALInference.c.

◆ LALInferenceBurninPTMCMC()

void LALInferenceBurninPTMCMC ( FILE *  filestream,
INT4  logl_idx,
INT4  nPar 
)

Determine burnin cycle from delta-logl criteria.

Find the cycle where the log(likelihood) is within nPar/2 of the maximum log(likelihood) sampled by the chain.

Parameters
filestreamThe PTMCMC input stream to be burned in.
[in]logl_idxThe column containing logl values.
nParUNDOCUMENTED

Definition at line 1137 of file LALInference.c.

◆ LALInferenceBurninStream()

void LALInferenceBurninStream ( FILE *  filestream,
INT4  burnin 
)

Burn-in a generic ASCII stream.

Reads past the desired number of lines of a filestream.

Parameters
filestreamThe filestream to be burned in.
[in]burninNumber of lines to read past in filestream.

Definition at line 1185 of file LALInference.c.

◆ LALInferenceReadAsciiHeader()

void LALInferenceReadAsciiHeader ( FILE *  input,
char  params[][VARNAME_MAX],
INT4 nCols 
)

Read column names from an ASCII file.

Reads the column names for an output file.

Parameters
[in]inputThe input files stream to parse.
[out]paramsThe column names found in the header.
[out]nColsNumber of columns found in the header.

Definition at line 1210 of file LALInference.c.

◆ LALInferenceSelectColsFromArray()

REAL8** LALInferenceSelectColsFromArray ( REAL8 **  inarray,
INT4  nRows,
INT4  nCols,
INT4  nSelCols,
INT4 selCols 
)

Utility for selecting columns from an array, in the specified order.

Selects a subset of columns from an existing array and create a new array with them.

Parameters
[in]inarrayThe array to select columns from.
[in]nRowsNumber of rows in inarray.
[in]nColsNumber of columns in inarray.
[in]nSelColsNumber of columns being extracted.
[in]selColsArray of column numbers to be extracted from inarray.
Returns
An array containing the requested columns in the order specified in selCols.

Definition at line 1239 of file LALInference.c.

◆ LALInferencePrintProposalStatsHeader()

int LALInferencePrintProposalStatsHeader ( FILE *  fp,
LALInferenceProposalCycle cycle 
)

Output proposal statistics header to file *fp.

Definition at line 1256 of file LALInference.c.

◆ LALInferencePrintProposalStats()

void LALInferencePrintProposalStats ( FILE *  fp,
LALInferenceProposalCycle cycle 
)

Output proposal statistics to file *fp.

Definition at line 1267 of file LALInference.c.

◆ LALInferenceProcessParamLine()

void LALInferenceProcessParamLine ( FILE *  inp,
char **  headers,
LALInferenceVariables vars 
)

Reads one line from the given file and stores the values there into the variable structure, using the given header array to name the columns.

Returns 0 on success.

Definition at line 2067 of file LALInference.c.

◆ LALInferenceSortVariablesByName()

void LALInferenceSortVariablesByName ( LALInferenceVariables vars)

Sorts the variable structure by name.

Definition at line 2205 of file LALInference.c.

◆ LALInferenceThinnedBufferToArray()

INT4 LALInferenceThinnedBufferToArray ( LALInferenceThreadState thread,
REAL8 **  DEarray,
INT4  step 
)

LALInferenceVariable buffer to array and vica versa.

Definition at line 1591 of file LALInference.c.

◆ LALInferenceBufferToArray()

INT4 LALInferenceBufferToArray ( LALInferenceThreadState thread,
REAL8 **  DEarray 
)

Definition at line 1613 of file LALInference.c.

◆ LALInferenceCopyVariablesToArray()

void LALInferenceCopyVariablesToArray ( LALInferenceVariables origin,
REAL8 target 
)

LALInference variables to an array, and vica versa.

Definition at line 1619 of file LALInference.c.

◆ LALInferenceCopyArrayToVariables()

void LALInferenceCopyArrayToVariables ( REAL8 origin,
LALInferenceVariables target 
)

Definition at line 1685 of file LALInference.c.

◆ LALInferenceLogSampleToFile()

void LALInferenceLogSampleToFile ( LALInferenceVariables algorithmParams,
LALInferenceVariables vars 
)

Append the sample to a file.

file pointer is stored in state->algorithmParams as a LALInferenceVariable called "outfile", as a void ptr. Caller is responsible for opening and closing file. Variables are alphabetically sorted before being written

Definition at line 2235 of file LALInference.c.

◆ LALInferenceLogSampleToArray()

void LALInferenceLogSampleToArray ( LALInferenceVariables algorithmParams,
LALInferenceVariables vars 
)

Append the sample to an array which can be later processed by the user.

Array is stored as a C array in a LALInferenceVariable in state->algorithmParams called "outputarray". Number of items in the array is stored as "N_outputarray". Will create the array and store it in this way if it does not exist. DOES NOT FREE ARRAY, user must clean up after use. Also outputs sample to disk if possible using LALInferenceLogSampleToFile()

Array is stored as a C array in a LALInferenceVariable in state->algorithmParams called "outputarray". Number of items in the array is stored as "N_outputarray". Will create the array and store it in this way if it does not exist. DOES NOT FREE ARRAY, user must clean up after use. Also outputs sample to disk if possible

Definition at line 2255 of file LALInference.c.

◆ LALInferenceMcEta2Masses()

void LALInferenceMcEta2Masses ( double  mc,
double  eta,
double *  m1,
double *  m2 
)

Convert from Mc, eta space to m1, m2 space (note m1 > m2).

Definition at line 2292 of file LALInference.c.

◆ LALInferenceMcQ2Masses()

void LALInferenceMcQ2Masses ( double  mc,
double  q,
double *  m1,
double *  m2 
)

Convert from Mc, q space to m1, m2 space (q = m2/m1, with m1 > m2).

Definition at line 2304 of file LALInference.c.

◆ LALInferenceQ2Eta()

void LALInferenceQ2Eta ( double  q,
double *  eta 
)

Convert from q to eta (q = m2/m1, with m1 > m2).

Definition at line 2314 of file LALInference.c.

◆ LALInferencedQuadMonSdQuadMonA()

void LALInferencedQuadMonSdQuadMonA ( REAL8  dQuadMonS,
REAL8  dQuadMonA,
REAL8 dQuadMon1,
REAL8 dQuadMon2 
)

Convert from dQuadMonS, dQuadMonA to dQuadMon1, dQuadMon2.

Definition at line 2322 of file LALInference.c.

◆ LALInferenceLambdaTsEta2Lambdas()

void LALInferenceLambdaTsEta2Lambdas ( REAL8  lambdaT,
REAL8  dLambdaT,
REAL8  eta,
REAL8 lambda1,
REAL8 lambda2 
)

Convert from lambdaT, dLambdaT, and eta to lambda1 and lambda2.

Definition at line 2329 of file LALInference.c.

◆ LALInferenceLogp1GammasMasses2Lambdas()

void LALInferenceLogp1GammasMasses2Lambdas ( REAL8  logp1,
REAL8  gamma1,
REAL8  gamma2,
REAL8  gamma3,
REAL8  mass1,
REAL8  mass2,
REAL8 lambda1,
REAL8 lambda2 
)

Calculate lambda1,2(m1,2|eos(logp1,gamma1,gamma2,gamma3))

Definition at line 2340 of file LALInference.c.

◆ LALInferenceSDGammasMasses2Lambdas()

void LALInferenceSDGammasMasses2Lambdas ( REAL8  gamma[],
REAL8  mass1,
REAL8  mass2,
REAL8 lambda1,
REAL8 lambda2,
int  size 
)

Convert from spectral parameters to lambda1, lambda2.

Definition at line 2370 of file LALInference.c.

◆ LALInferenceEOSPhysicalCheck()

int LALInferenceEOSPhysicalCheck ( LALInferenceVariables params,
ProcessParamsTable commandLine 
)

Check for causality violation and mass conflict given masses and eos.

Definition at line 2410 of file LALInference.c.

◆ AdiabaticIndex()

double AdiabaticIndex ( double  gamma[],
double  x,
int  size 
)

Specral decomposition of eos's adiabatic index.

Definition at line 2614 of file LALInference.c.

◆ LALInferenceSDGammaCheck()

int LALInferenceSDGammaCheck ( double  gamma[],
int  size 
)

Determine if the Adiabatic index is within 'prior' range.

Definition at line 2570 of file LALInference.c.

◆ LALInferenceKDTreeDelete()

void LALInferenceKDTreeDelete ( LALInferenceKDTree tree)

Delete a kD-tree.

Also deletes all contained cells, and points.

Definition at line 2668 of file LALInference.c.

◆ LALInferenceKDEmpty()

LALInferenceKDTree* LALInferenceKDEmpty ( REAL8 lowerLeft,
REAL8 upperRight,
size_t  ndim 
)

Constructs a fresh, empty kD tree.

The top-level cell will get the given bounds, which should enclose every point added by LALInferenceKDAddPoint().

Definition at line 2713 of file LALInference.c.

◆ LALInferenceKDAddPoint()

int LALInferenceKDAddPoint ( LALInferenceKDTree tree,
REAL8 pt 
)

Adds a point to the kD-tree, returns 0 on successful exit.

The memory for pt is owned by the tree, so should not be deallocated or modified except by LALInferenceKDTreeDelete().

Definition at line 3072 of file LALInference.c.

◆ LALInferenceKDFindCell()

LALInferenceKDTree* LALInferenceKDFindCell ( LALInferenceKDTree tree,
REAL8 pt,
size_t  Npts 
)

Returns the first cell that contains the given point that also contains fewer than Npts points, if possible.

If no cell containing the given point has fewer than Npts points, then returns the cell containing the fewest number of points and the given point. Non-positive Npts will give the fewest-point cell in the tree containing the given point. Returns NULL on error.

Definition at line 3112 of file LALInference.c.

◆ LALInferenceKDLogCellVolume()

double LALInferenceKDLogCellVolume ( LALInferenceKDTree cell)

Returns the log of the volume of the given cell, which is part of the given tree.

Definition at line 3116 of file LALInference.c.

◆ LALInferenceKDLogCellEigenVolume()

double LALInferenceKDLogCellEigenVolume ( LALInferenceKDTree cell)

Returns the log of the volume of the box aligned with the principal axes of the points in the given cell that tightly encloses those points.

Definition at line 3127 of file LALInference.c.

◆ LALInferenceKDVariablesToREAL8()

void LALInferenceKDVariablesToREAL8 ( LALInferenceVariables params,
REAL8 pt,
LALInferenceVariables templt 
)

Fills in the given REAL8 array with the parameter values from params; the ordering of the variables is taken from the order of the non-fixed variables in templt.

It is an error if pt does not point to enough storage to store all the non-fixed parameters from templt and params.

Definition at line 3143 of file LALInference.c.

◆ LALInferenceKDREAL8ToVariables()

void LALInferenceKDREAL8ToVariables ( LALInferenceVariables params,
REAL8 pt,
LALInferenceVariables templt 
)

Fills in the non-fixed variables in params from the given REAL8 array.

The ordering of variables is given by the order of the non-fixed variables in templt.

Definition at line 3155 of file LALInference.c.

◆ LALInferenceKDDrawEigenFrame()

void LALInferenceKDDrawEigenFrame ( gsl_rng *  rng,
LALInferenceKDTree tree,
REAL8 pt,
size_t  Npts 
)

Draws a pt uniformly from a randomly chosen cell of tree.

The chosen cell will be chosen to have (as nearly as possible) Npts in it.

Definition at line 3167 of file LALInference.c.

◆ LALInferenceKDLogProposalRatio()

REAL8 LALInferenceKDLogProposalRatio ( LALInferenceKDTree tree,
REAL8 current,
REAL8 proposed,
size_t  Npts 
)

Returns the log of the jump proposal probability ratio for the LALInferenceKDDrawEigenFrame() proposal to propose the point proposed given the current position current , where Npts is the parameter used to select the box to draw from in LALInferenceKDDrawEigenFrame().

Definition at line 3236 of file LALInference.c.

◆ LALInferenceCheckPositiveDefinite()

UINT4 LALInferenceCheckPositiveDefinite ( gsl_matrix *  matrix,
UINT4  dim 
)

Check matrix is positive definite.

dim is matrix dimensions

Definition at line 3281 of file LALInference.c.

◆ XLALMultiNormalDeviates()

void XLALMultiNormalDeviates ( REAL4Vector vector,
gsl_matrix *  matrix,
UINT4  dim,
RandomParams randParam 
)

Draw a random multivariate vector from Gaussian distr given covariance matrix.

Definition at line 3323 of file LALInference.c.

◆ XLALMultiStudentDeviates()

void XLALMultiStudentDeviates ( REAL4Vector vector,
gsl_matrix *  matrix,
UINT4  dim,
UINT4  n,
RandomParams randParam 
)

Draw a random multivariate vector from student-t distr given covariance matrix.

Definition at line 3374 of file LALInference.c.

◆ LALInferenceAngularDistance()

REAL8 LALInferenceAngularDistance ( REAL8  a1,
REAL8  a2 
)

Calculate shortest angular distance between a1 and a2 (modulo 2PI)

Definition at line 3425 of file LALInference.c.

◆ LALInferenceAngularVariance()

REAL8 LALInferenceAngularVariance ( LALInferenceVariables **  list,
const char pname,
int  N 
)

Calculate the variance of a distribution on an angle (modulo 2PI)

Definition at line 3431 of file LALInference.c.

◆ LALInferenceSanityCheck()

INT4 LALInferenceSanityCheck ( LALInferenceRunState state)

Sanity check the structures in the given state.

Will scan data for infinities and nans, and look for null pointers.

Definition at line 3450 of file LALInference.c.

◆ LALInferenceDumpWaveforms()

void LALInferenceDumpWaveforms ( LALInferenceModel model,
const char basefilename 
)

Dump all waveforms from the ifodata structure.

(currently: timeData, freqData) basefilename is optional text to append to file names

Definition at line 3590 of file LALInference.c.

◆ LALInferenceWriteVariablesBinary()

int LALInferenceWriteVariablesBinary ( FILE *  file,
LALInferenceVariables vars 
)

Write a LALInferenceVariables as binary to a given FILE pointer, returns the number of items written (should be the dimension of the variables) or -1 on error.

Definition at line 3695 of file LALInference.c.

◆ LALInferenceReadVariablesBinary()

LALInferenceVariables* LALInferenceReadVariablesBinary ( FILE *  stream)

Read from the given FILE * a LALInferenceVariables, which was previously written with LALInferenceWriteVariablesBinary() Returns a new LALInferenceVariables.

Definition at line 3778 of file LALInference.c.

◆ LALInferenceWriteVariablesArrayBinary()

int LALInferenceWriteVariablesArrayBinary ( FILE *  file,
LALInferenceVariables **  vars,
UINT4  N 
)

Write an array N of LALInferenceVariables to the given FILE * using LALInferenceWriteVariablesBinary().

Returns the number written (should be ==N)

Definition at line 3898 of file LALInference.c.

◆ LALInferenceReadVariablesArrayBinary()

int LALInferenceReadVariablesArrayBinary ( FILE *  file,
LALInferenceVariables **  vars,
UINT4  N 
)

Read N LALInferenceVariables from the binary FILE *file, previously written with LALInferenceWriteVariablesArrayBinary() returns the number read.

Definition at line 3910 of file LALInference.c.

◆ LALInferenceWriteRunStateBinary()

int LALInferenceWriteRunStateBinary ( FILE *  file,
LALInferenceRunState state 
)

Write the LALInferenceVariables contents of runState to a file in binary format.

Definition at line 3921 of file LALInference.c.

◆ LALInferenceReadRunStateBinary()

int LALInferenceReadRunStateBinary ( FILE *  file,
LALInferenceRunState state 
)

Reads the file and populates LALInferenceVariables in the runState that were saved with LALInferenceReadVariablesArrayBinary()

Definition at line 3931 of file LALInference.c.

◆ LALInferenceAddINT4Variable()

void LALInferenceAddINT4Variable ( LALInferenceVariables vars,
const char name,
INT4  value,
LALInferenceParamVaryType  vary 
)

Definition at line 3941 of file LALInference.c.

◆ LALInferenceGetINT4Variable()

INT4 LALInferenceGetINT4Variable ( LALInferenceVariables vars,
const char name 
)

Definition at line 3947 of file LALInference.c.

◆ LALInferenceSetINT4Variable()

void LALInferenceSetINT4Variable ( LALInferenceVariables vars,
const char name,
INT4  value 
)

Definition at line 3960 of file LALInference.c.

◆ LALInferenceAddINT8Variable()

void LALInferenceAddINT8Variable ( LALInferenceVariables vars,
const char name,
INT8  value,
LALInferenceParamVaryType  vary 
)

Definition at line 3964 of file LALInference.c.

◆ LALInferenceGetINT8Variable()

INT8 LALInferenceGetINT8Variable ( LALInferenceVariables vars,
const char name 
)

Definition at line 3970 of file LALInference.c.

◆ LALInferenceSetINT8Variable()

void LALInferenceSetINT8Variable ( LALInferenceVariables vars,
const char name,
INT8  value 
)

Definition at line 3983 of file LALInference.c.

◆ LALInferenceAddUINT4Variable()

void LALInferenceAddUINT4Variable ( LALInferenceVariables vars,
const char name,
UINT4  value,
LALInferenceParamVaryType  vary 
)

Definition at line 3987 of file LALInference.c.

◆ LALInferenceGetUINT4Variable()

UINT4 LALInferenceGetUINT4Variable ( LALInferenceVariables vars,
const char name 
)

Definition at line 3993 of file LALInference.c.

◆ LALInferenceSetUINT4Variable()

void LALInferenceSetUINT4Variable ( LALInferenceVariables vars,
const char name,
UINT4  value 
)

Definition at line 4006 of file LALInference.c.

◆ LALInferenceAddREAL4Variable()

void LALInferenceAddREAL4Variable ( LALInferenceVariables vars,
const char name,
REAL4  value,
LALInferenceParamVaryType  vary 
)

Definition at line 4010 of file LALInference.c.

◆ LALInferenceGetREAL4Variable()

REAL4 LALInferenceGetREAL4Variable ( LALInferenceVariables vars,
const char name 
)

Definition at line 4016 of file LALInference.c.

◆ LALInferenceSetREAL4Variable()

void LALInferenceSetREAL4Variable ( LALInferenceVariables vars,
const char name,
REAL4  value 
)

Definition at line 4029 of file LALInference.c.

◆ LALInferenceAddREAL8Variable()

void LALInferenceAddREAL8Variable ( LALInferenceVariables vars,
const char name,
REAL8  value,
LALInferenceParamVaryType  vary 
)

Definition at line 4033 of file LALInference.c.

◆ LALInferenceGetREAL8Variable()

REAL8 LALInferenceGetREAL8Variable ( LALInferenceVariables vars,
const char name 
)

Definition at line 4039 of file LALInference.c.

◆ LALInferenceSetREAL8Variable()

void LALInferenceSetREAL8Variable ( LALInferenceVariables vars,
const char name,
REAL8  value 
)

Definition at line 4052 of file LALInference.c.

◆ LALInferenceAddCOMPLEX8Variable()

void LALInferenceAddCOMPLEX8Variable ( LALInferenceVariables vars,
const char name,
COMPLEX8  value,
LALInferenceParamVaryType  vary 
)

Definition at line 4056 of file LALInference.c.

◆ LALInferenceGetCOMPLEX8Variable()

COMPLEX8 LALInferenceGetCOMPLEX8Variable ( LALInferenceVariables vars,
const char name 
)

Definition at line 4062 of file LALInference.c.

◆ LALInferenceSetCOMPLEX8Variable()

void LALInferenceSetCOMPLEX8Variable ( LALInferenceVariables vars,
const char name,
COMPLEX8  value 
)

Definition at line 4071 of file LALInference.c.

◆ LALInferenceAddCOMPLEX16Variable()

void LALInferenceAddCOMPLEX16Variable ( LALInferenceVariables vars,
const char name,
COMPLEX16  value,
LALInferenceParamVaryType  vary 
)

Definition at line 4075 of file LALInference.c.

◆ LALInferenceGetCOMPLEX16Variable()

COMPLEX16 LALInferenceGetCOMPLEX16Variable ( LALInferenceVariables vars,
const char name 
)

Definition at line 4081 of file LALInference.c.

◆ LALInferenceSetCOMPLEX16Variable()

void LALInferenceSetCOMPLEX16Variable ( LALInferenceVariables vars,
const char name,
COMPLEX16  value 
)

Definition at line 4090 of file LALInference.c.

◆ LALInferenceAddgslMatrixVariable()

void LALInferenceAddgslMatrixVariable ( LALInferenceVariables vars,
const char name,
gsl_matrix *  value,
LALInferenceParamVaryType  vary 
)

Definition at line 4094 of file LALInference.c.

◆ LALInferenceGetgslMatrixVariable()

gsl_matrix* LALInferenceGetgslMatrixVariable ( LALInferenceVariables vars,
const char name 
)

Definition at line 4100 of file LALInference.c.

◆ LALInferenceSetgslMatrixVariable()

void LALInferenceSetgslMatrixVariable ( LALInferenceVariables vars,
const char name,
gsl_matrix *  value 
)

Definition at line 4113 of file LALInference.c.

◆ LALInferenceAddREAL8VectorVariable()

void LALInferenceAddREAL8VectorVariable ( LALInferenceVariables vars,
const char name,
REAL8Vector value,
LALInferenceParamVaryType  vary 
)

Definition at line 4117 of file LALInference.c.

◆ LALInferenceGetREAL8VectorVariable()

REAL8Vector* LALInferenceGetREAL8VectorVariable ( LALInferenceVariables vars,
const char name 
)

Definition at line 4123 of file LALInference.c.

◆ LALInferenceSetREAL8VectorVariable()

void LALInferenceSetREAL8VectorVariable ( LALInferenceVariables vars,
const char name,
REAL8Vector value 
)

Definition at line 4136 of file LALInference.c.

◆ LALInferenceAddCOMPLEX16VectorVariable()

void LALInferenceAddCOMPLEX16VectorVariable ( LALInferenceVariables vars,
const char name,
COMPLEX16Vector value,
LALInferenceParamVaryType  vary 
)

Definition at line 4140 of file LALInference.c.

◆ LALInferenceGetCOMPLEX16VectorVariable()

COMPLEX16Vector* LALInferenceGetCOMPLEX16VectorVariable ( LALInferenceVariables vars,
const char name 
)

Definition at line 4146 of file LALInference.c.

◆ LALInferenceSetCOMPLEX16VectorVariable()

void LALInferenceSetCOMPLEX16VectorVariable ( LALInferenceVariables vars,
const char name,
COMPLEX16Vector value 
)

Definition at line 4159 of file LALInference.c.

◆ LALInferenceAddINT4VectorVariable()

void LALInferenceAddINT4VectorVariable ( LALInferenceVariables vars,
const char name,
INT4Vector value,
LALInferenceParamVaryType  vary 
)

Definition at line 4169 of file LALInference.c.

◆ LALInferenceAddUINT4VectorVariable()

void LALInferenceAddUINT4VectorVariable ( LALInferenceVariables vars,
const char name,
UINT4Vector value,
LALInferenceParamVaryType  vary 
)

Definition at line 4163 of file LALInference.c.

◆ LALInferenceGetINT4VectorVariable()

INT4Vector* LALInferenceGetINT4VectorVariable ( LALInferenceVariables vars,
const char name 
)

Definition at line 4188 of file LALInference.c.

◆ LALInferenceGetUINT4VectorVariable()

UINT4Vector* LALInferenceGetUINT4VectorVariable ( LALInferenceVariables vars,
const char name 
)

Definition at line 4175 of file LALInference.c.

◆ LALInferenceSetINT4VectorVariable()

void LALInferenceSetINT4VectorVariable ( LALInferenceVariables vars,
const char name,
INT4Vector value 
)

Definition at line 4205 of file LALInference.c.

◆ LALInferenceSetUINT4VectorVariable()

void LALInferenceSetUINT4VectorVariable ( LALInferenceVariables vars,
const char name,
UINT4Vector value 
)

Definition at line 4201 of file LALInference.c.

◆ LALInferenceAddMCMCrunphase_ptrVariable()

void LALInferenceAddMCMCrunphase_ptrVariable ( LALInferenceVariables vars,
const char name,
LALInferenceMCMCRunPhase value,
LALInferenceParamVaryType  vary 
)

Definition at line 4209 of file LALInference.c.

◆ LALInferenceGetMCMCrunphase_ptrVariable()

LALInferenceMCMCRunPhase* LALInferenceGetMCMCrunphase_ptrVariable ( LALInferenceVariables vars,
const char name 
)

Definition at line 4215 of file LALInference.c.

◆ LALInferenceSetMCMCrunphase_ptrVariable()

void LALInferenceSetMCMCrunphase_ptrVariable ( LALInferenceVariables vars,
const char name,
LALInferenceMCMCRunPhase value 
)

Definition at line 4228 of file LALInference.c.

◆ LALInferenceAddstringVariable()

void LALInferenceAddstringVariable ( LALInferenceVariables vars,
const char name,
const CHAR value,
LALInferenceParamVaryType  vary 
)

Definition at line 4232 of file LALInference.c.

◆ LALInferenceGetstringVariable()

const CHAR* LALInferenceGetstringVariable ( LALInferenceVariables vars,
const char name 
)

Definition at line 4238 of file LALInference.c.

◆ LALInferenceSetstringVariable()

void LALInferenceSetstringVariable ( LALInferenceVariables vars,
const char name,
const CHAR value 
)

Definition at line 4251 of file LALInference.c.

◆ LALInferenceFprintSplineCalibrationHeader()

void LALInferenceFprintSplineCalibrationHeader ( FILE *  output,
LALInferenceThreadState thread 
)

Print spline calibration parameter names as tab-separated ASCII.

Definition at line 4431 of file LALInference.c.

◆ LALInferenceBinaryLove()

void LALInferenceBinaryLove ( LALInferenceVariables vars,
REAL8 lambda1,
REAL8 lambda2 
)

Compute Tidal deformabilities following BinaryLove Universal relations.

Implelentation of the parametrization from Chatziioannou, Haster, Zimmerman (2018) https://doi.org/10.1103/PhysRevD.97.104036

Definition at line 4488 of file LALInference.c.

◆ LALInferenceDetFrameToEquatorial()

void LALInferenceDetFrameToEquatorial ( LALDetector det0,
LALDetector det1,
REAL8  t0,
REAL8  alpha,
REAL8  theta,
REAL8 tg,
REAL8 ra,
REAL8 dec 
)

Conversion routines between Equatorial (RA,DEC) and detector-based coordinate systems, where new "north pole" points along vector from det0 to det1.

theta - azimuth angle about vector joining det0 and det1 alpha - co-latitude (0,pi) relative to det0-det1 vector

Definition at line 18 of file DetectorFixedSkyCoords.c.

◆ LALInferenceEquatorialToDetFrame()

void LALInferenceEquatorialToDetFrame ( LALDetector det0,
LALDetector det1,
REAL8  tg,
REAL8  ra,
REAL8  dec,
REAL8 t0,
REAL8 alpha,
REAL8 theta 
)

Definition at line 101 of file DetectorFixedSkyCoords.c.

Typedef Documentation

◆ LALInferenceTemplateFunction

typedef void(* LALInferenceTemplateFunction) (struct tagLALInferenceModel *model)

Type declaration for template function, which operates on a LALInferenceIFOData structure *data.

Definition at line 377 of file LALInference.h.

◆ LALInferenceProposalFunction

typedef REAL8(* LALInferenceProposalFunction) (struct tagLALInferenceThreadState *thread, LALInferenceVariables *currentParams, LALInferenceVariables *proposedParams)

Jump proposal distribution Computes proposedParams based on currentParams and additional variables stored as proposalArgs inside runState, which could include correlation matrix, etc., as well as forward and reverse proposal probability.

The log of the Metropolis-Hasting proposal ratio is returned. A jump proposal distribution function could call other jump proposal distribution functions with various probabilities to allow for multiple jump proposal distributions

Definition at line 391 of file LALInference.h.

◆ LALInferencePriorFunction

typedef REAL8(* LALInferencePriorFunction) (struct tagLALInferenceRunState *runState, LALInferenceVariables *params, struct tagLALInferenceModel *model)

Type declaration for prior function which returns p(params) Can depend on runState ->priorArgs.

Definition at line 399 of file LALInference.h.

◆ LALInferenceCubeToPriorFunction

typedef UINT4(* LALInferenceCubeToPriorFunction) (struct tagLALInferenceRunState *runState, LALInferenceVariables *params, struct tagLALInferenceModel *model, double *cube, void *context)

Type declaration for CubeToPrior function which converts parameters in unit hypercube to their corresponding physical values according to the prior.

Can depend on runState ->priorArgs

Definition at line 407 of file LALInference.h.

◆ LALInferenceInitModelFunction

typedef LALInferenceModel*(* LALInferenceInitModelFunction) (struct tagLALInferenceRunState *runState)

Type declaration for variables init function, can be user-declared.

The function returns a pointer to a new LALInferenceVariables instance Reads runState->commandLine to get user config

Definition at line 475 of file LALInference.h.

◆ LALInferenceLikelihoodFunction

typedef REAL8(* LALInferenceLikelihoodFunction) (LALInferenceVariables *currentParams, struct tagLALInferenceIFOData *data, LALInferenceModel *model)

Type declaration for likelihood function Computes p(data | currentParams, templt ) templt is a LALInferenceTemplateFunction defined below.

Definition at line 487 of file LALInference.h.

◆ LALInferenceEvolveOneStepFunction

typedef INT4(* LALInferenceEvolveOneStepFunction) (struct tagLALInferenceRunState *runState)

Perform one step of an algorithm, replaces runState ->currentParams.

Definition at line 491 of file LALInference.h.

◆ LALInferenceSwapRoutine

typedef void(* LALInferenceSwapRoutine) (struct tagLALInferenceRunState *runState, FILE *)

Propose a swap between chain locations.

Definition at line 494 of file LALInference.h.

◆ LALInferenceAlgorithm

typedef void(* LALInferenceAlgorithm) (struct tagLALInferenceRunState *runState)

Type declaration for an algorithm function which is called by the driver code The user must initialise runState before use.

The Algorithm manipulates

Parameters
runStateto do its work

Definition at line 501 of file LALInference.h.

◆ LALInferenceLogFunction

typedef void(* LALInferenceLogFunction) (LALInferenceVariables *algorithmParams, LALInferenceVariables *vars)

Type declaration for output logging function, can be user-declared.

Definition at line 504 of file LALInference.h.

Enumeration Type Documentation

◆ LALInferenceVariableType

An enumerated type for denoting the type of a variable.

Several LAL types are supported as well as others.

Enumerator
LALINFERENCE_INT4_t 
LALINFERENCE_INT8_t 
LALINFERENCE_UINT4_t 
LALINFERENCE_REAL4_t 
LALINFERENCE_REAL8_t 
LALINFERENCE_COMPLEX8_t 
LALINFERENCE_COMPLEX16_t 
LALINFERENCE_gslMatrix_t 
LALINFERENCE_REAL8Vector_t 
LALINFERENCE_INT4Vector_t 
LALINFERENCE_UINT4Vector_t 
LALINFERENCE_COMPLEX16Vector_t 
LALINFERENCE_string_t 
LALINFERENCE_MCMCrunphase_ptr_t 
LALINFERENCE_void_ptr_t 

Definition at line 104 of file LALInference.h.

◆ LALInferenceParamVaryType

An enumerated type for denoting the topology of a parameter.

This information is used by the sampling routines when deciding what to vary in a proposal, etc.

Enumerator
LALINFERENCE_PARAM_LINEAR 
LALINFERENCE_PARAM_CIRCULAR 

A parameter that simply has a maximum and a minimum.

LALINFERENCE_PARAM_FIXED 

A parameter that is cyclic, such as an angle between 0 and 2pi.

LALINFERENCE_PARAM_OUTPUT 

A parameter that never changes, functions should respect this.

A parameter changed by an inner code and passed out

Definition at line 127 of file LALInference.h.

◆ LALInferenceMCMCRunPhase

Phase of MCMC run (depending on burn-in status, different actions are performed during the run, and this tag controls the activity).

Enumerator
LALINFERENCE_ONLY_PT 
LALINFERENCE_TEMP_PT 

Run only parallel tempers.

LALINFERENCE_ANNEALING 

In the parallel tempering phase of an annealed run.

LALINFERENCE_SINGLE_CHAIN 

In the annealing phase of an annealed run.

LALINFERENCE_LADDER_UPDATE 

In the single-chain sampling phase of an annealed run.

Move all temperature chains to cold chains location (helpful for burnin)

Definition at line 180 of file LALInference.h.

Macro Definition Documentation

◆ VARNAME_MAX

#define VARNAME_MAX   40

Definition at line 50 of file LALInference.h.

◆ VARVALSTRINGSIZE_MAX

#define VARVALSTRINGSIZE_MAX   128

Definition at line 51 of file LALInference.h.

◆ DETNAMELEN

#define DETNAMELEN   8

Definition at line 617 of file LALInference.h.

Variable Documentation

◆ LALInferenceTypeSize

size_t LALInferenceTypeSize[15]
extern

Definition at line 86 of file LALInference.c.