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
LALInference.h
Go to the documentation of this file.
1/*
2 *
3 * LALInference: Bayesian Followup
4 * include/LALInference.h: main header file
5 *
6 * Copyright (C) 2009,2012 Ilya Mandel, Vivien Raymond, Christian
7 * Roever, Marc van der Sluys, John Veitch, and Will M. Farr
8 *
9 *
10 * This program is free software; you can redistribute it and/or modify
11 * it under the terms of the GNU General Public License as published by
12 * the Free Software Foundation; either version 2 of the License, or
13 * (at your option) any later version.
14 *
15 * This program is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 * GNU General Public License for more details.
19 *
20 * You should have received a copy of the GNU General Public License
21 * along with with program; see the file COPYING. If not, write to the
22 * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
23 * MA 02110-1301 USA
24 */
25#ifndef LALInference_h
26#define LALInference_h
27
28/**
29 * \defgroup LALInference_h Header LALInference.h
30 * \ingroup lalinference_general
31 * \brief Main header file for LALInference common routines and structures
32 *
33 * LALInference is a Bayesian analysis toolkit for use with LAL. It contains
34 * common requirements for Bayesian codes such as Likelihood functions, data
35 * handling routines, MCMC and Nested Sampling algorithms and a template generation
36 * interface to the LALInspiral package.
37 *
38 * This file contains the basic structures for the algorithm state, interferometer
39 * data, manipulation of variables and type declarations for the standard function types.
40 *
41 */
42/** @{ */
43
44#include <math.h>
45#include <stdio.h>
46#include <stdlib.h>
47#include <gsl/gsl_errno.h>
48#include <gsl/gsl_spline.h>
49
50#define VARNAME_MAX 40
51#define VARVALSTRINGSIZE_MAX 128
52
53#include <lal/LALStdlib.h>
54#include <lal/LALConstants.h>
55#include <lal/SimulateCoherentGW.h>
56#include <lal/GeneratePPNInspiral.h>
57#include <lal/LIGOMetadataTables.h>
58#include <lal/LALDatatypes.h>
59#include <lal/FindChirp.h>
60#include <lal/Window.h>
61#include <lal/LALString.h>
62#include <lal/StringInput.h>
63#include <lal/StringVector.h>
64#include <lal/LALSimInspiral.h>
65#include <lal/LALSimInspiralWaveformCache.h>
66#include <lal/LALSimNeutronStar.h>
67#include <lal/LALHashTbl.h>
68#include <lal/RealFFT.h>
69#include <lal/LALDetectors.h>
70
71#include <gsl/gsl_linalg.h>
72#include <gsl/gsl_errno.h>
73#include <gsl/gsl_math.h>
74#include <gsl/gsl_min.h>
75#include <gsl/gsl_vector.h>
76#include <gsl/gsl_matrix.h>
77#include <gsl/gsl_blas.h>
78#include <gsl/gsl_linalg.h>
79#include <gsl/gsl_eigen.h>
80#include <gsl/gsl_rng.h>
81#include <gsl/gsl_randist.h>
82#include <gsl/gsl_statistics.h>
83#include <gsl/gsl_complex_math.h>
84#include <gsl/gsl_cdf.h>
85
86#include <sys/time.h>
87
88/*LIB imports*/
89#include <lal/LALInferenceBurstRoutines.h>
90
91//...other includes
92
93struct tagLALInferenceRunState;
94struct tagLALInferenceThreadState;
95struct tagLALInferenceIFOData;
96struct tagLALInferenceModel;
97
98/*Data storage type definitions*/
99
100/**
101 * An enumerated type for denoting the type of a variable. Several LAL
102 * types are supported as well as others.
103 */
104typedef enum tagLALInferenceVariableType {
121
122/**
123 * An enumerated type for denoting the topology of a parameter.
124 * This information is used by the sampling routines when deciding
125 * what to vary in a proposal, etc.
126 */
127typedef enum tagLALInferenceParamVaryType {
128 LALINFERENCE_PARAM_LINEAR, /** A parameter that simply has a maximum and a minimum */
129 LALINFERENCE_PARAM_CIRCULAR, /** A parameter that is cyclic, such as an angle between 0 and 2pi */
130 LALINFERENCE_PARAM_FIXED, /** A parameter that never changes, functions should respect this */
131 LALINFERENCE_PARAM_OUTPUT /** A parameter changed by an inner code and passed out */
133
134extern size_t LALInferenceTypeSize[15];
135
136/**
137 * The LALInferenceVariableItem list node structure
138 * This should only be accessed using the accessor functions below
139 * Implementation may change to hash table so please use only the
140 * accessor functions below.
141 *
142 * A note on memory ownership: each LALInferenceVariableItem() owns
143 * the memory for its contents. This means after calling
144 * LALInferenceAddVariable() or LALInferenceSetVariable(), you should
145 * not free the memory assoicated with the new variable's value. In
146 * fact, you should not access this memory through its original
147 * pointer again (access through pointers returned by
148 * LALInferenceGetVariable() is OK); the object may not be live in
149 * memory if its variable has been overwritten through another call to
150 * LALInferenceSetVariable().
151 */
152typedef struct
153tagVariableItem
154{
156 void *value;
159 struct tagVariableItem *next;
161
162
163/**
164 * The LALInferenceVariables structure to contain a set of parameters
165 * Implemented as a linked list of LALInferenceVariableItems.
166 * Should only be accessed using the accessor functions below
167 */
168typedef struct
169tagLALInferenceVariables
170{
173 LALHashTbl *hash_table;
175
176/**
177 * Phase of MCMC run (depending on burn-in status, different actions
178 * are performed during the run, and this tag controls the activity).
179 */
180typedef enum tagLALInferenceMCMCRunPhase {
181 LALINFERENCE_ONLY_PT, /** Run only parallel tempers. */
182 LALINFERENCE_TEMP_PT, /** In the parallel tempering phase of an annealed run */
183 LALINFERENCE_ANNEALING, /** In the annealing phase of an annealed run */
184 LALINFERENCE_SINGLE_CHAIN, /** In the single-chain sampling phase of an annealed run */
185 LALINFERENCE_LADDER_UPDATE /** Move all temperature chains to cold chains location (helpful for burnin) */
187
188/** Returns an array of header strings (terminated by NULL) from a common-format output file */
189char **LALInferenceGetHeaderLine(FILE *inp);
190
191/** Converts between internally used parameter names and those external (e.g. in SimInspiralTable?) */
192const char *LALInferenceTranslateInternalToExternalParamName(const char *inName);
193
194/** Converts between externally used parameter names and those internal */
195void LALInferenceTranslateExternalToInternalParamName(char *outName, const char *inName);
196
197/**
198 * Print the parameter names to a file as a tab-separated ASCII line
199 * \param out [in] pointer to output file
200 * \param params [in] LALInferenceVaraibles structure to print
201 */
203
204/**
205 * Print the parameters which do not vary to a file as a tab-separated ASCII line
206 * \param out [in] pointer to output file
207 * \param params [in] LALInferenceVaraibles structure to print
208 */
210
211/**
212 * Print the parameters which do not vary to a file as a tab-separated ASCII line, adding the given suffix
213 * \param out [in] pointer to output file
214 * \param params [in] LALInferenceVaraibles structure to print
215 * \param suffix [in] Suffix string to add to each parameter name
216 */
218
219/** Prints a variable item to a string. Print at most N characters. Returns the number of characters actually required
220 * to store the output (can be more or less than N) */
222
223/**
224 * Return a pointer to the memory the variable \c vars is stored in specified by \c name
225 * User must cast this pointer to the expected type before dereferencing
226 * it to get the value of the variable.
227 */
228void *LALInferenceGetVariable(const LALInferenceVariables * vars, const char * name);
229
230/** Get number of dimensions in variable \c vars */
232
233/** Get number of dimensions in \c vars which are not fixed to a certain value */
235
236/** Get number of dimensions in \c vars which are not fixed to a certain value,
237 * with a flag for skipping counting vectors */
239
240/**
241 * Get the LALInferenceVariableType of the \c idx -th item in the \c vars
242 * Indexing starts at 1
243 */
245
246/** Get the LALInferenceVariableType of the parameter named \c name in \c vars */
248
249/**
250 * Get the LALInferenceParamVaryType of the parameter named \c name in \c vars
251 * see the declaration of LALInferenceParamVaryType for possibilities
252 */
254
255/**
256 * Set the LALInferenceParamVaryType of the parameter named \c name in \c vars,
257 * see the declaration of LALInferenceParamVaryType for possibilities
258 */
260
261/**
262 * Get the name of the idx-th variable
263 * Indexing starts at 1
264 */
266
267/**
268 * Set a variable named \c name in \c vars with a value.
269 * Pass a void * in \c value to the value you wish to set,
270 * i.e. LALInferenceSetVariable(vars, "mu", (void *)&mu);
271 */
272void LALInferenceSetVariable(LALInferenceVariables * vars, const char * name, const void * value);
273
274/**
275 * Add a variable named \c name to \c vars with initial value referenced by \c value
276 * \param type is a LALInferenceVariableType (enumerated above)
277 * \param vary is a LALInferenceParamVaryType (enumerated above)
278 * \param vars UNDOCUMENTED
279 * \param name UNDOCUMENTED
280 * \param value UNDOCUMENTED
281 * If the variable already exists it will be over-written UNLESS IT HAS A CONFLICTING TYPE
282 */
283void LALInferenceAddVariable(LALInferenceVariables * vars, const char * name, const void * value,
285
286/**
287 * Remove \c name from \c vars
288 * Frees the memory for the \c name structure and its contents
289 */
290void LALInferenceRemoveVariable(LALInferenceVariables *vars,const char *name);
291
292/**
293 * Checks for \c name being present in \c vars
294 * returns 1(==true) or 0
295 */
296int LALInferenceCheckVariable(LALInferenceVariables *vars,const char *name);
297
298/**
299 * Checks for \c name being present in \c vars and having type LINEAR or CIRCULAR.
300 * returns 1 or 0
301 */
304/**
305 * Delete the variables in this structure.
306 * Does not free the LALInferenceVariables itself
307 * \param vars will have its dimension set to 0
308 */
310
311/** Deep copy the variables from one to another LALInferenceVariables structure */
313
314/* Copy REAL8s from "origin" to "target" if they weren't set on the command line */
316
317/** Print variables to stdout */
319
320/** Check for equality in two variables */
322
323/** Computes the factor relating the physical waveform to a measured
324 waveform for a spline-fit calibration model in amplitude and
325 phase. The spline points can be arbitrary frequencies, and the
326 values of the calibration offset at these points can be arbitary,
327 too. For amplitude offset, \f$\delta A\f$, and phase offset,
328 \f$\delta \phi\f$, the measured waveform is related to the
329 physical waveform via
330
331 \f[
332 h_\mathrm{meas} = h_\mathrm{phys} \left(1 + \delta A \right) \frac{2 + i \delta \phi}{2 - i \delta \phi}
333 \f]
334
335 The phase factor takes the form above rather than the more obvious
336 \f$\exp(i \delta \phi)\f$ or \f$1 + \delta \phi\f$ because it is
337 faster to compute than the former (but, for small \f$\phi\f$,
338 equivalent through second order in \f$\delta \phi\f$) and, unlike
339 the latter, it ensures that the complex amplitude of the
340 correction is always 1. (A similar technique is used when using
341 finite-difference approximations to solve the multi-dimensional
342 Schrodinger equation to ensure that the evolution remains
343 unitary.) Note that this implies that the phase calibration
344 parameter ranges over \f$\delta \phi \in [-\infty, \infty]\f$.
345
346 The values of \f$\delta A\f$ and \f$\delta \phi\f$ at arbitrary
347 frequency are obtained by a spline curve that passes through the
348 given values at the given frequencies.
349
350*/
352 REAL8Vector *deltaAmps,
353 REAL8Vector *deltaPhases,
354 COMPLEX16FrequencySeries *calFactor);
355
356 /** Modified version of LALInferenceSplineCalibrationFactor to compute the
357 * calibration factors for the specific frequency nodes used for
358 * Reduced Order Quadrature likelihoods.
359 */
360
362 REAL8Vector *deltaAmps,
363 REAL8Vector *deltaPhases,
364 REAL8Sequence *freqNodesLin,
365 COMPLEX16Sequence **calFactorROQLin,
366 REAL8Sequence *freqNodesQuad,
367 COMPLEX16Sequence **calFactorROQQuad);
368
369
370//Wrapper for template computation
371//(relies on LAL libraries for implementation) <- could be a #DEFINE ?
372//typedef void (LALTemplateFunction) (LALInferenceVariables *currentParams, struct tagLALInferenceIFOData *data); //Parameter Set is modelParams of LALInferenceIFOData
373/**
374 * Type declaration for template function, which operates on
375 * a LALInferenceIFOData structure \c *data
376 */
377typedef void (*LALInferenceTemplateFunction) (struct tagLALInferenceModel *model);
378
379
380/**
381 * Jump proposal distribution
382 * Computes \c proposedParams based on \c currentParams
383 * and additional variables stored as proposalArgs inside \c runState,
384 * which could include correlation matrix, etc.,
385 * as well as forward and reverse proposal probability. The log of the
386 * Metropolis-Hasting proposal ratio is returned.
387 * A jump proposal distribution function could call other jump proposal
388 * distribution functions with various probabilities to allow for multiple
389 * jump proposal distributions
390 */
391typedef REAL8 (*LALInferenceProposalFunction) (struct tagLALInferenceThreadState *thread,
393 LALInferenceVariables *proposedParams);
394
395/**
396 * Type declaration for prior function which returns p(\c params)
397 * Can depend on \c runState ->priorArgs
398 */
399typedef REAL8 (*LALInferencePriorFunction) (struct tagLALInferenceRunState *runState,
400 LALInferenceVariables *params, struct tagLALInferenceModel *model);
401
402/**
403 * Type declaration for CubeToPrior function which converts parameters in unit hypercube
404 * to their corresponding physical values according to the prior.
405 * Can depend on \c runState ->priorArgs
406 */
407typedef UINT4 (*LALInferenceCubeToPriorFunction) (struct tagLALInferenceRunState *runState,
408 LALInferenceVariables *params, struct tagLALInferenceModel *model, double *cube, void *context);
409
410/**
411 * Detector-dependent buffers and parameters
412 * A linked list meant for parameters and buffers that are separately specified for
413 * each detector. The list ordering should be the same as the IFO linked list.
414 */
415typedef struct
416tagLALInferenceIFOModel
417{
418 LALInferenceVariables *params; /** Parameters used when filling the buffers - template functions should copy to here */
419
420 LALDetector *detector; /** LALDetector structure for where this data came from */
421
422 void *extraData; /** Pointer to extra detector-dependent parameters, used by pulsar analyses */
423
424 REAL8TimeSeries *timehPlus, *timehCross; /** Time series model buffers */
425 COMPLEX16FrequencySeries *freqhPlus, *freqhCross; /** Freq series model buffers */
426 COMPLEX16TimeSeries *compTimeSignal; /** Complex time series signal buffer */
427 REAL8TimeSeries *timeData; /** Time series model buffer */
428
429 struct tagLALInferenceIFOModel *next; /** A pointer to the next set of parameters for linked list */
431
432/**
433 * Structure to constain a model and its parameters.
434 */
435typedef struct tagLALInferenceModel
436{
437 LALInferenceVariables *params; /** Parameters used when filling the buffers - template functions should copy to here */
438 LALInferenceIFOModel *ifo; /** IFO-dependent parameters and buffers */
439 LALSimulationDomain domain; /** Domain of model */
440 LALInferenceTemplateFunction templt; /** The template generation function */
441
442 REAL8 logprior; /** Prior value at *params* */
443 REAL8 loglikelihood; /** Likelihood value at *params* */
444 REAL8 SNR; /** Network SNR at *params* */
445 REAL8* ifo_loglikelihoods; /** Array of single-IFO likelihoods at *params* */
446 REAL8* ifo_SNRs; /** Array of single-IFO SNRs at *params* */
447
448 REAL8 fLow; /** Start frequency for waveform generation */
449 REAL8 fHigh; /** End frequency for waveform generation */
450 REAL8 deltaT, deltaF; /** Sampling rate information */
451 INT4 freqLength; /* Length of freq-domain buffer */
452
453 REAL8TimeSeries *timehPlus, *timehCross; /** Time series model buffers */
454 COMPLEX16FrequencySeries *freqhPlus, *freqhCross; /** Freq series model buffers */
455 COMPLEX16FrequencySeries **freqhs; /** Projected freq series model buffers */
456
457 LALDict *LALpars;
459 LALSimBurstWaveformCache *burstWaveformCache; /** Burst Waveform cache for LIB*/
460 REAL8FFTPlan *timeToFreqFFTPlan, *freqToTimeFFTPlan; /** Pre-calculated FFT plans for forward and reverse FFTs */
461 REAL8Window *window; /** A window */
462 REAL8 padding; /** The padding of the above window */
463 struct tagLALInferenceROQModel *roq; /** ROQ data */
464 int roq_flag; /** Is ROQ enabled */
465 LALSimNeutronStarFamily *eos_fam; /** Neutron Star equation of state family */
466
468
469
470/**
471 * Type declaration for variables init function, can be user-declared.
472 * The function returns a pointer to a new LALInferenceVariables instance
473 * Reads \c runState->commandLine to get user config
474 */
475typedef LALInferenceModel* (*LALInferenceInitModelFunction) (struct tagLALInferenceRunState *runState);
476
477
478//Likelihood calculator
479//Should take care to perform expensive evaluation of h+ and hx
480//only once if possible, unless necessary because different IFOs
481//have different data lengths or sampling rates
482/**
483 * Type declaration for likelihood function
484 * Computes p(\c data | \c currentParams, \c templt )
485 * templt is a LALInferenceTemplateFunction defined below
486 */
488 struct tagLALInferenceIFOData * data, LALInferenceModel *model);
489
490/** Perform one step of an algorithm, replaces \c runState ->currentParams */
491typedef INT4 (*LALInferenceEvolveOneStepFunction) (struct tagLALInferenceRunState *runState);
492
493/** Propose a swap between chain locations */
494typedef void (*LALInferenceSwapRoutine) (struct tagLALInferenceRunState *runState, FILE *);
495
496/**
497 * Type declaration for an algorithm function which is called by the driver code
498 * The user must initialise runState before use. The Algorithm manipulates
499 * \param runState to do its work
500 */
501typedef void (*LALInferenceAlgorithm) (struct tagLALInferenceRunState *runState);
502
503/** Type declaration for output logging function, can be user-declared */
505
506
507/**
508 * Structure for holding a LALInference proposal, along with name and stats.
509 */
510typedef struct
511tagLALInferenceProposal
512{
513 LALInferenceProposalFunction func; /* The actual proposal function */
514 char name[VARNAME_MAX]; /* The name of the proposal. This is used for printing stats */
515 INT4 weight; // Weight of proposal function in cycle
516 INT4 proposed; // Number of times proposal has been called
517 INT4 accepted; // Number of times a proposal from this function has been accepted
518 LALInferenceVariables *args; /** Local storage for arguments needed by the proposal (e.g. number of detectors) */
520
521/**
522 * Structure for holding a proposal cycle
523 */
524typedef struct
525tagLALInferenceProposalCycle
526{
527 LALInferenceProposal **proposals; /** Array of proposals (one per proposal function) */
528 INT4 *order; /* Array of proposal orders, with each element giving the index of the funcion in *proposals* */
529 INT4 length; /** Length of cycle */
530 INT4 nProposals; /* The number of unique proposals in the cycle */
531 INT4 counter; /** Counter for cycling through proposals */
532 char last_proposal_name[VARNAME_MAX]; /** Name of current proposal */
533 LALInferenceVariables *proposalArgs; /** Storage for arguments needed by proposal functions (e.g. number of detectors) */
535
536/**
537 * Structure containing chain-specific variables
538 */
539typedef struct
540tagLALInferenceThreadState
541{
542 INT4 id; /** Unique integer ID of this thread. Handy of I/O. */
544 INT4 step; /** Step counter for this thread. Negative numbers indicate burnin*/
545 INT4 effective_sample_size; /** Step counter for this thread. Negative numbers indicate burnin*/
546 LALInferenceProposalFunction proposal; /** The proposal cycle */
547 LALInferenceProposalCycle *cycle; /** Cycle of proposals to call */
548 LALInferenceModel *model; /** Stucture containing model buffers and parameters */
549 REAL8 currentPropDensity; /** Array containing multiple proposal densities */
551 LALInferenceVariables *proposalArgs, /** Arguments needed by proposals */
552 *algorithmParams, /** Stope things such as output arrays */
553 *priorArgs; /** Prior boundaries, etc. This is
554 stored at the thread level because proposals
555 often need to know about prior boundaries */
556 LALInferenceVariables *currentParams, /** The current parameters */
557 *preProposalParams, /** Current location going into jump proposal */
558 *proposedParams; /** Parameters proposed */
559 LALInferenceVariables **differentialPoints; /** Array of points for differential evolution */
560 size_t differentialPointsLength; /** Length of the current differential points stored in
561 differentialPoints. This should be removed can be given
562 as an algorithmParams entry */
563 size_t differentialPointsSize; /** Size of the differentialPoints memory block
564 (must be >= length of differential points).
565 Can also be removed. */
566 size_t differentialPointsSkip; /** When the DE buffer gets too long, start storing
567 only every n-th output point; this counter stores n */
568 REAL8 *currentIFOSNRs; /** Array storing single-IFO SNRs of current sample */
569 REAL8 *currentIFOLikelihoods; /** Array storing single-IFO likelihoods of current sample */
570 REAL8 currentSNR; /** Array storing network SNR of current sample */
572 REAL8 currentLikelihood; /** This should be removed, can be given as an algorithmParams or proposalParams entry */
573 REAL8 currentPrior; /** This should be removed, can be given as an algorithmParams entry */
576 gsl_rng *GSLrandom;
578 struct tagLALInferenceRunState *parent; /** Pointer to the parent RunState of the thread. e.g., Useful for getting data */
583
584
585/**
586 * Structure containing inference run state
587 * This includes pointers to the function types required to run
588 * the algorithm, and data structures as required
589 */
590typedef struct
591tagLALInferenceRunState
592{
593 ProcessParamsTable *commandLine; /** A ProcessParamsTable with command line arguments */
594 LALInferenceInitModelFunction initModel; /** A function that returns a new set of variables for the model */
595 LALInferenceAlgorithm algorithm; /** The algorithm function */
596 LALInferenceEvolveOneStepFunction evolve; /** The algorithm's single iteration function */
597 LALInferencePriorFunction prior; /** The prior for the parameters */
598 LALInferenceCubeToPriorFunction CubeToPrior; /** MultiNest prior for the parameters */
599 LALInferenceLikelihoodFunction likelihood; /** The likelihood function */
600 LALInferenceLogFunction logsample; /** Log sample, i.e. to disk */
601 struct tagLALInferenceIFOData *data; /** The data from the interferometers */
602 LALInferenceVariables *proposalArgs; /** Common arguments needed by proposals, to be copied to thread->cycle */
603 LALInferenceVariables *priorArgs, /** Any special arguments for the prior function */
604 *algorithmParams; /** Parameters which control the running of the algorithm*/
605 LALInferenceVariables **livePoints; /** Array of live points for Nested Sampling */
606 INT4 nthreads; /** Number of threads stored in ``threads``. */
608 gsl_rng *GSLrandom;
609 char *outFileName; /** Name for thread's output file */
610 char *resumeOutFileName; /** Name for thread's resume file */
611 char runID[VARNAME_MAX];
612 LALInferenceThreadState *threads; /** Array of chains for this run */
613
615
616
617#define DETNAMELEN 8
618
619/**
620 * Structure to contain IFO data.
621 * Some fields may be left empty if not needed
622 */
623typedef struct
624tagLALInferenceIFOData
625{
626 char name[DETNAMELEN]; /** Detector name */
627 REAL8TimeSeries *timeData, /** A time series from the detector */
628 *whiteTimeData, *windowedTimeData; /** white is not really white, but over-white. */
629 REAL8TimeSeries *varTimeData; /** A time series of the data noise variance */
630 /* Stores the log(L) for the model in presence of data. These were
631 added to allow for individual-detector log(L) output. The
632 convention is that loglikelihood always stores the log(L) for the
633 model in freqModel... or timeModel.... When a jump is accepted,
634 that value is copied into acceptedloglikelihood, which is the
635 quantity that is actually output in the output files. */
637 REAL8 fPlus, fCross; /** Detector responses */
638 REAL8 timeshift; /** What is this? */
639 COMPLEX16FrequencySeries *freqData, /** Buffer for frequency domain data */
640 *whiteFreqData; /* Over-white. */
641 COMPLEX16TimeSeries *compTimeData; /** Complex time series data buffers */
642 LALInferenceVariables *dataParams; /* Optional data parameters */
643 REAL8FrequencySeries *oneSidedNoisePowerSpectrum; /** one-sided Noise Power Spectrum */
644 REAL8FrequencySeries *noiseASD; /** (one-sided Noise Power Spectrum)^{-1/2} */
645// REAL8TimeSeries *timeDomainNoiseWeights; /** Roughly, InvFFT(1/Noise PSD). */
646 REAL8Window *window; /** A window */
647 REAL8 padding; /** Padding for the above window */
648 REAL8FFTPlan *timeToFreqFFTPlan, *freqToTimeFFTPlan; /** Pre-calculated FFT plans for forward and reverse FFTs */
649 REAL8FFTPlan *margFFTPlan; /** FFT plan needed for time/time-and-phase marginalisation */
650 REAL8 fLow, fHigh; /** integration limits for overlap integral in F-domain */
651 LALDetector *detector; /** LALDetector structure for where this data came from */
652 LIGOTimeGPS epoch; /** The epoch of this observation (the time of the first sample) */
653 REAL8 SNR; /** IF INJECTION ONLY, E(SNR) of the injection in the detector.*/
654 REAL8 STDOF; /** Degrees of freedom for IFO to be used in Student-T Likelihood. */
655 UINT4 likeli_counter; /** counts how many time the likelihood has been calculated */
656 UINT4 templa_counter; /** counts how many time the template has been calculated */
657 struct tagLALInferenceROQData *roq; /** ROQ data */
658
659 struct tagLALInferenceIFOData *next; /** A pointer to the next set of data for linked list */
661
662typedef struct
663tagLALInferenceROQData
664{
665 COMPLEX16 *weightsLinear; /** weights for <d|h>: NOTE: needs to be stored from data read from command line */
666 REAL8 *weightsQuadratic; /** weights for calculating <h|h>*/
672
673
674 struct tagLALInferenceROQSplineWeightsLinear *weights_linear;
675
676
677 /* Deprecated functions that should be removed at some point */
678 gsl_matrix_complex *weights; /** weights for the likelihood: NOTE: needs to be stored from data read from command line */
679 gsl_matrix_complex *mmweights; /** weights for calculating <h|h> if not using analytical formula */
680 double int_f_7_over_3; /** /int_{fmin}^{fmax} df f^(-7/3)/psd...for <h|h> part of the likelihood */
681 /* end deprecated function */
682
684
685/**
686 * * * Structure to contain spline of ROQ weights as a function of tc
687 * * */
688
689typedef struct
690tagLALInferenceROQSplineWeightsLinear
691{
692
693
696 gsl_interp_accel *acc_real_weight_linear;
697 gsl_interp_accel *acc_imag_weight_linear;
698
700/**
701 * * Structure to contain model-related Reduced Order Quadrature quantities
702 * */
703typedef struct
704tagLALInferenceROQModel
705{
710
712
714
715 REAL8Sequence * frequencyNodesLinear; /** empirical frequency nodes for the likelihood. NOTE: needs to be stored from data read from command line */
719
722
723 /* Deprecated functions that should be removed at some point */
724 gsl_vector_complex *hplus; /** waveform at frequency nodes. */
725 gsl_vector_complex *hcross;
726 gsl_vector_complex *hstrain;
727 gsl_vector *frequencyNodes; /** empirical frequency nodes for the likelihood. NOTE: needs to be stored from data read from command line */
729 /* end Deprecated functions */
730
732
733/**
734 * Structure to contain data-related Reduced Order Quadrature quantities
735 */
736/* Initialize an empty thread, saving a timestamp for benchmarking */
738
739/* Initialize a bunch of threads using LALInferenceInitThread */
741
742/** Returns the element of the process params table with "name" */
744
745/**
746 * parses a character string (passed as one of the options) and decomposes
747 * it into individual parameter character strings. \c input is of the form
748 * input : \"[one,two,three]\"
749 * and the resulting output \c strings is
750 * strings : {\"one\", \"two\", \"three\"}
751 * length of parameter names is for now limited to 512 characters.
752 * (should 'theoretically' (untested) be able to digest white space as well.
753 * Irrelevant for command line options, though.)
754 * \c n UNDOCUMENTED
755 */
756void LALInferenceParseCharacterOptionString(char *input, char **strings[], UINT4 *n);
757
758/** Return a ProcessParamsTable from the command line arguments */
760
761/** Return a ProcessParamsTrable from a string vector */
763
764/** Return a ProcessParamsTable from the command line arguments (SWIG version) */
766
767/** Output the command line based on the ProcessParamsTable \c procparams */
769
770/** Execute FFT for data in \c IFOdata */
772/** Execute Inverse FFT for data in \c IFOdata */
774
775/** Return the list node for "name" - do not rely on this */
777
778/**
779 * Return the list node for the idx-th item - do not rely on this
780 * Indexing starts at 1
781 */
783
784/**
785 * Pop the list node for "name". Returns a pointer to the node, which is removed from vars
786 */
788
789/** Output the sample to file *fp, in ASCII format */
791
792/** Output only non-fixed parameters */
794
795/** Output spline calibration parameters */
797
798/** Read in the non-fixed parameters from the given file (position in
799 the file must be arranged properly before calling this
800 function). */
802
803/** Utility for readling in delimited ASCII files. */
804REAL8 *LALInferenceParseDelimitedAscii(FILE *input, INT4 nCols, INT4 *wantedCols, INT4 *nLines);
805
806/* Parse a single line of delimited ASCII. */
807void parseLine(char *record, const char *delim, char arr[][VARNAME_MAX], INT4 *cnt);
808
809/* Discard the standard header of a PTMCMC chain file. */
810void LALInferenceDiscardPTMCMCHeader(FILE *filestream);
811
812/* Burn-in a PTMCMC output file. */
813void LALInferenceBurninPTMCMC(FILE *filestream, INT4 logl_idx, INT4 nPar);
814
815/* Burn-in a generic ASCII stream. */
816void LALInferenceBurninStream(FILE *filestream, INT4 burnin);
817
818/* Read column names from an ASCII file. */
819void LALInferenceReadAsciiHeader(FILE *input, char params[][VARNAME_MAX], INT4 *nCols);
820
821/* Utility for selecting columns from an array, in the specified order. */
822REAL8 **LALInferenceSelectColsFromArray(REAL8 **inarray, INT4 nRows, INT4 nCols, INT4 nSelCols, INT4 *selCols);
823
824/** Output proposal statistics header to file *fp */
826
827/** Output proposal statistics to file *fp */
829
830/**
831 * Reads one line from the given file and stores the values there into
832 * the variable structure, using the given header array to name the
833 * columns. Returns 0 on success.
834 */
835void LALInferenceProcessParamLine(FILE *inp, char **headers, LALInferenceVariables *vars);
836
837/** Sorts the variable structure by name */
839
840/** LALInferenceVariable buffer to array and vica versa */
843
844/** LALInference variables to an array, and vica versa */
846
848
849/**
850 * Append the sample to a file. file pointer is stored in state->algorithmParams as a
851 * LALInferenceVariable called "outfile", as a void ptr.
852 * Caller is responsible for opening and closing file.
853 * Variables are alphabetically sorted before being written
854 */
856
857/**
858 * Append the sample to an array which can be later processed by the user.
859 * Array is stored as a C array in a LALInferenceVariable in state->algorithmParams
860 * called "outputarray". Number of items in the array is stored as "N_outputarray".
861 * Will create the array and store it in this way if it does not exist.
862 * DOES NOT FREE ARRAY, user must clean up after use.
863 * Also outputs sample to disk if possible using LALInferenceLogSampleToFile()
864 */
866
867/** Convert from Mc, eta space to m1, m2 space (note m1 > m2).*/
868void LALInferenceMcEta2Masses(double mc, double eta, double *m1, double *m2);
869
870/** Convert from Mc, q space to m1, m2 space (q = m2/m1, with m1 > m2). */
871void LALInferenceMcQ2Masses(double mc, double q, double *m1, double *m2);
872
873/** Convert from q to eta (q = m2/m1, with m1 > m2). */
874void LALInferenceQ2Eta(double q, double *eta);
875/** Convert from dQuadMonS, dQuadMonA to dQuadMon1, dQuadMon2. */
876void LALInferencedQuadMonSdQuadMonA(REAL8 dQuadMonS, REAL8 dQuadMonA, REAL8 *dQuadMon1, REAL8 *dQuadMon2);
877/** Convert from lambdaT, dLambdaT, and eta to lambda1 and lambda2. */
879
880/** Calculate lambda1,2(m1,2|eos(logp1,gamma1,gamma2,gamma3)) */
881void LALInferenceLogp1GammasMasses2Lambdas(REAL8 logp1, REAL8 gamma1, REAL8 gamma2, REAL8 gamma3, REAL8 mass1, REAL8 mass2, REAL8 *lambda1, REAL8 *lambda2);
882
883/** Convert from spectral parameters to lambda1, lambda2 */
885
886/** Check for causality violation and mass conflict given masses and eos */
888
889/** Specral decomposition of eos's adiabatic index */
890double AdiabaticIndex(double gamma[],double x, int size);
891
892/** Determine if the Adiabatic index is within 'prior' range */
893int LALInferenceSDGammaCheck(double gamma[], int size);
894
895/**
896 * The kD trees in LALInference are composed of cells. Each cell
897 * represents a rectangular region in parameter space, defined by
898 * \f$\mathrm{lowerLeft} <= x <= \mathrm{upperRight}\f$. It also
899 * contains two sub-cells, each of which represents rectangular
900 * regions of half the size. The dimension along which a cell at a
901 * particular level in the tree is split in half is given by its
902 * level (mod the dimension of the space).
903 *
904 * Each cell contains some number (which may be zero) of points.
905 * Periodically, the cell will compute the mean and covariance of its
906 * points, and the eigenvectors of the covariance matrix (principal
907 * axes) of an ellipse fitting the the points. It will also compute
908 * the (tight) bounds of a box enclosing the points in a coordinate
909 * system aligned with the principal axes. When this has been done,
910 * the \c eigenFrameStale element will be set to zero. If points are
911 * subsequently added to the cell, the \c eigenFrameStale flag will
912 * be set to a non-zero value until the next re-computation of the
913 * principal axes.
914 */
915typedef struct tagLALInferenceKDTree {
916 size_t npts; /** Stores the number of tree points that lie in the cell. */
917 size_t ptsSize; /** Size of the pts buffer. */
918 size_t dim; /** Dimension of the system. */
920 REAL8 *ptsMean; /** Mean of pts. */
921 REAL8 *lowerLeft; /** Lower left (i.e. coordinate minimum) bound;
922 length is ndim from LALInferenceKDTree. */
923 REAL8 *upperRight; /** Upper right (i.e. coordinate maximum) bound. */
924 REAL8 **ptsCov; /** dim-by-dim covariance matrix. */
925 REAL8 **ptsCovEigenVects; /** Eigenvectors of the covariance matrix:
926 [i][j] is the jth component of the ith
927 eigenvector. */
928 REAL8 *eigenMin; /** Minimum coordinates of points in the eigen-frame. */
929 REAL8 *eigenMax; /** Maximum coordinates of points in the eigen-frame. */
930 int eigenFrameStale; /** == 1 when the mean, covariance, and
931 eigenvectors are out of date (i.e. more
932 points added). */
933 struct tagLALInferenceKDTree *left; /** Left (i.e. lower-coordinate)
934 sub-tree, may be NULL if
935 empty.*/
936 struct tagLALInferenceKDTree *right; /** Right
937 (i.e. upper-coordinate)
938 sub-tree, may be NULL if
939 empty. */
941
942/** Delete a kD-tree. Also deletes all contained cells, and points. */
944
945/**
946 * Constructs a fresh, empty kD tree. The top-level cell will get
947 * the given bounds, which should enclose every point added by
948 * LALInferenceKDAddPoint().
949 */
950LALInferenceKDTree *LALInferenceKDEmpty(REAL8 *lowerLeft, REAL8 *upperRight, size_t ndim);
951
952/**
953 * Adds a point to the kD-tree, returns 0 on successful exit. The
954 * memory for pt is owned by the tree, so should not be deallocated
955 * or modified except by LALInferenceKDTreeDelete().
956 */
958
959/**
960 * Returns the first cell that contains the given point that also
961 * contains fewer than Npts points, if possible. If no cell
962 * containing the given point has fewer than Npts points, then
963 * returns the cell containing the fewest number of points and the
964 * given point. Non-positive Npts will give the fewest-point cell in
965 * the tree containing the given point. Returns NULL on error.
966 */
968
969/**
970 * Returns the log of the volume of the given cell, which is part of
971 * the given tree.
972 */
974
975/**
976 * Returns the log of the volume of the box aligned with the
977 * principal axes of the points in the given cell that tightly
978 * encloses those points.
979 */
981
982/**
983 * Fills in the given REAL8 array with the parameter values from
984 * params; the ordering of the variables is taken from the order of
985 * the non-fixed variables in \c templt. It is an error if pt does
986 * not point to enough storage to store all the non-fixed parameters
987 * from \c templt and \c params.
988 */
990
991/**
992 * Fills in the non-fixed variables in params from the given REAL8
993 * array. The ordering of variables is given by the order of the
994 * non-fixed variables in \c templt.
995 */
997
998/**
999 * Draws a \c pt uniformly from a randomly chosen cell of \c
1000 * tree. The chosen cell will be chosen to have (as nearly as
1001 * possible) \c Npts in it.
1002 */
1003void LALInferenceKDDrawEigenFrame(gsl_rng * rng, LALInferenceKDTree *tree, REAL8 *pt, size_t Npts);
1004
1005/**
1006 * Returns the log of the jump proposal probability ratio for the
1007 * LALInferenceKDDrawEigenFrame() proposal to propose the point \c
1008 * proposed given the current position \c current , where \c Npts is
1009 * the parameter used to select the box to draw from in
1010 * LALInferenceKDDrawEigenFrame().
1011 */
1013 REAL8 *proposed, size_t Npts);
1014
1015/** Check matrix is positive definite. dim is matrix dimensions */
1017 gsl_matrix *matrix,
1018 UINT4 dim
1019 );
1020
1021/** Draw a random multivariate vector from Gaussian distr given covariance matrix */
1022void
1025 gsl_matrix *matrix,
1026 UINT4 dim,
1027 RandomParams *randParam
1028 );
1029/** Draw a random multivariate vector from student-t distr given covariance matrix */
1030void
1033 gsl_matrix *matrix,
1034 UINT4 dim,
1035 UINT4 n,
1036 RandomParams *randParam
1037 );
1038
1039
1040/** Calculate shortest angular distance between a1 and a2 (modulo 2PI) */
1042
1043/** Calculate the variance of a distribution on an angle (modulo 2PI) */
1044REAL8 LALInferenceAngularVariance(LALInferenceVariables **list,const char *pname, int N);
1045
1046/** Sanity check the structures in the given state. Will scan data for infinities and nans, and look for null pointers. */
1048
1049/**
1050 * Dump all waveforms from the ifodata structure. (currently: timeData, freqData)
1051 * basefilename is optional text to append to file names
1052 */
1053void LALInferenceDumpWaveforms(LALInferenceModel *model, const char *basefilename);
1054
1055/**
1056 * Write a LALInferenceVariables as binary to a given FILE pointer, returns the number
1057 * of items written (should be the dimension of the variables) or -1 on error
1058 */
1060
1061/**
1062 * Read from the given FILE * a LALInferenceVariables, which was previously
1063 * written with LALInferenceWriteVariablesBinary() Returns a new LALInferenceVariables
1064 */
1066
1067/**
1068 * Write an array N of LALInferenceVariables to the given FILE * using
1069 * LALInferenceWriteVariablesBinary(). Returns the number written (should be ==N)
1070 */
1072
1073/**
1074 * Read N LALInferenceVariables from the binary FILE *file, previously written with
1075 * LALInferenceWriteVariablesArrayBinary() returns the number read
1076 */
1078
1079/**
1080 * Write the LALInferenceVariables contents of runState to a file in binary format
1081 */
1083
1084/**
1085 * Reads the file and populates LALInferenceVariables in the runState that
1086 * were saved with LALInferenceReadVariablesArrayBinary()
1087 */
1089
1090
1091void LALInferenceAddINT4Variable(LALInferenceVariables * vars, const char * name, INT4 value, LALInferenceParamVaryType vary);
1092
1093INT4 LALInferenceGetINT4Variable(LALInferenceVariables * vars, const char * name);
1094
1095void LALInferenceSetINT4Variable(LALInferenceVariables* vars,const char* name,INT4 value);
1096
1097void LALInferenceAddINT8Variable(LALInferenceVariables * vars, const char * name, INT8 value, LALInferenceParamVaryType vary);
1098
1099INT8 LALInferenceGetINT8Variable(LALInferenceVariables * vars, const char * name);
1100
1101void LALInferenceSetINT8Variable(LALInferenceVariables* vars,const char* name,INT8 value);
1102
1103void LALInferenceAddUINT4Variable(LALInferenceVariables * vars, const char * name, UINT4 value, LALInferenceParamVaryType vary);
1104
1106
1107void LALInferenceSetUINT4Variable(LALInferenceVariables* vars,const char* name,UINT4 value);
1108
1109void LALInferenceAddREAL4Variable(LALInferenceVariables * vars, const char * name, REAL4 value, LALInferenceParamVaryType vary);
1110
1112
1113void LALInferenceSetREAL4Variable(LALInferenceVariables* vars,const char* name,REAL4 value);
1114
1115void LALInferenceAddREAL8Variable(LALInferenceVariables * vars, const char * name, REAL8 value, LALInferenceParamVaryType vary);
1116
1118
1119void LALInferenceSetREAL8Variable(LALInferenceVariables* vars,const char* name,REAL8 value);
1120
1122
1124
1125void LALInferenceSetCOMPLEX8Variable(LALInferenceVariables* vars,const char* name,COMPLEX8 value);
1126
1128
1130
1131void LALInferenceSetCOMPLEX16Variable(LALInferenceVariables* vars,const char* name,COMPLEX16 value);
1132
1133void LALInferenceAddgslMatrixVariable(LALInferenceVariables * vars, const char * name, gsl_matrix* value, LALInferenceParamVaryType vary);
1134
1135gsl_matrix* LALInferenceGetgslMatrixVariable(LALInferenceVariables * vars, const char * name);
1136
1137void LALInferenceSetgslMatrixVariable(LALInferenceVariables* vars,const char* name,gsl_matrix* value);
1138
1140
1142
1143void LALInferenceSetREAL8VectorVariable(LALInferenceVariables* vars,const char* name,REAL8Vector* value);
1144
1146
1148
1150
1152
1154
1156
1158
1159void LALInferenceSetINT4VectorVariable(LALInferenceVariables* vars, const char* name, INT4Vector* value);
1160
1161void LALInferenceSetUINT4VectorVariable(LALInferenceVariables* vars, const char* name, UINT4Vector* value);
1162
1164
1166
1168
1169#ifdef SWIG /* SWIG interface directives */
1170SWIGLAL(OWNS_THIS_ARG(const CHAR*, value));
1171#endif
1172
1173void LALInferenceAddstringVariable(LALInferenceVariables * vars, const char * name, const CHAR* value, LALInferenceParamVaryType vary);
1174
1175const CHAR* LALInferenceGetstringVariable(LALInferenceVariables * vars, const char * name);
1176
1177void LALInferenceSetstringVariable(LALInferenceVariables* vars,const char* name, const CHAR* value);
1178
1179#ifdef SWIG /* SWIG interface directives */
1180SWIGLAL_CLEAR(OWNS_THIS_ARG(const CHAR*, value));
1181#endif
1182
1183/**
1184 * Print spline calibration parameter names as tab-separated ASCII
1185 */
1187
1188/**
1189 * Compute Tidal deformabilities following BinaryLove Universal relations
1190 */
1192
1193
1194/**
1195 * Conversion routines between Equatorial (RA,DEC) and detector-based coordinate systems, where
1196 * new "north pole" points along vector from det0 to det1.
1197 * theta - azimuth angle about vector joining det0 and det1
1198 * alpha - co-latitude (0,pi) relative to det0-det1 vector
1199 */
1201 REAL8 t0, REAL8 alpha, REAL8 theta,
1202 REAL8 *tg, REAL8 *ra, REAL8 *dec);
1203
1205 REAL8 tg, REAL8 ra, REAL8 dec,
1206 REAL8 *t0, REAL8 *alpha, REAL8 *theta);
1207
1208/** @} */
1209
1210#endif
LALInferenceVariables currentParams
ProcessParamsTable * ptr
struct tagLALSimNeutronStarFamily LALSimNeutronStarFamily
const char *const name
sigmaKerr data[0]
double complex COMPLEX16
double REAL8
int64_t INT8
char CHAR
uint32_t UINT4
float complex COMPLEX8
int32_t INT4
float REAL4
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).
LALInferenceMCMCRunPhase * LALInferenceGetMCMCrunphase_ptrVariable(LALInferenceVariables *vars, const char *name)
REAL8 LALInferenceAngularDistance(REAL8 a1, REAL8 a2)
Calculate shortest angular distance between a1 and a2 (modulo 2PI)
INT4 LALInferenceBufferToArray(LALInferenceThreadState *thread, REAL8 **DEarray)
LALInferenceModel *(* LALInferenceInitModelFunction)(struct tagLALInferenceRunState *runState)
Type declaration for variables init function, can be user-declared.
Definition: LALInference.h:475
ProcessParamsTable * LALInferenceParseStringVector(LALStringVector *arglist)
Return a ProcessParamsTrable from a string vector.
int LALInferenceKDAddPoint(LALInferenceKDTree *tree, REAL8 *pt)
Adds a point to the kD-tree, returns 0 on successful exit.
INT4 LALInferenceThinnedBufferToArray(LALInferenceThreadState *thread, REAL8 **DEarray, INT4 step)
LALInferenceVariable buffer to array and vica versa.
void LALInferencePrintSampleNonFixed(FILE *fp, LALInferenceVariables *sample)
Output only non-fixed parameters.
Definition: LALInference.c:948
void LALInferenceReadAsciiHeader(FILE *input, char params[][VARNAME_MAX], INT4 *nCols)
Read column names from an ASCII file.
void LALInferenceAddINT8Variable(LALInferenceVariables *vars, const char *name, INT8 value, LALInferenceParamVaryType vary)
int LALInferenceReadVariablesArrayBinary(FILE *file, LALInferenceVariables **vars, UINT4 N)
Read N LALInferenceVariables from the binary FILE *file, previously written with LALInferenceWriteVar...
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 cou...
Definition: LALInference.c:260
void LALInferenceSetstringVariable(LALInferenceVariables *vars, const char *name, const CHAR *value)
void(* LALInferenceAlgorithm)(struct tagLALInferenceRunState *runState)
Type declaration for an algorithm function which is called by the driver code The user must initialis...
Definition: LALInference.h:501
void LALInferencedQuadMonSdQuadMonA(REAL8 dQuadMonS, REAL8 dQuadMonA, REAL8 *dQuadMon1, REAL8 *dQuadMon2)
Convert from dQuadMonS, dQuadMonA to dQuadMon1, dQuadMon2.
COMPLEX16Vector * LALInferenceGetCOMPLEX16VectorVariable(LALInferenceVariables *vars, const char *name)
LALInferenceVariableType LALInferenceGetVariableType(const LALInferenceVariables *vars, const char *name)
Get the LALInferenceVariableType of the parameter named name in vars.
Definition: LALInference.c:321
void LALInferenceDiscardPTMCMCHeader(FILE *filestream)
Discard the standard header of a PTMCMC chain file.
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: LALInference.h:399
void LALInferenceQ2Eta(double q, double *eta)
Convert from q to eta (q = m2/m1, with m1 > m2).
void(* LALInferenceLogFunction)(LALInferenceVariables *algorithmParams, LALInferenceVariables *vars)
Type declaration for output logging function, can be user-declared.
Definition: LALInference.h:504
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 corres...
Definition: LALInference.h:407
char * LALInferenceGetVariableName(LALInferenceVariables *vars, int idx)
Get the name of the idx-th variable Indexing starts at 1.
Definition: LALInference.c:339
#define VARNAME_MAX
Definition: LALInference.h:50
INT4 LALInferenceFprintParameterNonFixedHeaders(FILE *out, LALInferenceVariables *params)
Print the parameters which do not vary to a file as a tab-separated ASCII line.
UINT4 LALInferenceGetUINT4Variable(LALInferenceVariables *vars, const char *name)
void LALInferencePrintSample(FILE *fp, LALInferenceVariables *sample)
Output the sample to file *fp, in ASCII format.
Definition: LALInference.c:923
INT4 LALInferenceGetVariableDimensionNonFixed(LALInferenceVariables *vars)
Get number of dimensions in vars which are not fixed to a certain value.
Definition: LALInference.c:255
char ** LALInferenceGetHeaderLine(FILE *inp)
Returns an array of header strings (terminated by NULL) from a common-format output file.
int LALInferenceFprintParameterHeaders(FILE *out, LALInferenceVariables *params)
Print the parameter names to a file as a tab-separated ASCII line.
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,...
double AdiabaticIndex(double gamma[], double x, int size)
Specral decomposition of eos's adiabatic index.
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 sp...
void LALInferenceAddREAL8VectorVariable(LALInferenceVariables *vars, const char *name, REAL8Vector *value, LALInferenceParamVaryType vary)
void LALInferencePrintSplineCalibration(FILE *fp, LALInferenceThreadState *thread)
Output spline calibration parameters.
gsl_matrix * LALInferenceGetgslMatrixVariable(LALInferenceVariables *vars, const char *name)
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.
Definition: LALInference.c:395
REAL4 LALInferenceGetREAL4Variable(LALInferenceVariables *vars, const char *name)
double LALInferenceKDLogCellVolume(LALInferenceKDTree *cell)
Returns the log of the volume of the given cell, which is part of the given tree.
int LALInferenceWriteVariablesArrayBinary(FILE *file, LALInferenceVariables **vars, UINT4 N)
Write an array N of LALInferenceVariables to the given FILE * using LALInferenceWriteVariablesBinary(...
LALInferenceThreadState * LALInferenceInitThreads(INT4 nthreads)
Definition: LALInference.c:138
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 i...
void LALInferenceAddINT4VectorVariable(LALInferenceVariables *vars, const char *name, INT4Vector *value, LALInferenceParamVaryType vary)
int LALInferenceWriteRunStateBinary(FILE *file, LALInferenceRunState *state)
Write the LALInferenceVariables contents of runState to a file in binary format.
void LALInferenceLogSampleToFile(LALInferenceVariables *algorithmParams, LALInferenceVariables *vars)
Append the sample to a file.
UINT4 LALInferencePrintNVariableItem(char *out, UINT4 N, const LALInferenceVariableItem *const ptr)
Prints a variable item to a string.
Definition: LALInference.c:687
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))
REAL8 LALInferenceGetREAL8Variable(LALInferenceVariables *vars, const char *name)
INT4Vector * LALInferenceGetINT4VectorVariable(LALInferenceVariables *vars, const char *name)
void LALInferenceEquatorialToDetFrame(LALDetector *det0, LALDetector *det1, REAL8 tg, REAL8 ra, REAL8 dec, REAL8 *t0, REAL8 *alpha, REAL8 *theta)
const char * LALInferenceTranslateInternalToExternalParamName(const char *inName)
Converts between internally used parameter names and those external (e.g.
REAL8 LALInferenceKDLogProposalRatio(LALInferenceKDTree *tree, REAL8 *current, REAL8 *proposed, size_t Npts)
Returns the log of the jump proposal probability ratio for the LALInferenceKDDrawEigenFrame() proposa...
COMPLEX8 LALInferenceGetCOMPLEX8Variable(LALInferenceVariables *vars, const char *name)
void LALInferenceSetREAL8Variable(LALInferenceVariables *vars, const char *name, REAL8 value)
void LALInferenceExecuteFT(LALInferenceModel *model)
Execute FFT for data in IFOdata.
size_t LALInferenceTypeSize[15]
Definition: LALInference.c:86
double LALInferenceKDLogCellEigenVolume(LALInferenceKDTree *cell)
Returns the log of the volume of the box aligned with the principal axes of the points in the given c...
int LALInferenceCompareVariables(LALInferenceVariables *var1, LALInferenceVariables *var2)
Check for equality in two variables.
int LALInferenceEOSPhysicalCheck(LALInferenceVariables *params, ProcessParamsTable *commandLine)
Check for causality violation and mass conflict given masses and eos.
void LALInferenceRemoveVariable(LALInferenceVariables *vars, const char *name)
Remove name from vars Frees the memory for the name structure and its contents.
Definition: LALInference.c:450
void LALInferenceAddstringVariable(LALInferenceVariables *vars, const char *name, const CHAR *value, LALInferenceParamVaryType vary)
void(* LALInferenceSwapRoutine)(struct tagLALInferenceRunState *runState, FILE *)
Propose a swap between chain locations.
Definition: LALInference.h:494
char * LALInferencePrintCommandLine(ProcessParamsTable *procparams)
Output the command line based on the ProcessParamsTable procparams.
void LALInferenceLambdaTsEta2Lambdas(REAL8 lambdaT, REAL8 dLambdaT, REAL8 eta, REAL8 *lambda1, REAL8 *lambda2)
Convert from lambdaT, dLambdaT, and eta to lambda1 and lambda2.
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: LALInference.c:207
INT4 LALInferenceGetVariableTypeByIndex(LALInferenceVariables *vars, int idx)
Get the LALInferenceVariableType of the idx -th item in the vars Indexing starts at 1.
Definition: LALInference.c:326
int LALInferenceSDGammaCheck(double gamma[], int size)
Determine if the Adiabatic index is within 'prior' range.
LALInferenceParamVaryType
An enumerated type for denoting the topology of a parameter.
Definition: LALInference.h:127
void LALInferenceReadSampleNonFixed(FILE *fp, LALInferenceVariables *sample)
Read in the non-fixed parameters from the given file (position in the file must be arranged properly ...
Definition: LALInference.c:974
void LALInferenceSetCOMPLEX16Variable(LALInferenceVariables *vars, const char *name, COMPLEX16 value)
void LALInferenceSetUINT4VectorVariable(LALInferenceVariables *vars, const char *name, UINT4Vector *value)
void LALInferenceAddCOMPLEX16VectorVariable(LALInferenceVariables *vars, const char *name, COMPLEX16Vector *value, LALInferenceParamVaryType vary)
void LALInferenceDumpWaveforms(LALInferenceModel *model, const char *basefilename)
Dump all waveforms from the ifodata structure.
void LALInferenceSetCOMPLEX16VectorVariable(LALInferenceVariables *vars, const char *name, COMPLEX16Vector *value)
void LALInferenceCopyVariables(LALInferenceVariables *origin, LALInferenceVariables *target)
Deep copy the variables from one to another LALInferenceVariables structure.
Definition: LALInference.c:558
void LALInferenceBurninStream(FILE *filestream, INT4 burnin)
Burn-in a generic ASCII stream.
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 p...
Definition: LALInference.c:238
REAL8Vector * LALInferenceGetREAL8VectorVariable(LALInferenceVariables *vars, const char *name)
void LALInferenceAddREAL8Variable(LALInferenceVariables *vars, const char *name, REAL8 value, LALInferenceParamVaryType vary)
LALInferenceVariables * LALInferenceReadVariablesBinary(FILE *stream)
Read from the given FILE * a LALInferenceVariables, which was previously written with LALInferenceWri...
INT4 LALInferenceGetINT4Variable(LALInferenceVariables *vars, const char *name)
REAL8 ** LALInferenceSelectColsFromArray(REAL8 **inarray, INT4 nRows, INT4 nCols, INT4 nSelCols, INT4 *selCols)
Utility for selecting columns from an array, in the specified order.
void LALInferenceAddUINT4Variable(LALInferenceVariables *vars, const char *name, UINT4 value, LALInferenceParamVaryType vary)
REAL8 LALInferenceAngularVariance(LALInferenceVariables **list, const char *pname, int N)
Calculate the variance of a distribution on an angle (modulo 2PI)
ProcessParamsTable * LALInferenceGetProcParamVal(ProcessParamsTable *procparams, const char *name)
Returns the element of the process params table with "name".
int LALInferencePrintProposalStatsHeader(FILE *fp, LALInferenceProposalCycle *cycle)
Output proposal statistics header to file *fp.
void LALInferenceKDDrawEigenFrame(gsl_rng *rng, LALInferenceKDTree *tree, REAL8 *pt, size_t Npts)
Draws a pt uniformly from a randomly chosen cell of tree.
void LALInferencePrintProposalStats(FILE *fp, LALInferenceProposalCycle *cycle)
Output proposal statistics to file *fp.
void LALInferenceProcessParamLine(FILE *inp, char **headers, LALInferenceVariables *vars)
Reads one line from the given file and stores the values there into the variable structure,...
void LALInferenceAddCOMPLEX16Variable(LALInferenceVariables *vars, const char *name, COMPLEX16 value, LALInferenceParamVaryType vary)
ProcessParamsTable * LALInferenceParseCommandLine(int argc, char **argv)
Return a ProcessParamsTable from the command line arguments.
void LALInferenceAddINT4Variable(LALInferenceVariables *vars, const char *name, INT4 value, LALInferenceParamVaryType vary)
LALInferenceVariableItem * LALInferenceGetItem(const LALInferenceVariables *vars, const char *name)
Return the list node for "name" - do not rely on this.
Definition: LALInference.c:175
REAL8 * LALInferenceParseDelimitedAscii(FILE *input, INT4 nCols, INT4 *wantedCols, INT4 *nLines)
Utility for readling in delimited ASCII files.
void LALInferenceFprintSplineCalibrationHeader(FILE *output, LALInferenceThreadState *thread)
Print spline calibration parameter names as tab-separated ASCII.
LALInferenceMCMCRunPhase
Phase of MCMC run (depending on burn-in status, different actions are performed during the run,...
Definition: LALInference.h:180
void LALInferenceCopyArrayToVariables(REAL8 *origin, LALInferenceVariables *target)
void LALInferenceKDTreeDelete(LALInferenceKDTree *tree)
Delete a kD-tree.
void LALInferenceSetParamVaryType(LALInferenceVariables *vars, const char *name, LALInferenceParamVaryType vary)
Set the LALInferenceParamVaryType of the parameter named name in vars, see the declaration of LALInfe...
Definition: LALInference.c:231
void LALInferenceSetUINT4Variable(LALInferenceVariables *vars, const char *name, UINT4 value)
ProcessParamsTable * LALInferenceParseCommandLineStringVector(LALStringVector *args)
Return a ProcessParamsTable from the command line arguments (SWIG version)
void LALInferenceBurninPTMCMC(FILE *filestream, INT4 logl_idx, INT4 nPar)
Determine burnin cycle from delta-logl criteria.
INT4 LALInferenceGetVariableDimension(LALInferenceVariables *vars)
Get number of dimensions in variable vars.
Definition: LALInference.c:250
void LALInferenceAddMCMCrunphase_ptrVariable(LALInferenceVariables *vars, const char *name, LALInferenceMCMCRunPhase *value, LALInferenceParamVaryType vary)
void LALInferenceSetINT4Variable(LALInferenceVariables *vars, const char *name, INT4 value)
void(* LALInferenceTemplateFunction)(struct tagLALInferenceModel *model)
Type declaration for template function, which operates on a LALInferenceIFOData structure *data.
Definition: LALInference.h:377
int LALInferenceCheckVariableToPrint(LALInferenceVariables *vars, const char *name)
Definition: LALInference.c:514
LALInferenceKDTree * LALInferenceKDEmpty(REAL8 *lowerLeft, REAL8 *upperRight, size_t ndim)
Constructs a fresh, empty kD tree.
REAL8(* LALInferenceProposalFunction)(struct tagLALInferenceThreadState *thread, LALInferenceVariables *currentParams, LALInferenceVariables *proposedParams)
Jump proposal distribution Computes proposedParams based on currentParams and additional variables st...
Definition: LALInference.h:391
UINT4 LALInferenceCheckPositiveDefinite(gsl_matrix *matrix, UINT4 dim)
Check matrix is positive definite.
const CHAR * LALInferenceGetstringVariable(LALInferenceVariables *vars, const char *name)
void LALInferenceSDGammasMasses2Lambdas(REAL8 gamma[], REAL8 mass1, REAL8 mass2, REAL8 *lambda1, REAL8 *lambda2, int size)
Convert from spectral parameters to lambda1, lambda2.
UINT4Vector * LALInferenceGetUINT4VectorVariable(LALInferenceVariables *vars, const char *name)
void LALInferenceParseCharacterOptionString(char *input, char **strings[], UINT4 *n)
parses a character string (passed as one of the options) and decomposes it into individual parameter ...
void LALInferenceClearVariables(LALInferenceVariables *vars)
Delete the variables in this structure.
Definition: LALInference.c:532
#define DETNAMELEN
Definition: LALInference.h:617
void LALInferenceSortVariablesByName(LALInferenceVariables *vars)
Sorts the variable structure by name.
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.
void LALInferenceSetINT4VectorVariable(LALInferenceVariables *vars, const char *name, INT4Vector *value)
void LALInferenceAddgslMatrixVariable(LALInferenceVariables *vars, const char *name, gsl_matrix *value, LALInferenceParamVaryType vary)
void parseLine(char *record, const char *delim, char arr[][VARNAME_MAX], INT4 *cnt)
Parse a single line of delimited ASCII.
void LALInferenceBinaryLove(LALInferenceVariables *vars, REAL8 *lambda1, REAL8 *lambda2)
Compute Tidal deformabilities following BinaryLove Universal relations.
void LALInferenceCopyVariablesToArray(LALInferenceVariables *origin, REAL8 *target)
LALInference variables to an array, and vica versa.
void LALInferenceMcEta2Masses(double mc, double eta, double *m1, double *m2)
Convert from Mc, eta space to m1, m2 space (note m1 > m2).
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 suff...
int LALInferenceWriteVariablesBinary(FILE *file, LALInferenceVariables *vars)
Write a LALInferenceVariables as binary to a given FILE pointer, returns the number of items written ...
void XLALMultiNormalDeviates(REAL4Vector *vector, gsl_matrix *matrix, UINT4 dim, RandomParams *randParam)
Draw a random multivariate vector from Gaussian distr given covariance matrix.
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 calibratio...
INT4 LALInferenceSanityCheck(LALInferenceRunState *state)
Sanity check the structures in the given state.
INT4(* LALInferenceEvolveOneStepFunction)(struct tagLALInferenceRunState *runState)
Perform one step of an algorithm, replaces runState ->currentParams.
Definition: LALInference.h:491
void LALInferenceLogSampleToArray(LALInferenceVariables *algorithmParams, LALInferenceVariables *vars)
Append the sample to an array which can be later processed by the user.
void LALInferencePrintVariables(LALInferenceVariables *var)
Print variables to stdout.
Definition: LALInference.c:798
void LALInferenceExecuteInvFT(LALInferenceModel *model)
Execute Inverse FFT for data in IFOdata.
REAL8(* LALInferenceLikelihoodFunction)(LALInferenceVariables *currentParams, struct tagLALInferenceIFOData *data, LALInferenceModel *model)
Type declaration for likelihood function Computes p(data | currentParams, templt ) templt is a LALInf...
Definition: LALInference.h:487
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,...
void LALInferenceSetINT8Variable(LALInferenceVariables *vars, const char *name, INT8 value)
void LALInferenceSetREAL4Variable(LALInferenceVariables *vars, const char *name, REAL4 value)
int LALInferenceCheckVariable(LALInferenceVariables *vars, const char *name)
Checks for name being present in vars returns 1(==true) or 0.
Definition: LALInference.c:525
void LALInferenceSetgslMatrixVariable(LALInferenceVariables *vars, const char *name, gsl_matrix *value)
void LALInferenceCopyUnsetREAL8Variables(LALInferenceVariables *origin, LALInferenceVariables *target, ProcessParamsTable *commandLine)
Definition: LALInference.c:659
void LALInferenceSetMCMCrunphase_ptrVariable(LALInferenceVariables *vars, const char *name, LALInferenceMCMCRunPhase *value)
LALInferenceThreadState * LALInferenceInitThread(LALInferenceThreadState *thread)
Structure to contain data-related Reduced Order Quadrature quantities.
Definition: LALInference.c:105
void LALInferenceAddUINT4VectorVariable(LALInferenceVariables *vars, const char *name, UINT4Vector *value, LALInferenceParamVaryType vary)
void LALInferenceTranslateExternalToInternalParamName(char *outName, const char *inName)
Converts between externally used parameter names and those internal.
INT8 LALInferenceGetINT8Variable(LALInferenceVariables *vars, const char *name)
void LALInferenceSetCOMPLEX8Variable(LALInferenceVariables *vars, const char *name, COMPLEX8 value)
COMPLEX16 LALInferenceGetCOMPLEX16Variable(LALInferenceVariables *vars, const char *name)
LALInferenceVariableType
An enumerated type for denoting the type of a variable.
Definition: LALInference.h:104
LALInferenceParamVaryType LALInferenceGetVariableVaryType(LALInferenceVariables *vars, const char *name)
Get the LALInferenceParamVaryType of the parameter named name in vars see the declaration of LALInfer...
Definition: LALInference.c:226
void LALInferenceAddREAL4Variable(LALInferenceVariables *vars, const char *name, REAL4 value, LALInferenceParamVaryType vary)
void LALInferenceSetVariable(LALInferenceVariables *vars, const char *name, const void *value)
Set a variable named name in vars with a value.
Definition: LALInference.c:352
void LALInferenceSetREAL8VectorVariable(LALInferenceVariables *vars, const char *name, REAL8Vector *value)
int LALInferenceCheckVariableNonFixed(LALInferenceVariables *vars, const char *name)
Checks for name being present in vars and having type LINEAR or CIRCULAR.
Definition: LALInference.c:503
int LALInferenceReadRunStateBinary(FILE *file, LALInferenceRunState *state)
Reads the file and populates LALInferenceVariables in the runState that were saved with LALInferenceR...
void LALInferenceKDREAL8ToVariables(LALInferenceVariables *params, REAL8 *pt, LALInferenceVariables *templt)
Fills in the non-fixed variables in params from the given REAL8 array.
LALInferenceVariableItem * LALInferencePopVariableItem(LALInferenceVariables *vars, const char *name)
Pop the list node for "name".
void LALInferenceAddCOMPLEX8Variable(LALInferenceVariables *vars, const char *name, COMPLEX8 value, LALInferenceParamVaryType vary)
@ LALINFERENCE_PARAM_OUTPUT
A parameter that never changes, functions should respect this.
Definition: LALInference.h:131
@ LALINFERENCE_PARAM_FIXED
A parameter that is cyclic, such as an angle between 0 and 2pi.
Definition: LALInference.h:130
@ LALINFERENCE_PARAM_CIRCULAR
A parameter that simply has a maximum and a minimum.
Definition: LALInference.h:129
@ LALINFERENCE_PARAM_LINEAR
Definition: LALInference.h:128
@ LALINFERENCE_LADDER_UPDATE
In the single-chain sampling phase of an annealed run.
Definition: LALInference.h:185
@ LALINFERENCE_ANNEALING
In the parallel tempering phase of an annealed run.
Definition: LALInference.h:183
@ LALINFERENCE_SINGLE_CHAIN
In the annealing phase of an annealed run.
Definition: LALInference.h:184
@ LALINFERENCE_TEMP_PT
Run only parallel tempers.
Definition: LALInference.h:182
@ LALINFERENCE_ONLY_PT
Definition: LALInference.h:181
@ LALINFERENCE_string_t
Definition: LALInference.h:117
@ LALINFERENCE_COMPLEX16_t
Definition: LALInference.h:111
@ LALINFERENCE_MCMCrunphase_ptr_t
Definition: LALInference.h:118
@ LALINFERENCE_REAL8_t
Definition: LALInference.h:109
@ LALINFERENCE_INT4Vector_t
Definition: LALInference.h:114
@ LALINFERENCE_REAL4_t
Definition: LALInference.h:108
@ LALINFERENCE_INT4_t
Definition: LALInference.h:105
@ LALINFERENCE_REAL8Vector_t
Definition: LALInference.h:113
@ LALINFERENCE_void_ptr_t
Definition: LALInference.h:119
@ LALINFERENCE_UINT4_t
Definition: LALInference.h:107
@ LALINFERENCE_UINT4Vector_t
Definition: LALInference.h:115
@ LALINFERENCE_gslMatrix_t
Definition: LALInference.h:112
@ LALINFERENCE_INT8_t
Definition: LALInference.h:106
@ LALINFERENCE_COMPLEX8_t
Definition: LALInference.h:110
@ LALINFERENCE_COMPLEX16Vector_t
Definition: LALInference.h:116
LALSimulationDomain
static const INT4 q
args
type
Structure to contain IFO data.
Definition: LALInference.h:625
REAL8 padding
A window.
Definition: LALInference.h:647
REAL8TimeSeries * timeData
Detector name.
Definition: LALInference.h:627
REAL8TimeSeries * varTimeData
white is not really white, but over-white.
Definition: LALInference.h:629
LALDetector * detector
integration limits for overlap integral in F-domain
Definition: LALInference.h:651
REAL8FFTPlan * margFFTPlan
Pre-calculated FFT plans for forward and reverse FFTs.
Definition: LALInference.h:649
REAL8FrequencySeries * oneSidedNoisePowerSpectrum
Definition: LALInference.h:643
REAL8 timeshift
Detector responses.
Definition: LALInference.h:638
UINT4 templa_counter
counts how many time the likelihood has been calculated
Definition: LALInference.h:656
struct tagLALInferenceROQData * roq
counts how many time the template has been calculated
Definition: LALInference.h:657
REAL8FrequencySeries * noiseASD
one-sided Noise Power Spectrum
Definition: LALInference.h:644
REAL8Window * window
(one-sided Noise Power Spectrum)^{-1/2}
Definition: LALInference.h:646
COMPLEX16FrequencySeries * freqData
What is this?
Definition: LALInference.h:639
COMPLEX16FrequencySeries * whiteFreqData
Buffer for frequency domain data.
Definition: LALInference.h:640
UINT4 likeli_counter
Degrees of freedom for IFO to be used in Student-T Likelihood.
Definition: LALInference.h:655
REAL8 SNR
The epoch of this observation (the time of the first sample)
Definition: LALInference.h:653
REAL8 nullloglikelihood
A time series of the data noise variance.
Definition: LALInference.h:636
struct tagLALInferenceIFOData * next
ROQ data.
Definition: LALInference.h:659
LIGOTimeGPS epoch
LALDetector structure for where this data came from.
Definition: LALInference.h:652
REAL8FFTPlan * freqToTimeFFTPlan
Definition: LALInference.h:648
LALInferenceVariables * dataParams
Complex time series data buffers.
Definition: LALInference.h:642
REAL8 STDOF
IF INJECTION ONLY, E(SNR) of the injection in the detector.
Definition: LALInference.h:654
REAL8TimeSeries * whiteTimeData
A time series from the detector.
Definition: LALInference.h:628
COMPLEX16TimeSeries * compTimeData
Definition: LALInference.h:641
Detector-dependent buffers and parameters A linked list meant for parameters and buffers that are sep...
Definition: LALInference.h:417
COMPLEX16FrequencySeries * freqhCross
Definition: LALInference.h:425
REAL8TimeSeries * timeData
Complex time series signal buffer.
Definition: LALInference.h:427
struct tagLALInferenceIFOModel * next
Time series model buffer.
Definition: LALInference.h:429
REAL8TimeSeries * timehCross
Definition: LALInference.h:424
LALDetector * detector
Parameters used when filling the buffers - template functions should copy to here.
Definition: LALInference.h:420
void * extraData
LALDetector structure for where this data came from.
Definition: LALInference.h:422
LALInferenceVariables * params
Definition: LALInference.h:418
COMPLEX16TimeSeries * compTimeSignal
Freq series model buffers.
Definition: LALInference.h:426
The kD trees in LALInference are composed of cells.
Definition: LALInference.h:915
struct tagLALInferenceKDTree * right
Left (i.e.
Definition: LALInference.h:936
int eigenFrameStale
Maximum coordinates of points in the eigen-frame.
Definition: LALInference.h:930
REAL8 * eigenMax
Minimum coordinates of points in the eigen-frame.
Definition: LALInference.h:929
REAL8 * upperRight
Lower left (i.e.
Definition: LALInference.h:923
REAL8 ** ptsCovEigenVects
dim-by-dim covariance matrix.
Definition: LALInference.h:925
size_t ptsSize
Stores the number of tree points that lie in the cell.
Definition: LALInference.h:917
size_t dim
Size of the pts buffer.
Definition: LALInference.h:918
REAL8 * lowerLeft
Mean of pts.
Definition: LALInference.h:921
REAL8 ** pts
Dimension of the system.
Definition: LALInference.h:919
struct tagLALInferenceKDTree * left
== 1 when the mean, covariance, and eigenvectors are out of date (i.e.
Definition: LALInference.h:933
REAL8 * eigenMin
Eigenvectors of the covariance matrix: [i][j] is the jth component of the ith eigenvector.
Definition: LALInference.h:928
REAL8 ** ptsCov
Upper right (i.e.
Definition: LALInference.h:924
Structure to constain a model and its parameters.
Definition: LALInference.h:436
REAL8 * ifo_loglikelihoods
Network SNR at params
Definition: LALInference.h:445
REAL8 * ifo_SNRs
Array of single-IFO likelihoods at params
Definition: LALInference.h:446
LALSimNeutronStarFamily * eos_fam
Is ROQ enabled.
Definition: LALInference.h:465
REAL8 fHigh
Start frequency for waveform generation.
Definition: LALInference.h:449
REAL8TimeSeries * timehCross
Definition: LALInference.h:453
REAL8 logprior
The template generation function.
Definition: LALInference.h:442
LALSimBurstWaveformCache * burstWaveformCache
Waveform cache.
Definition: LALInference.h:459
LALSimInspiralWaveformCache * waveformCache
Definition: LALInference.h:458
REAL8 loglikelihood
Prior value at params
Definition: LALInference.h:443
REAL8 padding
A window.
Definition: LALInference.h:462
REAL8 SNR
Likelihood value at params
Definition: LALInference.h:444
int roq_flag
ROQ data.
Definition: LALInference.h:464
LALInferenceVariables * params
Definition: LALInference.h:437
COMPLEX16FrequencySeries ** freqhs
Freq series model buffers.
Definition: LALInference.h:455
COMPLEX16FrequencySeries * freqhCross
Definition: LALInference.h:454
LALDict * LALpars
Projected freq series model buffers.
Definition: LALInference.h:457
REAL8FFTPlan * freqToTimeFFTPlan
Definition: LALInference.h:460
INT4 freqLength
Sampling rate information.
Definition: LALInference.h:451
REAL8Window * window
Pre-calculated FFT plans for forward and reverse FFTs.
Definition: LALInference.h:461
REAL8 fLow
Array of single-IFO SNRs at params
Definition: LALInference.h:448
LALInferenceIFOModel * ifo
Parameters used when filling the buffers - template functions should copy to here.
Definition: LALInference.h:438
struct tagLALInferenceROQModel * roq
The padding of the above window.
Definition: LALInference.h:463
LALInferenceTemplateFunction templt
Domain of model.
Definition: LALInference.h:440
LALSimulationDomain domain
IFO-dependent parameters and buffers.
Definition: LALInference.h:439
Structure for holding a proposal cycle.
Definition: LALInference.h:526
LALInferenceProposal ** proposals
Definition: LALInference.h:527
INT4 * order
Array of proposals (one per proposal function)
Definition: LALInference.h:528
INT4 nProposals
Length of cycle.
Definition: LALInference.h:530
LALInferenceVariables * proposalArgs
Name of current proposal.
Definition: LALInference.h:533
Structure for holding a LALInference proposal, along with name and stats.
Definition: LALInference.h:512
LALInferenceVariables * args
Definition: LALInference.h:518
LALInferenceProposalFunction func
Definition: LALInference.h:513
REAL8 time_weights_width
weights for calculating <h|h>
Definition: LALInference.h:667
double int_f_7_over_3
weights for calculating <h|h> if not using analytical formula
Definition: LALInference.h:680
struct tagLALInferenceROQSplineWeightsLinear * weights_linear
Definition: LALInference.h:674
REAL8 * weightsQuadratic
weights for <d|h>: NOTE: needs to be stored from data read from command line
Definition: LALInference.h:666
gsl_matrix_complex * mmweights
weights for the likelihood: NOTE: needs to be stored from data read from command line
Definition: LALInference.h:679
gsl_matrix_complex * weights
Definition: LALInference.h:678
COMPLEX16 * weightsLinear
Definition: LALInference.h:665
COMPLEX16Sequence * calFactorLinear
Definition: LALInference.h:711
COMPLEX16FrequencySeries * hctildeQuadratic
Definition: LALInference.h:709
COMPLEX16FrequencySeries * hctildeLinear
Definition: LALInference.h:707
REAL8Sequence * frequencyNodesQuadratic
empirical frequency nodes for the likelihood.
Definition: LALInference.h:716
gsl_vector * frequencyNodes
Definition: LALInference.h:727
REAL8Sequence * frequencyNodesLinear
Definition: LALInference.h:715
gsl_vector_complex * hcross
waveform at frequency nodes.
Definition: LALInference.h:725
REAL8 * amp_squared
empirical frequency nodes for the likelihood.
Definition: LALInference.h:728
COMPLEX16Sequence * calFactorQuadratic
Definition: LALInference.h:713
gsl_vector_complex * hstrain
Definition: LALInference.h:726
gsl_vector_complex * hplus
Definition: LALInference.h:724
COMPLEX16FrequencySeries * hptildeQuadratic
Definition: LALInference.h:708
COMPLEX16FrequencySeries * hptildeLinear
Definition: LALInference.h:706
gsl_interp_accel * acc_real_weight_linear
Definition: LALInference.h:696
gsl_spline * spline_real_weight_linear
Definition: LALInference.h:694
gsl_interp_accel * acc_imag_weight_linear
Definition: LALInference.h:697
gsl_spline * spline_imag_weight_linear
Definition: LALInference.h:695
Structure containing inference run state This includes pointers to the function types required to run...
Definition: LALInference.h:592
ProcessParamsTable * commandLine
Definition: LALInference.h:593
LALInferenceVariables * proposalArgs
The data from the interferometers.
Definition: LALInference.h:602
LALInferenceVariables * algorithmParams
Any special arguments for the prior function.
Definition: LALInference.h:604
LALInferenceCubeToPriorFunction CubeToPrior
The prior for the parameters.
Definition: LALInference.h:598
INT4 nthreads
Array of live points for Nested Sampling.
Definition: LALInference.h:606
struct tagLALInferenceIFOData * data
Log sample, i.e.
Definition: LALInference.h:601
LALInferenceVariables * priorArgs
Common arguments needed by proposals, to be copied to thread->cycle.
Definition: LALInference.h:603
LALInferenceThreadState * threads
Definition: LALInference.h:612
LALInferenceAlgorithm algorithm
A function that returns a new set of variables for the model.
Definition: LALInference.h:595
LALInferenceEvolveOneStepFunction evolve
The algorithm function.
Definition: LALInference.h:596
LALInferencePriorFunction prior
The algorithm's single iteration function.
Definition: LALInference.h:597
LALInferenceLogFunction logsample
The likelihood function.
Definition: LALInference.h:600
char * resumeOutFileName
Name for thread's output file.
Definition: LALInference.h:610
LALInferenceInitModelFunction initModel
A ProcessParamsTable with command line arguments.
Definition: LALInference.h:594
LALInferenceVariables ** livePoints
Parameters which control the running of the algorithm.
Definition: LALInference.h:605
LALInferenceSwapRoutine parallelSwap
Number of threads stored in threads.
Definition: LALInference.h:607
LALInferenceLikelihoodFunction likelihood
MultiNest prior for the parameters.
Definition: LALInference.h:599
Structure containing chain-specific variables.
Definition: LALInference.h:541
LALInferenceVariables * proposedParams
Current location going into jump proposal.
Definition: LALInference.h:558
REAL8 currentSNR
Array storing single-IFO likelihoods of current sample.
Definition: LALInference.h:570
INT4 effective_sample_size
Step counter for this thread.
Definition: LALInference.h:545
LALInferenceVariables * algorithmParams
Arguments needed by proposals.
Definition: LALInference.h:552
REAL8 currentPrior
This should be removed, can be given as an algorithmParams or proposalParams entry.
Definition: LALInference.h:573
LALInferenceVariables * currentParams
Prior boundaries, etc.
Definition: LALInference.h:556
struct tagLALInferenceRunState * parent
Definition: LALInference.h:578
LALInferenceVariables * priorArgs
Stope things such as output arrays.
Definition: LALInference.h:553
REAL8 currentPropDensity
Stucture containing model buffers and parameters.
Definition: LALInference.h:549
LALInferenceModel * model
Cycle of proposals to call.
Definition: LALInference.h:548
REAL8 * currentIFOSNRs
When the DE buffer gets too long, start storing only every n-th output point; this counter stores n.
Definition: LALInference.h:568
INT4 accepted
This should be removed, can be given as an algorithmParams entry.
Definition: LALInference.h:574
size_t differentialPointsLength
Array of points for differential evolution.
Definition: LALInference.h:560
LALInferenceProposalFunction proposal
Step counter for this thread.
Definition: LALInference.h:546
LALInferenceProposalCycle * cycle
The proposal cycle.
Definition: LALInference.h:547
REAL8 * currentIFOLikelihoods
Array storing single-IFO SNRs of current sample.
Definition: LALInference.h:569
LALInferenceVariables * proposalArgs
Definition: LALInference.h:551
REAL8 temperature
Array containing multiple proposal densities.
Definition: LALInference.h:550
INT4 * temp_swap_accepts
Pointer to the parent RunState of the thread.
Definition: LALInference.h:579
size_t differentialPointsSize
Length of the current differential points stored in differentialPoints.
Definition: LALInference.h:563
size_t differentialPointsSkip
Size of the differentialPoints memory block (must be >= length of differential points).
Definition: LALInference.h:566
LALInferenceVariables ** differentialPoints
Parameters proposed.
Definition: LALInference.h:559
LALInferenceVariables * preProposalParams
The current parameters.
Definition: LALInference.h:557
REAL8 nullLikelihood
Array storing network SNR of current sample.
Definition: LALInference.h:571
The LALInferenceVariableItem list node structure This should only be accessed using the accessor func...
Definition: LALInference.h:154
LALInferenceVariableType type
Definition: LALInference.h:157
struct tagVariableItem * next
Definition: LALInference.h:159
LALInferenceParamVaryType vary
Definition: LALInference.h:158
The LALInferenceVariables structure to contain a set of parameters Implemented as a linked list of LA...
Definition: LALInference.h:170
LALHashTbl * hash_table
Definition: LALInference.h:173
LALInferenceVariableItem * head
Definition: LALInference.h:171
Stores previously-computed waveforms and parameters to take advantage of approximant- and parameter-s...