Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALPulsar 7.1.1.1-5e288d3
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
PulsarDataTypes.h
Go to the documentation of this file.
1/*
2 * Copyright (C) 2004, 2005 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 _PULSARDATATYPES_H /* Double-include protection. */
21#define _PULSARDATATYPES_H
22
23/* C++ protection. */
24#ifdef __cplusplus
25extern "C" {
26#endif
27
28/**
29 * \author Reinhard Prix
30 * \date 2005
31 * \defgroup PulsarDataTypes_h Header PulsarDataTypes.h
32 * \ingroup lalpulsar_general
33 * \brief Some common useful data-types for pulsar-searches.
34 *
35 * ### Synopsis ###
36 *
37 * \code
38 * #include <lal/PulsarDataTypes.h>
39 * \endcode
40 *
41 */
42/** @{ */
43
44#include <gsl/gsl_matrix.h>
45
46#include <lal/LALDatatypes.h>
47#include <lal/SkyCoordinates.h>
48
49// ---------- generic time/frequencies series vector types ----------
50
51/** A collection of (multi-IFO) REAL4 time-series */
52typedef struct tagMultiREAL4TimeSeries {
53#ifdef SWIG /* SWIG interface directives */
54 SWIGLAL( ARRAY_1D( MultiREAL4TimeSeries, REAL4TimeSeries *, data, UINT4, length ) );
55#endif /* SWIG */
56 UINT4 length; /**< Number of elements in array. */
57 REAL4TimeSeries **data; /**< Pointer to the data array. */
59
60/** A collection of (multi-IFO) REAL8 time-series */
61typedef struct tagMultiREAL8TimeSeries {
62#ifdef SWIG /* SWIG interface directives */
63 SWIGLAL( ARRAY_1D( MultiREAL8TimeSeries, REAL8TimeSeries *, data, UINT4, length ) );
64#endif /* SWIG */
65 UINT4 length; /**< Number of elements in array. */
66 REAL8TimeSeries **data; /**< Pointer to the data array. */
68
69/** A vector of COMPLEX8FrequencySeries */
70typedef struct tagCOMPLEX8FrequencySeriesVector {
71#ifdef SWIG /* SWIG interface directives */
72 SWIGLAL( ARRAY_1D( COMPLEX8FrequencySeriesVector, COMPLEX8FrequencySeries, data, UINT4, length ) );
73#endif /* SWIG */
74 UINT4 length; /**< Number of elements in array. */
75 COMPLEX8FrequencySeries *data; /**< Pointer to the data array. */
77
78/** A vector of COMPLEX16FrequencySeries */
79typedef struct tagCOMPLEX16FrequencySeriesVector {
80#ifdef SWIG /* SWIG interface directives */
81 SWIGLAL( ARRAY_1D( COMPLEX16FrequencySeriesVector, COMPLEX16FrequencySeries, data, UINT4, length ) );
82#endif /* SWIG */
83 UINT4 length; /**< Number of elements in array. */
84 COMPLEX16FrequencySeries *data; /**< Pointer to the data array. */
86
87/** A vector of REAL4FrequencySeries */
88typedef struct tagREAL4FrequencySeriesVector {
89#ifdef SWIG /* SWIG interface directives */
90 SWIGLAL( ARRAY_1D( REAL4FrequencySeriesVector, REAL4FrequencySeries, data, UINT4, length ) );
91#endif /* SWIG */
92 UINT4 length; /**< Number of elements in array. */
93 REAL4FrequencySeries *data; /**< Pointer to the data array. */
95
96/** A vector of REAL8FrequencySeries */
97typedef struct tagREAL8FrequencySeriesVector {
98#ifdef SWIG /* SWIG interface directives */
99 SWIGLAL( ARRAY_1D( REAL8FrequencySeriesVector, REAL8FrequencySeries, data, UINT4, length ) );
100#endif /* SWIG */
101 UINT4 length; /**< Number of elements in array. */
102 REAL8FrequencySeries *data; /**< Pointer to the data array. */
104
105// ---------- types describing the CW amplitude and phase parameters ----------
106
107/** maximal number of spin-parameters (Freq + spindowns) we can handle */
108#define PULSAR_MAX_SPINS 7
109
110/** maximal number of detectors we can handle (for static arrays of detector quantities) */
111#ifndef PULSAR_MAX_DETECTORS // allow this value to be overridden for e.g. E@H apps
112#define PULSAR_MAX_DETECTORS 10
113#endif
114
115/** Type containing the JKS 'amplitude parameters' {h0, cosi, phi0, psi} */
116typedef struct tagPulsarAmplitudeParams {
117 REAL8 psi; /**< polarization angle psi */
118 REAL8 phi0; /**< initial signal-phase (at some reference time) */
119 REAL8 aPlus; /**< Signal amplitude (plus) */
120 REAL8 aCross; /**< Signal amplitude (cross) */
122
123/** Struct for 'canonical' coordinates in amplitude-params space A^mu = {A1, A2, A3, A4} */
125
126/** Typedef for fixed-size array holding GW frequency and derivatives fk = d^k Freq/dt^k|(tau_ref) */
128
129/**
130 * Contains a "spin-range", ie spins \f$ f^{(k)} \f$ and corresponding bands \f$ \Delta f^{(k)} \f$
131 * at a given (SSB) reference GPS-time \f$ \tau \f$ .
132 * "Canonical" ordering refers to \f$ \Delta f^{(k)} >= 0 \f$ for all k.
133 */
134typedef struct tagPulsarSpinRange {
135 LIGOTimeGPS refTime; /**< SSB reference GPS-time at which spin-range is defined */
136 PulsarSpins fkdot; /**< Vector of spin-values \f$ f^{(k)} \f$ */
137 PulsarSpins fkdotBand; /**< Vector of spin-bands \f$ \Delta f^{(k)} \f$ , MUST be same length as fkdot */
139
140/** Type containing the 'Doppler-parameters' affecting the time-evolution of the phase */
141typedef struct tagPulsarDopplerParams {
142 LIGOTimeGPS refTime; /**< Reference time of pulsar parameters (in SSB!) */
143 REAL8 Alpha; /**< Sky position: RA (longitude) in equatorial coords and radians */
144 REAL8 Delta; /**< Sky position: DEC (latitude) in equatorial coords and radians */
145 PulsarSpins fkdot; /**< Intrinsic spins: [Freq, f1dot, f2dot, ... ] where fkdot = d^kFreq/dt^k */
146 REAL8 asini; /**< Binary: projected, normalized orbital semi-major axis (s).
147 If asini = 0 this is an \e isolated pulsar, and other binary parameters are ignored. */
148 REAL8 period; /**< Binary: orbital period (sec) */
149 REAL8 ecc; /**< Binary: orbital eccentricity */
150 LIGOTimeGPS tp; /**< Binary: time of observed periapsis passage (in SSB) */
151 REAL8 argp; /**< Binary: argument of periapsis (radians) */
153
154// ---------- transient-CW related types ----------
155
156/** Struct to define parameters of a 'transient window' to be applied to obtain transient signals */
157typedef enum tagtransientWindowType_t {
158 TRANSIENT_NONE = 0, /**< Note: in this case the window-parameters will be ignored, and treated as rect={data},
159 * i.e. a simple rectangular window covering all the data => this should always reproduce the
160 * standard F-statistic computation. */
161 TRANSIENT_RECTANGULAR = 1, /**< standard rectangular window covering [t0, t0+tau] */
162 TRANSIENT_EXPONENTIAL, /**< exponentially decaying window e^{-t0/tau} starting at t0.
163 * Note: we'll truncate this at some small (eg 3x) e-folding TRANSIENT_EXP_EFOLDING */
166
167/** Struct defining one transient window instance */
168typedef struct tagtransientWindow_t {
169 transientWindowType_t type; /**< window-type: none, rectangular, exponential, .... */
170 UINT4 t0; /**< GPS start-time 't0' */
171 UINT4 tau; /**< transient timescale tau in seconds */
173
174// ---------- 'integrated' types describing a complete CW signal ----------
175
176/** Type defining the parameters of a pulsar-source of CW Gravitational waves */
177typedef struct tagPulsarParams {
178 CHAR name[LALNameLength]; /**< 'name' for this sources, can be an empty string */
179 PulsarAmplitudeParams Amp; /**< 'Amplitude-parameters': h0, cosi, phi0, psi */
180 PulsarDopplerParams Doppler; /**< 'Phase-evolution parameters': {skypos, fkdot, orbital params } */
181 transientWindow_t Transient; /**< Transient window-parameters (start-time, duration, window-type) */
183
184/** Type containing a "candidate": parameter-space point with estimated errors and Fstat-value/significance */
185typedef struct tagPulsarCandidate {
186 PulsarAmplitudeParams Amp, dAmp; /**< amplitude-parameters and error-estimates */
187 PulsarDopplerParams Doppler, dDoppler;/**< Doppler-parameters and error-bars */
188 REAL8 significance; /**< a (user-chosen) measure of 'significance': Fstat, Hough-count,... */
189 gsl_matrix *AmpFisherMatrix; /**< Fisher-matrix of amplitude-subspace: has more info than dAmp! */
191
192/** @} */
193
194#ifdef __cplusplus
195}
196#endif
197/* C++ protection. */
198
199#endif /* Double-include protection. */
const char * name
Definition: SearchTiming.c:93
double REAL8
char CHAR
uint32_t UINT4
LALNameLength
REAL8 PulsarSpins[PULSAR_MAX_SPINS]
Typedef for fixed-size array holding GW frequency and derivatives fk = d^k Freq/dt^k|(tau_ref)
transientWindowType_t
Struct to define parameters of a 'transient window' to be applied to obtain transient signals.
#define PULSAR_MAX_SPINS
maximal number of spin-parameters (Freq + spindowns) we can handle
REAL8 PulsarAmplitudeVect[4]
Struct for 'canonical' coordinates in amplitude-params space A^mu = {A1, A2, A3, A4}.
@ TRANSIENT_RECTANGULAR
standard rectangular window covering [t0, t0+tau]
@ TRANSIENT_LAST
@ TRANSIENT_NONE
Note: in this case the window-parameters will be ignored, and treated as rect={data},...
@ TRANSIENT_EXPONENTIAL
exponentially decaying window e^{-t0/tau} starting at t0.
float data[BLOCKSIZE]
Definition: hwinject.c:360
A vector of COMPLEX16FrequencySeries.
COMPLEX16FrequencySeries * data
Pointer to the data array.
UINT4 length
Number of elements in array.
A vector of COMPLEX8FrequencySeries.
COMPLEX8FrequencySeries * data
Pointer to the data array.
UINT4 length
Number of elements in array.
A collection of (multi-IFO) REAL4 time-series.
UINT4 length
Number of elements in array.
REAL4TimeSeries ** data
Pointer to the data array.
A collection of (multi-IFO) REAL8 time-series.
REAL8TimeSeries ** data
Pointer to the data array.
UINT4 length
Number of elements in array.
Type containing the JKS 'amplitude parameters' {h0, cosi, phi0, psi}.
REAL8 aCross
Signal amplitude (cross)
REAL8 psi
polarization angle psi
REAL8 aPlus
Signal amplitude (plus)
REAL8 phi0
initial signal-phase (at some reference time)
Type containing a "candidate": parameter-space point with estimated errors and Fstat-value/significan...
REAL8 significance
a (user-chosen) measure of 'significance': Fstat, Hough-count,...
PulsarAmplitudeParams Amp
gsl_matrix * AmpFisherMatrix
Fisher-matrix of amplitude-subspace: has more info than dAmp!
PulsarDopplerParams dDoppler
Doppler-parameters and error-bars.
Type containing the 'Doppler-parameters' affecting the time-evolution of the phase.
REAL8 period
Binary: orbital period (sec)
LIGOTimeGPS tp
Binary: time of observed periapsis passage (in SSB)
PulsarSpins fkdot
Intrinsic spins: [Freq, f1dot, f2dot, ... ] where fkdot = d^kFreq/dt^k.
REAL8 Delta
Sky position: DEC (latitude) in equatorial coords and radians.
LIGOTimeGPS refTime
Reference time of pulsar parameters (in SSB!)
REAL8 Alpha
Sky position: RA (longitude) in equatorial coords and radians.
REAL8 ecc
Binary: orbital eccentricity.
REAL8 asini
Binary: projected, normalized orbital semi-major axis (s).
REAL8 argp
Binary: argument of periapsis (radians)
Type defining the parameters of a pulsar-source of CW Gravitational waves.
PulsarAmplitudeParams Amp
'Amplitude-parameters': h0, cosi, phi0, psi
transientWindow_t Transient
Transient window-parameters (start-time, duration, window-type)
PulsarDopplerParams Doppler
'Phase-evolution parameters': {skypos, fkdot, orbital params }
Contains a "spin-range", ie spins and corresponding bands at a given (SSB) reference GPS-time .
PulsarSpins fkdot
Vector of spin-values .
LIGOTimeGPS refTime
SSB reference GPS-time at which spin-range is defined.
PulsarSpins fkdotBand
Vector of spin-bands , MUST be same length as fkdot.
A vector of REAL4FrequencySeries.
REAL4FrequencySeries * data
Pointer to the data array.
UINT4 length
Number of elements in array.
A vector of REAL8FrequencySeries.
REAL8FrequencySeries * data
Pointer to the data array.
UINT4 length
Number of elements in array.
Struct defining one transient window instance.
UINT4 t0
GPS start-time 't0'.
UINT4 tau
transient timescale tau in seconds
transientWindowType_t type
window-type: none, rectangular, exponential, ....