LALPulsar  6.1.0.1-89842e6
SynthesizeCWDraws.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2010 Reinhard Prix
3  * Copyright (C) 2011, 2014 David Keitel
4  *
5  * This program is free software; you can redistribute it and/or modify
6  * it under the terms of the GNU General Public License as published by
7  * the Free Software Foundation; either version 2 of the License, or
8  * (at your option) any later version.
9  *
10  * This program is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13  * GNU General Public License for more details.
14  *
15  * You should have received a copy of the GNU General Public License
16  * along with with program; see the file COPYING. If not, write to the
17  * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
18  * MA 02110-1301 USA
19  */
20 
21 /**
22  * \defgroup SynthesizeCWDraws_h Header SynthesizeCWDraws.h
23  * \ingroup lalpulsar_general
24  * \author Reinhard Prix, David Keitel
25  *
26  * \brief
27  * Generate samples of various statistics (F-stat, F-atoms, B-stat,...) drawn from their
28  * respective distributions, assuming Gaussian noise, and drawing signal params from
29  * their (given) priors
30  *
31  * This is based on synthesizeBstat, and is mostly meant to be used for efficient
32  * Monte-Carlos studies, ROC curves etc
33  *
34  */
35 /** @{ */
36 
37 #ifndef _SYNTHESIZE_CW_DRAWS_H /* Double-include protection. */
38 #define _SYNTHESIZE_CW_DRAWS_H
39 
40 /* C++ protection. */
41 #ifdef __cplusplus
42 extern "C" {
43 #endif
44 
45 /*---------- INCLUDES ----------*/
46 
47 /* GSL includes */
48 #include <gsl/gsl_rng.h>
49 
50 /* LAL includes */
51 #include <lal/SkyCoordinates.h>
52 #include <lal/LALComputeAM.h>
53 
54 #include <lal/ProbabilityDensity.h>
55 #include <lal/TransientCW_utils.h>
56 
57 /*---------- exported types ----------*/
58 
59 /**
60  * Enumeration of allowed amplitude-prior types
61  */
62 typedef enum tagAmpPriorType_t {
63  AMP_PRIOR_TYPE_PHYSICAL = 0, /**< 'physical' priors: isotropic pdf{cosi,psi,phi0} AND flat pdf(h0) */
64  AMP_PRIOR_TYPE_CANONICAL, /**< 'canonical' priors: uniform in A^mu up to h_max */
67 
68 /**
69  * Signal (amplitude) parameter ranges
70  */
71 typedef struct tagAmplitudePrior_t {
72  pdf1D_t *pdf_h0Nat; /**< pdf for h0/sqrt{Sn} */
73  REAL8 fixedSNR; /**< alternative 1: adjust h0 to fix the optimal SNR of the signal */
74  BOOLEAN fixRhohMax; /**< alternative 2: draw h0 with fixed rhohMax = h0Max * (detM)^(1/8) <==> canonical Fstat prior */
75 
76  pdf1D_t *pdf_cosi; /**< pdf(cosi) */
77  pdf1D_t *pdf_psi; /**< pdf(psi) */
78  pdf1D_t *pdf_phi0; /**< pdf(phi0) */
80 
81 /**
82  * struct for buffering of AM-coeffs, if signal for same sky-position is injected
83  */
84 typedef struct tagmultiAMBuffer_t {
85  SkyPosition skypos; /**< sky-position for which we have AM-coeffs computed already */
86  MultiAMCoeffs *multiAM;; /**< pre-computed AM-coeffs for skypos */
88 
89 /**
90  * Hold all (generally) randomly drawn injection parameters: skypos, amplitude-params, M_mu_nu, transient-window, SNR
91  */
92 typedef struct tagInjParams_t {
99  REAL8 detM1o8; // (detMp)^(1/8): rescale param between h0, and rhoh = h0 * (detMp)^(1/8)
100 } InjParams_t;
101 
102 
103 /*---------- Global variables ----------*/
104 
105 /*---------- exported prototypes [API] ----------*/
106 int XLALDrawCorrelatedNoise( PulsarAmplitudeVect n_mu, const gsl_matrix *L, gsl_rng *rng );
107 
108 // ----- API to synthesize F-stat atoms for transient CW searches
109 FstatAtomVector *XLALGenerateFstatAtomVector( const DetectorStateSeries *detStates, const AMCoeffs *amcoeffs );
111 
112 int XLALAddNoiseToFstatAtomVector( FstatAtomVector *atoms, gsl_rng *rng );
113 int XLALAddNoiseToMultiFstatAtomVector( MultiFstatAtomVector *multiAtoms, gsl_rng *rng );
114 
117 
119 int write_InjParams_to_fp( FILE *fp, const InjParams_t *par, const UINT4 dataStartGPS, const BOOLEAN outputMmunuX, const UINT4 numDetectors );
120 
123  SkyPosition skypos,
124  AmplitudePrior_t AmpPrior,
125  transientWindowRange_t transientInjectRange,
126  const MultiDetectorStateSeries *multiDetStates,
127  BOOLEAN SignalOnly,
128  multiAMBuffer_t *multiAMBuffer,
129  gsl_rng *rng,
130  INT4 lineX,
131  const MultiNoiseWeights *multiNoiseWeights
132  );
133 
134 /** @} */
135 
136 #ifdef __cplusplus
137 }
138 #endif
139 /* C++ protection. */
140 
141 #endif /* Double-include protection. */
#define L(i, j)
unsigned char BOOLEAN
double REAL8
uint32_t UINT4
int32_t INT4
REAL8 PulsarAmplitudeVect[4]
Struct for 'canonical' coordinates in amplitude-params space A^mu = {A1, A2, A3, A4}.
FstatAtomVector * XLALGenerateFstatAtomVector(const DetectorStateSeries *detStates, const AMCoeffs *amcoeffs)
Generate an FstatAtomVector for given antenna-pattern functions.
int XLALAddNoiseToMultiFstatAtomVector(MultiFstatAtomVector *multiAtoms, gsl_rng *rng)
Add Gaussian-noise components to given multi-FstatAtomVector.
REAL8 XLALAddSignalToMultiFstatAtomVector(MultiFstatAtomVector *multiAtoms, AntennaPatternMatrix *M_mu_nu, const PulsarAmplitudeVect A_Mu, transientWindow_t transientWindow, INT4 lineX)
Add given signal s_mu = M_mu_nu A^nu within the given transient-window to multi-IFO noise-atoms.
int write_InjParams_to_fp(FILE *fp, const InjParams_t *par, const UINT4 dataStartGPS, const BOOLEAN outputMmunuX, const UINT4 numDetectors)
Write an injection-parameters structure to the given file-pointer, adding one line with the injection...
MultiFstatAtomVector * XLALGenerateMultiFstatAtomVector(const MultiDetectorStateSeries *multiDetStates, const MultiAMCoeffs *multiAM)
Generate a multi-FstatAtomVector for given antenna-pattern functions.
MultiFstatAtomVector * XLALSynthesizeTransientAtoms(InjParams_t *injParamsOut, SkyPosition skypos, AmplitudePrior_t AmpPrior, transientWindowRange_t transientInjectRange, const MultiDetectorStateSeries *multiDetStates, BOOLEAN SignalOnly, multiAMBuffer_t *multiAMBuffer, gsl_rng *rng, INT4 lineX, const MultiNoiseWeights *multiNoiseWeights)
Generates a multi-Fstat atoms vector for given parameters, drawing random parameters wherever require...
int XLALRescaleMultiFstatAtomVector(MultiFstatAtomVector *multiAtoms, REAL8 rescale)
Rescale a given multi-Fstat atoms vector {Fa,Fb} by given scalar factor.
int XLALAddNoiseToFstatAtomVector(FstatAtomVector *atoms, gsl_rng *rng)
Add Gaussian-noise components to given FstatAtomVector.
AmpPriorType_t
Enumeration of allowed amplitude-prior types.
int XLALDrawCorrelatedNoise(PulsarAmplitudeVect n_mu, const gsl_matrix *L, gsl_rng *rng)
Generate 4 random-noise draws n_mu = {n_1, n_2, n_3, n_4} with correlations according to the matrix M...
REAL8 XLALAddSignalToFstatAtomVector(FstatAtomVector *atoms, AntennaPatternMatrix *M_mu_nu, const PulsarAmplitudeVect A_Mu, transientWindow_t transientWindow)
Add signal s_mu = M_mu_nu A^nu within the given transient-window to given atoms.
@ AMP_PRIOR_TYPE_LAST
@ AMP_PRIOR_TYPE_CANONICAL
'canonical' priors: uniform in A^mu up to h_max
@ AMP_PRIOR_TYPE_PHYSICAL
'physical' priors: isotropic pdf{cosi,psi,phi0} AND flat pdf(h0)
This structure contains the per-SFT (weighted) antenna-pattern functions , with the SFT-index,...
Definition: LALComputeAM.h:63
Signal (amplitude) parameter ranges.
pdf1D_t * pdf_h0Nat
pdf for h0/sqrt{Sn}
REAL8 fixedSNR
alternative 1: adjust h0 to fix the optimal SNR of the signal
BOOLEAN fixRhohMax
alternative 2: draw h0 with fixed rhohMax = h0Max * (detM)^(1/8) <==> canonical Fstat prior
pdf1D_t * pdf_phi0
pdf(phi0)
pdf1D_t * pdf_cosi
pdf(cosi)
pdf1D_t * pdf_psi
pdf(psi)
Struct holding the "antenna-pattern" matrix , in terms of the multi-detector scalar product.
Definition: LALComputeAM.h:127
Timeseries of DetectorState's, representing the detector-info at different timestamps.
A vector of -statistic 'atoms', i.e.
Definition: ComputeFstat.h:174
Hold all (generally) randomly drawn injection parameters: skypos, amplitude-params,...
SkyPosition skypos
transientWindow_t transientWindow
PulsarAmplitudeParams ampParams
MultiAMCoeffs multiAM
PulsarAmplitudeVect ampVect
Multi-IFO container for antenna-pattern coefficients and atenna-pattern matrix .
Definition: LALComputeAM.h:137
Multi-IFO time-series of DetectorStates.
A multi-detector vector of FstatAtomVector.
Definition: ComputeFstat.h:186
One noise-weight (number) per SFT (therefore indexed over IFOs and SFTs.
Definition: PSDutils.h:71
Type containing the JKS 'amplitude parameters' {h0, cosi, phi0, psi}.
struct for buffering of AM-coeffs, if signal for same sky-position is injected
MultiAMCoeffs * multiAM
SkyPosition skypos
sky-position for which we have AM-coeffs computed already
Struct defining one transient window instance.
Struct defining a range of transient windows.