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
DopplerFullScan.h
Go to the documentation of this file.
1/*
2 * Copyright (C) 2007, 2008, 2012 Karl Wette
3 * Copyright (C) 2006 Reinhard Prix
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 * \author Reinhard Prix, Karl Wette
23 * \date 2006
24 * \file
25 * \brief Header file defining the API for DopplerFullScan.
26 *
27 */
28
29#ifndef _DOPPLERFULLSCAN_H /* Double-include protection. */
30#define _DOPPLERFULLSCAN_H
31
32/* C++ protection. */
33#ifdef __cplusplus
34extern "C" {
35#endif
36
37/*---------- INCLUDES ----------*/
38#include <lal/LALDatatypes.h>
39#include <lal/SkyCoordinates.h>
40#include <lal/PtoleMetric.h>
41#include <lal/LALBarycenter.h>
42#include <lal/PulsarDataTypes.h>
43#include <lal/ComputeFstat.h>
44#include <lal/DopplerScan.h>
45#include <lal/LatticeTiling.h>
46
47/*---------- DEFINES ----------*/
48
49/*---------- external types ----------*/
50
51/* ==================== FULL MULTIDIMENSIONAL-GRID types ==================== */
52/**
53 * Structure describing a region in paramter-space (a,d,f,f1dot,..).
54 * Currently this is simply a direct product of skyRegion x FreqBand x f1dotBand.
55 */
56
57/** initialization struct for full InitDopplerScan() [UNDER CONSTRUCTION] */
58#ifdef SWIG /* SWIG interface directives */
59SWIGLAL( IMMUTABLE_MEMBERS( tagDopplerFullScanInit, Detector, ephemeris, gridFile ) );
60#endif /* SWIG */
61typedef struct tagDopplerFullScanInit {
62 DopplerRegion searchRegion; /**< Doppler-space region to be covered + scanned */
63 DopplerGridType gridType; /**< which type of grid to generate */
64 LALPulsarMetricType metricType; /**< which metric to use if GRID_METRIC */
65 BOOLEAN projectMetric; /**< project metric on f=const subspace */
66 PulsarDopplerParams stepSizes; /**< user-settings for stepsizes if GRID_FLAT */
67 REAL8 metricMismatch; /**< for GRID_METRIC and GRID_ISOTROPIC */
68 LIGOTimeGPS startTime; /**< start-time of the observation */
69 REAL8 Tspan; /**< total time spanned by observation */
70 const LALDetector *Detector; /**< Current detector */
71 const EphemerisData *ephemeris; /**< ephemeris-data for numerical metrics */
72 const CHAR *gridFile; /**< filename for sky-grid or full-grid if GRID_FILE_SKYGRID or GRID_FILE_FULLGRID */
73 REAL8 extraArgs[3]; /**< extra grid-specific setup arguments */
75
76/** opaque type to reflects the current state of a full multi-dimensional DopplerScan */
77typedef struct tagDopplerFullScanState DopplerFullScanState; /* opaque type */
78
79
80/*---------- Global variables ----------*/
81/* some empty structs for initializations */
82
83/*---------- external prototypes [API] ----------*/
84
85/* ----- full-fledged multi-dimensional Doppler scanner functions ----- */
86DopplerFullScanState *XLALInitDopplerFullScan( const DopplerFullScanInit *init );
87
88int XLALNextDopplerPos( PulsarDopplerParams *pos, DopplerFullScanState *scan );
89REAL8 XLALNumDopplerTemplates( DopplerFullScanState *scan );
90UINT8 XLALNumDopplerPointsAtDimension( DopplerFullScanState *scan, const size_t dim );
91int XLALGetDopplerSpinRange( PulsarSpinRange *spinRange, const DopplerFullScanState *scan );
92void XLALDestroyDopplerFullScan( DopplerFullScanState *scan );
93
94/* ----- variout utility functions ----- */
95
96///
97/// Set a first spindown bound derived from spindown age and braking indices
98///
100 LatticeTiling *tiling, ///< [in] Lattice tiling
101 const size_t freq_dimension, ///< [in] Frequency dimension
102 const size_t f1dot_dimension, ///< [in] First spindown dimension
103 const double age, ///< [in] Spindown age
104 const double min_braking, ///< [in] Minimum braking index
105 const double max_braking ///< [in] Maximum braking index
106);
107
108///
109/// Set a second spindown bound derived from braking indices
110///
112 LatticeTiling *tiling, ///< [in] Lattice tiling
113 const size_t freq_dimension, ///< [in] Frequency dimension
114 const size_t f1dot_dimension, ///< [in] First spindown dimension
115 const size_t f2dot_dimension, ///< [in] Second spindown dimension
116 const double min_braking, ///< [in] Minimum braking index
117 const double max_braking ///< [in] Maximum braking index
118);
119
120///
121/// Return the step size of the spindown lattice tiling in a given dimension, or 0 for non-tiled dimensions.
122///
124 DopplerFullScanState *scan, ///< [in] Doppler scan state object
125 const size_t dim ///< [in] Dimension of which to return step size
126);
127
128// ========== deprecated LAL functions ==========
129void InitDopplerFullScan( LALStatus *, DopplerFullScanState **scanState, const DopplerFullScanInit *init );
130
131#ifdef __cplusplus
132}
133#endif
134/* C++ protection. */
135
136#endif /* Double-include protection. */
void XLALDestroyDopplerFullScan(DopplerFullScanState *scan)
Destroy the a full DopplerFullScanState structure.
REAL8 XLALNumDopplerTemplates(DopplerFullScanState *scan)
Return (and compute if not done so yet) the total number of Doppler templates of the DopplerScan scan...
void InitDopplerFullScan(LALStatus *, DopplerFullScanState **scanState, const DopplerFullScanInit *init)
int XLALSetLatticeTilingF2DotBrakingBound(LatticeTiling *tiling, const size_t freq_dimension, const size_t f1dot_dimension, const size_t f2dot_dimension, const double min_braking, const double max_braking)
Set a second spindown bound derived from braking indices.
UINT8 XLALNumDopplerPointsAtDimension(DopplerFullScanState *scan, const size_t dim)
Compute and return the total number of frequency and spindown points up to and along a certain dimens...
DopplerFullScanState * XLALInitDopplerFullScan(const DopplerFullScanInit *init)
Set up a full multi-dimensional grid-scan.
int XLALGetDopplerSpinRange(PulsarSpinRange *spinRange, const DopplerFullScanState *scan)
Return the spin-range spanned by the given template bank stored in the opaque DopplerFullScanState.
REAL8 XLALGetDopplerLatticeTilingStepSize(DopplerFullScanState *scan, const size_t dim)
Return the step size of the spindown lattice tiling in a given dimension, or 0 for non-tiled dimensio...
int XLALNextDopplerPos(PulsarDopplerParams *pos, DopplerFullScanState *scan)
Function to step through the full template grid point by point.
int XLALSetLatticeTilingF1DotAgeBrakingBound(LatticeTiling *tiling, const size_t freq_dimension, const size_t f1dot_dimension, const double age, const double min_braking, const double max_braking)
Set a first spindown bound derived from spindown age and braking indices.
DopplerGridType
different types of grids:
Definition: DopplerScan.h:89
unsigned char BOOLEAN
uint64_t UINT8
double REAL8
char CHAR
LALPulsarMetricType
Constants defining different types of pulsar-metrics.
Definition: PtoleMetric.h:95
Structure describing a region in paramter-space (a,d,f,f1dot,..).
LALPulsarMetricType metricType
which metric to use if GRID_METRIC
const CHAR * gridFile
filename for sky-grid or full-grid if GRID_FILE_SKYGRID or GRID_FILE_FULLGRID
REAL8 metricMismatch
for GRID_METRIC and GRID_ISOTROPIC
const LALDetector * Detector
Current detector.
LIGOTimeGPS startTime
start-time of the observation
REAL8 Tspan
total time spanned by observation
DopplerRegion searchRegion
Doppler-space region to be covered + scanned.
BOOLEAN projectMetric
project metric on f=const subspace
PulsarDopplerParams stepSizes
user-settings for stepsizes if GRID_FLAT
const EphemerisData * ephemeris
ephemeris-data for numerical metrics
DopplerGridType gridType
which type of grid to generate
This structure contains all information about the center-of-mass positions of the Earth and Sun,...
Type containing the 'Doppler-parameters' affecting the time-evolution of the phase.
Contains a "spin-range", ie spins and corresponding bands at a given (SSB) reference GPS-time .
--— internal [opaque] type to store the state of a FULL multidimensional grid-scan --—
PulsarSpinRange spinRange
spin-range covered by template bank