LALPulsar  6.1.0.1-89842e6
CWMakeFakeData.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2013 Reinhard Prix
3  *
4  * This program is free software; you can redistribute it and/or modify
5  * it under the terms of the GNU General Public License as published by
6  * the Free Software Foundation; either version 2 of the License, or
7  * (at your option) any later version.
8  *
9  * This program is distributed in the hope that it will be useful,
10  * but WITHOUT ANY WARRANTY; without even the implied warranty of
11  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12  * GNU General Public License for more details.
13  *
14  * You should have received a copy of the GNU General Public License
15  * along with with program; see the file COPYING. If not, write to the
16  * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
17  * MA 02110-1301 USA
18  */
19 
20 #ifndef _CWMAKEFAKEDATA_H // Double-include protection
21 #define _CWMAKEFAKEDATA_H
22 
23 /* C++ protection. */
24 #ifdef __cplusplus
25 extern "C" {
26 #endif
27 
28 /**
29  * \defgroup CWMakeFakeData_h Header CWMakeFakeData.h
30  * \ingroup lalpulsar_inject
31  * \author Reinhard Prix, Karl Wette
32  *
33  * \brief Module for generating 'fake' data containing CW signals and/or Gaussian noise.
34  * This basically presents a high-level wrapper API to the lower-level CW signal-generation
35  * functions in lalsuite.
36  */
37 /** @{ */
38 
39 // ---------- exported INCLUDES ----------
40 #include <math.h>
41 
42 // gsl includes
43 
44 // lal includes
45 #include <lal/LALDatatypes.h>
46 #include <lal/SFTfileIO.h>
47 #include <lal/PulsarDataTypes.h>
48 #include <lal/ConfigFile.h>
49 #include <lal/DetectorStates.h>
50 #include <lal/FITSFileIO.h>
51 
52 // ---------- exported TYPES ----------
53 /**
54  * Straightforward vector type of N PulsarParams structs
55  */
56 typedef struct tagPulsarParamsVector {
57 #ifdef SWIG /* SWIG interface directives */
58  SWIGLAL( ARRAY_1D( PulsarParamsVector, PulsarParams, data, UINT4, length ) );
59 #endif /* SWIG */
60  UINT4 length; //!< number of pulsar-signals
61  PulsarParams *data; //!< array of pulsar-signal parameters
63 
64 /**
65  * Struct controlling all the aspects of the fake data (time-series + SFTs)
66  * to be produced by XLALCWMakeFakeData() and XLALCWMakeFakeMultiData()
67  */
68 #ifdef SWIG /* SWIG interface directives */
69 SWIGLAL( IMMUTABLE_MEMBERS( tagCWMFDataParams, SFTWindowType ) );
70 #endif /* SWIG */
71 typedef struct tagCWMFDataParams {
72  REAL8 fMin; //!< smallest frequency guaranteed to be generated [returned fMin can be smaller]
73  REAL8 Band; //!< smallest frequency band guaranteed to be generated [returned Band can be larger]
74  MultiLALDetector multiIFO; //!< detectors to generate data for
75  MultiNoiseFloor multiNoiseFloor; //!< ... and corresponding noise-floors to generate Gaussian white noise for
76  MultiLIGOTimeGPSVector *multiTimestamps; //!< timestamps to generate SFTs for
77  const char *SFTWindowType; //!< window to apply to the SFT timeseries
78  REAL8 SFTWindowParam; //!< parameter required for *some* windows [otherwise must be 0]
79  UINT4 randSeed; //!< seed value for random-number generator
80  MultiREAL8TimeSeries *inputMultiTS; //!< [optional] input time-series for signals+noise to be added to
81  REAL8 sourceDeltaT; //!< [optional] source-frame sampling period. '0' means to use the previous internal defaults
83 
84 // ---------- Global variables ----------
85 
86 extern const char *const InjectionSourcesHelpString;
87 
88 // ---------- exported prototypes [API] ----------
89 
90 #ifdef SWIG // SWIG interface directives
91 SWIGLAL( INOUT_STRUCTS( MultiSFTVector **, multiSFTs ) );
92 SWIGLAL( INOUT_STRUCTS( MultiREAL8TimeSeries **, multiTseries ) );
93 SWIGLAL( INOUT_STRUCTS( SFTVector **, SFTVect ) );
94 SWIGLAL( INOUT_STRUCTS( REAL8TimeSeries **, Tseries ) );
95 #endif
96 
100 int XLALCWMakeFakeData( SFTVector **SFTVect, REAL8TimeSeries **Tseries,
101  const PulsarParamsVector *injectionSources, const CWMFDataParams *dataParams, UINT4 detectorIndex, const EphemerisData *edat );
102 
104 XLALGenerateCWSignalTS( const PulsarParams *pulsarParams, const LALDetector *site, LIGOTimeGPS startTime, REAL8 duration, REAL8 fSamp, REAL8 fHet, const EphemerisData *edat, REAL8 sourceDeltaT );
105 SFTVector *
106 XLALMakeSFTsFromREAL8TimeSeries( const REAL8TimeSeries *timeseries, const LIGOTimeGPSVector *timestamps, const char *windowType, REAL8 windowParam );
107 
108 int XLALReadPulsarParams( PulsarParams *pulsarParams, LALParsedDataFile *cfgdata, const CHAR *secName, const LIGOTimeGPS *refTimeDef );
109 PulsarParamsVector *XLALPulsarParamsFromFile( const char *fname, const LIGOTimeGPS *refTimeDef );
110 PulsarParamsVector *XLALPulsarParamsFromUserInput( const LALStringVector *UserInput, const LIGOTimeGPS *refTimeDef );
111 
114 
116 
117 int XLALFITSWritePulsarParamsVector( FITSFile *file, const CHAR *tableName, const PulsarParamsVector *list );
118 
120 
121 /** @} */
122 
123 #ifdef __cplusplus
124 }
125 #endif
126 /* C++ protection. */
127 
128 #endif /* Double-include protection. */
LIGOTimeGPSVector * timestamps
int XLALFindSmallestValidSamplingRate(UINT4 *n1, UINT4 n0, const LIGOTimeGPSVector *timestamps)
Find the smallest sampling rate of the form fsamp = n / Tsft, with n>=n0, such that all gap sizes Dt_...
int XLALCWMakeFakeData(SFTVector **SFTvect, REAL8TimeSeries **Tseries, const PulsarParamsVector *injectionSources, const CWMFDataParams *dataParams, UINT4 detectorIndex, const EphemerisData *edat)
Single-IFO version of XLALCWMakeFakeMultiData(), handling the actual work, but same input API.
PulsarParamsVector * XLALPulsarParamsVectorAppend(PulsarParamsVector *list, const PulsarParamsVector *add)
Append the given PulsarParamsVector 'add' to the vector 'list' ( which can be NULL),...
PulsarParamsVector * XLALCreatePulsarParamsVector(UINT4 numPulsars)
Create zero-initialized PulsarParamsVector for numPulsars.
int XLALReadPulsarParams(PulsarParams *pulsarParams, LALParsedDataFile *cfgdata, const CHAR *secName, const LIGOTimeGPS *refTimeDef)
Function to parse a config-file-type string (or section thereof) into a PulsarParams struct.
SFTVector * XLALMakeSFTsFromREAL8TimeSeries(const REAL8TimeSeries *timeseries, const LIGOTimeGPSVector *timestamps, const char *windowType, REAL8 windowParam)
Make SFTs from given REAL8TimeSeries at given timestamps, potentially applying a time-domain window o...
PulsarParamsVector * XLALPulsarParamsFromUserInput(const LALStringVector *UserInput, const LIGOTimeGPS *refTimeDef)
Function to determine the PulsarParamsVector input from a user-input defining CW sources.
const char *const InjectionSourcesHelpString
void XLALDestroyPulsarParamsVector(PulsarParamsVector *ppvect)
Destructor for PulsarParamsVector type.
int XLALFITSWritePulsarParamsVector(FITSFile *file, const CHAR *tableName, const PulsarParamsVector *list)
Write a PulsarParamsVector to a FITS file.
void XLALDestroyCWMFDataParams(CWMFDataParams *params)
Destructor for a CWMFDataParams type.
REAL4TimeSeries * XLALGenerateCWSignalTS(const PulsarParams *pulsarParams, const LALDetector *site, LIGOTimeGPS startTime, REAL8 duration, REAL8 fSamp, REAL8 fHet, const EphemerisData *edat, REAL8 sourceDeltaT)
Generate a (heterodyned) REAL4 timeseries of a CW signal for given pulsarParams, site,...
PulsarParamsVector * XLALPulsarParamsFromFile(const char *fname, const LIGOTimeGPS *refTimeDef)
Parse a given 'CWsources' config file for PulsarParams, return vector of all pulsar definitions found...
int XLALCWMakeFakeMultiData(MultiSFTVector **multiSFTs, MultiREAL8TimeSeries **multiTseries, const PulsarParamsVector *injectionSources, const CWMFDataParams *dataParams, const EphemerisData *edat)
Generate fake 'CW' data, returned either as SFTVector or REAL4Timeseries or both, for given CW-signal...
struct tagFITSFile FITSFile
Representation of a FITS file.
Definition: FITSFileIO.h:54
double REAL8
char CHAR
uint32_t UINT4
A vector of COMPLEX8FrequencySeries.
Struct controlling all the aspects of the fake data (time-series + SFTs) to be produced by XLALCWMake...
REAL8 sourceDeltaT
[optional] source-frame sampling period. '0' means to use the previous internal defaults
REAL8 fMin
smallest frequency guaranteed to be generated [returned fMin can be smaller]
MultiLALDetector multiIFO
detectors to generate data for
MultiLIGOTimeGPSVector * multiTimestamps
timestamps to generate SFTs for
REAL8 Band
smallest frequency band guaranteed to be generated [returned Band can be larger]
MultiREAL8TimeSeries * inputMultiTS
[optional] input time-series for signals+noise to be added to
MultiNoiseFloor multiNoiseFloor
... and corresponding noise-floors to generate Gaussian white noise for
UINT4 randSeed
seed value for random-number generator
const char * SFTWindowType
window to apply to the SFT timeseries
REAL8 SFTWindowParam
parameter required for some windows [otherwise must be 0]
This structure contains all information about the center-of-mass positions of the Earth and Sun,...
A vector of 'timestamps' of type LIGOTimeGPS.
Definition: SFTfileIO.h:188
array of detectors definitions 'LALDetector'
A collection of (multi-IFO) LIGOTimeGPSVector time-stamps vectors.
Definition: SFTfileIO.h:198
array of detector-specific 'noise floors' (ie PSD values), assumed constant over the frequency-band o...
A collection of (multi-IFO) REAL8 time-series.
A collection of SFT vectors – one for each IFO in a multi-IFO search.
Definition: SFTfileIO.h:179
Type defining the parameters of a pulsar-source of CW Gravitational waves.
Straightforward vector type of N PulsarParams structs.
PulsarParams * data
array of pulsar-signal parameters
UINT4 length
number of pulsar-signals