LALPulsar  6.1.0.1-c9a8ef6
DopplerFullScan.c
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2007, 2008, 2012 Karl Wette
3  * Copyright (C) 2006, 2007 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, 2007, 2008
24  * \file
25  * \brief Functions for handling "full" multi-dimensional search-grids for CFS.
26  * as opposed to the "factored" grids implemented in DopplerScan.[ch]
27  */
28 
29 /*---------- INCLUDES ----------*/
30 #include <math.h>
31 
32 #include <lal/LALStdlib.h>
33 #include <lal/FileIO.h>
34 #include <lal/DetectorSite.h>
35 #include <lal/LALError.h>
36 #include <lal/ConfigFile.h>
37 #include <lal/LogPrintf.h>
38 
39 #include <lal/DopplerFullScan.h>
40 
41 /*---------- DEFINES ----------*/
42 #define MIN(x,y) (x < y ? x : y)
43 #define MAX(x,y) (x > y ? x : y)
44 
45 #define TRUE (1==1)
46 #define FALSE (1==0)
47 
48 #ifdef __GNUC__
49 #define UNUSED __attribute__ ((unused))
50 #else
51 #define UNUSED
52 #endif
53 
54 /*---------- internal types ----------*/
55 typedef struct {
56  PulsarDopplerParams thisPoint; /**< current doppler-position of the scan */
57  DopplerSkyScanState skyScan; /**< keep track of sky-grid stepping */
59 
60 /** doubly linked list of REAL8-vectors (physical vectors) */
61 typedef struct tagREAL8VectorList {
63  struct tagREAL8VectorList *next;
64  struct tagREAL8VectorList *prev;
66 
67 /** ----- internal [opaque] type to store the state of a FULL multidimensional grid-scan ----- */
69  INT2 state; /**< idle, ready or finished */
70  DopplerGridType gridType; /**< what type of grid are we dealing with */
71  REAL8 numTemplates; /**< total number of templates in the grid */
72  PulsarSpinRange spinRange; /**< spin-range covered by template bank */
73  SkyRegion skyRegion; /**< sky-range covered by template bank */
74 
75  /* ----- full multi-dim parameter-space grid stuff ----- */
76  REAL8VectorList *covering; /**< multi-dimensional covering */
77  REAL8VectorList *thisGridPoint; /**< pointer to current grid-point */
78 
79  /* ----- spindown lattice tiling ----- */
80  LatticeTiling *spindownTiling; /**< spindown lattice tiling */
81  LatticeTilingIterator *spindownTilingItr; /**< iterator over spindown lattice tiling */
82  gsl_vector *spindownTilingPoint; /**< current point in spindown lattice tiling */
83 
84  /* ----- emulate old-style factored grids */
85  factoredGridScan_t *factoredScan; /**< only used to emulate FACTORED grids sky x Freq x f1dot */
86 
87 } /* struct DopplerFullScanState */;
88 
89 /*---------- internal prototypes ----------*/
90 int XLALInitFactoredGrid( DopplerFullScanState *scan, const DopplerFullScanInit *init );
91 int nextPointInFactoredGrid( PulsarDopplerParams *pos, DopplerFullScanState *scan );
92 int XLALLoadFullGridFile( DopplerFullScanState *scan, const DopplerFullScanInit *init );
93 
94 /*==================== FUNCTION DEFINITIONS ====================*/
95 
96 /**
97  * Set up a full multi-dimensional grid-scan.
98  * Currently this only emulates a 'factored' grid-scan with 'sky x Freq x f1dot ...' , but
99  * keeps all details within the DopplerScan module for future extension to real multidimensional
100  * grids.
101  *
102  * NOTE: Use 'XLALNextDopplerPos()' to step through this template grid.
103  *
104  */
105 DopplerFullScanState *
106 XLALInitDopplerFullScan( const DopplerFullScanInit *init /**< [in] initialization parameters */
107  )
108 {
109  XLAL_CHECK_NULL( init != NULL, XLAL_EINVAL );
110 
111  DopplerFullScanState *thisScan;
112  XLAL_CHECK_NULL( ( thisScan = LALCalloc( 1, sizeof( *thisScan ) ) ) != NULL, XLAL_ENOMEM );
113 
114  thisScan->gridType = init->gridType;
115 
116  /* store the user-input spinRange (includes refTime) in DopplerFullScanState */
117  thisScan->spinRange.refTime = init->searchRegion.refTime;
118  memcpy( thisScan->spinRange.fkdot, init->searchRegion.fkdot, sizeof( PulsarSpins ) );
119  memcpy( thisScan->spinRange.fkdotBand, init->searchRegion.fkdotBand, sizeof( PulsarSpins ) );
120 
121  // check that some old metric-codes aren't used with refTime!=startTime, which they don't handle correctly
122  switch ( thisScan->gridType ) {
123  case GRID_METRIC:
124  case GRID_METRIC_SKYFILE:
125  case GRID_SPINDOWN_SQUARE: /* square parameter space */
126  case GRID_SPINDOWN_AGEBRK: /* age-braking index parameter space */
127 
129  "NOTE: gridType={metric,4,spin-square,spin-age-brk} only work for refTime (%f) == startTime (%f)!\n",
130  XLALGPSGetREAL8( &( init->searchRegion.refTime ) ), XLALGPSGetREAL8( &( init->startTime ) ) );;
131 
132  break;
133  default:
134  break;
135  }
136 
137  /* which "class" of template grid to generate?: factored, or full-multidim ? */
138  switch ( thisScan->gridType ) {
139  /* emulate old 'factored' grids 'sky x f0dot x f1dot x f2dot x f3dot': */
140  case GRID_FLAT:
141  case GRID_ISOTROPIC:
142  case GRID_METRIC:
143  case GRID_FILE_SKYGRID:
144  case GRID_METRIC_SKYFILE:
145  /* backwards-compatibility mode */
147  break;
148 
149  /* ----- multi-dimensional covering of full parameter space ----- */
150  case GRID_FILE_FULLGRID:
152  break;
153 
154  case GRID_SPINDOWN_SQUARE: /* square parameter space */
155  case GRID_SPINDOWN_AGEBRK: { /* age-braking index parameter space */
156  const size_t n = 2 + PULSAR_MAX_SPINS;
157 
158  /* Check that the reference time is the same as the start time */
159  XLAL_CHECK_NULL( XLALGPSCmp( &thisScan->spinRange.refTime, &init->startTime ) == 0, XLAL_EINVAL,
160  "\nGRID_SPINDOWN_{SQUARE,AGEBRK}: This option currently restricts the reference time to be the same as the start time.\n" );
161 
162  /* Create a vector to hold lattice tiling parameter-space points */
163  XLAL_CHECK_NULL( ( thisScan->spindownTilingPoint = gsl_vector_alloc( n ) ) != NULL, XLAL_ENOMEM,
164  "\nGRID_SPINDOWN_{SQUARE,AGEBRK}: gsl_vector_alloc failed\n" );
165 
166  /* Create a lattice tiling */
167  XLAL_CHECK_NULL( ( thisScan->spindownTiling = XLALCreateLatticeTiling( n ) ) != NULL, XLAL_EFUNC );
168 
169  /* Parse the sky region string and check that it consists of only one point, and set bounds on it */
170  SkyRegion XLAL_INIT_DECL( sky );
172  XLAL_CHECK_NULL( sky.numVertices == 1, XLAL_EINVAL, "\nGRID_SPINDOWN_{SQUARE,AGEBRK}: This option can only handle a single sky position.\n" );
173  XLAL_CHECK_NULL( sky.vertices[0].system == COORDINATESYSTEM_EQUATORIAL, XLAL_EINVAL, "\nGRID_SPINDOWN_{SQUARE,AGEBRK}: This option only understands COORDINATESYSTEM_EQUATORIAL\n" );
174 
175  XLAL_CHECK_NULL( XLALSetLatticeTilingConstantBound( thisScan->spindownTiling, 0, sky.vertices[0].longitude, sky.vertices[0].longitude ) == XLAL_SUCCESS, XLAL_EFUNC );
176 
177  XLAL_CHECK_NULL( XLALSetLatticeTilingConstantBound( thisScan->spindownTiling, 1, sky.vertices[0].latitude, sky.vertices[0].latitude ) == XLAL_SUCCESS, XLAL_EFUNC );
178  if ( sky.vertices ) {
179  XLALFree( sky.vertices );
180  }
181 
182  /* Set up parameter space */
183  if ( thisScan->gridType == GRID_SPINDOWN_SQUARE ) { /* square parameter space */
184 
185  /* Set square bounds on the frequency and spindowns */
186  for ( size_t i = 0; i < PULSAR_MAX_SPINS; ++i ) {
187  XLAL_CHECK_NULL( XLALSetLatticeTilingConstantBound( thisScan->spindownTiling, 2 + i, init->searchRegion.fkdot[i], init->searchRegion.fkdot[i] + init->searchRegion.fkdotBand[i] ) == XLAL_SUCCESS, XLAL_EFUNC );
188  }
189  /* if requested by user, suppress padding in the fkdot (k>0) dimensions;
190  * if not, default padding is kept the same as in XLALCreateLatticeTiling()
191  */
192  if ( !init->extraArgs[0] ) {
193  for ( size_t i = 1; i < PULSAR_MAX_SPINS; ++i ) {
194  XLAL_CHECK_NULL( XLALSetLatticeTilingPadding( thisScan->spindownTiling, 2 + i, 0, 0, 0, 0, -1 ) == XLAL_SUCCESS, XLAL_EFUNC );
195  }
196  }
197 
198  } else if ( thisScan->gridType == GRID_SPINDOWN_AGEBRK ) { /* age-braking index parameter space */
199 
200  /* Get age and braking index from extra arguments */
201  const REAL8 spindownAge = init->extraArgs[0];
202  const REAL8 minBraking = init->extraArgs[1];
203  const REAL8 maxBraking = init->extraArgs[2];
204 
205  /* Set age-braking index parameter space */
206  XLAL_CHECK_NULL( XLAL_SUCCESS == XLALSetLatticeTilingConstantBound( thisScan->spindownTiling, 2, init->searchRegion.fkdot[0], init->searchRegion.fkdot[0] + init->searchRegion.fkdotBand[0] ), XLAL_EFUNC );
207  XLAL_CHECK_NULL( XLAL_SUCCESS == XLALSetLatticeTilingF1DotAgeBrakingBound( thisScan->spindownTiling, 2, 3, spindownAge, minBraking, maxBraking ), XLAL_EFUNC );
208  XLAL_CHECK_NULL( XLAL_SUCCESS == XLALSetLatticeTilingF2DotBrakingBound( thisScan->spindownTiling, 2, 3, 4, minBraking, maxBraking ), XLAL_EFUNC );
209 
210  /* This current only goes up to second spindown, so bound higher dimensions */
211  for ( size_t i = 3; i < PULSAR_MAX_SPINS; ++i ) {
212  XLAL_CHECK_NULL( XLAL_SUCCESS == XLALSetLatticeTilingConstantBound( thisScan->spindownTiling, 2 + i, init->searchRegion.fkdot[i], init->searchRegion.fkdot[i] + init->searchRegion.fkdotBand[i] ), XLAL_EFUNC );
213  }
214 
215  }
216 
217  /* Create a lattice tiling with Anstar lattice and spindown metric */
218  gsl_matrix *metric;
219  XLAL_CHECK_NULL( ( metric = gsl_matrix_alloc( n, n ) ) != NULL, XLAL_ENOMEM );
220  gsl_matrix_set_identity( metric );
221  gsl_matrix_view spin_metric = gsl_matrix_submatrix( metric, 2, 2, PULSAR_MAX_SPINS, PULSAR_MAX_SPINS );
222  XLAL_CHECK_NULL( XLALSpindownMetric( &spin_metric.matrix, init->Tspan ) == XLAL_SUCCESS, XLAL_EFUNC );
224 
225  /* Create iterator over flat lattice tiling */
226  XLAL_CHECK_NULL( ( thisScan->spindownTilingItr = XLALCreateLatticeTilingIterator( thisScan->spindownTiling, n ) ) != NULL, XLAL_EFUNC );
227 
228  /* Cleanup */
229  gsl_matrix_free( metric );
230 
231  }
232 
233  break;
234 
235  default:
236  XLAL_ERROR_NULL( XLAL_EINVAL, "\nInvalid grid type '%d'\n\n", init->gridType );
237  break;
238  } /* switch gridType */
239 
240  /* we're ready */
241  thisScan->state = STATE_READY;
242 
243  /* return result */
244  return thisScan;
245 
246 } // XLALInitDopplerFullScan()
247 
248 
249 /** \deprecated Use XLALInitDopplerFullScan() instead.
250  */
251 void
252 InitDopplerFullScan( LALStatus *status, /**< pointer to LALStatus structure */
253  DopplerFullScanState **scan, /**< [out] initialized Doppler scan state */
254  const DopplerFullScanInit *init /**< [in] initialization parameters */
255  )
256 {
257  INITSTATUS( status );
258 
259  XLAL_PRINT_DEPRECATION_WARNING( "XLALInitDopplerFullScan" );
260 
264 
265  if ( ( ( *scan ) = XLALInitDopplerFullScan( init ) ) == NULL ) {
267  }
268 
269  RETURN( status );
270 
271 } /* InitDopplerFullScan() */
272 
273 /**
274  * Return the spin-range spanned by the given template bank stored in the
275  * *opaque* DopplerFullScanState.
276  *
277  * \note The user cannot directly access any internal fields of that opaque type,
278  * which is why we need this API.
279  */
280 int
281 XLALGetDopplerSpinRange( PulsarSpinRange *spinRange, const DopplerFullScanState *scan )
282 {
283  if ( spinRange == NULL ) {
284  XLAL_ERROR( XLAL_EINVAL, "\nPulsarSpinRange pointer 'spinRange' is NULL\n" );
285  }
286 
287  if ( scan == NULL ) {
288  XLAL_ERROR( XLAL_EINVAL, "\nDopplerFullScanState pointer 'scan' is NULL\n" );
289  }
290 
291  ( *spinRange ) = scan->spinRange; // simple struct-copy is all that's needed
292 
293  return XLAL_SUCCESS;
294 
295 } /* XLALGetDopplerSpinRange() */
296 
297 /**
298  * Initialize Doppler-scanner to emulate an old-style factored template grid: 'sky x f0dot x f1dot x f2dot x f3dot'.
299  * This is a compatiblity-mode with the previous implementation currently also used in ComputeFstatistic.c.
300  */
301 int
302 XLALInitFactoredGrid( DopplerFullScanState *scan, /**< [bout] initialized Doppler scan state */
303  const DopplerFullScanInit *init /**< [in] initialization parameters */
304  )
305 {
306  XLAL_CHECK( scan, XLAL_EINVAL );
307  XLAL_CHECK( scan->state == STATE_IDLE, XLAL_EINVAL );
308  XLAL_CHECK( init, XLAL_EINVAL );
309 
310  DopplerSkyScanInit XLAL_INIT_DECL( skyScanInit );
311  SkyPosition skypos;
312  factoredGridScan_t *fscan = NULL;
313 
314  /* prepare initialization of DopplerSkyScanner to step through paramter space */
315  skyScanInit.dAlpha = init->stepSizes.Alpha;
316  skyScanInit.dDelta = init->stepSizes.Delta;
317  skyScanInit.gridType = init->gridType;
318  skyScanInit.metricType = init->metricType;
319  skyScanInit.metricMismatch = init->metricMismatch;
320  skyScanInit.projectMetric = init->projectMetric;
321  skyScanInit.obsBegin = init->startTime;
322  skyScanInit.obsDuration = init->Tspan;
323 
324  skyScanInit.Detector = init->Detector;
325  skyScanInit.ephemeris = init->ephemeris; /* used only by Ephemeris-based metric */
326  skyScanInit.skyGridFile = init->gridFile;
327  skyScanInit.skyRegionString = init->searchRegion.skyRegionString;
328  skyScanInit.Freq = init->searchRegion.fkdot[0] + init->searchRegion.fkdotBand[0];
329 
330  XLAL_CHECK( ( fscan = LALCalloc( 1, sizeof( *fscan ) ) ) != NULL, XLAL_ENOMEM );
331  scan->factoredScan = fscan;
332  XLAL_CHECK( XLALInitDopplerSkyScan( &( fscan->skyScan ), &skyScanInit ) == XLAL_SUCCESS, XLAL_EFUNC );
333 
334  /* overload spin step-sizes with user-settings if given */
335  for ( UINT4 i = 0; i < PULSAR_MAX_SPINS; i ++ ) {
336  if ( init->stepSizes.fkdot[i] ) {
337  fscan->skyScan.dfkdot[i] = init->stepSizes.fkdot[i];
338  }
339  }
340 
341  /* ----- set Doppler-scanner to start-point ----- */
342  fscan->thisPoint.refTime = init->searchRegion.refTime; /* set proper reference time for spins */
343 
344  /* normalize skyposition: correctly map into [0,2pi]x[-pi/2,pi/2] */
345  skypos.longitude = fscan->skyScan.skyNode->Alpha;
346  skypos.latitude = fscan->skyScan.skyNode->Delta;
348  XLALNormalizeSkyPosition( &skypos.longitude, &skypos.latitude );
349  fscan->thisPoint.Alpha = skypos.longitude;
350  fscan->thisPoint.Delta = skypos.latitude;
351  /* set spins to start */
352  for ( UINT4 i = 0; i < PULSAR_MAX_SPINS; i ++ ) {
353  fscan->thisPoint.fkdot[i] = scan->spinRange.fkdot[i];
354  }
355 
356  { /* count total number of templates */
357  REAL8 nSky, nTot;
358  REAL8 nSpins[PULSAR_MAX_SPINS];
359  nSky = fscan->skyScan.numSkyGridPoints;
360  for ( UINT4 i = 0; i < PULSAR_MAX_SPINS; i ++ ) {
361  nSpins[i] = floor( scan->spinRange.fkdotBand[i] / fscan->skyScan.dfkdot[i] ) + 1.0;
362  }
363  nTot = nSky;
364  for ( UINT4 i = 0; i < PULSAR_MAX_SPINS; i ++ ) {
365  nTot *= nSpins[i];
366  }
367  scan->numTemplates = nTot;
368  XLALPrintInfo( "Template grid: nSky x nFreq x nf1dot = %.0f x %.0f x %.0f = %.0f \n", nSky, nSpins[0], nSpins[1], nTot );
369  }
370 
371  return XLAL_SUCCESS;
372 
373 } // XLALInitFactoredGrid()
374 
375 
376 /**
377  * Return (and compute if not done so yet) the total number of Doppler templates
378  * of the DopplerScan \a scan
379  */
380 REAL8
381 XLALNumDopplerTemplates( DopplerFullScanState *scan )
382 {
383  if ( ! scan->numTemplates ) { /* not pre-computed already ? */
384  switch ( scan->gridType ) {
385  /* case GRID_METRIC_LATTICE: */
386  /* LogPrintf ( LOG_DEBUG, "Now counting number of templates in lattice ... "); */
387  /* scan->numTemplates = XLALCountLatticeTemplates ( scan->latticeScan ); */
388  /* LogPrintfVerbatim( LOG_DEBUG, " done. (%.0f)\n", scan->numTemplates ); */
389  /* break; */
390 
391  case GRID_SPINDOWN_SQUARE: /* square parameter space */
392  case GRID_SPINDOWN_AGEBRK: /* age-braking index parameter space */
393  LogPrintf( LOG_DEBUG, "Counting spindown lattice templates ... " );
394  scan->numTemplates = ( REAL8 )XLALTotalLatticeTilingPoints( scan->spindownTilingItr );
395  LogPrintfVerbatim( LOG_DEBUG, "%0.0f\n", scan->numTemplates );
396  break;
397 
398  default:
399  /* FIXME: not implemented yet */
400  LogPrintf( LOG_NORMAL, "template counting not implemented yet for gridType=%d!\n", scan->gridType );
401  return -1;
402  break;
403  } /* switch() */
404  } /* ! numTemplates */
405 
406  return scan->numTemplates;
407 
408 } /* XLALNumDopplerTemplates() */
409 
410 /**
411  * Function to step through the full template grid point by point.
412  * Normal return = 0,
413  * errors return -1,
414  * end of scan is signalled by return = 1
415  */
416 int
417 XLALNextDopplerPos( PulsarDopplerParams *pos, DopplerFullScanState *scan )
418 {
419 
420  /* This traps coding errors in the calling routine. */
421  if ( pos == NULL || scan == NULL ) {
423  }
424  if ( scan->state == STATE_IDLE ) {
425  XLALPrintError( "\nCalled XLALNextDopplerPos() on un-initialized DopplerFullScanState !\n\n" );
427  }
428 
429  /* is this search finished? then return '1' */
430  if ( scan->state == STATE_FINISHED ) {
431  return 1;
432  }
433 
434  // set refTime in returned template to the refTime of the grid
435  pos->refTime = scan->spinRange.refTime;
436 
437  /* ----- step foward one template in full grid ----- */
438  /* Which "class" of template grid are we dealing with: factored, or full-multidim ? */
439  switch ( scan->gridType ) {
440  /* emulate old 'factored' grids 'sky x f0dot x f1dot x f2dot x f3dot': */
441  case GRID_FLAT:
442  case GRID_ISOTROPIC:
443  case GRID_METRIC:
444  case GRID_FILE_SKYGRID:
445  case GRID_METRIC_SKYFILE:
446  /* backwards-compatibility mode */
447  nextPointInFactoredGrid( pos, scan );
448  break;
449 
450  case GRID_FILE_FULLGRID:
451  XLAL_INIT_MEM( pos->fkdot );
452  pos->fkdot[0] = scan->thisGridPoint->entry.data[0];
453  pos->Alpha = scan->thisGridPoint->entry.data[1];
454  pos->Delta = scan->thisGridPoint->entry.data[2];
455  pos->fkdot[1] = scan->thisGridPoint->entry.data[3];
456  pos->fkdot[2] = scan->thisGridPoint->entry.data[4];
457  pos->fkdot[3] = scan->thisGridPoint->entry.data[5];
458  pos->asini = 0; // isolated pulsar
459  /* advance to next grid point */
460  if ( ( scan->thisGridPoint = scan->thisGridPoint->next ) == NULL ) {
461  scan->state = STATE_FINISHED;
462  }
463 
464  break;
465 
466  /* case GRID_METRIC_LATTICE: */
467  /* if ( XLALgetCurrentDopplerPos ( pos, scan->latticeScan, COORDINATESYSTEM_EQUATORIAL ) ) { */
468  /* XLAL_ERROR ( XLAL_EFUNC ); */
469  /* } */
470  /* /\* advance to next point *\/ */
471  /* ret = XLALadvanceLatticeIndex ( scan->latticeScan ); */
472  /* if ( ret < 0 ) { */
473  /* XLAL_ERROR ( XLAL_EFUNC ); */
474  /* } */
475  /* else if ( ret == 1 ) */
476  /* { */
477  /* XLALPrintError ( "\n\nXLALadvanceLatticeIndex(): this was the last lattice points!\n\n"); */
478  /* scan->state = STATE_FINISHED; */
479  /* } */
480  /* #if 0 */
481  /* { /\* debugging *\/ */
482  /* gsl_vector_int *lal_index = NULL; */
483  /* XLALgetCurrentLatticeIndex ( &lal)index, scan->latticeScan ); */
484  /* XLALfprintfGSLvector_int ( stderr, "%d", lal_index ); */
485  /* gsl_vector_int_free ( lal_index ); */
486  /* } */
487  /* #endif */
488 
489  break;
490 
491  case GRID_SPINDOWN_SQUARE: /* square parameter space */
492  case GRID_SPINDOWN_AGEBRK: { /* age-braking index parameter space */
493 
494  /* Advance to next tile */
495  int retn = XLALNextLatticeTilingPoint( scan->spindownTilingItr, scan->spindownTilingPoint );
496  if ( retn < 0 ) {
497  XLALPrintError( "\nGRID_SPINDOWN_{SQUARE,AGEBRK}: XLALNextLatticeTilingPoint() failed\n" );
498  return -1;
499  }
500 
501  if ( retn > 0 ) {
502 
503  /* Found a point */
504  pos->Alpha = gsl_vector_get( scan->spindownTilingPoint, 0 );
505  pos->Delta = gsl_vector_get( scan->spindownTilingPoint, 1 );
506  pos->fkdot[0] = gsl_vector_get( scan->spindownTilingPoint, 2 );
507  for ( size_t i = 1; i < PULSAR_MAX_SPINS; ++i ) {
508  pos->fkdot[i] = gsl_vector_get( scan->spindownTilingPoint, i + 2 );
509  }
510 
511  return 0;
512 
513  } else {
514 
515  /* No more points */
516  scan->state = STATE_FINISHED;
517  return 1;
518 
519  }
520 
521  }
522 
523  break;
524 
525  default:
526  XLALPrintError( "\nInvalid grid type '%d'\n\n", scan->gridType );
528  return -1;
529  break;
530  } /* switch gridType */
531 
532  return 0;
533 
534 } /* XLALNextDopplerPos() */
535 
536 /**
537  * return current grid-point and step forward one template in 'factored' grids (sky x f0dot x f1dot ... )
538  * return 0 = OK, -1 = ERROR
539  */
540 int
541 nextPointInFactoredGrid( PulsarDopplerParams *pos, DopplerFullScanState *scan )
542 {
544  factoredGridScan_t *fscan;
546  PulsarSpins fkdotMax;
547  SkyPosition skypos;
548 
549  if ( pos == NULL || scan == NULL ) {
550  return -1;
551  }
552 
553  if ( ( fscan = scan->factoredScan ) == NULL ) {
554  return -1;
555  }
556 
557  range = &( scan->spinRange ); /* shortcut */
558 
559  ( *pos ) = fscan->thisPoint; /* RETURN current Doppler-point (struct-copy) */
560 
561  nextPos = fscan->thisPoint; /* struct-copy: start from current point to get next one */
562 
563  /* shortcuts: spin boundaries */
564  fkdotMax[3] = range->fkdot[3] + range->fkdotBand[3];
565  fkdotMax[2] = range->fkdot[2] + range->fkdotBand[2];
566  fkdotMax[1] = range->fkdot[1] + range->fkdotBand[1];
567  fkdotMax[0] = range->fkdot[0] + range->fkdotBand[0];
568 
569  /* try to advance to next template */
570  nextPos.fkdot[0] += fscan->skyScan.dfkdot[0]; /* f0dot one forward */
571  if ( nextPos.fkdot[0] > fkdotMax[0] ) { /* too far? */
572  nextPos.fkdot[0] = range->fkdot[0]; /* f0dot return to start */
573  nextPos.fkdot[1] += fscan->skyScan.dfkdot[1]; /* f1dot one step forward */
574  if ( nextPos.fkdot[1] > fkdotMax[1] ) {
575  nextPos.fkdot[1] = range->fkdot[1]; /* f1dot return to start */
576  nextPos.fkdot[2] += fscan->skyScan.dfkdot[2]; /* f2dot one forward */
577  if ( nextPos.fkdot[2] > fkdotMax[2] ) {
578  nextPos.fkdot[2] = range->fkdot[2]; /* f2dot return to start */
579  nextPos.fkdot[3] += fscan->skyScan.dfkdot[3]; /* f3dot one forward */
580  if ( nextPos.fkdot[3] > fkdotMax[3] ) {
581  nextPos.fkdot[3] = range->fkdot[3]; /* f3dot return to start */
582  /* skygrid one forward */
583  if ( ( fscan->skyScan.skyNode = fscan->skyScan.skyNode->next ) == NULL ) { /* no more sky-points ?*/
584  fscan->skyScan.state = STATE_FINISHED; /* avoid warning when freeing */
585  scan->state = STATE_FINISHED; /* we're done */
586  } else {
587  /* normalize next skyposition: correctly map into [0,2pi]x[-pi/2,pi/2] */
588  skypos.longitude = fscan->skyScan.skyNode->Alpha;
589  skypos.latitude = fscan->skyScan.skyNode->Delta;
591  XLALNormalizeSkyPosition( &skypos.longitude, &skypos.latitude );
592  nextPos.Alpha = skypos.longitude;
593  nextPos.Delta = skypos.latitude;
594  }
595  } /* f3dot */
596  } /* f2dot */
597  } /* f1dot */
598  } /* f0dot */
599 
600  /* prepare next step */
601  fscan->thisPoint = nextPos;
602 
603  return 0;
604 
605 } /* nextPointInFactoredGrid() */
606 
607 static void
609 {
610  REAL8VectorList *ptr, *next;
611 
612  if ( !head ) {
613  return;
614  }
615 
616  next = head;
617 
618  do {
619  /* step to next element */
620  ptr = next;
621  /* remember pointer to next element */
622  next = ptr->next;
623  /* free current element */
624  LALFree( ptr->entry.data );
625  LALFree( ptr );
626 
627  } while ( ( ptr = next ) != NULL );
628 
629  return;
630 } /* XLALREAL8VectorListDestroy() */
631 
632 /**
633  * Destroy the a full DopplerFullScanState structure
634  */
635 void
636 XLALDestroyDopplerFullScan( DopplerFullScanState *scan )
637 {
638  if ( scan == NULL ) {
639  return;
640  }
641 
642  if ( scan->factoredScan ) {
643  XLALDestroyDopplerSkyScan( &( scan->factoredScan->skyScan ) );
644  XLALFree( scan->factoredScan );
645  }
646 
647  if ( scan->covering ) {
648  XLALREAL8VectorListDestroy( scan->covering );
649  }
650 
651  if ( scan->spindownTiling ) {
652  XLALDestroyLatticeTiling( scan->spindownTiling );
653  }
654  if ( scan->spindownTilingItr ) {
655  XLALDestroyLatticeTilingIterator( scan->spindownTilingItr );
656  }
657  if ( scan->spindownTilingPoint ) {
658  gsl_vector_free( scan->spindownTilingPoint );
659  }
660 
661  if ( scan->skyRegion.vertices ) {
662  XLALFree( scan->skyRegion.vertices );
663  }
664 
665  XLALFree( scan );
666 
667  return;
668 
669 } // XLALDestroyDopplerFullScan()
670 
671 
672 static REAL8VectorList *
674 {
675  UINT4 dim;
676  REAL8VectorList *ptr = NULL; /* running list-pointer */
677  REAL8VectorList *newElement = NULL; /* new list-element */
678  /* check illegal input */
679  if ( ( head == NULL ) || ( entry == NULL ) ) {
680  return NULL;
681  }
682 
683  /* find tail of list */
684  ptr = head;
685  while ( ptr->next ) {
686  ptr = ptr->next;
687  }
688 
689  /* construct new list-element */
690  dim = entry->length;
691  if ( ( newElement = LALCalloc( 1, sizeof( *newElement ) ) ) == NULL ) {
692  return NULL;
693  }
694  if ( ( newElement->entry.data = LALCalloc( dim, sizeof( entry->data[0] ) ) ) == NULL ) {
695  LALFree( newElement );
696  return NULL;
697  }
698  newElement->entry.length = dim;
699  memcpy( newElement->entry.data, entry->data, dim * sizeof( entry->data[0] ) );
700 
701  /* link this to the tail of list */
702  ptr->next = newElement;
703  newElement->prev = ptr;
704 
705  return newElement;
706 
707 } /* XLALREAL8VectorListAddEntry() */
708 
709 /**
710  * load a full multi-dim template grid from the file init->gridFile,
711  * the file-format is: lines of 6 columns, which are:
712  *
713  * Freq Alpha Delta f1dot f2dot f3dot
714  *
715  * \note
716  * *) this function returns the effective spinRange covered by the read-in template bank
717  * by storing it in scan->spinRange, potentially overwriting any previous user-input values in there.
718  *
719  * *) a possible future extension should probably *clip* the template-bank to the user-specified ranges,
720  * then return the effective ranges spanned by the resultant template bank.
721  *
722  * *) in order to avoid surprises until such a feature is implemented, we currently return an error if
723  * any of the input spinRanges are non-zero
724  *
725  */
726 int
727 XLALLoadFullGridFile( DopplerFullScanState *scan,
728  const DopplerFullScanInit *init
729  )
730 {
731  XLAL_CHECK( ( scan != NULL ) && ( init != NULL ), XLAL_EINVAL );
732  XLAL_CHECK( init->gridFile != NULL, XLAL_EINVAL );
733  XLAL_CHECK( scan->state == STATE_IDLE, XLAL_EINVAL );
734 
736  REAL8VectorList *tail = NULL;
737  REAL8Vector *entry = NULL;
738  UINT4 numTemplates;
739  FILE *fp;
740 
741 
742  /* Check that all user-input spin- and sky-ranges are zero, otherwise fail!
743  *
744  * NOTE: In the future we should allow combining the user-input ranges with
745  * those found in the grid-file by forming the intersection, ie *clipping*
746  * of the read-in grids to the user-input ranges.
747  * Right now we require empty ranges input, and report back the ranges from the grid-file
748  *
749  */
750  if ( init->searchRegion.skyRegionString != NULL ) {
751  XLAL_ERROR( XLAL_EINVAL, "\nnon-NULL skyRegion input currently not supported! skyRegion = '%s'\n\n", init->searchRegion.skyRegionString );
752  }
753  for ( UINT4 s = 0; s < PULSAR_MAX_SPINS; s ++ ) {
754  if ( ( init->searchRegion.fkdot[s] != 0 ) || ( init->searchRegion.fkdotBand[s] != 0 ) ) {
755  XLAL_ERROR( XLAL_EINVAL, "\nnon-zero input spinRanges currently not supported! fkdot[%d] = %g, fkdotBand[%d] = %g\n\n",
756  s, init->searchRegion.fkdot[s], s, init->searchRegion.fkdotBand[s] );
757  }
758  } /* for s < max_spins */
759 
760  /* open input data file */
761  XLAL_CHECK( ( fp = fopen( init->gridFile, "r" ) ) != NULL, XLAL_ESYS, "Could not open data-file: `%s`\n\n", init->gridFile );
762 
763  /* prepare grid-entry buffer */
764  XLAL_CHECK( ( entry = XLALCreateREAL8Vector( 6 ) ) != NULL, XLAL_EFUNC );
765 
766  /* keep track of the sky- and spinRanges spanned by the template bank */
767  REAL8 FreqMax = - LAL_REAL4_MAX, FreqMin = LAL_REAL4_MAX; // only using REAL4 ranges to avoid over/under flows, and should be enough
768  REAL8 f1dotMax = - LAL_REAL4_MAX, f1dotMin = LAL_REAL4_MAX;
769  REAL8 f2dotMax = - LAL_REAL4_MAX, f2dotMin = LAL_REAL4_MAX;
770  REAL8 f3dotMax = - LAL_REAL4_MAX, f3dotMin = LAL_REAL4_MAX;
771 
772  REAL8 alphaMax = - LAL_REAL4_MAX, alphaMin = LAL_REAL4_MAX;
773  REAL8 deltaMax = - LAL_REAL4_MAX, deltaMin = LAL_REAL4_MAX;
774 
775  /* parse this list of lines into a full grid */
776  numTemplates = 0;
777  tail = &head; /* head will remain empty! */
778  CHAR line[2048];
779  while ( ! feof( fp ) && ! ferror( fp ) && fgets( line, sizeof( line ), fp ) != NULL ) {
780 
781  // Skip over any comment lines
782  if ( line[0] == '#' || line[0] == '%' ) {
783  continue;
784  }
785 
786  // File format expects lines containing 6 columns: Freq Alpha Delta f1dot f2dot f3dot
787  REAL8 Freq, Alpha, Delta, f1dot, f2dot, f3dot;
788  if ( 6 != sscanf( line, "%" LAL_REAL8_FORMAT " %" LAL_REAL8_FORMAT " %" LAL_REAL8_FORMAT " %" LAL_REAL8_FORMAT " %" LAL_REAL8_FORMAT " %" LAL_REAL8_FORMAT "\n",
789  &Freq, &Alpha, &Delta, &f1dot, &f2dot, &f3dot ) ) {
790  XLALPrintError( "ERROR: Failed to parse 6 REAL8's from line %d in grid-file '%s'\n\n", numTemplates + 1, init->gridFile );
791  if ( head.next ) {
793  }
795  }
796 
797  /* keep track of maximal spans */
798  alphaMin = fmin( Alpha, alphaMin );
799  deltaMin = fmin( Delta, deltaMin );
800  FreqMin = fmin( Freq, FreqMin );
801  f1dotMin = fmin( f1dot, f1dotMin );
802  f2dotMin = fmin( f2dot, f2dotMin );
803  f3dotMin = fmin( f3dot, f3dotMin );
804 
805  alphaMax = fmax( Alpha, alphaMax );
806  deltaMax = fmax( Delta, deltaMax );
807  FreqMax = fmax( Freq, FreqMax );
808  f1dotMax = fmax( f1dot, f1dotMax );
809  f2dotMax = fmax( f2dot, f2dotMax );
810  f3dotMax = fmax( f3dot, f3dotMax );
811 
812  /* add this entry to template-bank list */
813  entry->data[0] = Freq;
814  entry->data[1] = Alpha;
815  entry->data[2] = Delta;
816  entry->data[3] = f1dot;
817  entry->data[4] = f2dot;
818  entry->data[5] = f3dot;
819 
820  if ( ( tail = XLALREAL8VectorListAddEntry( tail, entry ) ) == NULL ) {
821  if ( head.next ) {
823  }
825  }
826 
827  numTemplates ++ ;
828 
829  } /* while !feof(fp) && ... */
830  if ( ferror( fp ) ) {
831  XLAL_ERROR( XLAL_EIO );
832  }
833 
834  XLALDestroyREAL8Vector( entry );
835 
836  /* ---------- update scan-state ---------- */
837 
838  // ----- report back ranges actually spanned by grid-file
839 
840  CHAR *skyRegionString = NULL;
842  XLAL_CHECK( ( skyRegionString = XLALSkySquare2String( alphaMin, deltaMin, ( alphaMax - alphaMin ) + eps, ( deltaMax - deltaMin ) + eps ) ) != NULL, XLAL_EFUNC );
843 
844  // note: we slight expanded the enclosing sky-square by eps to avoid complaints when a grid-file contains
845  // only points in a line, which is perfectly valid here.
846  XLAL_CHECK( XLALParseSkyRegionString( &scan->skyRegion, skyRegionString ) == XLAL_SUCCESS, XLAL_EFUNC );
847  XLALFree( skyRegionString );
848 
849  scan->spinRange.fkdot[0] = FreqMin;
850  scan->spinRange.fkdotBand[0] = FreqMax - FreqMin;
851 
852  scan->spinRange.fkdot[1] = f1dotMin;
853  scan->spinRange.fkdotBand[1] = f1dotMax - f1dotMin;
854 
855  scan->spinRange.fkdot[2] = f2dotMin;
856  scan->spinRange.fkdotBand[2] = f2dotMax - f2dotMin;
857 
858  scan->spinRange.fkdot[3] = f3dotMin;
859  scan->spinRange.fkdotBand[3] = f3dotMax - f3dotMin;
860 
861  scan->numTemplates = numTemplates;
862  scan->covering = head.next; /* pass result (without head!) */
863  scan->thisGridPoint = scan->covering; /* init to start */
864 
865  XLALPrintInfo( "Template grid: nTot = %.0f\n", 1.0 * numTemplates );
866  XLALPrintInfo( "Spanned ranges: Freq in [%g, %g], f1dot in [%g, %g], f2dot in [%g, %g], f3dot in [%g, %g]\n",
867  FreqMin, FreqMax, f1dotMin, f1dotMax, f2dotMin, f2dotMax, f3dotMin, f3dotMax );
868 
869  return XLAL_SUCCESS;
870 
871 } /* XLALLoadFullGridFile() */
872 
873 typedef struct {
874  size_t freq_dim;
875  double scale;
877 
878 static double F1DotAgeBrakingBound(
879  const void *data,
880  const size_t dim UNUSED,
881  const gsl_matrix *cache UNUSED,
882  const gsl_vector *point
883 )
884 {
885 
886  // Get bounds info
888 
889  // Get current value of frequency
890  const double freq = gsl_vector_get( point, info->freq_dim );
891 
892  // Return first spindown bounds
893  return info->scale * freq;
894 
895 }
896 
898  LatticeTiling *tiling,
899  const size_t freq_dimension,
900  const size_t f1dot_dimension,
901  const double age,
902  const double min_braking,
903  const double max_braking
904 )
905 {
906 
907  // Check input
908  XLAL_CHECK( tiling != NULL, XLAL_EFAULT );
909  XLAL_CHECK( freq_dimension < f1dot_dimension, XLAL_EINVAL );
910  XLAL_CHECK( age > 0.0, XLAL_EINVAL );
911  XLAL_CHECK( min_braking > 1.0, XLAL_EINVAL );
912  XLAL_CHECK( max_braking > 1.0, XLAL_EINVAL );
913  XLAL_CHECK( min_braking <= max_braking, XLAL_EINVAL );
914 
915  // Set the parameter-space bound
918  info_lower.freq_dim = info_upper.freq_dim = freq_dimension;
919  info_lower.scale = -1.0 / ( ( min_braking - 1.0 ) * age );
920  info_upper.scale = -1.0 / ( ( max_braking - 1.0 ) * age );
921  XLAL_CHECK( XLALSetLatticeTilingBound( tiling, f1dot_dimension, F1DotAgeBrakingBound, sizeof( info_lower ), &info_lower, &info_upper ) == XLAL_SUCCESS, XLAL_EFAILED );
922 
923  return XLAL_SUCCESS;
924 
925 }
926 
927 typedef struct {
928  size_t freq_dim;
929  size_t f1dot_dim;
930  double scale;
932 
933 static double F2DotBrakingBound(
934  const void *data,
935  const size_t dim UNUSED,
936  const gsl_matrix *cache UNUSED,
937  const gsl_vector *point
938 )
939 {
940 
941  // Get bounds info
943 
944  // Get current values of frequency and first spindown
945  const double freq = gsl_vector_get( point, info->freq_dim );
946  const double f1dot = gsl_vector_get( point, info->f1dot_dim );
947 
948  // Return second spindown bounds
949  return info->scale * f1dot * f1dot / freq;
950 
951 }
952 
954  LatticeTiling *tiling,
955  const size_t freq_dimension,
956  const size_t f1dot_dimension,
957  const size_t f2dot_dimension,
958  const double min_braking,
959  const double max_braking
960 )
961 {
962 
963  // Check input
964  XLAL_CHECK( tiling != NULL, XLAL_EFAULT );
965  XLAL_CHECK( freq_dimension < f1dot_dimension, XLAL_EINVAL );
966  XLAL_CHECK( f1dot_dimension < f2dot_dimension, XLAL_EINVAL );
967  XLAL_CHECK( min_braking > 0.0, XLAL_EINVAL );
968  XLAL_CHECK( max_braking > 0.0, XLAL_EINVAL );
969  XLAL_CHECK( min_braking <= max_braking, XLAL_EINVAL );
970 
971  // Set the parameter-space bound
974  info_lower.freq_dim = info_upper.freq_dim = freq_dimension;
975  info_lower.f1dot_dim = info_upper.f1dot_dim = f1dot_dimension;
976  info_lower.scale = min_braking;
977  info_upper.scale = max_braking;
978  XLAL_CHECK( XLALSetLatticeTilingBound( tiling, f2dot_dimension, F2DotBrakingBound, sizeof( info_lower ), &info_lower, &info_upper ) == XLAL_SUCCESS, XLAL_EFAILED );
979 
980  return XLAL_SUCCESS;
981 
982 }
const REAL8 eps
static double F2DotBrakingBound(const void *data, const size_t dim UNUSED, const gsl_matrix *cache UNUSED, const gsl_vector *point)
int nextPointInFactoredGrid(PulsarDopplerParams *pos, DopplerFullScanState *scan)
return current grid-point and step forward one template in 'factored' grids (sky x f0dot x f1dot ....
int XLALInitFactoredGrid(DopplerFullScanState *scan, const DopplerFullScanInit *init)
Initialize Doppler-scanner to emulate an old-style factored template grid: 'sky x f0dot x f1dot x f2d...
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...
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.
static void XLALREAL8VectorListDestroy(REAL8VectorList *head)
static REAL8VectorList * XLALREAL8VectorListAddEntry(REAL8VectorList *head, const REAL8Vector *entry)
void InitDopplerFullScan(LALStatus *status, DopplerFullScanState **scan, const DopplerFullScanInit *init)
int XLALGetDopplerSpinRange(PulsarSpinRange *spinRange, const DopplerFullScanState *scan)
Return the spin-range spanned by the given template bank stored in the opaque DopplerFullScanState.
int XLALLoadFullGridFile(DopplerFullScanState *scan, const DopplerFullScanInit *init)
load a full multi-dim template grid from the file init->gridFile, the file-format is: lines of 6 colu...
int XLALNextDopplerPos(PulsarDopplerParams *pos, DopplerFullScanState *scan)
Function to step through the full template grid point by point.
DopplerFullScanState * XLALInitDopplerFullScan(const DopplerFullScanInit *init)
Set up a full multi-dimensional grid-scan.
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.
static double F1DotAgeBrakingBound(const void *data, const size_t dim UNUSED, const gsl_matrix *cache UNUSED, const gsl_vector *point)
CHAR * XLALSkySquare2String(REAL8 Alpha, REAL8 Delta, REAL8 AlphaBand, REAL8 DeltaBand)
parse a 'classical' sky-square (Alpha, Delta, AlphaBand, DeltaBand) into a "SkyRegion"-string of the ...
Definition: DopplerScan.c:1295
int XLALParseSkyRegionString(SkyRegion *region, const CHAR *input)
parse a skyRegion-string into a SkyRegion structure: the expected string-format is " (longitude1,...
Definition: DopplerScan.c:1203
int XLALInitDopplerSkyScan(DopplerSkyScanState *skyScan, const DopplerSkyScanInit *init)
Initialize the Doppler sky-scanner.
Definition: DopplerScan.c:170
void XLALDestroyDopplerSkyScan(DopplerSkyScanState *skyScan)
Destroy the DopplerSkyScanState structure.
Definition: DopplerScan.c:294
#define DOPPLERSCANH_MSGENONULL
Definition: DopplerScan.h:71
#define DOPPLERSCANH_ENONULL
Definition: DopplerScan.h:55
#define DOPPLERSCANH_MSGEXLAL
Definition: DopplerScan.h:76
#define DOPPLERSCANH_EINPUT
Definition: DopplerScan.h:58
#define DOPPLERSCANH_MSGEINPUT
Definition: DopplerScan.h:74
#define DOPPLERSCANH_EXLAL
Definition: DopplerScan.h:60
@ STATE_READY
initialized and ready
Definition: DopplerScan.h:83
@ STATE_FINISHED
all templates have been read
Definition: DopplerScan.h:84
@ STATE_IDLE
not initialized yet
Definition: DopplerScan.h:82
DopplerGridType
different types of grids:
Definition: DopplerScan.h:89
@ GRID_METRIC
generate grid using a 2D sky-metric
Definition: DopplerScan.h:93
@ GRID_SPINDOWN_SQUARE
spindown tiling for a single sky position and square parameter space
Definition: DopplerScan.h:100
@ GRID_FLAT
"flat" sky-grid: fixed step-size (dAlpha,dDelta)
Definition: DopplerScan.h:91
@ GRID_FILE_SKYGRID
read skygrid from a file
Definition: DopplerScan.h:94
@ GRID_FILE_FULLGRID
load the full D-dim grid from a file
Definition: DopplerScan.h:98
@ GRID_METRIC_SKYFILE
'hybrid' grid-construction: use skygrid from file, metric for others
Definition: DopplerScan.h:95
@ GRID_ISOTROPIC
approximately isotropic sky-grid
Definition: DopplerScan.h:92
@ GRID_SPINDOWN_AGEBRK
spindown tiling for a single sky position and non-square parameter space defined by spindown age and ...
Definition: DopplerScan.h:101
ProcessParamsTable * ptr
#define LALCalloc(m, n)
#define LALFree(p)
#define ABORT(statusptr, code, mesg)
#define ASSERT(assertion, statusptr, code, mesg)
#define INITSTATUS(statusptr)
#define RETURN(statusptr)
int s
#define LAL_REAL4_MAX
#define LAL_REAL8_EPS
#define XLAL_INIT_MEM(x)
double REAL8
#define XLAL_INIT_DECL(var,...)
int16_t INT2
char CHAR
uint32_t UINT4
void XLALFree(void *p)
#define LAL_REAL8_FORMAT
int XLALSetLatticeTilingBound(LatticeTiling *tiling, const size_t dim, const LatticeTilingBound func, const size_t data_len, const void *data_lower, const void *data_upper)
Set a parameter-space bound on a dimension of the lattice tiling.
UINT8 XLALTotalLatticeTilingPoints(const LatticeTilingIterator *itr)
Return the total number of points covered by the lattice tiling iterator.
int XLALSetLatticeTilingConstantBound(LatticeTiling *tiling, const size_t dim, const double bound1, const double bound2)
Set a constant lattice tiling parameter-space bound, given by the minimum and maximum of the two supp...
int XLALSetTilingLatticeAndMetric(LatticeTiling *tiling, const TilingLattice lattice, const gsl_matrix *metric, const double max_mismatch)
Set the tiling lattice, parameter-space metric, and maximum prescribed mismatch.
LatticeTilingIterator * XLALCreateLatticeTilingIterator(const LatticeTiling *tiling, const size_t itr_ndim)
Create a new lattice tiling iterator.
void XLALDestroyLatticeTilingIterator(LatticeTilingIterator *itr)
Destroy a lattice tiling iterator.
void XLALDestroyLatticeTiling(LatticeTiling *tiling)
Destroy a lattice tiling.
int XLALNextLatticeTilingPoint(LatticeTilingIterator *itr, gsl_vector *point)
Advance lattice tiling iterator, and optionally return the next point in point.
LatticeTiling * XLALCreateLatticeTiling(const size_t ndim)
Create a new lattice tiling.
int XLALSetLatticeTilingPadding(LatticeTiling *tiling, const size_t dim, const double lower_bbox_pad, const double upper_bbox_pad, const int lower_intp_pad, const int upper_intp_pad, const int find_bound_extrema)
Control the padding of lattice tiling parameter-space bounds in the given dimension.
@ TILING_LATTICE_ANSTAR
An-star ( ) lattice.
Definition: LatticeTiling.h:65
void void LogPrintfVerbatim(LogLevel_t, const char *format,...) _LAL_GCC_PRINTF_FORMAT_(2
void LogPrintf(LogLevel_t, const char *format,...) _LAL_GCC_PRINTF_FORMAT_(2
LOG_DEBUG
LOG_NORMAL
int XLALSpindownMetric(gsl_matrix *metric, double Tspan)
Frequency and frequency derivative components of the metric, suitable for a directed search with only...
Definition: PtoleMetric.c:828
REAL8 PulsarSpins[PULSAR_MAX_SPINS]
Typedef for fixed-size array holding GW frequency and derivatives fk = d^k Freq/dt^k|(tau_ref)
#define PULSAR_MAX_SPINS
maximal number of spin-parameters (Freq + spindowns) we can handle
void XLALNormalizeSkyPosition(double *restrict longitude, double *restrict latitude)
COORDINATESYSTEM_EQUATORIAL
REAL8Vector * XLALCreateREAL8Vector(UINT4 length)
void XLALDestroyREAL8Vector(REAL8Vector *vector)
#define XLAL_ERROR_NULL(...)
#define xlalErrno
int int int XLALPrintInfo(const char *fmt,...) _LAL_GCC_PRINTF_FORMAT_(1
#define XLAL_ERROR(...)
#define XLAL_CHECK(assertion,...)
int XLALPrintError(const char *fmt,...) _LAL_GCC_PRINTF_FORMAT_(1
#define XLAL_CHECK_NULL(assertion,...)
#define XLAL_PRINT_DEPRECATION_WARNING(replacement)
XLAL_ENOMEM
XLAL_SUCCESS
XLAL_EFAULT
XLAL_EFUNC
XLAL_EIO
XLAL_ESYS
XLAL_EINVAL
XLAL_EFAILED
int XLALGPSCmp(const LIGOTimeGPS *t0, const LIGOTimeGPS *t1)
REAL8 XLALGPSGetREAL8(const LIGOTimeGPS *epoch)
REAL8 XLALGPSDiff(const LIGOTimeGPS *t1, const LIGOTimeGPS *t0)
head
pos
n
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
REAL8 extraArgs[3]
extra grid-specific setup arguments
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
PulsarSpins fkdotBand
spin-intervals
Definition: DopplerScan.h:119
CHAR * skyRegionString
sky-region string '(a1,d1), (a2,d2), ..'
Definition: DopplerScan.h:116
LIGOTimeGPS refTime
Definition: DopplerScan.h:117
PulsarSpins fkdot
reference time of definition of Doppler parameters
Definition: DopplerScan.h:118
struct tagDopplerSkyGrid * next
Definition: DopplerScan.h:127
initialization-structure passed to InitDopplerSkyScan()
Definition: DopplerScan.h:134
this structure reflects the current state of a DopplerSkyScan
Definition: DopplerScan.h:152
UINT4 numSkyGridPoints
how many skygrid-points
Definition: DopplerScan.h:155
PulsarSpins dfkdot
fixed-size steps in spins
Definition: DopplerScan.h:156
DopplerSkyGrid * skyNode
pointer to current grid-node in skygrid
Definition: DopplerScan.h:158
scan_state_t state
idle, ready or finished
Definition: DopplerScan.h:153
Type containing the 'Doppler-parameters' affecting the time-evolution of the phase.
PulsarSpins fkdot
Intrinsic spins: [Freq, f1dot, f2dot, ...
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.
Contains a "spin-range", ie spins and corresponding bands at a given (SSB) reference GPS-time .
REAL8 * data
doubly linked list of REAL8-vectors (physical vectors)
REAL8Vector entry
struct tagREAL8VectorList * next
struct tagREAL8VectorList * prev
REAL8 longitude
REAL8 latitude
CoordinateSystem system
structure describing a polygon-region in the sky
Definition: DopplerScan.h:108
PulsarDopplerParams thisPoint
current doppler-position of the scan
DopplerSkyScanState skyScan
keep track of sky-grid stepping
--— internal [opaque] type to store the state of a FULL multidimensional grid-scan --—
REAL8VectorList * covering
multi-dimensional covering
SkyRegion skyRegion
sky-range covered by template bank
factoredGridScan_t * factoredScan
only used to emulate FACTORED grids sky x Freq x f1dot
REAL8VectorList * thisGridPoint
pointer to current grid-point
INT2 state
idle, ready or finished
LatticeTiling * spindownTiling
spindown lattice tiling
DopplerGridType gridType
what type of grid are we dealing with
PulsarSpinRange spinRange
spin-range covered by template bank
REAL8 numTemplates
total number of templates in the grid
LatticeTilingIterator * spindownTilingItr
iterator over spindown lattice tiling
gsl_vector * spindownTilingPoint
current point in spindown lattice tiling