Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALApps 10.1.0.1-5e288d3
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
tmpltbank.c
Go to the documentation of this file.
1/*
2* Copyright (C) 2007 Alexander Dietz, Duncan Brown, Eirini Messaritaki, Gareth Jones, Benjamin Owen, Patrick Brady, Robert Adam Mercer, Stephen Fairhurst, Craig Robinson , Thomas Cokelaer, Evan Ochsner
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/*-----------------------------------------------------------------------
21 *
22 * File Name: tmpltbank.c
23 *
24 * Author: Brown, D. A.
25 *
26 *
27 *-----------------------------------------------------------------------
28 */
29
30/**
31 * \file
32 * \ingroup lalapps_inspiral
33 *
34 * <dl>
35 *
36 * <dt>Name</dt><dd>
37 * \c lalapps_tmpltbank --- program to generate inspiral template banks.</dd>
38 *
39 * <dt>Synopsis</dt><dd>
40 * <tt>lalapps_tmpltbank</tt>
41 * [<tt>--help</tt>]
42 * [<tt>--verbose</tt>]
43 * [<tt>--version</tt>]
44 * [<tt>--user-tag</tt> <i>usertag</i>]
45 * [<tt>--comment</tt> <i>comment</i>]
46 * <tt>--gps-start-time</tt> <i>gps_start</i>
47 * <tt>--gps-end-time</tt> <i>gps_end</i>
48 * [<tt>--pad-data</tt> <i>time_pad</i>]
49 * [<tt>--glob-frame-data</tt>]
50 * [<tt>--frame-type</tt> <i>type</i>]
51 * [<tt>--frame-cache</tt> <i>cache_file</i>]
52 * <tt>--calibration-cache</tt> <i>cal_file</i>
53 * <tt>--glob-calibration-data</tt>
54 * <tt>--channel-name</tt> <i>channel</i>
55 * [<tt>--calibrated-data</tt> <i>cal_type</i>]
56 * [<tt>--geo-high-pass-freq</tt> <i>geo_freq</i>]
57 * [<tt>--geo-high-pass-order</tt> <i>geo_order</i>]
58 * [<tt>--geo-high-pass-atten</tt> <i>geo_atten</i>]
59 * <tt>--sample-rate</tt> <i>sample_freq</i>
60 * <tt>--resample-filter</tt> <i>filter_type</i>
61 * [<tt>--disable-high-pass</tt>]
62 * [<tt>--enable-high-pass</tt> <i>high_freq</i>]
63 * [<tt>--high-pass-order</tt> <i>high_order</i>]
64 * [<tt>--high-pass-attenuation</tt> <i>high_atten</i>]
65 * <tt>--spectrum-type</tt> <i>spectype</i>
66 * [<tt>--dynamic-range-exponent</tt> <i>exp</i>]
67 * <tt>--segment-length</tt> <i>seglen</i>
68 * [<tt>--number-of-segments</tt> <i>segnum</i>]
69 * [<tt>--standard-candle</tt>]
70 * [<tt>--candle-snr</tt> <i>candle_snr</i>]
71 * [<tt>--candle-mass1</tt> <i>candle_mass1</i>]
72 * [<tt>--candle-mass2</tt> <i>candle_mass2</i>]
73 * <tt>--low-frequency-cutoff</tt> <i>cutlow</i>
74 * <tt>--high-frequency-cutoff</tt> <i>cuthigh</i>
75 * [<tt>--minimum-mass</tt> <i>minmass</i>]
76 * [<tt>--maximum-mass</tt> <i>maxmass</i>]
77 * [<tt>--minimum-psi0</tt> <i>psi0min</i>]
78 * [<tt>--maximum-psi0</tt> <i>psi0max</i>]
79 * [<tt>--minimum-psi3</tt> <i>psi3min</i>]
80 * [<tt>--maximum-psi3</tt> <i>psi3max</i>]
81 * [<tt>--maximum-fcut-tmplts</tt> <i>maxTemp</i>]
82 * [<tt>--alpha</tt> <i>alpha</i>]
83 * <tt>--minimal-match</tt> <i>match</i>
84 * <tt>--order</tt> <i>order</i>
85 * <tt>--approximant</tt> <i>approx</i>
86 * <tt>--space</tt> <i>space</i>
87 * [<tt>--write-raw-data</tt>]
88 * [<tt>--write-response</tt>]
89 * [<tt>--write-spectrum</tt>]
90 * [<tt>--write-strain-spectrum</tt>]</dd>
91 *
92 * <dt>Options</dt><dd>
93 * The following command line arguments are available when running tmpltbank.c
94 * \</dd>
95 *
96 * <dt><tt>--alpha</tt> <i> alpha</i></dt><dd>
97 * Set BCV amplitude correction to <i>alpha</i>.</dd>
98 *
99 * <dt><tt>--approximant</tt> <i> approx</i></dt><dd>
100 * Sets the approximant of the waveform to <i>approx</i>. <tt>TaylorT2</tt> is the standard stationary phase frequency domain chirp used in the BNS search. Available parameters: <tt>TaylorT1</tt>, <tt>TaylorT2</tt>, <tt>TaylorT3</tt>, <tt>TaylorF1</tt>, <tt>TaylorF2</tt>, <tt>PadeT1</tt>, <tt>PadeF1</tt>, <tt>EOB</tt>, <tt>BCV</tt>, <tt>SpinTaylorT3</tt>, <tt>BCVSpin</tt>.</dd>
101 *
102 * <dt><tt>--calibrated-data</tt> <i>type</i></dt><dd>
103 * Calibrated data of <i>type</i> <tt>real_4</tt> or <tt>real_8</tt>.</dd>
104 *
105 * <dt><tt>--calibration-cache</tt> <i>cal_file</i></dt><dd>
106 * Obtain calibration from LAL frame cache <i>cal_file</i>.</dd>
107 *
108 * <dt><tt>--candle-mass1</tt> <i>candle_mass1</i></dt><dd>
109 * Mass <i>candle_mass1</i> of first component in candle binary. Must be specified is the option <tt>--standard-candle</tt> is set.</dd>
110 *
111 * <dt><tt>--candle-mass2</tt> <i>candle_mass2</i></dt><dd>
112 * Mass <i>candle_mass2</i> of second component in candle binary. Must be specified is the option <tt>--standard-candle</tt> is set.</dd>
113 *
114 * <dt><tt>--candle-snr</tt> <i>candle_snr</i></dt><dd>
115 * Set the signal-to-noise ratio of standard candle to <i>candle_snr</i>. Must be specified is the option <tt>--standard-candle</tt> is set.</dd>
116 *
117 * <dt><tt>--channel-name</tt> <i>channel</i></dt><dd>
118 * Read data from interferometer channel <i>channel</i>.</dd>
119 *
120 * <dt><tt>--comment</tt> <i>comment</i></dt><dd>
121 * Set the process table comment to <i>comment</i>.</dd>
122 *
123 * <dt><tt>--disable-high-pass</tt></dt><dd>
124 * Turn off the IIR highpass filter. This is an optimistc option. Someday the data will be so good we won't need high pass filtering! </dd>
125 *
126 * <dt><tt>--dynamic-range-exponent</tt> <i>exp</i></dt><dd>
127 * Set dynamic range scaling to \f${2}^{exp}\f$.</dd>
128 *
129 * <dt><tt>--enable-high-pass</tt> <i>high_freq</i></dt><dd>
130 * High pass data above <i>high_freq</i> Hz using an IIR filter.</dd>
131 *
132 * <dt><tt>--frame-cache</tt> <i>cache_file</i></dt><dd>
133 * This option is used instead of <tt>--glob-frame-data</tt> to read frame data from a frame cache file <i>cache_file</i>. </dd>
134 *
135 * <dt><tt>--frame-type</tt> <i>type</i></dt><dd>
136 * This option specified the type of frames containing the input data. This option must be specified with the <tt>--glob-frame-data</tt> option.???????????</dd>
137 *
138 * <dt><tt>--geo-high-pass-atten</tt> <i>geo_atten</i></dt><dd>
139 * Set the attenuation of the high pass filter to <i>geo_atten</i>. Only if <tt>--calibrated-data</tt> is set to <tt>real_8</tt>.</dd>
140 *
141 * <dt><tt>--geo-high-pass-freq</tt> <i>geo_freq</i></dt><dd>
142 * This sets the high pass filter frequency for GEO data above <i>geo_freq</i> Hz using an IIR filter. Only if <tt>--calibrated-data</tt> is set to <tt>real_8</tt>.</dd>
143 *
144 * <dt><tt>--geo-high-pass-order</tt> <i>geo_order</i></dt><dd>
145 * Set the order of the GEO high pass filter to <i>geo_order</i>. Only if <tt>--calibrated-data</tt> is set to <tt>real_8</tt>.</dd>
146 *
147 * <dt><tt>--glob-calibration-data</tt></dt><dd>
148 * Is this option is specified, the calibration is obtained by globbing in the working directory.?????????</dd>
149 *
150 * <dt><tt>--glob-frame-data</tt></dt><dd>
151 * This option along with <tt>--frame-type</tt>
152 * can be used instead of <tt>--frame-cache</tt> to read data stored locally in
153 * the working directory. It finds files of the specified frame type with a *.gwf
154 * extension. </dd>
155 *
156 * <dt><tt>--gps-end-time</tt> <i>gps_end</i></dt><dd>
157 * Set the integer part of the GPS time <i>gps_end</i> you want
158 * to stop reading data. </dd>
159 *
160 * <dt><tt>--gps-start-time</tt> <i>gps_start</i></dt><dd>
161 * Set the integer part of the GPS time <i>gps_start</i> from which you wish to begin reading data.</dd>
162 *
163 * <dt><tt>--help</tt></dt><dd> display the help message which gives brief explanations
164 * of the command arguments. </dd>
165 *
166 * <dt><tt>--high-frequency-cutoff</tt> <i>cuthigh</i></dt><dd>
167 * Do not filter above <i>cuthigh</i> Hz.</dd>
168 *
169 * <dt><tt>--high-pass-attenuation</tt> <i>high_atten</i></dt><dd>
170 * Set the attenuation of the high pass filter to <i>high_atten</i>.</dd>
171 *
172 * <dt><tt>--high-pass-order</tt> <i>high_order</i></dt><dd>
173 * Set the order of the high pass filter to <i>high_order</i>.</dd>
174 *
175 * <dt><tt>--low-frequency-cutoff</tt> <i>cutlow</i></dt><dd>
176 * Do not filter below <i>cutlow</i> Hz.</dd>
177 *
178 * <dt><tt>-maximum-fcut-tmplts</tt> <i> maxTemp</i></dt><dd>
179 * Set the maximum number of templates in fcut direction to <i>maxTemp</i>.</dd>
180 *
181 * <dt><tt>--maximum-mass</tt> <i>maxmass</i></dt><dd>
182 * Set maximum component mass of bank to <i>maxmass</i>.</dd>
183 *
184 * <dt><tt>--maximum-psi0</tt> <i>psi0max</i></dt><dd>
185 * Set maximum range of BCV parameter psi0 to <i> psi0max</i>.</dd>
186 *
187 * <dt><tt>--maximum-psi3</tt> <i>psi3max</i></dt><dd>
188 * Set maximum range of BCV parameter psi3 to <i> psi3max</i>.</dd>
189 *
190 * <dt><tt>--minimal-match</tt> <i>match</i></dt><dd>
191 * Specifies the minimal match <i>match</i> between templates in the
192 * bank and all possible signals in the parameter space.</dd>
193 *
194 * <dt><tt>--minimum-mass</tt> <i>minmass</i></dt><dd>
195 * Set minimum component mass of bank to <i>minmass</i>.</dd>
196 *
197 * <dt><tt>--minimum-psi0</tt> <i>psi0min</i></dt><dd>
198 * Set minimum range of BCV parameter psi0 to <i> psi0min</i>.</dd>
199 *
200 * <dt><tt>--minimum-psi3</tt> <i>psi3min</i></dt><dd>
201 * Set minimum range of BCV parameter psi3 to <i> psi3min</i>.</dd>
202 *
203 * <dt><tt>--number-of-segments</tt> <i>segnum</i></dt><dd>
204 * Set number of data segments to <i>segnum</i>.</dd>
205 *
206 * <dt><tt>--order</tt> <i>order</i></dt><dd>
207 * This sets the order of the waveform to <i>order</i>. Usually it is set to <tt>twoPN</tt> (second order post newtonian). Available parameters: <tt>newtonian</tt>, <tt>oneHalfPN</tt>, <tt>onePN</tt>, <tt>onePointFivePN</tt>, <tt>twoPN</tt>, <tt>twoPointFivePN</tt>, <tt>threePN</tt>, <tt>threePointFivePN</tt>.</dd>
208 *
209 * <dt><tt>--pad-data</tt> <i>time_pad</i></dt><dd>
210 * This flag specifies an amount of time <i>time_pad</i> to add to
211 * the beginning and end of the input time series data. Padding the data is
212 * necessary because resampling and filtering corrupts these portions.
213 * 8 seconds is the accepted choice for this paramenter. See LAL documentation
214 * for a description of resampling and high pass filtering. </dd>
215 *
216 * <dt><tt>--resample-filter</tt> <i>filter_type</i></dt><dd>
217 * Set resample filter <i>filter_type</i> to <tt>ldas</tt> or <tt>butterworth</tt>. In the normal case the <i>ldas</i> filter is used.</dd>
218 *
219 * <dt><tt>--sample-rate</tt> <i>sample_freq</i></dt><dd>
220 * Specifies the sampling frequency <i>sample_freq</i> at which you
221 * want to filter the data downsampling if necessary.</dd>
222 *
223 * <dt><tt>--segment-length</tt> <i>seglen</i></dt><dd>
224 * Set data segment length to <i>seglen</i> points.</dd>
225 *
226 * <dt><tt>--space</tt> <i> space</i></dt><dd>
227 * In order to make the template bank coordinates nice and friendly these parameters are used instead of masses.
228 * Usually \c Tau0Tau3 is used. Available parameters:
229 * <tt>Tau0Tau2</tt>, <tt>Tau0Tau3</tt>, <tt>Psi0Psi3</tt>.</dd>
230 *
231 * <dt><tt>--spectrum-type</tt> <i>spec_type</i></dt><dd>
232 * Use PSD estimator <i>spec_type</i> <tt>mean</tt> or <tt>median</tt> to choose how the average is calculated. Since the median average is less affected by a loud glitch <tt>median</tt> is used generally.</dd>
233 *
234 * <dt><tt>--standard-candle</tt></dt><dd>
235 * Compute a standard candle from the PSD. In that case the arguments <tt>candle-mass1</tt>, <tt>candle-mass2</tt> and <tt>candle-snr</tt> must also be specified.</dd>
236 *
237 * <dt><tt>--verbose</tt></dt><dd> print progress information as the code executes.</dd>
238 *
239 * <dt><tt>--version</tt></dt><dd> print version information and exit without running
240 * the tmpltbank code. </dd>
241 *
242 * <dt><tt>--user-tag</tt> <i>usertag</i></dt><dd>
243 * Set the user tag to the string <i>usertag</i>.
244 * This string must not contain spaces or dashes ("-"). This string will appear
245 * in the name of the file to which output information is written, and is recorded
246 * in the various XML tables within the file.</dd>
247 *
248 * <dt><tt>--write-raw-data</tt></dt><dd>
249 * Write raw data to a frame file.</dd>
250 *
251 * <dt><tt>--write-response</tt></dt><dd>
252 * Write the computed response function to a frame.</dd>
253 *
254 * <dt><tt>--write-spectrum</tt></dt><dd>
255 * Write the uncalibrated psd to a frame.</dd>
256 *
257 * <dt><tt>--write-strain-spectrum</tt></dt><dd>
258 * Write the calibrated strain psd to a text file.
259 *
260 * </dd>
261 *
262 * <dt>Description</dt><dd>
263 * \c lalapps_tmpltbank is a stand alone code for generating inspiral
264 * template banks for LIGO or GEO data with the LAL bank package. The code
265 * generates a calibrated power spectrum at the specified time for the
266 * requested channel and uses this to compute the template bank.
267 * The number of templates and the
268 * values of the bank parameters in the bank also depend on the minimal
269 * match, the
270 * minimum and maximum values of mass1 and mass2 (for the BNS search) or the
271 * minimum and maximum values of psi0, psi3, the bank-alpha and the number of
272 * fcut values (for the BCV search), which are all command-line arguments.
273 * Other necessary pieces of information are the approximant and its order and
274 * the space that the template bank will be laid on. The output of the code is
275 * an xml file and the bank is contained in a \c sngl_inspiral table. The code has
276 * also the capability of outputing the raw data, the response function and the
277 * calibrated and unclibrated power spectra to frame files.
278 * See the LAL bank package
279 * documentation for detailed information on the algorithms used to generate the
280 * template banks.</dd>
281 *
282 * <dt>Example</dt><dd>
283 *
284 * \code
285 * lalapps_tmpltbank \
286 * --gps-start-time 734357353 --gps-end-time 734358377 \
287 * --frame-cache cache/L-734357345-734361107.cache \
288 * --segment-length 1048576 --number-of-segments 7 \
289 * --pad-data 7 --sample-rate 4096 --resample-filter ldas \
290 * --enable-high-pass 5.000000e+01 --spectrum-type median
291 * --low-frequency-cutoff 7.000000e+01 --high-frequency-cutoff 2.048000e+03 \
292 * --minimum-mass 1.000000e+00 --maximum-mass 3.000000e+00 \
293 * --minimal-match 9.700000e-01 --calibration-cache \
294 * /ldas_outgoing/calibration/cache_files/L1-CAL-V03-729273600-734367600.cache \
295 * --space Tau0Tau3 --approximant TaylorT1 --order twoPN \
296 * --channel-name L1:LSC-AS_Q
297 *
298 * \endcode</dd>
299 *
300 * <dt>Author</dt><dd>
301 * Duncan Brown and Alexander Dietz</dd>
302 * </dl>
303 */
304
305#include <config.h>
306#include <stdio.h>
307#include <stdlib.h>
308#include <string.h>
309#include <sys/types.h>
310#include <sys/stat.h>
311#include <fcntl.h>
312#include <regex.h>
313#include <time.h>
314
315#include <series.h>
316#include <lalappsfrutils.h>
317
318#ifdef LALAPPS_CUDA_ENABLED
319#include <cuda_runtime_api.h>
320#endif
321
322#include <lal/LALConfig.h>
323#include <lal/LALgetopt.h>
324#include <lal/LALStdio.h>
325#include <lal/LALStdlib.h>
326#include <lal/LALError.h>
327#include <lal/LALDatatypes.h>
328#include <lal/AVFactories.h>
329#include <lal/LALConstants.h>
330#include <lal/PrintFTSeries.h>
331#include <lal/LALFrStream.h>
332#include <lal/Window.h>
333#include <lal/TimeFreqFFT.h>
334#include <lal/IIRFilter.h>
335#include <lal/ResampleTimeSeries.h>
336#include <lal/BandPassTimeSeries.h>
337#include <lal/LIGOLwXML.h>
338#include <lal/LIGOLwXMLRead.h>
339#include <lal/LIGOMetadataTables.h>
340#include <lal/LIGOMetadataUtils.h>
341#include <lal/LIGOMetadataInspiralUtils.h>
342#include <lal/Date.h>
343#include <lal/Units.h>
344#include <lal/LALInspiral.h>
345#include <lal/LALInspiralBank.h>
346#include <lal/LALSimNoise.h>
347
348#include <LALAppsVCSInfo.h>
349
350#include "inspiral.h"
351
352#define CVS_ID_STRING "$Id$"
353#define CVS_NAME_STRING "$Name$"
354#define CVS_REVISION "$Revision$"
355#define CVS_SOURCE "$Source$"
356#define CVS_DATE "$Date$"
357#define PROGRAM_NAME "tmpltbank"
358
359int arg_parse_check( int argc, char *argv[], ProcessParamsTable *procparams );
360
361/* type of data to analyze */
362enum
363{
366 real_8
368
369/*
370 * which type of PSD to use - 'simulated' refers to LALSimNoise PSD functions
371 *
372 */
373enum
374{
380
381/*
382 *
383 * variables that control program behaviour
384 *
385 */
386
387
388/* debugging */
389extern int vrbflg; /* verbosity of lal function */
390
391/* parameters used to generate calibrated power spectrum */
392LIGOTimeGPS gpsStartTime = { 0, 0 }; /* input data GPS start time */
393LIGOTimeGPS gpsEndTime = { 0, 0 }; /* input data GPS end time */
394INT4 padData = 0; /* saftety margin on input data */
395CHAR *fqChanName = NULL; /* name of data channel */
396INT4 globFrameData = 0; /* glob *.gwf to get frame data */
397CHAR *frInCacheName = NULL; /* cache file containing frames */
398CHAR *frInType = NULL; /* type of data frames */
399INT4 numPoints = -1; /* points in a segment */
400INT4 numSegments = -1; /* number of segments */
401CHAR ifo[3]; /* two character ifo code */
402CHAR *channelName = NULL; /* channel string */
403INT4 inputDataLength = 0; /* number of points in input */
404INT4 resampFiltType = -1; /* low pass filter used for res */
405INT4 sampleRate = -1; /* sample rate of filter data */
406INT4 highPass = -1; /* enable high pass on raw data */
407REAL4 highPassFreq = 0; /* high pass frequency */
408INT4 highPassOrder = -1; /* order of the td iir filter */
409REAL4 highPassAtten = -1; /* attenuation of the td filter */
410REAL4 fLow = -1; /* low frequency cutoff */
411CHAR *calCacheName = NULL; /* location of calibration data */
412INT4 globCalData = 0; /* glob for calibration frames */
413INT4 pointCal = 0; /* don't average cal over chunk */
414REAL4 dynRangeExponent = 0; /* exponent of dynamic range */
415REAL4 strainHighPassFreq = -1; /* h(t) high pass frequency */
416INT4 strainHighPassOrder = -1; /* h(t) high pass filter order */
417REAL4 strainHighPassAtten = -1; /* h(t) high pass attenuation */
418REAL8 (*specFunc)(REAL8) = NULL; /* pointer to simPSD functions */
419
420/* template bank generation parameters */
421REAL4 minMass = -1; /* minimum component mass */
422REAL4 maxMass = -1; /* maximum component mass */
423REAL4 minTotalMass = -1; /* minimum total mass */
424REAL4 maxTotalMass = -1; /* maximum total mass */
425REAL4 chirpMassCutoff = -1; /* maximum chirp mass to keep */
426REAL4 etaMinCutoff = -1; /* minimum eta to keep */
427REAL4 etaMaxCutoff = -1; /* maximum eta to keep */
428REAL4 psi0Min = 0; /* minimum value of psi0 */
429REAL4 psi0Max = 0; /* maximum value of psi0 */
430REAL4 psi3Min = 0; /* minimum value of psi3 */
431REAL4 psi3Max = 0; /* maximum value of psi3 */
432REAL4 alpha = 0; /* BCV amplitude correction */
433REAL4 betaMin = 0; /* minimum BCV spin parameter */
434REAL4 betaMax = 0; /* maximum BCV spin parameter */
435INT4 maxFcutTmplts = -1; /* num tmplts in fcut direction */
436REAL4 minMatch = -1; /* minimum requested match */
437REAL4 fUpper = -1; /* upper frequency cutoff */
438REAL4 chiMin = 0.0; /* minimum value of chi for PTF */
439REAL4 chiMax = 1.0; /* maximum value of chi for PTF */
440REAL4 kappaMin = -1.0; /* minimum value of kappa for PTF */
441REAL4 kappaMax = 1.0; /* maximum value of kappa for PTF */
442INT4 nPointsChi = 3; /* PTF template bank density */
443INT4 nPointsKappa = 5; /* PTF templated bank density */
444LALPNOrder order; /* post-Newtonian order */
445Approximant approximant; /* approximation method */
446CoordinateSpace space; /* coordinate space used */
447INT4 haveGridSpacing = 0; /* flag to indicate gridspacing */
449FreqCut maxFreqCut; /* Max. upper frequency cutoff */
450FreqCut minFreqCut; /* Min. upper frequency cutoff */
451INT4 numFreqCut = 0; /* # of upper freq. cuts to use */
452
453GridSpacing gridSpacing = SquareNotOriented; /* grid spacing (square or hexa)*/
454int polygonFit = 1; /* fit a polygon around BCV bank */
455
456/* standard candle parameters */
457INT4 computeCandle = 0; /* should we compute a candle? */
458REAL4 candleSnr = -1; /* candle signal to noise ratio */
459REAL4 candleMinMass = -1; /* standard candle mass (solar) */
460REAL4 candleMaxMass = 50; /* standard candle mass (solar) */
461
462/* TD follow up filenames */
465
466/* output parameters */
467CHAR *userTag = NULL;
468CHAR *ifoTag = NULL;
469int writeRawData = 0; /* write the raw data to a file */
470int writeResponse = 0; /* write response function used */
471int writeSpectrum = 0; /* write computed psd to file */
472int writeStrainSpec = 0; /* write computed strain spec */
474
475/* other command line args */
476CHAR comment[LIGOMETA_COMMENT_MAX]; /* process param comment */
477
478int main ( int argc, char *argv[] )
479{
480 /* lal function variables */
482
483 /* frame input data */
484 LALCache *frInCache = NULL;
485 LALCache *frGlobCache = NULL;
486 LALFrStream *frStream = NULL;
487 FrChanIn frChan;
488
489 /* frame output data */
490 struct FrFile *frOutFile = NULL;
491 struct FrameH *outFrame = NULL;
492
493 /* raw input data storage */
494 REAL4TimeSeries chan;
495 REAL8TimeSeries strainChan;
498
499 /* structures for preconditioning */
500 ResampleTSParams resampleParams;
501
502 /* templates */
504 SnglInspiralTable *tmplt = NULL;
505 INT4 numCoarse = 0;
506
507 /* output data */
508 SnglInspiralTable *templateBank = NULL;
509 ProcessTable *proctable;
510 ProcessParamsTable *procparams;
512 LIGOLwXMLStream *results;
513
514
515 /* counters and other variables */
516 UINT4 cut, i, j, k;
517 const LALUnit strainPerCount = {0,{0,0,0,0,0,1,-1},{0,0,0,0,0,0,0}};
518 CHAR fname[256];
519 REAL8 respRe, respIm;
520 REAL8 shf;
522 UINT4 numInputPoints;
523 const REAL8 epsilon = 1.0e-8;
524 UINT4 resampleChan = 0;
525 REAL8 tsLength;
526 REAL8 dynRange = 0;
527
528 /* TD follow-up variables */
529 SnglInspiralTable *tdFollowUp = NULL;
530
531 /*
532 *
533 * initialization
534 *
535 */
536
537
538 /* set up inital debugging values */
540
541 /* create the process and process params tables */
542 proctable = (ProcessTable *) calloc( 1, sizeof(ProcessTable) );
543 XLALGPSTimeNow(&(proctable->start_time));
546 procparams = (ProcessParamsTable *) calloc( 1, sizeof(ProcessParamsTable) );
547 memset( comment, 0, LIGOMETA_COMMENT_MAX * sizeof(CHAR) );
548
549 /* create the search summary table */
550 searchsumm = (SearchSummaryTable *) calloc( 1, sizeof(SearchSummaryTable) );
551
552 /* call the argument parse and check function */
553 arg_parse_check( argc, argv, procparams );
554
555 /* can use LALMalloc() / LALCalloc() from here */
556
557 /* fill the comment, if a user has specified on, or leave it blank */
558 if ( ! *comment )
559 {
560 snprintf( proctable->comment, LIGOMETA_COMMENT_MAX, " " );
561 snprintf( searchsumm->comment, LIGOMETA_COMMENT_MAX,
562 " " );
563 }
564 else
565 {
566 snprintf( proctable->comment, LIGOMETA_COMMENT_MAX,
567 "%s", comment );
568 snprintf( searchsumm->comment, LIGOMETA_COMMENT_MAX,
569 "%s", comment );
570 }
571
572 /* the number of nodes for a standalone job is always 1 */
573 searchsumm->nnodes = 1;
574
575
576 /* If we're doing a td follow-up we dont want to generate a bank
577 * if there was no BCV trigger in this time
578 */
579 if ( tdFileNames )
580 {
581 SnglInspiralTable **thisTdFollow = &tdFollowUp;
582
583 if ( vrbflg )
584 {
585 fprintf( stdout, "We are doing a TD follow-up\n" );
586 fprintf( stdout, "Following up %d files...\n", numTDFiles );
587 }
588
589 for (i = 0; i < (UINT4)numTDFiles; i++ )
590 {
592 if ( !*thisTdFollow )
593 {
594 fprintf( stderr, "Error reading file %s\n", tdFileNames[i] );
595 exit( 1 );
596 }
597
598 while ( (*thisTdFollow)->next )
599 {
600 thisTdFollow = &(*thisTdFollow)->next;
601 }
602 }
603
604 tdFollowUp = XLALIfoCutSingleInspiral( &tdFollowUp, ifo );
605 if ( tdFollowUp )
606 tdFollowUp = XLALTimeCutSingleInspiral( tdFollowUp, &gpsStartTime,
607 &gpsEndTime );
608
609 /* If there are no events to follow up, we just exit */
610 if ( !tdFollowUp ) goto cleanExit;
611
612 /* Free the follow-up events */
613 XLALDestroySnglInspiralTable( tdFollowUp );
614 }
615
616
617 if ( dynRangeExponent )
618 {
619 /* compute the dynamic range scaling for the psd computation */
620 dynRange = (REAL8) pow( 2.0, dynRangeExponent );
621 }
622 else
623 {
624 dynRange = 1.0;
625 }
626 if ( vrbflg )
627 fprintf( stdout, "using dynamic range scaling %e\n", dynRange );
628
629
630
632 {
633
634 /*
635 *
636 * read in the input data channel
637 *
638 */
639
640 /* set the time series parameters of the input data and resample params */
641 memset( &resampleParams, 0, sizeof(ResampleTSParams) );
642 resampleParams.deltaT = 1.0 / (REAL8) sampleRate;
643
644 /* set the params of the input data time series */
645 memset( &chan, 0, sizeof(REAL4TimeSeries) );
646 memset( &strainChan, 0, sizeof(REAL8TimeSeries) );
647 chan.epoch = gpsStartTime;
648 chan.epoch.gpsSeconds -= padData; /* subtract pad seconds from start */
649
650 /* copy the start time into the REAL8 h(t) time series */
651 strainChan.epoch = chan.epoch;
652
653 if ( globFrameData )
654 {
655 CHAR ifoRegExPattern[6];
656
657 if ( vrbflg ) fprintf( stdout, "globbing for *.gwf frame files from %c "
658 "of type %s in current directory\n", fqChanName[0], frInType );
659
660 frGlobCache = NULL;
661
662 /* create a frame cache by globbing all *.gwf files in the pwd */
663 frGlobCache = XLALCacheGlob(NULL, NULL);
664
665 /* check we globbed at least one frame file */
666 if ( ! frGlobCache->length )
667 {
668 fprintf( stderr, "error: no frame files of type %s found\n",
669 frInType );
670 exit( 1 );
671 }
672
673 /* sieve out the requested data type */
674 snprintf( ifoRegExPattern,
675 XLAL_NUM_ELEM(ifoRegExPattern), ".*%c.*",
676 fqChanName[0] );
677 frInCache = XLALCacheDuplicate(frGlobCache);
678 XLALCacheSieve(frInCache, 0, 0, ifoRegExPattern, frInType, NULL);
679
680 /* check we got at least one frame file back after the sieve */
681 if ( ! frInCache->length )
682 {
683 fprintf( stderr, "error: no frame files of type %s globbed as input\n",
684 frInType );
685 exit( 1 );
686 }
687
688 XLALDestroyCache( frGlobCache );
689 }
690 else
691 {
692 if ( vrbflg ) fprintf( stdout,
693 "reading frame file locations from cache file: %s\n", frInCacheName );
694
695 /* read a frame cache from the specified file */
696 frInCache = XLALCacheImport(frInCacheName);
697 }
698
699 /* open the input data frame stream from the frame cache */
700 LAL_CALL( LALFrCacheOpen( &status, &frStream, frInCache ), &status );
701
702 /* set the mode of the frame stream to fail on gaps or time errors */
704
705 /* enable frame-file checksum checking */
707
708 /* seek to required epoch and set chan name */
709 LAL_CALL( LALFrSeek( &status, &(chan.epoch), frStream ), &status );
710 frChan.name = fqChanName;
711
712 if ( calData == real_8 )
713 {
714 /* determine the sample rate of the raw data */
715 LAL_CALL( LALFrGetREAL8TimeSeries( &status, &strainChan, &frChan,
716 frStream ), &status );
717
718 /* copy the data parameters from the h(t) channel to input data channel */
719 snprintf( chan.name, LALNameLength, "%s", strainChan.name );
720 chan.epoch = strainChan.epoch;
721 chan.deltaT = strainChan.deltaT;
722 chan.f0 = strainChan.f0;
723 chan.sampleUnits = strainChan.sampleUnits;
724 }
725 else
726 {
727 /* determine the sample rate of the raw data and allocate enough memory */
728 LAL_CALL( LALFrGetREAL4TimeSeries( &status, &chan, &frChan, frStream ),
729 &status );
730 }
731
732 /* determine if we need to resample the channel */
733 if ( vrbflg )
734 {
735 fprintf( stdout, "resampleParams.deltaT = %e\n", resampleParams.deltaT );
736 fprintf( stdout, "chan.deltaT = %e\n", chan.deltaT );
737 }
738 if ( ! ( fabs( resampleParams.deltaT - chan.deltaT ) < epsilon ) )
739 {
740 resampleChan = 1;
741 if ( vrbflg )
742 fprintf( stdout, "input channel will be resampled\n" );
743
744 if ( resampFiltType == 0 )
745 {
746 resampleParams.filterType = LDASfirLP;
747 }
748 else if ( resampFiltType == 1 )
749 {
750 resampleParams.filterType = defaultButterworth;
751 }
752 }
753
754 /* determine the number of points to get and create storage for the data */
755 inputLengthNS = (REAL8) ( LAL_INT8_C(1000000000) *
757 chan.deltaT *= 1.0e9;
758 numInputPoints = (UINT4) floor( inputLengthNS / chan.deltaT + 0.5 );
759 if ( calData == real_8 )
760 {
761 /* create storage for the REAL8 h(t) input data */
762 LAL_CALL( LALDCreateVector( &status, &(strainChan.data), numInputPoints ),
763 &status );
764 }
765 LAL_CALL( LALSCreateVector( &status, &(chan.data), numInputPoints ),
766 &status );
767
768 if ( vrbflg ) fprintf( stdout, "input channel %s has sample interval "
769 "(deltaT) = %e\nreading %d points from frame stream\n", fqChanName,
770 chan.deltaT / 1.0e9, numInputPoints );
771
772 if ( calData == real_8 )
773 {
774 /* read in the REAL8 h(t) data here */
775 PassBandParamStruc strainHighpassParam;
776
777 /* read the REAL8 h(t) data from the time series into strainChan */
778 /* which already has the correct amount of memory allocated */
779 if ( vrbflg ) fprintf( stdout, "reading REAL8 h(t) data from frames... " );
780
781 LAL_CALL( LALFrGetREAL8TimeSeries( &status, &strainChan, &frChan,
782 frStream ), &status);
783
784 if ( vrbflg ) fprintf( stdout, "done\n" );
785
786 /* high pass the h(t) data using the parameters specified on the cmd line*/
787 strainHighpassParam.nMax = strainHighPassOrder;
788 strainHighpassParam.f1 = -1.0;
789 strainHighpassParam.f2 = (REAL8) strainHighPassFreq;
790 strainHighpassParam.a1 = -1.0;
791 strainHighpassParam.a2 = (REAL8)(1.0 - strainHighPassAtten);
792 if ( vrbflg ) fprintf( stdout,
793 "applying %d order high pass to REAL8 h(t) data: "
794 "%3.2f of signal passes at %4.2f Hz\n",
795 strainHighpassParam.nMax, strainHighpassParam.a2,
796 strainHighpassParam.f2 );
797
799 &strainHighpassParam ), &status );
800
801 /* cast the REAL8 h(t) data to REAL4 in the chan time series */
802 /* which already has the correct amount of memory allocated */
803 for ( j = 0 ; j < numInputPoints ; ++j )
804 {
805 chan.data->data[j] = (REAL4) ( strainChan.data->data[j] * dynRange );
806 }
807
808 /* re-copy the data parameters from h(t) channel to input data channel */
809 snprintf( chan.name, LALNameLength, "%s", strainChan.name );
810 chan.epoch = strainChan.epoch;
811 chan.deltaT = strainChan.deltaT;
812 chan.f0 = strainChan.f0;
813 chan.sampleUnits = strainChan.sampleUnits;
814
815 /* free the REAL8 h(t) input data */
816 LAL_CALL( LALDDestroyVector( &status, &(strainChan.data) ), &status );
817 strainChan.data = NULL;
818 }
819 else
820 {
821 /* read the data channel time series from frames */
822 LAL_CALL( LALFrGetREAL4TimeSeries( &status, &chan, &frChan, frStream ),
823 &status );
824
825 if ( calData == real_4 )
826 {
827 /* multiply the input data by dynRange */
828 for ( j = 0 ; j < numInputPoints ; ++j )
829 {
830 chan.data->data[j] *= dynRange;
831 }
832 }
833 }
834 memcpy( &(chan.sampleUnits), &lalADCCountUnit, sizeof(LALUnit) );
835
836 /* store the start and end time of the raw channel in the search summary */
837 /* FIXME: loss of precision; consider
838 searchsumm->in_start_time = searchsumm->in_end_time = chan.epoch;
839 XLALGPSAdd(&searchsumm->in_end_time, chan.deltaT * (REAL8) chan.data->length);
840 */
841 searchsumm->in_start_time = chan.epoch;
842 tsLength = XLALGPSGetREAL8(&(chan.epoch) );
843 tsLength += chan.deltaT * (REAL8) chan.data->length;
844 XLALGPSSetREAL8( &(searchsumm->in_end_time), tsLength );
845
846 /* close the frame file stream and destroy the cache */
847 LAL_CALL( LALFrClose( &status, &frStream ), &status );
848 XLALDestroyCache( frInCache );
849
850 /* write the raw channel data as read in from the frame files */
851 if ( writeRawData ) outFrame = fr_add_proc_REAL4TimeSeries( outFrame,
852 &chan, "ct", "RAW" );
853
854 if ( vrbflg ) fprintf( stdout, "read channel %s from frame stream\n"
855 "got %d points with deltaT %e\nstarting at GPS time %d sec %d ns\n",
856 chan.name, chan.data->length, chan.deltaT,
858
859 /* resample the input data */
860 if ( resampleChan )
861 {
862 if (vrbflg) fprintf( stdout, "resampling input data from %e to %e\n",
863 chan.deltaT, resampleParams.deltaT );
864
865 LAL_CALL( LALResampleREAL4TimeSeries( &status, &chan, &resampleParams ),
866 &status );
867
868 if ( vrbflg ) fprintf( stdout, "channel %s resampled:\n"
869 "%d points with deltaT %e\nstarting at GPS time %d sec %d ns\n",
870 chan.name, chan.data->length, chan.deltaT,
872
873 /* write the resampled channel data as read in from the frame files */
874 if ( writeRawData ) outFrame = fr_add_proc_REAL4TimeSeries( outFrame,
875 &chan, "ct", "RAW_RESAMP" );
876 }
877 }
878
879 /*
880 *
881 * compute a calibrated strain spectrum
882 *
883 */
884
885 /* create storage for the response and spectrum */
886 memset( &spec, 0, sizeof(REAL4FrequencySeries) );
887 LAL_CALL( LALSCreateVector( &status, &(spec.data), numPoints / 2 + 1 ),
888 &status );
889 memset( &resp, 0, sizeof(COMPLEX8FrequencySeries) );
890 LAL_CALL( LALCCreateVector( &status, &(resp.data), numPoints / 2 + 1 ),
891 &status );
892 resp.epoch = spec.epoch = gpsStartTime;
893
895 {
896 /* iir filter to remove low frequencies from data channel */
897 if ( highPass )
898 {
899 PassBandParamStruc highpassParam;
900 highpassParam.nMax = highPassOrder;
901 highpassParam.f1 = -1.0;
902 highpassParam.f2 = (REAL8) highPassFreq;
903 highpassParam.a1 = -1.0;
904 highpassParam.a2 = (REAL8)(1.0 - highPassAtten); /* a2 is not attenuation */
905
906 if ( vrbflg ) fprintf( stdout, "applying %d order high pass: "
907 "%3.2f of signal passes at %4.2f Hz\n",
908 highpassParam.nMax, highpassParam.a2, highpassParam.f2 );
909
910 LAL_CALL( LALDButterworthREAL4TimeSeries( &status, &chan, &highpassParam ),
911 &status );
912 }
913
914 /* remove pad from requested data from start and end of time series */
915 memmove( chan.data->data, chan.data->data + padData * sampleRate,
916 (chan.data->length - 2 * padData * sampleRate) * sizeof(REAL4) );
917 XLALRealloc( chan.data->data,
918 (chan.data->length - 2 * padData * sampleRate) * sizeof(REAL4) );
919 chan.data->length -= 2 * padData * sampleRate;
920 chan.epoch.gpsSeconds += padData;
921
922 if ( vrbflg ) fprintf( stdout, "after removal of %d second padding at "
923 "start and end:\ndata channel sample interval (deltaT) = %e\n"
924 "data channel length = %d\nstarting at %d sec %d ns\n",
925 padData , chan.deltaT , chan.data->length,
927
928 /* store the start and end time of the filter channel in the search summ */
929 /* FIXME: loss of precision; consider
930 searchsumm->out_start_time = chan.epoch;
931 XLALGPSAdd(&searchsumm->out_start_time, chan.deltaT * (REAL8) chan.data->length);
932 */
933 searchsumm->out_start_time = chan.epoch;
934 tsLength = XLALGPSGetREAL8( &(chan.epoch) );
935 tsLength += chan.deltaT * (REAL8) chan.data->length;
936 XLALGPSSetREAL8( &(searchsumm->out_end_time), tsLength );
937
938 /* compute the windowed power spectrum for the data channel */
939 { /* Kipp: ported this block from LALStatus 20190806 */
941 REAL4FFTPlan *plan = XLALCreateForwardREAL4FFTPlan( numPoints, 0 );
942 switch ( specType )
943 {
944 case specType_mean:
945 if ( vrbflg ) fprintf( stdout, "computing mean psd with overlap %d\n", numPoints / 2 );
946 XLALREAL4AverageSpectrumWelch( &spec, &chan, numPoints, numPoints / 2, window, plan );
947 break;
948 case specType_median:
949 if ( vrbflg ) fprintf( stdout, "computing median psd with overlap %d\n", numPoints / 2 );
950 XLALREAL4AverageSpectrumMedian( &spec, &chan, numPoints, numPoints / 2, window, plan );
951 break;
953 fprintf( stderr, "This should never happen! Exiting..." );
954 exit( 1 );
955 default:
956 fprintf( stderr, "unknown spectrum type %d\n", specType );
957 exit( 1 );
958 }
959
960 XLALDestroyREAL4Window( window );
962 }
963 LAL_CALL( LALSDestroyVector( &status, &(chan.data) ), &status );
964 }
965
967 /* initialize spectrum frequency bins */
968 {
969 spec.f0 = 0.0;
970 /* frequency bin (Hz) is 1/duration of segment */
971 /* = sample rate/number of points */
972 spec.deltaF = ( (REAL8) sampleRate ) / ( (REAL8) numPoints );
973
974 /* evaluate the simulated PSD */
975 for ( k = 0; k < spec.data->length; ++k )
976 {
977 REAL8 sim_psd_freq = (REAL8) k * spec.deltaF;
978
979 /* PSD can be a very small number, rescale it before casting to REAL4 */
980 spec.data->data[k] = (REAL4) (dynRange * dynRange * specFunc( sim_psd_freq ));
981 }
982 }
983
984 /* write the spectrum data to a file */
985 if ( writeSpectrum )
986 {
987 strcpy( spec.name, chan.name );
988 outFrame = fr_add_proc_REAL4FrequencySeries( outFrame,
989 &spec, "ct/sqrtHz", "PSD" );
990 }
991
992 /* set the parameters of the response to match the data and spectrum */
993 resp.deltaF = spec.deltaF;
994 resp.f0 = spec.f0;
996
998 {
999 /* if we are using calibrated data or a design PSD set the response to unity */
1000 /* i.e. 1 over the dynamic range scaling factor */
1001 if ( vrbflg ) fprintf( stdout, "generating unity response function\n" );
1002 for( k = 0; k < resp.data->length; ++k )
1003 {
1004 resp.data->data[k] = (REAL4) (1.0 / dynRange);
1005 }
1006 }
1007 else
1008 {
1009 fprintf( stderr, "uncalibrated data no longer supported" );
1010 exit( 1 );
1011 }
1012
1013 /* write the calibration data to a file */
1014 if ( writeResponse )
1015 {
1016 strcpy( resp.name, chan.name );
1017 outFrame = fr_add_proc_COMPLEX8FrequencySeries( outFrame,
1018 &resp, "strain/ct", "RESPONSE" );
1019 }
1020
1021 /* set low frequency cutoff of power spectrum */
1022 cut = fLow / spec.deltaF > 1 ? fLow / spec.deltaF : 1;
1023
1024 /* compute a calibrated strain power spectrum */
1025 bankIn.shf.epoch = spec.epoch;
1026 memcpy( bankIn.shf.name, spec.name, LALNameLength * sizeof(CHAR) );
1027 bankIn.shf.deltaF = spec.deltaF;
1028 bankIn.shf.f0 = spec.f0;
1029 bankIn.shf.data = NULL;
1030 if (XLALUnitMultiply( &(bankIn.shf.sampleUnits), &(spec.sampleUnits), &(resp.sampleUnits) ) == NULL) {
1031 return LAL_EXLAL;
1032 }
1033 LAL_CALL( LALDCreateVector( &status, &(bankIn.shf.data), spec.data->length ),
1034 &status );
1035 memset( bankIn.shf.data->data, 0,
1036 bankIn.shf.data->length * sizeof(REAL8) );
1037
1038 shf = spec.data->data[cut] *
1039 ( crealf(resp.data->data[cut]) * crealf(resp.data->data[cut]) +
1040 cimagf(resp.data->data[cut]) * cimagf(resp.data->data[cut]) );
1041 for ( k = 1; k < cut ; ++k )
1042 {
1043 bankIn.shf.data->data[k] = shf;
1044 }
1045 for ( k = cut; k < bankIn.shf.data->length; ++k )
1046 {
1047 respRe = (REAL8) crealf(resp.data->data[k]);
1048 respIm = (REAL8) cimagf(resp.data->data[k]);
1049 bankIn.shf.data->data[k] = (REAL8) spec.data->data[k] *
1050 ( respRe * respRe + respIm * respIm );
1051 }
1052
1053 /* write the scaled strain spectrum data to a file */
1054 if ( writeStrainSpec )
1055 {
1056#if 0
1057 strcpy( spec.name, chan.name );
1058 outFrame = fr_add_proc_REAL8FrequencySeries( outFrame,
1059 &(bankIn.shf), "strain/sqrtHz", "STRAIN_PSD" );
1060#endif
1061 snprintf( fname, sizeof(fname), "%s-TMPLTBANK-%d-%d.strainspec.txt",
1064 LALDPrintFrequencySeries( &(bankIn.shf), fname );
1065 }
1066
1067
1068 /*
1069 *
1070 * compute the standard candle distance
1071 *
1072 */
1073
1074 if ( computeCandle )
1075 {
1076 while ( candleMinMass <= candleMaxMass )
1077 {
1078 if ( vrbflg ) fprintf( stdout, "maximum distance for (%3.2f,%3.2f) "
1079 "at signal-to-noise %3.2f = ", candleMinMass, candleMinMass, candleSnr );
1080
1082 }
1083 }
1084
1085
1086 /*
1087 *
1088 * compute the template bank
1089 *
1090 */
1091
1092
1093 /* bank generation parameters */
1094 bankIn.mMin = (REAL8) minMass;
1095 bankIn.mMax = (REAL8) maxMass;
1096 if (maxTotalMass > 0)
1097 {
1098 bankIn.MMax = (REAL8) maxTotalMass;
1099 bankIn.mMax = bankIn.MMax - bankIn.mMin;
1100 if (minTotalMass > 0)
1101 {
1103 bankIn.MMin = (REAL8) minTotalMass;
1104 }
1105 else
1106 {
1108 }
1109 }
1110 else
1111 {
1113 bankIn.MMax = bankIn.mMax * 2.0;
1114 }
1115 bankIn.psi0Min = (REAL8) psi0Min;
1116 bankIn.psi0Max = (REAL8) psi0Max;
1117 bankIn.psi3Min = (REAL8) psi3Min;
1118 bankIn.psi3Max = (REAL8) psi3Max;
1120 bankIn.alpha = (REAL8) alpha;
1121 bankIn.betaMin = (REAL8) betaMin;
1122 bankIn.betaMax = (REAL8) betaMax;
1123 bankIn.chiMin = (REAL8) chiMin;
1124 bankIn.chiMax = (REAL8) chiMax;
1125 bankIn.kappaMin = (REAL8) kappaMin;
1126 bankIn.kappaMax = (REAL8) kappaMax;
1127 bankIn.nPointsChi = nPointsChi;
1128 bankIn.nPointsKappa = nPointsKappa;
1129 bankIn.mmCoarse = (REAL8) minMatch;
1130 bankIn.mmFine = 0.99; /* doesn't matter since no fine bank yet */
1131 bankIn.fLower = (REAL8) fLow;
1132 bankIn.fUpper = (REAL8) fUpper;
1133 bankIn.iflso = 0; /* currently not implemented */
1134 bankIn.tSampling = (REAL8) sampleRate;
1135 bankIn.order = order;
1136 bankIn.approximant = approximant;
1137 bankIn.gridSpacing = gridSpacing;
1138 bankIn.space = space;
1139 bankIn.insidePolygon = polygonFit; /*by default it uses a polygon fitting. */
1140 bankIn.etamin = bankIn.mMin * ( bankIn.MMax - bankIn.mMin) /
1141 ( bankIn.MMax * bankIn.MMax );
1142 bankIn.LowGM = -4.;
1143 bankIn.HighGM = 6.;
1144 bankIn.computeMoments = computeMoments; /* by default, gammas/moments are recomputed */
1145 bankIn.maxFreqCut = maxFreqCut;
1146 bankIn.minFreqCut = minFreqCut;
1147 bankIn.numFreqCut = numFreqCut;
1148
1149 /* generate the template bank */
1150 if ( vrbflg )
1151 {
1152 fprintf( stdout, "generating template bank parameters... " );
1153 fflush( stdout );
1154 }
1155 LAL_CALL( LALInspiralBankGeneration( &status, &bankIn, &tmplt, &numCoarse),
1156 &status );
1157
1158 /* Do chirp mass cut */
1159 if ( chirpMassCutoff > 0 )
1160 {
1161 tmplt = XLALMassCut( tmplt, "mchirp", 0, chirpMassCutoff, -1, -1 );
1162 /* count the remaining tmplts */
1163 numCoarse = XLALCountSnglInspiral( tmplt );
1164 }
1165
1166 /* Do eta cut */
1167 if ( etaMinCutoff >= 0 || etaMaxCutoff > 0 )
1168 {
1169 tmplt = XLALMassCut( tmplt, "eta", etaMinCutoff, etaMaxCutoff, -1, -1 );
1170 /* count the remaining tmplts */
1171 numCoarse = XLALCountSnglInspiral( tmplt );
1172 }
1173
1174 if ( vrbflg )
1175 {
1176 fprintf( stdout, "done. Got %d templates\n", numCoarse );
1177 fflush( stdout );
1178 }
1179
1180 if ( numCoarse )
1181 {
1182 templateBank = tmplt;
1183 snprintf( tmplt->ifo, LIGOMETA_IFO_MAX, "%s", ifo );
1184 snprintf( tmplt->search, LIGOMETA_SEARCH_MAX, "tmpltbank" );
1185 snprintf( tmplt->channel, LIGOMETA_CHANNEL_MAX,
1186 "%s", channelName );
1187 while( (tmplt = tmplt->next) )
1188 {
1189 snprintf( tmplt->ifo, LIGOMETA_IFO_MAX, "%s", ifo );
1190 snprintf( tmplt->search, LIGOMETA_SEARCH_MAX, "tmpltbank" );
1191 snprintf( tmplt->channel, LIGOMETA_CHANNEL_MAX,
1192 "%s", channelName );
1193 }
1194 }
1195
1196 /* save the number of templates in the search summary table */
1197 searchsumm->nevents = numCoarse;
1198
1199
1200 /*
1201 *
1202 * free the data storage
1203 *
1204 */
1205
1206
1207 LAL_CALL( LALDDestroyVector( &status, &(bankIn.shf.data) ), &status );
1208 LAL_CALL( LALSDestroyVector( &status, &(spec.data) ), &status );
1209 LAL_CALL( LALCDestroyVector( &status, &(resp.data) ), &status );
1210
1211
1212 /*
1213 *
1214 * write the results to disk
1215 *
1216 */
1217
1218
1219 /* write the output frame */
1221 {
1222 snprintf( fname, sizeof(fname), "%s-TMPLTBANK-%d-%d.gwf",
1225 frOutFile = FrFileONew( fname, 0 );
1226 FrameWrite( outFrame, frOutFile );
1227 FrFileOEnd( frOutFile );
1228 }
1229
1230cleanExit:
1231 /* open the output xml file */
1232 if ( userTag && ifoTag && !outCompress)
1233 {
1234 snprintf( fname, sizeof(fname), "%s-TMPLTBANK_%s_%s-%d-%d.xml",
1237 }
1238 else if (userTag && !ifoTag && !outCompress)
1239 {
1240 snprintf( fname, sizeof(fname), "%s-TMPLTBANK_%s-%d-%d.xml",
1243 }
1244 else if (!userTag && ifoTag && !outCompress)
1245 {
1246 snprintf( fname, sizeof(fname), "%s-TMPLTBANK_%s-%d-%d.xml",
1249 }
1250 else if ( userTag && ifoTag && outCompress)
1251 {
1252 snprintf( fname, sizeof(fname), "%s-TMPLTBANK_%s_%s-%d-%d.xml.gz",
1255 }
1256 else if (userTag && !ifoTag && outCompress)
1257 {
1258 snprintf( fname, sizeof(fname), "%s-TMPLTBANK_%s-%d-%d.xml.gz",
1261 }
1262 else if (!userTag && ifoTag && outCompress)
1263 {
1264 snprintf( fname, sizeof(fname), "%s-TMPLTBANK_%s-%d-%d.xml.gz",
1267 }
1268 else if (!userTag && !ifoTag && outCompress)
1269 {
1270 snprintf( fname, sizeof(fname), "%s-TMPLTBANK-%d-%d.xml.gz",
1273 }
1274 else
1275 {
1276 snprintf( fname, sizeof(fname), "%s-TMPLTBANK-%d-%d.xml",
1279 }
1280 results = XLALOpenLIGOLwXMLFile( fname );
1281
1282 /* write the process table */
1283 snprintf( proctable->ifos, LIGOMETA_IFO_MAX, "%s", ifo );
1284 XLALGPSTimeNow(&(proctable->end_time));
1285 XLALWriteLIGOLwXMLProcessTable( results, proctable );
1286 free( proctable );
1287
1288 /* erase the first empty process params entry */
1289 {
1290 ProcessParamsTable *emptyPPtable = procparams;
1291 procparams = procparams->next;
1292 free( emptyPPtable );
1293 }
1294
1295 /* write the process params table */
1296 XLALWriteLIGOLwXMLProcessParamsTable( results, procparams );
1297 XLALDestroyProcessParamsTable( procparams );
1298
1299 /* write the search summary table */
1301
1302 /* write the template bank to the file */
1303 if ( templateBank )
1304 XLALWriteLIGOLwXMLSnglInspiralTable( results, templateBank );
1305 while ( templateBank )
1306 {
1307 tmplt = templateBank;
1308 templateBank = templateBank->next;
1309 LALFree( tmplt );
1310 }
1311
1312 /* close the output xml file */
1313 XLALCloseLIGOLwXMLFile( results );
1314
1315 /* free the rest of the memory, check for memory leaks and exit */
1316 if ( tdFileNames ) free( tdFileNames );
1317 if ( calCacheName ) free( calCacheName );
1318 if ( frInCacheName ) free( frInCacheName );
1319 if ( frInType ) free( frInType );
1320 if ( channelName ) free( channelName );
1321 if ( fqChanName ) free( fqChanName );
1323
1324#ifdef LALAPPS_CUDA_ENABLED
1325 cudaDeviceReset();
1326#endif
1327
1328 exit( 0 );
1329}
1330
1331/* ------------------------------------------------------------------------- */
1332
1333#define ADD_PROCESS_PARAM( pptype, format, ppvalue ) \
1334this_proc_param = this_proc_param->next = (ProcessParamsTable *) \
1335 calloc( 1, sizeof(ProcessParamsTable) );\
1336 snprintf( this_proc_param->program, LIGOMETA_PROGRAM_MAX, "%s", \
1337 PROGRAM_NAME );\
1338 snprintf( this_proc_param->param, LIGOMETA_PARAM_MAX, "--%s", \
1339 long_options[option_index].name );\
1340 snprintf( this_proc_param->type, LIGOMETA_TYPE_MAX, "%s", pptype );\
1341 snprintf( this_proc_param->value, LIGOMETA_VALUE_MAX, format, ppvalue );
1342
1343#define USAGE( a ) \
1344fprintf(a, " --help display this message\n");\
1345fprintf(a, " --verbose print progress information\n");\
1346fprintf(a, " --version print version information and exit\n");\
1347fprintf(a, " --user-tag STRING set the process_params usertag to STRING\n");\
1348fprintf(a, " --ifo-tag STRING set the ifotag to STRING - for file naming\n");\
1349fprintf(a, " --comment STRING set the process table comment to STRING\n");\
1350fprintf(a, " --write-compress write a compressed xml file\n");\
1351fprintf(a, "\n");\
1352fprintf(a, " NOTE: Data-related options not required when using a simulated PSD are labelled @\n");\
1353fprintf(a, "\n");\
1354fprintf(a, " --gps-start-time SEC GPS second of data start time\n");\
1355fprintf(a, " --gps-end-time SEC GPS second of data end time\n");\
1356fprintf(a, "@ --pad-data T pad the data start and end time by T seconds\n");\
1357fprintf(a, "\n");\
1358fprintf(a, "@ --glob-frame-data glob *.gwf files in the pwd to obtain frame data\n");\
1359fprintf(a, "@ --frame-type TAG input data is contained in frames of type TAG\n");\
1360fprintf(a, "@ --frame-cache obtain frame data from LAL frame cache FILE\n");\
1361fprintf(a, "@ --calibration-cache FILE obtain calibration from LAL frame cache FILE\n");\
1362fprintf(a, "@ --glob-calibration-data obtain calibration by globbing in working dir\n");\
1363fprintf(a, "\n");\
1364fprintf(a, "@ --channel-name CHAN read data from interferometer channel CHAN\n");\
1365fprintf(a, "@ --calibrated-data TYPE calibrated data of TYPE real_4 or real_8\n");\
1366fprintf(a, "@ --strain-high-pass-freq F high pass REAL8 h(t) data above F Hz\n");\
1367fprintf(a, "@ --strain-high-pass-order O set the order of the h(t) high pass filter to O\n");\
1368fprintf(a, "@ --strain-high-pass-atten A set the attenuation of the high pass filter to A\n");\
1369fprintf(a, "@ --point-calibration use the first point in the chunk to calibrate\n");\
1370fprintf(a, "\n");\
1371fprintf(a, " --sample-rate F filter data at F Hz, downsampling if necessary\n");\
1372fprintf(a, "@ --resample-filter TYPE set resample filter to TYPE [ldas|butterworth]\n");\
1373fprintf(a, "\n");\
1374fprintf(a, " --disable-high-pass turn off the IIR highpass filter\n");\
1375fprintf(a, "@ --enable-high-pass F high pass data above F Hz using an IIR filter\n");\
1376fprintf(a, "@ --high-pass-order O set the order of the high pass filter to O\n");\
1377fprintf(a, "@ --high-pass-attenuation A set the attenuation of the high pass filter to A\n");\
1378fprintf(a, " --spectrum-type TYPE use PSD estimator TYPE \n");\
1379fprintf(a, " (mean|median|iLIGOSRD|eLIGOModel|GEOModel|\n");\
1380fprintf(a, " |aLIGONoSRMLoP|aLIGONoSRMHiP|aLIGOZDLoP|aLIGOZDHiP|\n");\
1381fprintf(a, " |iVirgoModel|aVirgoModel|KAGRAModel)\n");\
1382fprintf(a, " --dynamic-range-exponent X set dynamic range scaling to 2^X (eg X=69.0)\n");\
1383fprintf(a, "\n");\
1384fprintf(a, " --segment-length N set data segment length to N points\n");\
1385fprintf(a, "@ --number-of-segments N set number of data segments to N\n");\
1386fprintf(a, "\n");\
1387fprintf(a, " --td-follow-up FILE follow up BCV events contained in FILE\n");\
1388fprintf(a, "\n");\
1389fprintf(a, " --standard-candle compute a standard candle from the PSD\n");\
1390fprintf(a, " --candle-snr SNR signal-to-noise ratio of standard candle\n");\
1391fprintf(a, " --candle-minmass M minimum component mass for (equal-mass) candle binary\n");\
1392fprintf(a, " --candle-maxmass M maximum component mass for candle binary, default=50\n");\
1393fprintf(a, "\n");\
1394fprintf(a, " --low-frequency-cutoff F do not filter below F Hz\n");\
1395fprintf(a, " --high-frequency-cutoff F upper frequency cutoff in Hz\n");\
1396fprintf(a, " --disable-compute-moments do not recompute the moments stored in the template bank. \n");\
1397fprintf(a, "\n");\
1398fprintf(a, " --minimum-mass MASS set minimum component mass of bank to MASS: required\n");\
1399fprintf(a, " --maximum-mass MASS set maximum component mass of bank to MASS\n");\
1400fprintf(a, " --max-total-mass MASS set maximum total mass of the bank to MASS. Will override --maximum-mass option\n");\
1401fprintf(a, " --min-total-mass MASS set minimum total mass of the bank to MASS: --max-total-mass must also be given\n");\
1402fprintf(a, " --chirp-mass-cutoff MASS set chirp mass cutoff to MASS\n");\
1403fprintf(a, " --max-eta ETA set maximum symmetric mass ratio of the bank to ETA\n");\
1404fprintf(a, " --min-eta ETA set minimum symmetric mass ratio of the bank to ETA\n");\
1405fprintf(a, "\n");\
1406fprintf(a, " --minimum-psi0 PSI0 set minimum range of BCV parameter psi0 to PSI0\n");\
1407fprintf(a, " --maximum-psi0 PSI0 set maximum range of BCV parameter psi0 to PSI0\n");\
1408fprintf(a, " --minimum-psi3 PSI3 set minimum range of BCV parameter psi3 to PSI3\n");\
1409fprintf(a, " --maximum-psi3 PSI3 set maximum range of BCV parameter psi3 to PSI3\n");\
1410fprintf(a, " --maximum-fcut-tmplts N maximum number of tmplts in fcut direction is N\n");\
1411fprintf(a, " --disable-polygon-fit disable the polygon fitting for BCV bank\n");\
1412fprintf(a, " --alpha ALPHA set alpha for the BCV bank generation\n");\
1413fprintf(a, " --minimum-beta BETA set minimum BCV spin parameter beta to BETA\n");\
1414fprintf(a, " --maximum-beta BETA set maximum BCV spin parameter beta to BETA\n");\
1415fprintf(a, "\n");\
1416fprintf(a, " --minimum-spin1 SPIN1_MIN set minimum value of chi for PTF to SPIN1_MIN (0.0)\n");\
1417fprintf(a, " --maximum-spin1 SPIN1_MAX set maximum value of chi for PTF to SPIN1_MAX (1.0)\n");\
1418fprintf(a, " --minimum-kappa1 KAPPA1_MIN set minimum value of kappa for PTF to KAPPA1_MIN (-1.0)\n");\
1419fprintf(a, " --maximum-kappa1 KAPPA1_MAX set maximum value of kappa for PTF to KAPPA1_MAX (1.0)\n");\
1420fprintf(a, " --npoints-chi N-CHI set number of points in the Chi direction for PTF template bank to N-CHI (3)\n");\
1421fprintf(a, " --npoints-kappa N-KAPPA set number of points in the Kappa direction for PTF template bank to N-KAPPA (5)\n");\
1422fprintf(a, "\n");\
1423fprintf(a, " --minimal-match M generate bank with minimal match M\n");\
1424fprintf(a, "\n");\
1425fprintf(a, " --order ORDER set post-Newtonian order of the waveform to ORDER\n");\
1426fprintf(a, " (newtonian|oneHalfPN|onePN|onePointFivePN|\n");\
1427fprintf(a, " twoPN|twoPointFive|threePN|threePointFivePN)\n");\
1428fprintf(a, " --approximant APPROX set approximant of the waveform to APPROX\n");\
1429fprintf(a, " (TaylorT1|TaylorT2|TaylorT3|TaylorF1|TaylorF2|\n");\
1430fprintf(a, " PadeT1|PadeT2|EOB|EOBNR|BCV|SpinTaylorT3|BCVSpin)\n");\
1431fprintf(a, " --num-freq-cutoffs Ncut create a template bank with Ncut different upper \n");\
1432fprintf(a, " frequency cutoffs (must be a positive integer) \n");\
1433fprintf(a, " --max-high-freq-cutoff MAX formula to compute the largest high freq. cutoff\n");\
1434fprintf(a, " possible choices in ascending order: (SchwarzISCO|BKLISCO|LightRing|FRD|ERD|LRD)\n");\
1435fprintf(a, " --min-high-freq-cutoff MIN formula to compute the smallest high freq. cutoff\n");\
1436fprintf(a, " possible choices in ascending order: (SchwarzISCO|BKLISCO|LightRing|FRD|ERD|LRD)\n");\
1437fprintf(a, " --space SPACE grid up template bank with mass parameters SPACE\n");\
1438fprintf(a, " (Tau0Tau2|Tau0Tau3|Psi0Psi3)\n");\
1439fprintf(a, " --grid-spacing GRIDSPACING grid up template bank with GRIDSPACING\n");\
1440fprintf(a, " (Hexagonal|SquareNotOriented)\n");\
1441fprintf(a, "\n");\
1442fprintf(a, " --write-response write the computed response function to a frame\n");\
1443fprintf(a, " --write-spectrum write the uncalibrated psd to a frame\n");\
1444fprintf(a, " --write-strain-spectrum write the calibrated strain psd to a text file\n");
1445
1446
1447int arg_parse_check( int argc, char *argv[], ProcessParamsTable *procparams )
1448{
1449 /* LALgetopt arguments */
1450 struct LALoption long_options[] =
1451 {
1452 /* these options set a flag */
1453 {"verbose", no_argument, &vrbflg, 1 },
1454 {"write-compress", no_argument, &outCompress, 1 },
1455 {"disable-high-pass", no_argument, &highPass, 0 },
1456 {"standard-candle", no_argument, &computeCandle, 1 },
1457 {"glob-frame-data", no_argument, &globFrameData, 1 },
1458 {"glob-calibration-data", no_argument, &globCalData, 1 },
1459 {"point-calibration", no_argument, &pointCal, 1 },
1460 /* parameters used to generate calibrated power spectrum */
1461 {"gps-start-time", required_argument, 0, 'a'},
1462 {"gps-end-time", required_argument, 0, 'b'},
1463 {"channel-name", required_argument, 0, 'c'},
1464 {"segment-length", required_argument, 0, 'd'},
1465 {"number-of-segments", required_argument, 0, 'e'},
1466 {"sample-rate", required_argument, 0, 'g'},
1467#ifdef LALAPPS_CUDA_ENABLED
1468 {"gpu-device-id", required_argument, 0, '+'},
1469#endif
1470 {"calibrated-data", required_argument, 0, 'M'},
1471 {"strain-high-pass-freq", required_argument, 0, 'J'},
1472 {"strain-high-pass-order", required_argument, 0, 'K'},
1473 {"strain-high-pass-atten", required_argument, 0, 'L'},
1474 {"help", no_argument, 0, 'h'},
1475 {"low-frequency-cutoff", required_argument, 0, 'i'},
1476 {"spectrum-type", required_argument, 0, 'j'},
1477 {"dynamic-range-exponent", required_argument, 0, 'f'},
1478 {"calibration-cache", required_argument, 0, 'p'},
1479 {"comment", required_argument, 0, 's'},
1480 {"enable-high-pass", required_argument, 0, 't'},
1481 {"high-pass-order", required_argument, 0, 'H'},
1482 {"high-pass-attenuation", required_argument, 0, 'I'},
1483 {"frame-cache", required_argument, 0, 'u'},
1484 {"frame-type", required_argument, 0, 'n'},
1485 {"pad-data", required_argument, 0, 'x'},
1486 {"user-tag", required_argument, 0, 'Z'},
1487 {"ifo-tag", required_argument, 0, 'Y'},
1488 {"version", no_argument, 0, 'V'},
1489 {"resample-filter", required_argument, 0, 'r'},
1490 /* template bank generation parameters */
1491 {"minimum-mass", required_argument, 0, 'A'},
1492 {"maximum-mass", required_argument, 0, 'B'},
1493 {"minimum-psi0", required_argument, 0, 'P'},
1494 {"maximum-psi0", required_argument, 0, 'Q'},
1495 {"minimum-psi3", required_argument, 0, 'R'},
1496 {"maximum-psi3", required_argument, 0, 'S'},
1497 {"maximum-fcut-tmplts", required_argument, 0, 'U'},
1498 {"minimum-beta", required_argument, 0, 'o'},
1499 {"maximum-beta", required_argument, 0, 'O'},
1500 {"alpha", required_argument, 0, 'T'},
1501 {"minimum-spin1", required_argument, 0, '4'},
1502 {"maximum-spin1", required_argument, 0, '5'},
1503 {"minimum-kappa1", required_argument, 0, '6'},
1504 {"maximum-kappa1", required_argument, 0, '7'},
1505 {"npoints-chi", required_argument, 0, '8'},
1506 {"npoints-kappa", required_argument, 0, '9'},
1507 {"minimal-match", required_argument, 0, 'C'},
1508 {"high-frequency-cutoff", required_argument, 0, 'D'},
1509 {"order", required_argument, 0, 'E'},
1510 {"approximant", required_argument, 0, 'F'},
1511 {"num-freq-cutoffs", required_argument, 0, '1'},
1512 {"max-high-freq-cutoff", required_argument, 0, '2'},
1513 {"min-high-freq-cutoff", required_argument, 0, '3'},
1514 {"space", required_argument, 0, 'G'},
1515 {"grid-spacing", required_argument, 0, 'v'},
1516 {"max-total-mass", required_argument, 0, 'y'},
1517 {"min-total-mass", required_argument, 0, 'W'},
1518 {"chirp-mass-cutoff", required_argument, 0, 'q'},
1519 {"max-eta", required_argument, 0, '0'},
1520 {"min-eta", required_argument, 0, 'X'},
1521 {"disable-polygon-fit", no_argument, &polygonFit, 0 },
1522 {"disable-compute-moments", no_argument, &computeMoments, 0 },
1523 /* standard candle parameters */
1524 {"candle-snr", required_argument, 0, 'k'},
1525 {"candle-minmass", required_argument, 0, 'l'},
1526 {"candle-maxmass", required_argument, 0, 'm'},
1527 /* frame writing options */
1528 {"write-raw-data", no_argument, &writeRawData, 1 },
1529 {"write-response", no_argument, &writeResponse, 1 },
1530 {"write-spectrum", no_argument, &writeSpectrum, 1 },
1531 {"write-strain-spectrum", no_argument, &writeStrainSpec, 1 },
1532 /* td follow up file */
1533 {"td-follow-up", required_argument, 0, 'w'},
1534 {0, 0, 0, 0}
1535 };
1536
1537 int c;
1538#ifdef LALAPPS_CUDA_ENABLED
1539 INT4 gpuDeviceID = 0;
1540 cudaError_t cudaError = cudaSuccess;
1541#endif
1542 ProcessParamsTable *this_proc_param = procparams;
1543 UINT4 haveOrder = 0;
1544 UINT4 haveApprox = 0;
1545 UINT4 haveSpace = 0;
1546 UINT4 havePsi0Min = 0;
1547 UINT4 havePsi0Max = 0;
1548 UINT4 havePsi3Min = 0;
1549 UINT4 havePsi3Max = 0;
1550 UINT4 haveAlpha = 0;
1551 UINT4 haveNumFcut = 0;
1552 UINT4 haveMaxFcut = 0;
1553 UINT4 haveMinFcut = 0;
1554
1555 /*
1556 *
1557 * parse command line arguments
1558 *
1559 */
1560
1561 while ( 1 )
1562 {
1563 /* LALgetopt_long stores long option here */
1564 int option_index = 0;
1565 size_t LALoptarg_len;
1566
1567 c = LALgetopt_long_only( argc, argv,
1568#ifdef LALAPPS_CUDA_ENABLED
1569 "a:b:c:d:e:f:g:hi:j:k:l:m:n:o:p:r:s:t:u:v:x:yX:0:"
1570 "A:B:C:D:E:F:G:H:I:J:K:L:M:O:P:Q:R:S:T:U:VZ:1:2:3:4:5:6:7:8:9:+:",
1571#else
1572 "a:b:c:d:e:f:g:hi:j:k:l:m:n:o:p:r:s:t:u:v:x:yX:0:"
1573 "A:B:C:D:E:F:G:H:I:J:K:L:M:O:P:Q:R:S:T:U:VZ:1:2:3:4:5:6:7:8:9:",
1574#endif
1575 long_options, &option_index );
1576
1577 /* detect the end of the options */
1578 if ( c == - 1 )
1579 {
1580 break;
1581 }
1582
1583 switch ( c )
1584 {
1585 case 0:
1586 /* if this option set a flag, do nothing else now */
1587 if ( long_options[option_index].flag != 0 )
1588 {
1589 break;
1590 }
1591 else
1592 {
1593 fprintf( stderr, "error parsing option %s with argument %s\n",
1594 long_options[option_index].name, LALoptarg );
1595 exit( 1 );
1596 }
1597 break;
1598
1599 case 'a':
1600 {
1601 long int gstartt = atol( LALoptarg );
1602 if ( gstartt < 441417609 )
1603 {
1604 fprintf( stderr, "invalid argument to --%s:\n"
1605 "GPS start time is prior to "
1606 "Jan 01, 1994 00:00:00 UTC:\n"
1607 "(%ld specified)\n",
1608 long_options[option_index].name, gstartt );
1609 exit( 1 );
1610 }
1611 gpsStartTime.gpsSeconds = (INT4) gstartt;
1613 ADD_PROCESS_PARAM( "int", "%ld", gstartt );
1614 }
1615 break;
1616
1617 case 'b':
1618 {
1619 long int gendt = atol( LALoptarg );
1620 if ( gendt < 441417609 )
1621 {
1622 fprintf( stderr, "invalid argument to --%s:\n"
1623 "GPS end time is prior to "
1624 "Jan 01, 1994 00:00:00 UTC:\n"
1625 "(%ld specified)\n",
1626 long_options[option_index].name, gendt );
1627 exit( 1 );
1628 }
1629 gpsEndTime.gpsSeconds = (INT4) gendt;
1631 ADD_PROCESS_PARAM( "int", "%ld", gendt );
1632 }
1633 break;
1634
1635 case 'c':
1636 {
1637 /* create storage for the channel name and copy it */
1638 char *channamptr = NULL;
1639 LALoptarg_len = strlen( LALoptarg ) + 1;
1640 fqChanName = (CHAR *) calloc( LALoptarg_len, sizeof(CHAR) );
1641 memcpy( fqChanName, LALoptarg, LALoptarg_len );
1642 ADD_PROCESS_PARAM( "string", "%s", LALoptarg );
1643
1644 /* check that we have a proper channel name */
1645 if ( ! (channamptr = strstr( fqChanName, ":" ) ) )
1646 {
1647 fprintf( stderr, "invalid argument to --%s:\n"
1648 "channel name must be a full LIGO channel name "
1649 "e.g. L1:LSC-AS_Q\n(%s specified)\n",
1650 long_options[option_index].name, LALoptarg );
1651 exit( 1 );
1652 }
1653 LALoptarg_len = strlen( ++channamptr ) + 1;
1654 channelName = (CHAR *) calloc( LALoptarg_len, sizeof(CHAR) );
1655 memcpy( channelName, channamptr, LALoptarg_len );
1656
1657 /* copy the first two characters to ifo */
1658 memset( ifo, 0, sizeof(ifo) );
1659 memcpy( ifo, LALoptarg, sizeof(ifo) - 1 );
1660 }
1661 break;
1662
1663 case 'd':
1664 numPoints = (INT4) atoi( LALoptarg );
1665 if ( numPoints < 2 || numPoints % 2 )
1666 {
1667 fprintf( stderr, "invalid argument to --%s:\n"
1668 "number of points must be a non-zero power of 2: "
1669 "(%d specified) \n",
1670 long_options[option_index].name, numPoints );
1671 exit( 1 );
1672 }
1673 ADD_PROCESS_PARAM( "int", "%d", numPoints );
1674 break;
1675
1676 case 'e':
1677 numSegments = (INT4) atoi( LALoptarg );
1678 if ( numSegments < 1 )
1679 {
1680 fprintf( stderr, "invalid argument to --%s:\n"
1681 "number of data segment must be greater than 0: "
1682 "(%d specified)\n",
1683 long_options[option_index].name, numSegments );
1684 exit( 1 );
1685 }
1686 ADD_PROCESS_PARAM( "int", "%d", numSegments );
1687 break;
1688
1689 case 'g':
1690 sampleRate = (INT4) atoi( LALoptarg );
1691 if ( sampleRate < 2 || sampleRate > 16384 || sampleRate % 2 )
1692 {
1693 fprintf( stderr, "invalid argument to --%s:\n"
1694 "rate must be power of 2 between 2 and 16384 inclusive: "
1695 "(%d specified)\n",
1696 long_options[option_index].name, sampleRate );
1697 exit( 1 );
1698 }
1699 ADD_PROCESS_PARAM( "int", "%d", sampleRate );
1700 break;
1701
1702#ifdef LALAPPS_CUDA_ENABLED
1703 case '+':
1704 gpuDeviceID = (INT4) atoi( LALoptarg );
1705 cudaError = cudaSetDevice( gpuDeviceID );
1706 if ( cudaError != cudaSuccess )
1707 {
1708 fprintf( stderr, "invalid argument to --%s:\n"
1709 "could not associate thread to GPU %d\n"
1710 "CudaError: %s\n",
1711 long_options[option_index].name, gpuDeviceID,
1712 cudaGetErrorString(cudaError));
1713 exit( 1 );
1714 }
1715 ADD_PROCESS_PARAM( "int", "%d", gpuDeviceID );
1716 break;
1717#endif
1718
1719 case 'M':
1720 /* specify which type of calibrated data */
1721 {
1722 if ( ! strcmp( "real_4", LALoptarg ) )
1723 {
1724 calData = real_4;
1725 }
1726 else if ( ! strcmp( "real_8", LALoptarg ) )
1727 {
1728 calData = real_8;
1729 }
1730 else
1731 {
1732 fprintf( stderr, "invalid argument to --%s:\n"
1733 "unknown data type specified;\n"
1734 "%s (must be one of: real_4, real_8)\n",
1735 long_options[option_index].name, LALoptarg);
1736 }
1737 ADD_PROCESS_PARAM( "string", "%s", LALoptarg );
1738 }
1739 break;
1740
1741 case 'J':
1743 if ( strainHighPassFreq <= 0 )
1744 {
1745 fprintf( stderr, "invalid argument to --%s:\n"
1746 "REAL8 h(t) high pass filter frequency must be greater than 0 Hz:"
1747 "(%f Hz specified)\n",
1748 long_options[option_index].name, strainHighPassFreq );
1749 exit( 1 );
1750 }
1751 ADD_PROCESS_PARAM( "float", "%e", strainHighPassFreq );
1752 break;
1753
1754 case 'K':
1756 if ( strainHighPassOrder <= 0 )
1757 {
1758 fprintf( stderr, "invalid argument to --%s:\n"
1759 "REAL8 h(t) high pass filter order must be greater than 0: "
1760 "(%d specified)\n",
1761 long_options[option_index].name, strainHighPassOrder );
1762 exit( 1 );
1763 }
1764 ADD_PROCESS_PARAM( "int", "%d", strainHighPassOrder );
1765 break;
1766
1767 case 'L':
1769 if ( strainHighPassAtten < 0.0 || strainHighPassAtten > 1.0 )
1770 {
1771 fprintf( stderr, "invalid argument to --%s:\n"
1772 "REAL8 h(t) high pass attenuation must be in the range [0:1]: "
1773 "(%f specified)\n",
1774 long_options[option_index].name, strainHighPassAtten );
1775 exit( 1 );
1776 }
1777 ADD_PROCESS_PARAM( "float", "%e", strainHighPassAtten );
1778 break;
1779
1780 case 'h':
1781 USAGE( stdout );
1782 exit( 0 );
1783 break;
1784
1785 case 'i':
1786 fLow = (REAL4) atof( LALoptarg );
1787 if ( fLow < 0 )
1788 {
1789 fprintf( stdout, "invalid argument to --%s:\n"
1790 "low frequency cutoff is less than 0 Hz: "
1791 "(%f Hz specified)\n",
1792 long_options[option_index].name, fLow );
1793 exit( 1 );
1794 }
1795 ADD_PROCESS_PARAM( "float", "%e", fLow );
1796 break;
1797
1798 case 'j':
1799 if ( ! strcmp( "mean", LALoptarg ) )
1800 {
1802 }
1803 else if ( ! strcmp( "median", LALoptarg ) )
1804 {
1806 }
1807 else if ( ! strcmp( "iLIGOSRD", LALoptarg ) )
1810 else if ( ! strcmp( "eLIGOModel", LALoptarg ) )
1813 else if ( ! strcmp( "GEOModel", LALoptarg ) )
1816 else if ( ! strcmp( "aLIGONoSRMLoP", LALoptarg ) )
1819 else if ( ! strcmp( "aLIGONoSRMHiP", LALoptarg ) )
1822 else if ( ! strcmp( "aLIGOZDLoP", LALoptarg ) )
1825 else if ( ! strcmp( "aLIGOZDHiP", LALoptarg ) )
1828 else if ( ! strcmp( "iVirgoModel", LALoptarg ) )
1831 else if ( ! strcmp( "aVirgoModel", LALoptarg ) )
1834 else if ( ! strcmp( "KAGRAModel", LALoptarg ) )
1837
1838 else
1839 {
1840 fprintf( stderr, "invalid argument to --%s:\n"
1841 "unknown power spectrum type: %s\n",
1842 long_options[option_index].name, LALoptarg );
1843 exit( 1 );
1844 }
1845 ADD_PROCESS_PARAM( "string", "%s", LALoptarg );
1846 break;
1847
1848 case 'f':
1849 dynRangeExponent = (REAL4) atof( LALoptarg );
1850 ADD_PROCESS_PARAM( "float", "%e", dynRangeExponent );
1851 break;
1852
1853 case 'p':
1854 /* create storage for the calibration frame cache name */
1855 LALoptarg_len = strlen( LALoptarg ) + 1;
1856 calCacheName = (CHAR *) calloc( LALoptarg_len, sizeof(CHAR));
1857 memcpy( calCacheName, LALoptarg, LALoptarg_len );
1858 ADD_PROCESS_PARAM( "string", "%s", LALoptarg );
1859 break;
1860
1861 case 'r':
1862 if ( ! strcmp( "ldas", LALoptarg ) )
1863 {
1864 resampFiltType = 0;
1865 }
1866 else if ( ! strcmp( "butterworth", LALoptarg ) )
1867 {
1868 resampFiltType = 1;
1869 }
1870 else
1871 {
1872 fprintf( stderr, "invalid argument to --%s:\n"
1873 "unknown resampling filter type: "
1874 "%s (must be ldas or butterworth)\n",
1875 long_options[option_index].name, LALoptarg );
1876 exit( 1 );
1877 }
1878 ADD_PROCESS_PARAM( "string", "%s", LALoptarg );
1879 break;
1880
1881 case 's':
1882 if ( strlen( LALoptarg ) > LIGOMETA_COMMENT_MAX - 1 )
1883 {
1884 fprintf( stderr, "invalid argument to --%s:\n"
1885 "comment must be less than %d characters\n",
1886 long_options[option_index].name, LIGOMETA_COMMENT_MAX );
1887 exit( 1 );
1888 }
1889 else
1890 {
1891 snprintf( comment, LIGOMETA_COMMENT_MAX, "%s", LALoptarg);
1892 }
1893 break;
1894
1895 case 't':
1896 highPass = 1;
1897 highPassFreq = (REAL4) atof( LALoptarg );
1898 if ( highPassFreq < 0 )
1899 {
1900 fprintf( stdout, "invalid argument to --%s:\n"
1901 "low frequency cutoff is less than 0 Hz: "
1902 "(%f Hz specified)\n",
1903 long_options[option_index].name, highPassFreq );
1904 exit( 1 );
1905 }
1906 ADD_PROCESS_PARAM( "float", "%e", highPassFreq );
1907 break;
1908
1909 case 'H':
1910 highPassOrder = (INT4) atoi( LALoptarg );
1911 if ( highPassOrder <= 0 )
1912 {
1913 fprintf( stdout, "invalid argument to --%s:\n"
1914 "high pass filter order must be greater than 0: "
1915 "(%d specified)\n",
1916 long_options[option_index].name, highPassOrder );
1917 exit( 1 );
1918 }
1919 ADD_PROCESS_PARAM( "int", "%d", highPassOrder );
1920 break;
1921
1922 case 'I':
1923 highPassAtten = (REAL4) atof( LALoptarg );
1924 if ( highPassAtten < 0.0 || highPassAtten > 1.0 )
1925 {
1926 fprintf( stdout, "invalid argument to --%s:\n"
1927 "high pass attenuation must be in the range [0:1]: "
1928 "(%f specified)\n",
1929 long_options[option_index].name, highPassAtten );
1930 exit( 1 );
1931 }
1932 ADD_PROCESS_PARAM( "float", "%e", highPassAtten );
1933 break;
1934
1935 case 'u':
1936 LALoptarg_len = strlen( LALoptarg ) + 1;
1937 frInCacheName = (CHAR *) calloc( LALoptarg_len, sizeof(CHAR) );
1938 memcpy( frInCacheName, LALoptarg, LALoptarg_len );
1939 ADD_PROCESS_PARAM( "string", "%s", LALoptarg );
1940 break;
1941
1942 case 'n':
1943 LALoptarg_len = strlen( LALoptarg ) + 1;
1944 frInType = (CHAR *) calloc( LALoptarg_len, sizeof(CHAR) );
1945 memcpy( frInType, LALoptarg, LALoptarg_len );
1946 ADD_PROCESS_PARAM( "string", "%s", LALoptarg );
1947 break;
1948
1949 case 'x':
1950 padData = (UINT4) atoi( LALoptarg );
1951 if ( padData < 0 )
1952 {
1953 fprintf( stderr, "invalid argument to --%s:\n"
1954 "number of seconds to pad from input data"
1955 "must be greater than 0: (%d specified)\n",
1956 long_options[option_index].name, padData );
1957 exit( 1 );
1958 }
1959 ADD_PROCESS_PARAM( "int", "%d", padData );
1960 break;
1961
1962 case 'Z':
1963 /* create storage for the usertag */
1964 LALoptarg_len = strlen( LALoptarg ) + 1;
1965 userTag = (CHAR *) calloc( LALoptarg_len, sizeof(CHAR) );
1966 memcpy( userTag, LALoptarg, LALoptarg_len );
1967
1968 this_proc_param = this_proc_param->next = (ProcessParamsTable *)
1969 calloc( 1, sizeof(ProcessParamsTable) );
1970 snprintf( this_proc_param->program, LIGOMETA_PROGRAM_MAX, "%s",
1971 PROGRAM_NAME );
1972 snprintf( this_proc_param->param, LIGOMETA_PARAM_MAX, "--user-tag" );
1973 snprintf( this_proc_param->type, LIGOMETA_TYPE_MAX, "string" );
1974 snprintf( this_proc_param->value, LIGOMETA_VALUE_MAX, "%s",
1975 LALoptarg );
1976 break;
1977
1978 case 'A':
1979 minMass = (REAL4) atof( LALoptarg );
1980 if ( minMass <= 0 )
1981 {
1982 fprintf( stdout, "invalid argument to --%s:\n"
1983 "minimum component mass must be > 0: "
1984 "(%f solar masses specified)\n",
1985 long_options[option_index].name, minMass );
1986 exit( 1 );
1987 }
1988 ADD_PROCESS_PARAM( "float", "%e", minMass );
1989 break;
1990
1991 case 'B':
1992 maxMass = (REAL4) atof( LALoptarg );
1993 if ( maxMass <= 0 )
1994 {
1995 fprintf( stdout, "invalid argument to --%s:\n"
1996 "maximum component mass must be > 0: "
1997 "(%f solar masses specified)\n",
1998 long_options[option_index].name, maxMass );
1999 exit( 1 );
2000 }
2001 ADD_PROCESS_PARAM( "float", "%e", maxMass );
2002 break;
2003
2004 case 'P':
2005 psi0Min = (REAL4) atof( LALoptarg );
2006 if ( psi0Min <= 0 )
2007 {
2008 fprintf( stdout, "invalid argument to --%s:\n"
2009 "minimum value of psi0 must be > 0: "
2010 "(%f specified)\n",
2011 long_options[option_index].name, psi0Min );
2012 exit( 1 );
2013 }
2014 ADD_PROCESS_PARAM( "float", "%e", psi0Min );
2015 havePsi0Min = 1;
2016 break;
2017
2018 case 'Q':
2019 psi0Max = (REAL4) atof( LALoptarg );
2020 if ( psi0Max <= 0 )
2021 {
2022 fprintf( stdout, "invalid argument to --%s:\n"
2023 "maximum value of psi0 must be > 0: "
2024 "(%f specified)\n",
2025 long_options[option_index].name, psi0Max );
2026 exit( 1 );
2027 }
2028 ADD_PROCESS_PARAM( "float", "%e", psi0Max );
2029 havePsi0Max = 1;
2030 break;
2031
2032 case 'R':
2033 psi3Min = (REAL4) atof( LALoptarg );
2034 if ( psi3Min >= 0 )
2035 {
2036 fprintf( stdout, "invalid argument to --%s:\n"
2037 "miniumum value of psi3 must be < 0: "
2038 "(%f specified)\n",
2039 long_options[option_index].name, psi3Min );
2040 exit( 1 );
2041 }
2042 ADD_PROCESS_PARAM( "float", "%e", psi3Min );
2043 havePsi3Min = 1;
2044 break;
2045
2046 case 'S':
2047 psi3Max = (REAL4) atof( LALoptarg );
2048 ADD_PROCESS_PARAM( "float", "%e", psi3Max );
2049 havePsi3Max = 1;
2050 break;
2051
2052 case 'o':
2053 betaMin = (REAL4) atof( LALoptarg );
2054 ADD_PROCESS_PARAM( "float", "%e", betaMin );
2055 break;
2056
2057 case 'O':
2058 betaMax = (REAL4) atof( LALoptarg );
2059 ADD_PROCESS_PARAM( "float", "%e", betaMax );
2060 break;
2061
2062 case 'U':
2063 maxFcutTmplts = (INT4) atof( LALoptarg );
2064 if ( maxFcutTmplts < 0 )
2065 {
2066 fprintf( stdout, "invalid argument to --%s:\n"
2067 "number of templates in f_final direction must be >= 0"
2068 "(%d specified)\n",
2069 long_options[option_index].name, maxFcutTmplts );
2070 exit( 1 );
2071 }
2072 ADD_PROCESS_PARAM( "int", "%d", maxFcutTmplts );
2073 break;
2074
2075 case 'T':
2076 alpha = (REAL4) atof( LALoptarg );
2077 if ( alpha < -1 || alpha > 1 )
2078 {
2079 fprintf( stdout, "invalid argument to --%s:\n"
2080 "value of alpha must be the range [0:1]"
2081 "(%f specified)\n",
2082 long_options[option_index].name, alpha );
2083 exit( 1 );
2084 }
2085 ADD_PROCESS_PARAM( "float", "%e", alpha );
2086 haveAlpha = 1;
2087 break;
2088
2089 case 'C':
2090 minMatch = (REAL4) atof( LALoptarg );
2091 if ( minMatch <= 0 )
2092 {
2093 fprintf( stdout, "invalid argument to --%s:\n"
2094 "minimum match of bank must be > 0: "
2095 "(%f specified)\n",
2096 long_options[option_index].name, minMatch );
2097 exit( 1 );
2098 }
2099 ADD_PROCESS_PARAM( "float", "%e", minMatch );
2100 break;
2101
2102 case 'D':
2103 fUpper = (REAL4) atof( LALoptarg );
2104 if ( fUpper <= 0 )
2105 {
2106 fprintf( stdout, "invalid argument to --%s:\n"
2107 "miniumu component mass must be > 0: "
2108 "(%f specified)\n",
2109 long_options[option_index].name, fUpper );
2110 exit( 1 );
2111 }
2112 ADD_PROCESS_PARAM( "float", "%e", fUpper );
2113 break;
2114
2115 case 'E':
2116 if ( ! strcmp( "newtonian", LALoptarg ) )
2117 {
2119 }
2120 else if ( ! strcmp( "oneHalfPN", LALoptarg ) )
2121 {
2123 }
2124 else if ( ! strcmp( "onePN", LALoptarg ) )
2125 {
2127 }
2128 else if ( ! strcmp( "onePointFivePN", LALoptarg ) )
2129 {
2131 }
2132 else if ( ! strcmp( "twoPN", LALoptarg ) )
2133 {
2135 }
2136 else if ( ! strcmp( "twoPointFive", LALoptarg ) )
2137 {
2139 }
2140 else if ( ! strcmp( "threePN", LALoptarg ) )
2141 {
2143 }
2144 else if ( ! strcmp( "threePointFivePN", LALoptarg ) )
2145 {
2147 }
2148 else
2149 {
2150 fprintf( stderr, "invalid argument to --%s:\n"
2151 "unknown order specified: "
2152 "%s (must be one of: newtonian, oneHalfPN, onePN,\n"
2153 "onePointFivePN, twoPN, twoPointFivePN, threePN or\n"
2154 "threePointFivePN)\n",
2155 long_options[option_index].name, LALoptarg );
2156 exit( 1 );
2157 }
2158 haveOrder = 1;
2159 ADD_PROCESS_PARAM( "string", "%s", LALoptarg );
2160 break;
2161
2162 case 'F':
2163 if ( ! strcmp( "TaylorT1", LALoptarg ) )
2164 {
2166 }
2167 else if ( ! strcmp( "TaylorT2", LALoptarg ) )
2168 {
2170 }
2171 else if ( ! strcmp( "TaylorT3", LALoptarg ) )
2172 {
2174 }
2175 else if ( ! strcmp( "TaylorF1", LALoptarg ) )
2176 {
2178 }
2179 else if ( ! strcmp( "TaylorF2", LALoptarg ) )
2180 {
2182 }
2183 else if ( ! strcmp( "PadeT1", LALoptarg ) )
2184 {
2186 }
2187 else if ( ! strcmp( "PadeF1", LALoptarg ) )
2188 {
2190 }
2191 else if ( ! strcmp( "EOB", LALoptarg ) )
2192 {
2193 approximant = EOB;
2194 }
2195 else if ( ! strcmp( "EOBNR", LALoptarg ) )
2196 {
2198 }
2199 else if ( ! strcmp( "EOBNRv2", LALoptarg ) )
2200 {
2202 }
2203 else if ( ! strcmp( "IMRPhenomA", LALoptarg ) )
2204 {
2206 }
2207 else if ( ! strcmp( "IMRPhenomB", LALoptarg ) )
2208 {
2210 }
2211 else if ( ! strcmp( "BCV", LALoptarg ) )
2212 {
2213 approximant = BCV;
2214 }
2215 else if ( ! strcmp( "SpinTaylorT3", LALoptarg ) )
2216 {
2218 }
2219 else if ( ! strcmp( "BCVSpin", LALoptarg ) )
2220 {
2222 }
2223 else if ( ! strcmp( "FindChirpPTF", LALoptarg ) )
2224 {
2226 }
2227 else
2228 {
2229 fprintf( stderr, "invalid argument to --%s:\n"
2230 "unknown order specified: "
2231 "%s (must be one of: TaylorT1, TaylorT2, TaylorT3, TaylorF1,\n"
2232 "TaylorF2, PadeT1, PadeF1, EOB, EOBNR, EOBNRv2, IMRPhenomA,\n"
2233 "IMRPhenomB, BCV, SpinTaylorT3, BCVSpin\n"
2234 "or FindChirpPTF)\n", long_options[option_index].name, LALoptarg );
2235 exit( 1 );
2236 }
2237 haveApprox = 1;
2238 ADD_PROCESS_PARAM( "string", "%s", LALoptarg );
2239 break;
2240
2241 case 'G':
2242 if ( ! strcmp( "Tau0Tau2", LALoptarg ) )
2243 {
2244 space = Tau0Tau2;
2245 }
2246 else if ( ! strcmp( "Tau0Tau3", LALoptarg ) )
2247 {
2248 space = Tau0Tau3;
2249 }
2250 else if ( ! strcmp( "Psi0Psi3", LALoptarg ) )
2251 {
2252 space = Psi0Psi3;
2253 }
2254 else
2255 {
2256 fprintf( stderr, "invalid argument to --%s:\n"
2257 "unknown space specified: "
2258 "%s (must be one of: Tau0Tau2, Tau0Tau3 or Psi0Psi3)\n",
2259 long_options[option_index].name, LALoptarg );
2260 exit( 1 );
2261 }
2262 haveSpace = 1;
2263 ADD_PROCESS_PARAM( "string", "%s", LALoptarg );
2264 break;
2265
2266 case 'v':
2267 if ( ! strcmp( "Hexagonal", LALoptarg) )
2268 {
2269 haveGridSpacing = 1;
2271 }
2272 else if ( ! strcmp( "SquareNotOriented", LALoptarg) )
2273 {
2274 haveGridSpacing = 1;
2276 }
2277 else
2278 {
2279 fprintf(stderr, "invalid argument to --%s:\n"
2280 "unknown grid spacing specified: "
2281 "%s (must be one of Hexagonal, SquareNotOriented )\n",
2282 long_options[option_index].name, LALoptarg );
2283 exit(1);
2284 }
2285 ADD_PROCESS_PARAM( "string", "%s", LALoptarg );
2286 break;
2287
2288 case 'y':
2289 maxTotalMass = (REAL4) atof( LALoptarg );
2290 if ( maxTotalMass <= 0 )
2291 {
2292 fprintf( stdout, "invalid argument to --%s:\n"
2293 "maximum total mass must be > 0: "
2294 "(%f solar masses specified)\n",
2295 long_options[option_index].name, maxTotalMass );
2296 exit( 1 );
2297 }
2298 ADD_PROCESS_PARAM( "float", "%e", maxTotalMass );
2299 break;
2300
2301 case 'W':
2302 minTotalMass = (REAL4) atof( LALoptarg );
2303 if ( minTotalMass <= 0 )
2304 {
2305 fprintf( stdout, "invalid argument to --%s:\n"
2306 "minimum total mass must be > 0: "
2307 "(%f solar masses specified)\n",
2308 long_options[option_index].name, minTotalMass );
2309 exit( 1 );
2310 }
2311 ADD_PROCESS_PARAM( "float", "%e", minTotalMass );
2312 break;
2313
2314 case 'q':
2315 chirpMassCutoff = (REAL4) atof( LALoptarg );
2316 if ( chirpMassCutoff <= 0 )
2317 {
2318 fprintf( stdout, "invalid argument to --%s:\n"
2319 "chirp mass cutoff must be > 0: "
2320 "(%f solar masses specified)\n",
2321 long_options[option_index].name, chirpMassCutoff );
2322 exit( 1 );
2323 }
2324 ADD_PROCESS_PARAM( "float", "%e", chirpMassCutoff );
2325 break;
2326
2327 case 'X':
2328 etaMinCutoff = (REAL4) atof( LALoptarg );
2329 if ( etaMinCutoff < 0 || etaMinCutoff >= 0.25 )
2330 {
2331 fprintf( stdout, "invalid argument to --%s:\n"
2332 "minimum eta must be >= 0 and < 0.25: "
2333 "(%f specified)\n",
2334 long_options[option_index].name, etaMinCutoff);
2335 exit( 1 );
2336 }
2337 ADD_PROCESS_PARAM( "float", "%e", etaMinCutoff);
2338 break;
2339
2340 case '0':
2341 etaMaxCutoff = (REAL4) atof( LALoptarg );
2342 if ( etaMaxCutoff <= 0 || etaMaxCutoff > 0.25 )
2343 {
2344 fprintf( stdout, "invalid argument to --%s:\n"
2345 "maximum eta must be > 0 and <= 0.25: "
2346 "(%f specified)\n",
2347 long_options[option_index].name, etaMaxCutoff );
2348 exit( 1 );
2349 }
2350 ADD_PROCESS_PARAM( "float", "%e", etaMaxCutoff );
2351 break;
2352
2353 case 'k':
2354 candleSnr = (REAL4) atof( LALoptarg );
2355 if ( candleSnr <= 0 )
2356 {
2357 fprintf( stdout, "invalid argument to --%s:\n"
2358 "standard candle signal-to-noise ratio must be > 0: "
2359 "(%f specified)\n",
2360 long_options[option_index].name, candleSnr );
2361 exit( 1 );
2362 }
2363 ADD_PROCESS_PARAM( "float", "%e", candleSnr );
2364 break;
2365
2366 case 'l':
2367 candleMinMass = (REAL4) atof( LALoptarg );
2368 if ( candleMinMass <= 0 )
2369 {
2370 fprintf( stdout, "invalid argument to --%s:\n"
2371 "standard candle minimum component mass must be > 0: "
2372 "(%f specified)\n",
2373 long_options[option_index].name, candleMinMass );
2374 exit( 1 );
2375 }
2376 ADD_PROCESS_PARAM( "float", "%e", candleMinMass );
2377 break;
2378
2379 case 'm':
2380 candleMaxMass = (REAL4) atof( LALoptarg );
2381 if ( ( candleMaxMass <= 0 ) || ( candleMaxMass < candleMinMass ) )
2382 {
2383 fprintf( stdout, "invalid argument to --%s:\n"
2384 "standard candle maximum component mass must be > 0 and greater than min mass: "
2385 "(%f specified)\n",
2386 long_options[option_index].name, candleMaxMass );
2387 exit( 1 );
2388 }
2389 ADD_PROCESS_PARAM( "float", "%e", candleMaxMass );
2390 break;
2391
2392 case 'w':
2393 numTDFiles = 1;
2394 /* Count the number of thinca files to follow up */
2395 while ( !strstr( argv[LALoptind], "--" ) )
2396 {
2397 numTDFiles++;
2398 LALoptind++;
2399 }
2401
2402 /* Set pointers to the relevant filenames */
2403 tdFileNames = (CHAR **) calloc( numTDFiles, sizeof(CHAR *));
2404 numTDFiles = 0;
2405
2406 while ( !strstr( argv[LALoptind], "--" ) )
2407 {
2408 tdFileNames[numTDFiles++] = argv[LALoptind];
2409 ADD_PROCESS_PARAM( "string", "%s", argv[LALoptind] );
2410 LALoptind++;
2411
2412 }
2413 break;
2414
2415 case 'Y':
2416 /* create storage for the ifo-tag */
2417 LALoptarg_len = strlen( LALoptarg ) + 1;
2418 ifoTag = (CHAR *) calloc( LALoptarg_len, sizeof(CHAR) );
2419 memcpy( ifoTag, LALoptarg, LALoptarg_len );
2420 ADD_PROCESS_PARAM( "string", "%s", LALoptarg );
2421 break;
2422
2423 case 'V':
2424 /* print version information and exit */
2425 fprintf( stdout, "LIGO/LSC Standalone Inspiral Template Bank Code\n"
2426 "Duncan Brown <duncan@gravity.phys.uwm.edu>\n");
2427 XLALOutputVCSInfo(stderr, lalAppsVCSInfoList, 0, "%% ");
2428 exit( 0 );
2429 break;
2430
2431 case '?':
2432 USAGE( stderr );
2433 exit( 1 );
2434 break;
2435
2436 case '1':
2437 numFreqCut = (INT4) atof( LALoptarg );
2438 if( numFreqCut < 1 )
2439 {
2440 fprintf( stdout, "invalid argument to --%s:\n"
2441 "Value must be a positive integer "
2442 "(%d specified)\n",
2443 long_options[option_index].name, numFreqCut );
2444 exit( 1 );
2445 }
2446 ADD_PROCESS_PARAM( "int", "%d", numFreqCut );
2447 haveNumFcut = 1;
2448 break;
2449
2450 case '2':
2451 if ( ! strcmp( "SchwarzISCO", LALoptarg ) )
2452 {
2454 }
2455 else if( ! strcmp( "BKLISCO", LALoptarg ) )
2456 {
2458 }
2459 else if ( ! strcmp( "LightRing", LALoptarg ) )
2460 {
2462 }
2463 else if ( ! strcmp( "FRD", LALoptarg ) )
2464 {
2466 }
2467 else if ( ! strcmp( "ERD", LALoptarg ) )
2468 {
2470 }
2471 else if ( ! strcmp( "LRD", LALoptarg ) )
2472 {
2474 }
2475 else
2476 {
2477 fprintf( stderr, "invalid argument to --%s:\n"
2478 "unknown cutoff frequency specified: "
2479 "%s (must be one of: SchwarzISCO, BKLISCO, LightRing, FRD, ERD or LRD)\n",
2480 long_options[option_index].name, LALoptarg );
2481 exit( 1 );
2482 }
2483 ADD_PROCESS_PARAM( "string", "%s", LALoptarg );
2484 haveMaxFcut = 1;
2485 break;
2486
2487 case '3':
2488 if ( ! strcmp( "SchwarzISCO", LALoptarg ) )
2489 {
2491 }
2492 else if ( ! strcmp( "BKLISCO", LALoptarg ) )
2493 {
2495 }
2496 else if ( ! strcmp( "LightRing", LALoptarg ) )
2497 {
2499 }
2500 else if ( ! strcmp( "FRD", LALoptarg ) )
2501 {
2503 }
2504 else if ( ! strcmp( "ERD", LALoptarg ) )
2505 {
2507 }
2508 else if ( ! strcmp( "LRD", LALoptarg ) )
2509 {
2511 }
2512 else
2513 {
2514 fprintf( stderr, "invalid argument to --%s:\n"
2515 "unknown cutoff frequency specified: "
2516 "%s (must be one of: SchwarzISCO, BKLISCO, LightRing, FRD, ERD, or LRD)\n",
2517 long_options[option_index].name, LALoptarg );
2518 exit( 1 );
2519 }
2520 ADD_PROCESS_PARAM( "string", "%s", LALoptarg );
2521 haveMinFcut = 1;
2522 break;
2523
2524 case '4':
2525 chiMin = atof( LALoptarg );
2526 if ( chiMin < 0. )
2527 {
2528 fprintf( stdout, "invalid argument to --%s:\n"
2529 "Spin magnitude can only take values between 0 and 1. : "
2530 "(%f specified)\n",
2531 long_options[option_index].name, chiMin );
2532 exit( 1 );
2533 }
2534 ADD_PROCESS_PARAM( "float", "%e", chiMin );
2535 break;
2536
2537 case '5':
2538 chiMax = atof( LALoptarg );
2539 if ( chiMax > 1. )
2540 {
2541 fprintf( stdout, "invalid argument to --%s:\n"
2542 "Spin magnitude can only take values between 0 and 1. : "
2543 "(%f specified)\n",
2544 long_options[option_index].name, chiMax );
2545 exit( 1 );
2546 }
2547 ADD_PROCESS_PARAM( "float", "%e", chiMax );
2548 break;
2549
2550 case '6':
2551 kappaMin = atof( LALoptarg );
2552 if ( kappaMin < -1. )
2553 {
2554 fprintf( stdout, "invalid argument to --%s:\n"
2555 "Kappa can only take values between -1. and 1. : "
2556 "(%f specified)\n",
2557 long_options[option_index].name, kappaMin );
2558 exit( 1 );
2559 }
2560 ADD_PROCESS_PARAM( "float", "%e", kappaMin );
2561 break;
2562
2563 case '7':
2564 kappaMax = atof( LALoptarg );
2565 if ( kappaMax > 1. )
2566 {
2567 fprintf( stdout, "invalid argument to --%s:\n"
2568 "Kappa can only take values between -1. and 1. : "
2569 "(%f specified)\n",
2570 long_options[option_index].name, kappaMax );
2571 exit( 1 );
2572 }
2573 ADD_PROCESS_PARAM( "float", "%e", kappaMax );
2574 break;
2575
2576 case '8':
2577 nPointsChi = atof( LALoptarg );
2578 if ( nPointsChi < 1 )
2579 {
2580 fprintf( stdout, "invalid argument to --%s:\n"
2581 "Number of points in the Chi direction must be greater than 0 : "
2582 "(%d specified)\n",
2583 long_options[option_index].name, nPointsChi );
2584 exit( 1 );
2585 }
2586 ADD_PROCESS_PARAM( "int", "%d", nPointsChi );
2587 break;
2588
2589 case '9':
2590 nPointsKappa = atof( LALoptarg );
2591 if ( nPointsKappa < 1 )
2592 {
2593 fprintf( stdout, "invalid argument to --%s:\n"
2594 "Number of points in the Kappa direction must be greater than 0 : "
2595 "(%d specified)\n",
2596 long_options[option_index].name, nPointsKappa );
2597 exit( 1 );
2598 }
2599 ADD_PROCESS_PARAM( "int", "%d", nPointsKappa );
2600 break;
2601
2602 default:
2603 fprintf( stderr, "unknown error while parsing options\n" );
2604 USAGE( stderr );
2605 exit( 1 );
2606 }
2607 }
2608
2609
2610
2611
2612 if ( LALoptind < argc )
2613 {
2614 fprintf( stderr, "extraneous command line arguments:\n" );
2615 while ( LALoptind < argc )
2616 {
2617 fprintf ( stderr, "%s\n", argv[LALoptind++] );
2618 }
2619 exit( 1 );
2620 }
2621
2622 /* add option without arguments into the process param table */
2623 if (vrbflg==1)
2624 {
2625 this_proc_param = this_proc_param->next = (ProcessParamsTable *)
2626 calloc( 1, sizeof(ProcessParamsTable) );
2627 snprintf( this_proc_param->program, LIGOMETA_PROGRAM_MAX,
2628 "%s", PROGRAM_NAME );
2629 snprintf( this_proc_param->param, LIGOMETA_PARAM_MAX,
2630 "--verbose" );
2631 snprintf( this_proc_param->type, LIGOMETA_TYPE_MAX, "string" );
2632 snprintf( this_proc_param->value, LIGOMETA_TYPE_MAX, " " );
2633 }
2634 if (outCompress==1)
2635 {
2636 this_proc_param = this_proc_param->next = (ProcessParamsTable *)
2637 calloc( 1, sizeof(ProcessParamsTable) );
2638 snprintf( this_proc_param->program, LIGOMETA_PROGRAM_MAX,
2639 "%s", PROGRAM_NAME );
2640 snprintf( this_proc_param->param, LIGOMETA_PARAM_MAX,
2641 "--write-compress" );
2642 snprintf( this_proc_param->type, LIGOMETA_TYPE_MAX, "string" );
2643 snprintf( this_proc_param->value, LIGOMETA_TYPE_MAX, " " );
2644 }
2645 if (computeMoments==0)
2646 {
2647 this_proc_param = this_proc_param->next = (ProcessParamsTable *)
2648 calloc( 1, sizeof(ProcessParamsTable) );
2649 snprintf( this_proc_param->program, LIGOMETA_PROGRAM_MAX,
2650 "%s", PROGRAM_NAME );
2651 snprintf( this_proc_param->param, LIGOMETA_PARAM_MAX,
2652 "--disable-compute-moments" );
2653 snprintf( this_proc_param->type, LIGOMETA_TYPE_MAX, "string" );
2654 snprintf( this_proc_param->value, LIGOMETA_TYPE_MAX, " " );
2655 }
2656 if (polygonFit==0)
2657 {
2658 this_proc_param = this_proc_param->next = (ProcessParamsTable *)
2659 calloc( 1, sizeof(ProcessParamsTable) );
2660 snprintf( this_proc_param->program, LIGOMETA_PROGRAM_MAX,
2661 "%s", PROGRAM_NAME );
2662 snprintf( this_proc_param->param, LIGOMETA_PARAM_MAX,
2663 "--disable-polygon-fit" );
2664 snprintf( this_proc_param->type, LIGOMETA_TYPE_MAX, "string" );
2665 snprintf( this_proc_param->value, LIGOMETA_TYPE_MAX, " " );
2666 }
2667 if (globFrameData==1)
2668 {
2669 this_proc_param = this_proc_param->next = (ProcessParamsTable *)
2670 calloc( 1, sizeof(ProcessParamsTable) );
2671 snprintf( this_proc_param->program, LIGOMETA_PROGRAM_MAX,
2672 "%s", PROGRAM_NAME );
2673 snprintf( this_proc_param->param, LIGOMETA_PARAM_MAX,
2674 "--glob-frame-data" );
2675 snprintf( this_proc_param->type, LIGOMETA_TYPE_MAX, "string" );
2676 snprintf( this_proc_param->value, LIGOMETA_TYPE_MAX, " " );
2677 }
2678 if (globCalData==1)
2679 {
2680 this_proc_param = this_proc_param->next = (ProcessParamsTable *)
2681 calloc( 1, sizeof(ProcessParamsTable) );
2682 snprintf( this_proc_param->program, LIGOMETA_PROGRAM_MAX,
2683 "%s", PROGRAM_NAME );
2684 snprintf( this_proc_param->param, LIGOMETA_PARAM_MAX,
2685 "--glob-calibration-data" );
2686 snprintf( this_proc_param->type, LIGOMETA_TYPE_MAX, "string" );
2687 snprintf( this_proc_param->value, LIGOMETA_TYPE_MAX, " " );
2688 }
2689 if (pointCal==1)
2690 {
2691 this_proc_param = this_proc_param->next = (ProcessParamsTable *)
2692 calloc( 1, sizeof(ProcessParamsTable) );
2693 snprintf( this_proc_param->program, LIGOMETA_PROGRAM_MAX,
2694 "%s", PROGRAM_NAME );
2695 snprintf( this_proc_param->param, LIGOMETA_PARAM_MAX,
2696 "--point-calibration" );
2697 snprintf( this_proc_param->type, LIGOMETA_TYPE_MAX, "string" );
2698 snprintf( this_proc_param->value, LIGOMETA_TYPE_MAX, " " );
2699 }
2700
2701 /*
2702 *
2703 * check validity of arguments
2704 *
2705 */
2706
2707 /* if using a simulated PSD, use the first two characters of */
2708 /* ifoTag as an ifo prefix */
2710 {
2711 if ( ! ifoTag )
2712 {
2713 fprintf( stderr, "--ifo-tag must be specified if using a simulated PSD\n" );
2714 exit( 1 );
2715 }
2716 else
2717 {
2718 memset( ifo, 0, sizeof(ifo) );
2719 memcpy( ifo, ifoTag, sizeof(ifo) - 1 );
2720 }
2721 }
2722
2723 /* check validity of input data time */
2724 if ( ! gpsStartTime.gpsSeconds )
2725 {
2726 fprintf( stderr, "--gps-start-time must be specified\n" );
2727 exit( 1 );
2728 }
2729 if ( ! gpsEndTime.gpsSeconds )
2730 {
2731 fprintf( stderr, "--gps-end-time must be specified\n" );
2732 exit( 1 );
2733 }
2735 {
2736 fprintf( stderr, "invalid gps time range: "
2737 "start time: %d, end time %d\n",
2739 exit( 1 );
2740 }
2741
2742 /* check validity of data length parameters */
2743 if ( numPoints < 0 )
2744 {
2745 fprintf( stderr, "--segment-length must be specified\n" );
2746 exit( 1 );
2747 }
2749 {
2750 fprintf( stderr, "--number-of-segments must be specified if using frama data\n" );
2751 exit( 1 );
2752 }
2753
2754 /* check sample rate has been given */
2755 if ( sampleRate < 0 )
2756 {
2757 fprintf( stderr, "--sample-rate must be specified\n" );
2758 exit( 1 );
2759 }
2760
2761 /* check high pass option has been given */
2762 if ( highPass < 0 )
2763 {
2764 fprintf( stderr, "--disable-high-pass or --enable-high-pass (freq)"
2765 " must be specified\n" );
2766 exit( 1 );
2767 }
2768 else if ( ! highPass )
2769 {
2770 this_proc_param = this_proc_param->next = (ProcessParamsTable *)
2771 calloc( 1, sizeof(ProcessParamsTable) );
2772 snprintf( this_proc_param->program, LIGOMETA_PROGRAM_MAX,
2773 "%s", PROGRAM_NAME );
2774 snprintf( this_proc_param->param, LIGOMETA_PARAM_MAX,
2775 "--disable-high-pass" );
2776 snprintf( this_proc_param->type, LIGOMETA_TYPE_MAX, "string" );
2777 snprintf( this_proc_param->value, LIGOMETA_TYPE_MAX, " " );
2778 }
2779 else
2780 {
2781 /* check that all the high pass parameters have been specified */
2782 if ( highPassOrder < 0 )
2783 {
2784 fprintf( stderr, "--high-pass-order must be specified\n" );
2785 exit( 1 );
2786 }
2787 if ( highPassAtten < 0 )
2788 {
2789 fprintf( stderr, "--high-pass-attenuation must be specified\n" );
2790 exit( 1 );
2791 }
2792 }
2793
2795 {
2796 /* check that strain high pass parameters have been specified */
2797 if ( strainHighPassFreq < 0 )
2798 {
2799 fprintf( stderr,
2800 "--strain-high-pass-freq must be specified for REAL8 h(t) data\n" );
2801 exit( 1 );
2802 }
2803 if ( strainHighPassOrder < 0 )
2804 {
2805 fprintf( stderr,
2806 "--strain-high-pass-order must be specified for REAL8 h(t) data\n");
2807 exit( 1 );
2808 }
2809 if ( strainHighPassAtten < 0 )
2810 {
2811 fprintf( stderr,
2812 "--strain-high-pass-atten must be specified for REAL8 h(t) data\n");
2813 exit( 1 );
2814 }
2815 }
2816
2817 /* check validity of input data length */
2819 {
2821 (numPoints / 2);
2822 {
2823 UINT8 gpsChanIntervalNS = gpsEndTime.gpsSeconds * LAL_INT8_C(1000000000) -
2824 gpsStartTime.gpsSeconds * LAL_INT8_C(1000000000);
2825 UINT8 inputDataLengthNS = (UINT8) inputDataLength * LAL_INT8_C(1000000000)
2826 / (UINT8) sampleRate;
2827
2828 if ( inputDataLengthNS != gpsChanIntervalNS )
2829 {
2830 fprintf( stderr, "length of input data and data chunk do not match\n" );
2831 fprintf( stderr, "start time: %d, end time %d\n",
2833 fprintf( stderr, "gps channel time interval: %" LAL_UINT8_FORMAT " ns\n"
2834 "computed input data length: %" LAL_UINT8_FORMAT " ns\n",
2835 gpsChanIntervalNS, inputDataLengthNS );
2836 exit( 1 );
2837 }
2838 }
2839 }
2840
2841 /* check standard candle arguments */
2842 if ( computeCandle )
2843 {
2844 this_proc_param = this_proc_param->next = (ProcessParamsTable *)
2845 calloc( 1, sizeof(ProcessParamsTable) );
2846 snprintf( this_proc_param->program, LIGOMETA_PROGRAM_MAX,
2847 "%s", PROGRAM_NAME );
2848 snprintf( this_proc_param->param, LIGOMETA_PARAM_MAX,
2849 "--standard-candle" );
2850 snprintf( this_proc_param->type, LIGOMETA_TYPE_MAX, "string" );
2851 snprintf( this_proc_param->value, LIGOMETA_TYPE_MAX, " " );
2852
2853 if ( candleSnr < 0 )
2854 {
2855 fprintf( stderr,
2856 "--candle-snr must be specified if --standard-candle is given\n" );
2857 exit( 1 );
2858 }
2859 if ( candleMinMass < 0 )
2860 {
2861 fprintf( stderr,
2862 "--candle-minmass must be specified if --standard-candle is given\n" );
2863 exit( 1 );
2864 }
2865 if ( candleMaxMass < 0 )
2866 {
2867 fprintf( stderr,
2868 "--candle-maxmass must be specified if --standard-candle is given\n" );
2869 exit( 1 );
2870 }
2871 if ( ( specType == specType_simulated ) && numPoints < 0 )
2872 {
2873 fprintf( stderr,
2874 "--num-points must be specified if --standard-candle is given\n"
2875 "when using a simulated PSD\n" );
2876 exit( 1 );
2877 }
2878 }
2879
2880 /* check that the spectrum generation parameters have been given */
2881 if ( fLow < 0 )
2882 {
2883 fprintf( stderr, "--low-frequency-cutoff must be specified\n" );
2884 exit( 1 );
2885 }
2887 {
2888 fprintf( stderr, "--resample-filter must be specified for frame data\n" );
2889 exit( 1 );
2890 }
2892 {
2893 fprintf( stderr, "--spectrum-type must be specified\n" );
2894 exit( 1 );
2895 }
2896
2897 /* check for potential underflows in the simulated spectrum */
2898 if ( ( specType == specType_simulated ) && dynRangeExponent < 10 )
2899 {
2900 fprintf( stderr, "If using a simulated PSD, a suitable dynamic \n"
2901 "range exponent must be given, eg 69.0. Exiting...\n" );
2902 exit( 1 );
2903 }
2904
2905 /* check that a channel has been requested and fill the ifo */
2907 {
2908 fprintf( stderr, "--channel-name must be specified for frame data\n" );
2909 exit( 1 );
2910 }
2911
2912 /* check that we can correctly obtain the input frame data */
2914 {
2915 if ( globFrameData )
2916 {
2917 if ( frInCacheName )
2918 {
2919 fprintf( stderr,
2920 "--frame-cache must not be specified when globbing frame data\n" );
2921 exit( 1 );
2922 }
2923
2924 if ( ! frInType )
2925 {
2926 fprintf( stderr,
2927 "--frame-type must be specified when globbing frame data\n" );
2928 exit( 1 );
2929 }
2930 }
2931 else
2932 {
2933 if ( ! frInCacheName )
2934 {
2935 fprintf( stderr,
2936 "--frame-cache must be specified when not globbing frame data\n" );
2937 exit( 1 );
2938 }
2939
2940 if ( frInType )
2941 {
2942 fprintf( stderr, "--frame-type must not be specified when obtaining "
2943 "frame data from a cache file\n" );
2944 exit( 1 );
2945 }
2946 }
2947 }
2948
2949 /* record the glob frame data option in the process params */
2950 if ( globFrameData )
2951 {
2952 this_proc_param = this_proc_param->next = (ProcessParamsTable *)
2953 calloc( 1, sizeof(ProcessParamsTable) );
2954 snprintf( this_proc_param->program, LIGOMETA_PROGRAM_MAX,
2955 "%s", PROGRAM_NAME );
2956 snprintf( this_proc_param->param, LIGOMETA_PARAM_MAX,
2957 "--glob-frame-data" );
2958 snprintf( this_proc_param->type, LIGOMETA_TYPE_MAX, "string" );
2959 snprintf( this_proc_param->value, LIGOMETA_TYPE_MAX, " " );
2960 }
2961
2962 /* store point calibration option */
2963 if ( pointCal )
2964 {
2965 this_proc_param = this_proc_param->next = (ProcessParamsTable *)
2966 calloc( 1, sizeof(ProcessParamsTable) );
2967 snprintf( this_proc_param->program, LIGOMETA_PROGRAM_MAX,
2968 "%s", PROGRAM_NAME );
2969 snprintf( this_proc_param->param, LIGOMETA_PARAM_MAX,
2970 "--point-calibration" );
2971 snprintf( this_proc_param->type, LIGOMETA_TYPE_MAX, "string" );
2972 snprintf( this_proc_param->value, LIGOMETA_TYPE_MAX, " " );
2973 }
2974
2975 /* check we can calibrate the data if it's not h(t) */
2977 {
2978 if ( ! ( calCacheName || globCalData ) )
2979 {
2980 fprintf( stderr, "either --calibration-cache or "
2981 "--glob-calibration-data must be specified\n" );
2982 exit( 1 );
2983 }
2984 else if ( calCacheName && globCalData )
2985 {
2986 fprintf( stderr, "only one of --calibration-cache or "
2987 "--glob-calibration-data can be specified\n" );
2988 exit( 1 );
2989 }
2990 }
2991 else
2992 {
2993 if ( calCacheName || globCalData )
2994 {
2995 fprintf( stderr, "neither --calibration-cache nor --glob-calibration-data\n"
2996 "should be given when using calibrated data or a simulated PSD\n" );
2997 exit( 1 );
2998 }
2999 }
3000
3001 /* record the glob calibration data option in the process params */
3002 if ( globCalData )
3003 {
3004 this_proc_param = this_proc_param->next = (ProcessParamsTable *)
3005 calloc( 1, sizeof(ProcessParamsTable) );
3006 snprintf( this_proc_param->program, LIGOMETA_PROGRAM_MAX,
3007 "%s", PROGRAM_NAME );
3008 snprintf( this_proc_param->param, LIGOMETA_PARAM_MAX,
3009 "--glob-calibration-data" );
3010 snprintf( this_proc_param->type, LIGOMETA_TYPE_MAX, "string" );
3011 snprintf( this_proc_param->value, LIGOMETA_TYPE_MAX, " " );
3012 }
3013
3014 /* check that the bank type has been specified */
3015 if ( ! haveOrder )
3016 {
3017 fprintf( stderr, "--order must be specified\n" );
3018 exit( 1 );
3019 }
3020 if ( ! haveApprox )
3021 {
3022 fprintf( stderr, "--approximant must be specified\n" );
3023 exit( 1 );
3024 }
3025 if ( ! haveSpace )
3026 {
3027 fprintf( stderr, "--space must be specified\n" );
3028 exit( 1 );
3029 }
3030
3031 /* check validity of grid spacing with respect to approximant */
3032 if ( ! haveGridSpacing )
3033 {
3034 fprintf( stderr, "--grid-spacing must be specified\n" );
3035 exit( 1 );
3036 }
3038 {
3039 fprintf( stderr, "--grid-spacing must be either SquareNotOriented or Hexagonal\n" );
3040 exit( 1 );
3041 }
3042
3043 /* check that the correct range parameters have been given for the bank */
3044 if ( approximant == BCV )
3045 {
3046 if ( ! havePsi0Min )
3047 {
3048 fprintf( stderr, "--minimum-psi0 must be specified\n" );
3049 exit( 1 );
3050 }
3051 if ( ! havePsi0Max )
3052 {
3053 fprintf( stderr, "--maximum-psi0 must be specified\n" );
3054 exit( 1 );
3055 }
3056 if ( ! havePsi3Min )
3057 {
3058 fprintf( stderr, "--minimum-psi3 must be specified\n" );
3059 exit( 1 );
3060 }
3061 if ( ! havePsi3Max )
3062 {
3063 fprintf( stderr, "--maximum-psi3 must be specified\n" );
3064 exit( 1 );
3065 }
3066 if ( ! haveAlpha )
3067 {
3068 fprintf( stderr, "--alpha must be specified\n" );
3069 exit( 1 );
3070 }
3071 if ( maxFcutTmplts < 0 )
3072 {
3073 fprintf( stderr, "--maximum-fcut-tmplts must be specified\n" );
3074 exit( 1 );
3075 }
3076
3077 if ( psi3Max <= psi3Min )
3078 {
3079 fprintf( stdout, "invalid argument to --maximum-psi3:\n"
3080 "maximum value of psi3 must be greater than minimum value of psi3: "
3081 "(%f specified)\n",
3082 psi3Max );
3083 exit( 1 );
3084 }
3085
3086 minMass = maxMass = 0;
3087 }
3088 else if ( approximant == BCVSpin )
3089 {
3090 if ( minMass < 0 && ( ! havePsi0Min || ! havePsi3Min ) )
3091 {
3092 fprintf( stderr, "--minimum-mass or --minimum-psi0 and --minimum-psi3 "
3093 "must be specified\n" );
3094 exit( 1 );
3095 }
3096 if ( maxMass < 0 && ( ! havePsi0Max || ! havePsi3Max ) )
3097 {
3098 fprintf( stderr, "--maximum-mass or --maximum-psi0 and --maximum-psi3 "
3099 "must be specified\n" );
3100 exit( 1 );
3101 }
3102 }
3103 else
3104 {
3105 if ( minMass < 0 )
3106 {
3107 fprintf( stderr, "--minimum-mass must be specified\n" );
3108 exit( 1 );
3109 }
3110 if ( maxMass < 0 && maxTotalMass < 0 )
3111 {
3112 fprintf( stderr, "Either --maximum-mass or --max-total-mass must be specified\n" );
3113 exit( 1 );
3114 }
3115 if ( minTotalMass > 0 && maxTotalMass < 0 )
3116 {
3117 fprintf( stderr, "--max-total-mass must be specified with --min-total-mass\n" );
3118 exit( 1 );
3119 }
3120 }
3121
3122 /* check that the bank parameters have been specified */
3123 if ( minMatch < 0 )
3124 {
3125 fprintf( stderr, "--minimal-match must be specified\n" );
3126 exit( 1 );
3127 }
3128 if ( fUpper < 0 )
3129 {
3130 fprintf( stderr, "--high-frequency-cutoff must be specified\n" );
3131 exit( 1 );
3132 }
3133
3134 /* Check that multiple cutoff freq. options specified */
3135 if ( ! haveNumFcut )
3136 {
3137 fprintf( stderr, "must specify --num-freq-cutoffs\n" );
3138 exit( 1 );
3139 }
3140 if ( ! haveMaxFcut )
3141 {
3142 fprintf( stderr, "must specify --max-high-freq-cutoff\n" );
3143 exit( 1 );
3144 }
3145 if ( ! haveMinFcut )
3146 {
3147 fprintf( stderr, "must specify --min-high-freq-cutoff\n" );
3148 exit( 1 );
3149 }
3150 /* Check Min and Max upper freq. cuts are the same if NumFreqCut = 1 */
3151 if ( numFreqCut == 1 )
3152 {
3153 if( maxFreqCut < minFreqCut || maxFreqCut > minFreqCut )
3154 {
3155 fprintf(stderr, "--max-high-freq-cutoff must equal --min-high-freq-cutoff when --num-freq-cutoffs = 1\n" );
3156 exit( 1 );
3157 }
3158 }
3159
3161 {
3162 /* check max and mins are the correct way around */
3163 if (chiMin > chiMax )
3164 {
3165 fprintf( stderr,
3166 "Error: argument to --minimum-spin1 must be less than --maximum-spin1 .\n" );
3167 exit( 1 );
3168 }
3169
3170 /* check that kappa min-max are set correctly */
3171 if (kappaMin > kappaMax)
3172 {
3173 fprintf( stderr,
3174 "Error: argument to --minimum-kappa1 must be less than --maximum-kappa1 .\n" );
3175 exit( 1 );
3176 }
3177 }
3178
3179 if( etaMaxCutoff > 0 || etaMinCutoff >= 0 )
3180 {
3181 if( etaMaxCutoff <= 0 )
3182 {
3183 fprintf( stderr,
3184 "Error: argument --max-eta must be given if --min-eta is given\n");
3185 exit(1);
3186 }
3187
3188 if( etaMinCutoff < 0 )
3189 {
3190 fprintf( stderr,
3191 "Error: argument --min-eta must be given if --max-eta is given\n");
3192 exit(1);
3193 }
3194
3196 {
3197 fprintf( stderr,
3198 "Error: value for --max-eta must be greater than or equal to value for --min-eta\n");
3199 exit(1);
3200 }
3201 }
3202
3203 return 0;
3204}
3205
3206#undef ADD_PROCESS_PARAM
const LALVCSInfoList lalAppsVCSInfoList
NULL-terminated list of VCS and build information for LALApps and its dependencies
const LALVCSInfo lalAppsVCSIdentInfo
Identable VCS and build information for LALApps.
lal_errhandler_t lal_errhandler
int LAL_ERR_EXIT(LALStatus *stat, const char *func, const char *file, const int line, volatile const char *id)
#define LAL_CALL(function, statusptr)
void LALFrSeek(LALStatus *status, const LIGOTimeGPS *epoch, LALFrStream *stream)
void LALFrGetREAL8TimeSeries(LALStatus *status, REAL8TimeSeries *series, FrChanIn *chanin, LALFrStream *stream)
void LALFrCacheOpen(LALStatus *status, LALFrStream **output, LALCache *cache)
void LALFrGetREAL4TimeSeries(LALStatus *status, REAL4TimeSeries *series, FrChanIn *chanin, LALFrStream *stream)
void LALFrClose(LALStatus *status, LALFrStream **stream)
int j
int k
void LALCheckMemoryLeaks(void)
#define LALFree(p)
#define LAL_EXLAL
int LALgetopt_long_only(int argc, char *const *argv, const char *options, const struct LALoption *long_options, int *opt_index)
int LALoptind
char * LALoptarg
#define no_argument
#define required_argument
int XLALCloseLIGOLwXMLFile(LIGOLwXMLStream *xml)
LIGOLwXMLStream * XLALOpenLIGOLwXMLFile(const char *path)
int XLALWriteLIGOLwXMLProcessTable(LIGOLwXMLStream *, const ProcessTable *)
int XLALWriteLIGOLwXMLProcessParamsTable(LIGOLwXMLStream *, const ProcessParamsTable *)
int XLALWriteLIGOLwXMLSnglInspiralTable(LIGOLwXMLStream *xml, const SnglInspiralTable *sngl_inspiral)
int XLALWriteLIGOLwXMLSearchSummaryTable(LIGOLwXMLStream *, const SearchSummaryTable *)
SnglInspiralTable * XLALSnglInspiralTableFromLIGOLw(const char *fileName)
SnglInspiralTable * XLALTimeCutSingleInspiral(SnglInspiralTable *eventList, LIGOTimeGPS *startTime, LIGOTimeGPS *endTime)
SnglInspiralTable * XLALIfoCutSingleInspiral(SnglInspiralTable **eventHead, char *ifo)
INT4 XLALCountSnglInspiral(SnglInspiralTable *head)
SnglInspiralTable * XLALMassCut(SnglInspiralTable *eventHead, const char *massCut, REAL4 massRangeLow, REAL4 massRangeHigh, REAL4 mass2RangeLow, REAL4 mass2RangeHigh)
#define LIGOMETA_SEARCH_MAX
#define LIGOMETA_COMMENT_MAX
#define LIGOMETA_TYPE_MAX
#define LIGOMETA_CHANNEL_MAX
#define LIGOMETA_PROGRAM_MAX
#define LIGOMETA_VALUE_MAX
#define LIGOMETA_PARAM_MAX
#define LIGOMETA_IFO_MAX
void XLALDestroyProcessParamsTable(ProcessParamsTable *)
int XLALPopulateProcessTable(ProcessTable *ptable, const char *program_name, const char *cvs_revision, const char *cvs_source, const char *cvs_date, long process_id)
void XLALDestroySnglInspiralTable(SnglInspiralTable *head)
#define fprintf
REAL8 dynRange
Definition: blindinj.c:122
static LALUnit strainPerCount
Definition: getresp.c:41
void LALButterworthREAL8TimeSeries(LALStatus *stat, REAL8TimeSeries *series, PassBandParamStruc *params)
void LALDButterworthREAL4TimeSeries(LALStatus *stat, REAL4TimeSeries *series, PassBandParamStruc *params)
LIGOTimeGPS * XLALGPSTimeNow(LIGOTimeGPS *gpstime)
void XLALDestroyCache(LALCache *cache)
LALCache * XLALCacheGlob(const char *dirstr, const char *fnptrn)
LALCache * XLALCacheImport(const char *fname)
LALCache * XLALCacheDuplicate(const LALCache *cache)
int XLALCacheSieve(LALCache *cache, INT4 t0, INT4 t1, const char *srcregex, const char *dscregex, const char *urlregex)
#define XLAL_NUM_ELEM(x)
uint64_t UINT8
#define LAL_INT8_C(c)
double REAL8
#define XLAL_INIT_DECL(var,...)
char CHAR
uint32_t UINT4
int32_t INT4
float REAL4
LALNameLength
int XLALFrStreamSetMode(LALFrStream *stream, int mode)
LAL_FR_STREAM_CHECKSUM_MODE
LAL_FR_STREAM_VERBOSE_MODE
GridSpacing
void LALInspiralBankGeneration(LALStatus *status, InspiralCoarseBankIn *in, SnglInspiralTable **out, INT4 *count)
FreqCut
CoordinateSpace
MinMaxComponentTotalMass
MinMaxComponentMass
MinComponentMassMaxTotalMass
SquareNotOriented
Hexagonal
FreqCut_SchwarzISCO
FreqCut_FRD
FreqCut_LightRing
FreqCut_ERD
FreqCut_BKLISCO
FreqCut_LRD
Tau0Tau2
Psi0Psi3
Tau0Tau3
void * XLALRealloc(void *p, size_t n)
Approximant
LALPNOrder
EOB
PadeT1
BCVSpin
EOBNRv2
EOBNR
PadeF1
IMRPhenomA
TaylorF1
TaylorT3
FindChirpPTF
TaylorF2
IMRPhenomB
BCV
TaylorT1
SpinTaylorT3
TaylorT2
LAL_PNORDER_TWO_POINT_FIVE
LAL_PNORDER_THREE
LAL_PNORDER_TWO
LAL_PNORDER_ONE
LAL_PNORDER_THREE_POINT_FIVE
LAL_PNORDER_HALF
LAL_PNORDER_ONE_POINT_FIVE
LAL_PNORDER_NEWTONIAN
double XLALSimNoisePSDGEO(double f)
double XLALSimNoisePSDaLIGOZeroDetHighPower(double f)
double XLALSimNoisePSDaLIGONoSRMHighPower(double f)
double XLALSimNoisePSDeLIGOModel(double f)
double XLALSimNoisePSDaLIGONoSRMLowPower(double f)
double XLALSimNoisePSDAdvVirgo(double f)
double XLALSimNoisePSDaLIGOZeroDetLowPower(double f)
double XLALSimNoisePSDiLIGOSRD(double f)
double XLALSimNoisePSDKAGRA(double f)
double XLALSimNoisePSDVirgo(double f)
#define LAL_UINT8_FORMAT
int XLALOutputVCSInfo(FILE *fp, const LALVCSInfoList vcs_list, const int verbose, const char *prefix)
void LALDPrintFrequencySeries(REAL8FrequencySeries *series, const CHAR *filename)
REAL4FFTPlan * XLALCreateForwardREAL4FFTPlan(UINT4 size, int measurelvl)
void XLALDestroyREAL4FFTPlan(REAL4FFTPlan *plan)
void LALResampleREAL4TimeSeries(LALStatus *status, REAL4TimeSeries *ts, ResampleTSParams *params)
LDASfirLP
defaultButterworth
int XLALREAL4AverageSpectrumMedian(REAL4FrequencySeries *spectrum, const REAL4TimeSeries *tseries, UINT4 seglen, UINT4 stride, const REAL4Window *window, const REAL4FFTPlan *plan)
int XLALREAL4AverageSpectrumWelch(REAL4FrequencySeries *spectrum, const REAL4TimeSeries *tseries, UINT4 seglen, UINT4 stride, const REAL4Window *window, const REAL4FFTPlan *plan)
const LALUnit lalADCCountUnit
LALUnit * XLALUnitMultiply(LALUnit *output, const LALUnit *unit1, const LALUnit *unit2)
void LALCCreateVector(LALStatus *, COMPLEX8Vector **, UINT4)
void LALDCreateVector(LALStatus *, REAL8Vector **, UINT4)
void LALCDestroyVector(LALStatus *, COMPLEX8Vector **)
void LALDDestroyVector(LALStatus *, REAL8Vector **)
void LALSDestroyVector(LALStatus *, REAL4Vector **)
void LALSCreateVector(LALStatus *, REAL4Vector **, UINT4)
void XLALDestroyREAL4Window(REAL4Window *window)
REAL4Window * XLALCreateHannREAL4Window(UINT4 length)
LIGOTimeGPS * XLALGPSSetREAL8(LIGOTimeGPS *epoch, REAL8 t)
REAL8 XLALGPSGetREAL8(const LIGOTimeGPS *epoch)
INT8 inputLengthNS
Definition: inspfrinj.c:317
char name[LIGOMETA_SOURCE_MAX]
Definition: inspinj.c:561
static LALStatus status
Definition: inspinj.c:552
int i
Definition: inspinj.c:596
FrameH * fr_add_proc_REAL8FrequencySeries(FrameH *frame, REAL8FrequencySeries *chan, const char *unit, const char *suffix)
FrameH * fr_add_proc_COMPLEX8FrequencySeries(FrameH *frame, COMPLEX8FrequencySeries *chan, const char *unit, const char *suffix)
FrameH * fr_add_proc_REAL4TimeSeries(FrameH *frame, REAL4TimeSeries *chan, const char *unit, const char *suffix)
FrameH * fr_add_proc_REAL4FrequencySeries(FrameH *frame, REAL4FrequencySeries *chan, const char *unit, const char *suffix)
searchsumm
c
float atol
CHAR fname[256]
Definition: spininj.c:211
CHAR name[LALNameLength]
COMPLEX8Sequence * data
COMPLEX8 * data
const CHAR * name
InspiralBankMassRange massRange
GridSpacing gridSpacing
ComputeMoments computeMoments
REAL8FrequencySeries shf
InsidePolygonEnum insidePolygon
Approximant approximant
CoordinateSpace space
UINT4 length
const char *const vcsDate
const char *const vcsStatus
const char *const vcsId
int * flag
INT4 gpsNanoSeconds
CHAR type[LIGOMETA_TYPE_MAX]
CHAR param[LIGOMETA_PARAM_MAX]
CHAR value[LIGOMETA_VALUE_MAX]
struct tagProcessParamsTable * next
CHAR program[LIGOMETA_PROGRAM_MAX]
LIGOTimeGPS start_time
CHAR ifos[LIGOMETA_IFOS_MAX]
LIGOTimeGPS end_time
CHAR comment[LIGOMETA_COMMENT_MAX]
REAL4Sequence * data
CHAR name[LALNameLength]
CHAR name[LALNameLength]
REAL4Sequence * data
LALUnit sampleUnits
LIGOTimeGPS epoch
REAL4 * data
REAL8Sequence * data
CHAR name[LALNameLength]
REAL8Sequence * data
LALUnit sampleUnits
LIGOTimeGPS epoch
CHAR name[LALNameLength]
REAL8 * data
ResampleTSFilter filterType
CHAR ifo[LIGOMETA_IFO_MAX]
struct tagSnglInspiralTable * next
CHAR channel[LIGOMETA_CHANNEL_MAX]
CHAR search[LIGOMETA_SEARCH_MAX]
REAL4 kappaMin
Definition: tmpltbank.c:440
@ real_4
Definition: tmpltbank.c:365
@ real_8
Definition: tmpltbank.c:366
@ undefined
Definition: tmpltbank.c:364
REAL4 dynRangeExponent
Definition: tmpltbank.c:414
int main(int argc, char *argv[])
Definition: tmpltbank.c:478
REAL4 fLow
Definition: tmpltbank.c:410
REAL4 alpha
Definition: tmpltbank.c:432
int writeStrainSpec
Definition: tmpltbank.c:472
INT4 numSegments
Definition: tmpltbank.c:400
INT4 nPointsChi
Definition: tmpltbank.c:442
Approximant approximant
Definition: tmpltbank.c:445
INT4 maxFcutTmplts
Definition: tmpltbank.c:435
int writeRawData
Definition: tmpltbank.c:469
INT4 numTDFiles
Definition: tmpltbank.c:464
REAL4 fUpper
Definition: tmpltbank.c:437
REAL4 candleSnr
Definition: tmpltbank.c:458
GridSpacing gridSpacing
Definition: tmpltbank.c:453
LALPNOrder order
Definition: tmpltbank.c:444
INT4 strainHighPassOrder
Definition: tmpltbank.c:416
INT4 computeCandle
Definition: tmpltbank.c:457
REAL4 psi3Min
Definition: tmpltbank.c:430
REAL4 candleMaxMass
Definition: tmpltbank.c:460
#define PROGRAM_NAME
Definition: tmpltbank.c:357
REAL4 minMatch
Definition: tmpltbank.c:436
CHAR * calCacheName
Definition: tmpltbank.c:411
REAL4 strainHighPassFreq
Definition: tmpltbank.c:415
#define USAGE(a)
Definition: tmpltbank.c:1343
CHAR * channelName
Definition: tmpltbank.c:402
REAL4 chirpMassCutoff
Definition: tmpltbank.c:425
INT4 computeMoments
Definition: tmpltbank.c:448
REAL8(* specFunc)(REAL8)
Definition: tmpltbank.c:418
CHAR * ifoTag
Definition: tmpltbank.c:468
INT4 globFrameData
Definition: tmpltbank.c:396
REAL4 psi0Max
Definition: tmpltbank.c:429
INT4 padData
Definition: tmpltbank.c:394
int polygonFit
Definition: tmpltbank.c:454
REAL4 psi0Min
Definition: tmpltbank.c:428
REAL4 highPassAtten
Definition: tmpltbank.c:409
INT4 outCompress
Definition: tmpltbank.c:473
INT4 numFreqCut
Definition: tmpltbank.c:451
INT4 sampleRate
Definition: tmpltbank.c:405
REAL4 highPassFreq
Definition: tmpltbank.c:407
#define ADD_PROCESS_PARAM(pptype, format, ppvalue)
Definition: tmpltbank.c:1333
REAL4 minMass
Definition: tmpltbank.c:421
REAL4 kappaMax
Definition: tmpltbank.c:441
LIGOTimeGPS gpsStartTime
Definition: tmpltbank.c:392
REAL4 betaMax
Definition: tmpltbank.c:434
INT4 nPointsKappa
Definition: tmpltbank.c:443
INT4 resampFiltType
Definition: tmpltbank.c:404
int vrbflg
defined in lal/lib/std/LALError.c
CHAR comment[LIGOMETA_COMMENT_MAX]
Definition: tmpltbank.c:476
REAL4 maxTotalMass
Definition: tmpltbank.c:424
CHAR ifo[3]
Definition: tmpltbank.c:401
INT4 globCalData
Definition: tmpltbank.c:412
CHAR ** tdFileNames
Definition: tmpltbank.c:463
int writeResponse
Definition: tmpltbank.c:470
CHAR * fqChanName
Definition: tmpltbank.c:395
INT4 pointCal
Definition: tmpltbank.c:413
@ specType_median
Definition: tmpltbank.c:376
@ specType_mean
Definition: tmpltbank.c:375
@ specType_undefined
Definition: tmpltbank.c:378
@ specType_simulated
Definition: tmpltbank.c:377
INT4 highPassOrder
Definition: tmpltbank.c:408
REAL4 psi3Max
Definition: tmpltbank.c:431
REAL4 betaMin
Definition: tmpltbank.c:433
INT4 haveGridSpacing
Definition: tmpltbank.c:447
enum @7 calData
enum @8 specType
REAL4 maxMass
Definition: tmpltbank.c:422
INT4 numPoints
Definition: tmpltbank.c:399
CHAR * frInCacheName
Definition: tmpltbank.c:397
REAL4 candleMinMass
Definition: tmpltbank.c:459
REAL4 chiMin
Definition: tmpltbank.c:438
int writeSpectrum
Definition: tmpltbank.c:471
REAL4 etaMaxCutoff
Definition: tmpltbank.c:427
CoordinateSpace space
Definition: tmpltbank.c:446
LIGOTimeGPS gpsEndTime
Definition: tmpltbank.c:393
REAL4 etaMinCutoff
Definition: tmpltbank.c:426
CHAR * frInType
Definition: tmpltbank.c:398
FreqCut maxFreqCut
Definition: tmpltbank.c:449
INT4 highPass
Definition: tmpltbank.c:406
INT4 inputDataLength
Definition: tmpltbank.c:403
CHAR * userTag
Definition: tmpltbank.c:467
FreqCut minFreqCut
Definition: tmpltbank.c:450
REAL4 chiMax
Definition: tmpltbank.c:439
REAL4 minTotalMass
Definition: tmpltbank.c:423
REAL4 strainHighPassAtten
Definition: tmpltbank.c:417
int arg_parse_check(int argc, char *argv[], ProcessParamsTable *procparams)
Definition: tmpltbank.c:1447