binj.c

Go to the documentation of this file.
00001 /*
00002  * Copyright (C) 2007 Kipp Cannon, Lisa M. Goggin, Patrick Brady, Saikat
00003  * Ray-Majumder, Xavier Siemens
00004  *
00005  * This program is free software; you can redistribute it and/or modify it
00006  * under the terms of the GNU General Public License as published by the
00007  * Free Software Foundation; either version 2 of the License, or (at your
00008  * option) any later version.
00009  *
00010  * This program is distributed in the hope that it will be useful, but
00011  * WITHOUT ANY WARRANTY; without even the implied warranty of
00012  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00013  * General Public License for more details.
00014  *
00015  * You should have received a copy of the GNU General Public License along
00016  * with with program; see the file COPYING. If not, write to the Free
00017  * Software Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA
00018  * 02111-1307  USA
00019  */
00020 
00021 
00022 /*
00023  * ============================================================================
00024  *
00025  *                                  Preamble
00026  *
00027  * ============================================================================
00028  */
00029 
00030 
00031 #include <ctype.h>
00032 #include <stdio.h>
00033 #include <stdlib.h>
00034 #include <string.h>
00035 #include <getopt.h>
00036 #include <unistd.h>
00037 #include <time.h>
00038 #include <math.h>
00039 #include <limits.h>
00040 #include <gsl/gsl_rng.h>
00041 #include <gsl/gsl_randist.h>
00042 
00043 
00044 #include <lal/Date.h>
00045 #include <lal/GenerateBurst.h>
00046 #include <lal/LALConstants.h>
00047 #include <lal/LALSimBurst.h>
00048 #include <lal/LALStdio.h>
00049 #include <lal/LALStdlib.h>
00050 #include <lal/LIGOLwXML.h>
00051 #include <lal/LIGOMetadataTables.h>
00052 #include <lal/LIGOMetadataUtils.h>
00053 #include <lal/TimeDelay.h>
00054 #include <lal/TimeSeries.h>
00055 #include <lal/XLALError.h>
00056 
00057 
00058 #include <lalapps.h>
00059 #include <processtable.h>
00060 
00061 
00062 int snprintf(char *str, size_t size, const char *format, ...);
00063 
00064 
00065 RCSID("$Id: binj.c,v 1.69 2008/04/26 02:27:33 kipp Exp $");
00066 
00067 
00068 #define CVS_REVISION "$Revision: 1.69 $"
00069 #define CVS_SOURCE "$Source: /usr/local/cvs/lscsoft/lalapps/src/power/binj.c,v $"
00070 #define CVS_DATE "$Date: 2008/04/26 02:27:33 $"
00071 #define PROGRAM_NAME "lalapps_binj"
00072 
00073 
00074 /*
00075  * ============================================================================
00076  *
00077  *                                Command Line
00078  *
00079  * ============================================================================
00080  */
00081 
00082 
00083 enum population {
00084         POPULATION_TARGETED,
00085         POPULATION_ALL_SKY_SINEGAUSSIAN,
00086         POPULATION_STRING_CUSP
00087 };
00088 
00089 
00090 struct options {
00091         INT8 gps_start_time;
00092         INT8 gps_end_time;
00093         enum population population;
00094         double ra;
00095         double dec;
00096         double maxA;
00097         double minA;
00098         double maxbandwidth;
00099         double minbandwidth;
00100         double maxduration;
00101         double minduration;
00102         double maxEoverr2;
00103         double minEoverr2;
00104         double maxf;
00105         double minf;
00106         double maxhrss;
00107         double minhrss;
00108         double q;
00109         unsigned long seed;
00110         double time_step;
00111         char *user_tag;
00112 };
00113 
00114 
00115 static struct options options_defaults(void)
00116 {
00117         struct options defaults;
00118         
00119         defaults.gps_start_time = 0;
00120         defaults.gps_end_time = 0;
00121         defaults.population = -1;
00122         defaults.ra = XLAL_REAL8_FAIL_NAN;
00123         defaults.dec = XLAL_REAL8_FAIL_NAN;
00124         defaults.maxbandwidth = XLAL_REAL8_FAIL_NAN;
00125         defaults.minbandwidth = XLAL_REAL8_FAIL_NAN;
00126         defaults.maxduration = XLAL_REAL8_FAIL_NAN;
00127         defaults.minduration = XLAL_REAL8_FAIL_NAN;
00128         defaults.maxEoverr2 = XLAL_REAL8_FAIL_NAN;
00129         defaults.minEoverr2 = XLAL_REAL8_FAIL_NAN;
00130         defaults.maxf = XLAL_REAL8_FAIL_NAN;
00131         defaults.minf = XLAL_REAL8_FAIL_NAN;
00132         defaults.maxhrss = XLAL_REAL8_FAIL_NAN;
00133         defaults.minhrss = XLAL_REAL8_FAIL_NAN;
00134         defaults.minA = XLAL_REAL8_FAIL_NAN;
00135         defaults.maxA = XLAL_REAL8_FAIL_NAN;
00136         defaults.q = XLAL_REAL8_FAIL_NAN;
00137         defaults.seed = 0;
00138         defaults.time_step = 210.0 / LAL_PI;
00139         defaults.user_tag = NULL;
00140 
00141         return defaults;
00142 }
00143 
00144 
00145 static void print_usage(void)
00146 {
00147         fprintf(stderr, 
00148 "lalapps_binj [options]\n" 00149 "\n" 00150 "Options:\n" 00151 "\n" 00152 "--gps-end-time seconds\n" 00153 "--gps-start-time seconds\n" 00154 "       Bounds of interval in which to synthesize injections.\n" 00155 "\n" 00156 "--help\n" 00157 "       display this message\n" 00158 "\n" 00159 "--max-amplitude value\n" 00160 "--min-amplitude value\n" 00161 "       Set the bounds of the injection ampltiudes.  These only affect\n" 00162 "       string cusp injections.\n" 00163 "\n" 00164 "--max-bandwidth hertz\n" 00165 "--min-bandwidth hertz\n" 00166 "       Set the bounds of the injection bandwidthds.  These only affect\n" 00167 "       btlwnb waveforms.\n" 00168 "\n" 00169         ); fprintf(stderr, 
00170 "--max-duration seconds\n" 00171 "--min-duration seconds\n" 00172 "       Set the bounds of the injection durations.  These only affect\n" 00173 "       btlwnb waveforms.\n" 00174 "\n" 00175 "--max-e-over-r2 value\n" 00176 "--min-e-over-r2 value\n" 00177 "       Set the bounds of the range of equivalent isotropic radiated\n" 00178 "       energies of btlwnb waveforms.  The units are M_{sun} / pc^{2} (solar\n" 00179 "       masses per parsec^{2}).\n" 00180 "\n" 00181         ); fprintf(stderr, 
00182 "--max-frequency hertz\n" 00183 "--min-frequency hertz\n" 00184 "       Set the bounds of the injection frequencies.  These are the centre\n" 00185 "       frequencies of sine-Gaussians and btlwnb waveforms, and the\n" 00186 "       high-frequency cut-offs of string cusp waveforms.\n" 00187 "\n" 00188 "--max-hrss value\n" 00189 "--min-hrss value\n" 00190         ); fprintf(stderr, 
00191 "       Set the bounds of the injection h_{rss} values.  These only affect\n" 00192 "       sine-Gaussian injections.  (Actually, these set the bounds of the\n" 00193 "       product of the waveform's hrss and its duration, which makes the\n" 00194 "       injections lie along the 50 efficiency curve better.  To convert to\n" 00195 "       real hrss multiply by sqrt(2) pi f/Q.) \n" 00196 "\n" 00197         ); fprintf(stderr, 
00198 "--population name\n" 00199 "       Select the injection population to synthesize.  Allowed values are\n" 00200 "       \"targeted\", \"string_cusp\", and \"all_sky_sinegaussian\".\n" 00201 "\n" 00202 "--q value\n" 00203 "       Set the Q for sine-Gaussian injections.\n" 00204 "\n" 00205 "--ra-dec ra,dec\n" 00206 "       Set the right-ascension and declination of the sky location from which\n" 00207 "       injections should originate when generating a targeted population.\n" 00208 "       Co-ordinates are in radians.\n" 00209 "\n" 00210         ); fprintf(stderr, 
00211 "--seed value\n" 00212 "       Set the random number generator's seed (0 = off, default = 0).\n" 00213 "\n" 00214 "--time-step value\n" 00215 "       Set the time betwen injections to value/pi seconds (default = 210).\n" 00216 "\n" 00217 "--user-tag string\n" 00218 "       Set the user tag in the process and search summary tables to this.\n"
00219         );
00220 }
00221 
00222 
00223 static ProcessParamsTable **add_process_param(ProcessParamsTable **proc_param, const ProcessTable *process, const char *type, const char *param, const char *value)
00224 {
00225         *proc_param = XLALCreateProcessParamsTableRow(process);
00226         snprintf((*proc_param)->program, sizeof((*proc_param)->program), PROGRAM_NAME);
00227         snprintf((*proc_param)->type, sizeof((*proc_param)->type), type);
00228         snprintf((*proc_param)->param, sizeof((*proc_param)->param), "--%s", param);
00229         snprintf((*proc_param)->value, sizeof((*proc_param)->value), value);
00230 
00231         return(&(*proc_param)->next);
00232 }
00233 
00234         
00235 #define ADD_PROCESS_PARAM(process, type) \
00236         do { paramaddpoint = add_process_param(paramaddpoint, process, type, long_options[option_index].name, optarg); } while(0)
00237 
00238 
00239 static struct options parse_command_line(int *argc, char **argv[], const ProcessTable *process, ProcessParamsTable **paramaddpoint)
00240 {
00241         struct options options = options_defaults();
00242         int c;
00243         int option_index;
00244         struct option long_options[] = {
00245                 {"gps-end-time", required_argument, NULL, 'A'},
00246                 {"gps-start-time", required_argument, NULL, 'B'},
00247                 {"help", no_argument, NULL, 'C'},
00248                 {"max-amplitude", required_argument, NULL, 'D'},
00249                 {"min-amplitude", required_argument, NULL, 'E'},
00250                 {"max-bandwidth", required_argument, NULL, 'F'},
00251                 {"min-bandwidth", required_argument, NULL, 'G'},
00252                 {"max-duration", required_argument, NULL, 'H'},
00253                 {"min-duration", required_argument, NULL, 'I'},
00254                 {"max-e-over-r2", required_argument, NULL, 'S'},
00255                 {"min-e-over-r2", required_argument, NULL, 'T'},
00256                 {"max-frequency", required_argument, NULL, 'J'},
00257                 {"min-frequency", required_argument, NULL, 'K'},
00258                 {"max-hrss", required_argument, NULL, 'L'},
00259                 {"min-hrss", required_argument, NULL, 'M'},
00260                 {"population", required_argument, NULL, 'N'},
00261                 {"q", required_argument, NULL, 'O'},
00262                 {"ra-dec", required_argument, NULL, 'U'},
00263                 {"seed", required_argument, NULL, 'P'},
00264                 {"time-step", required_argument, NULL, 'Q'},
00265                 {"user-tag", required_argument, NULL, 'R'},
00266                 {NULL, 0, NULL, 0}
00267         };
00268 
00269         do switch(c = getopt_long(*argc, *argv, "", long_options, &option_index)) {
00270         case 'A':
00271                 XLALClearErrno();
00272                 {
00273                         LIGOTimeGPS tmp;
00274                         XLALStrToGPS(&tmp, optarg, NULL);
00275                         options.gps_end_time = XLALGPSToINT8NS(&tmp);
00276                 }
00277                 if(xlalErrno) {
00278                         fprintf(stderr, "invalid --%s (%s specified)\n", long_options[option_index].name, optarg);
00279                         exit(1);
00280                 }
00281                 ADD_PROCESS_PARAM(process, "lstring");
00282                 break;
00283 
00284         case 'B':
00285                 XLALClearErrno();
00286                 {
00287                         LIGOTimeGPS tmp;
00288                         XLALStrToGPS(&tmp, optarg, NULL);
00289                         options.gps_start_time = XLALGPSToINT8NS(&tmp);
00290                 }
00291                 if(xlalErrno) {
00292                         fprintf(stderr, "invalid --%s (%s specified)\n", long_options[option_index].name, optarg);
00293                         exit(1);
00294                 }
00295                 ADD_PROCESS_PARAM(process, "lstring");
00296                 break;
00297 
00298         case 'C':
00299                 print_usage();
00300                 exit(0);
00301 
00302         case 'D':
00303                 options.maxA = atof(optarg);
00304                 ADD_PROCESS_PARAM(process, "real_8");
00305                 break;
00306 
00307         case 'E':
00308                 options.minA = atof(optarg);
00309                 ADD_PROCESS_PARAM(process, "real_8");
00310                 break;
00311 
00312         case 'F':
00313                 options.maxbandwidth = atof(optarg);
00314                 ADD_PROCESS_PARAM(process, "real_8");
00315                 break;
00316 
00317         case 'G':
00318                 options.minbandwidth = atof(optarg);
00319                 ADD_PROCESS_PARAM(process, "real_8");
00320                 break;
00321 
00322         case 'H':
00323                 options.maxduration = atof(optarg);
00324                 ADD_PROCESS_PARAM(process, "real_8");
00325                 break;
00326 
00327         case 'I':
00328                 options.minduration = atof(optarg);
00329                 ADD_PROCESS_PARAM(process, "real_8");
00330                 break;
00331 
00332         case 'J':
00333                 options.maxf = atof(optarg);
00334                 ADD_PROCESS_PARAM(process, "real_8");
00335                 break;
00336 
00337         case 'K':
00338                 options.minf = atof(optarg);
00339                 ADD_PROCESS_PARAM(process, "real_8");
00340                 break;
00341 
00342         case 'L':
00343                 options.maxhrss = atof(optarg);
00344                 ADD_PROCESS_PARAM(process, "real_8");
00345                 break;
00346 
00347         case 'M':
00348                 options.minhrss = atof(optarg);
00349                 ADD_PROCESS_PARAM(process, "real_8");
00350                 break;
00351 
00352         case 'N':
00353                 if(!strcmp(optarg, "targeted"))
00354                         options.population = POPULATION_TARGETED;
00355                 else if(!strcmp(optarg, "string_cusp"))
00356                         options.population = POPULATION_STRING_CUSP;
00357                 else if(!strcmp(optarg, "all_sky_sinegaussian"))
00358                         options.population = POPULATION_ALL_SKY_SINEGAUSSIAN;
00359                 else {
00360                         fprintf(stderr, "error: unrecognized population \"%s\"", optarg);
00361                         exit(1);
00362                 }
00363                 ADD_PROCESS_PARAM(process, "lstring");
00364                 break;
00365 
00366         case 'O':
00367                 options.q = atof(optarg);
00368                 ADD_PROCESS_PARAM(process, "real_8");
00369                 break;
00370 
00371         case 'P':
00372                 options.seed = atol(optarg);
00373                 ADD_PROCESS_PARAM(process, "int_8u");
00374                 break;
00375 
00376         case 'Q':
00377                 options.time_step = atof(optarg) / LAL_PI;
00378                 ADD_PROCESS_PARAM(process, "real_8");
00379                 break;
00380 
00381         case 'R':
00382                 options.user_tag = optarg;
00383                 ADD_PROCESS_PARAM(process, "lstring");
00384                 break;
00385 
00386         case 'S':
00387                 options.maxEoverr2 = atof(optarg);
00388                 ADD_PROCESS_PARAM(process, "real_8");
00389                 break;
00390 
00391         case 'T':
00392                 options.minEoverr2 = atof(optarg);
00393                 ADD_PROCESS_PARAM(process, "real_8");
00394                 break;
00395 
00396         case 'U':
00397                 {
00398                         char *end;
00399                         options.ra = strtod(optarg, &end);
00400                         while(isspace(*end))
00401                                 end++;
00402                         if(*end != ',') {
00403                                 fprintf(stderr, "error: cannot parse --ra-dec \"%s\"\n", optarg);
00404                                 exit(1);
00405                         }
00406                         options.dec = strtod(end + 1, &end);
00407                         while(isspace(*end))
00408                                 end++;
00409                         if(*end != '\0') {
00410                                 fprintf(stderr, "error: cannot parse --ra-dec \"%s\"\n", optarg);
00411                                 exit(1);
00412                         }
00413                 }
00414                 ADD_PROCESS_PARAM(process, "lstring");
00415                 break;
00416 
00417         case 0:
00418                 /* option sets a flag */
00419                 break;
00420 
00421         case -1:
00422                 /* end of arguments */
00423                 break;
00424 
00425         case '?':
00426                 /* unrecognized option */
00427                 print_usage();
00428                 exit(1);
00429 
00430         case ':':
00431                 /* missing argument for an option */
00432                 print_usage();
00433                 exit(1);
00434         } while(c != -1);
00435 
00436         /* check some of the input parameters for consistency */
00437         if(options.maxA < options.minA) {
00438                 fprintf(stderr, "error: --max-amplitude < --min-amplitude\n");
00439                 exit(1);
00440         }
00441         if(options.maxbandwidth < options.minbandwidth) {
00442                 fprintf(stderr, "error: --max-bandwidth < --min-bandwidth\n");
00443                 exit(1);
00444         }
00445         if(options.maxduration < options.minduration) {
00446                 fprintf(stderr, "error: --max-duration < --min-duration\n");
00447                 exit(1);
00448         }
00449         if(options.maxf < options.minf) {
00450                 fprintf(stderr, "error: --max-frequency < --min-frequency\n");
00451                 exit(1);
00452         }
00453         if(options.maxhrss < options.minhrss) {
00454                 fprintf(stderr, "error: --max-hrss < --min-hrss\n");
00455                 exit(1);
00456         }
00457 
00458         if(!options.gps_start_time || !options.gps_end_time) {
00459                 fprintf(stderr, "--gps-start-time and --gps-end-time are both required\n");
00460                 exit(1);
00461         }
00462         if(options.gps_end_time < options.gps_start_time) {
00463                 fprintf(stderr, "error: --gps-end-time < --gps-start-time\n");
00464                 exit(1);
00465         }
00466 
00467         return options;
00468 }
00469 
00470 
00471 /* 
00472  * ============================================================================
00473  *
00474  *                                 Sequences
00475  *
00476  * ============================================================================
00477  */
00478 
00479 
00480 /*
00481  * Repeating arithmetic sequence.
00482  */
00483 
00484 
00485 static double sequence_arithmetic_next(double low, double high, double delta)
00486 {
00487         static int i = 0;
00488         double x = low + delta * i++;
00489 
00490         if(high < low || high - low < delta)
00491                 return XLAL_REAL8_FAIL_NAN;
00492 
00493         if(x > high) {
00494                 i = 1;
00495                 return low;
00496         }
00497 
00498         return x;
00499 }
00500 
00501 
00502 /*
00503  * Repeating geometric sequence.
00504  */
00505 
00506 
00507 static double sequence_geometric_next(double low, double high, double ratio)
00508 {
00509         static unsigned i = 0;
00510         double x = low * pow(ratio, i++);
00511 
00512         if(high < low || high / low < ratio)
00513                 return XLAL_REAL8_FAIL_NAN;
00514 
00515         if(x > high) {
00516                 i = 1;
00517                 return low;
00518         }
00519 
00520         return x;
00521 }
00522 
00523 
00524 /*
00525  * Random amplitude presets.
00526  */
00527 
00528 
00529 static double sequence_preset_next(gsl_rng *rng)
00530 {
00531         static const double presets[] = {
00532                 1e-22,
00533                 2e-22,
00534                 5e-22,
00535                 1e-21,
00536                 2e-21,
00537                 5e-21,
00538                 1e-20
00539         };
00540 
00541         return presets[gsl_rng_uniform_int(rng, sizeof(presets)/sizeof(*presets))];
00542 }
00543 
00544 
00545 /* 
00546  * ============================================================================
00547  *
00548  *                Logarithmically-Distributed Random Variable
00549  *
00550  * ============================================================================
00551  */
00552 
00553 
00554 /*
00555  * Return a random number in the range [a, b), whose logarithm is uniformly
00556  * distributed
00557  */
00558 
00559 
00560 static double ran_flat_log(gsl_rng *rng, double a, double b)
00561 {
00562         return exp(gsl_ran_flat(rng, log(a), log(b)));
00563 }
00564 
00565 
00566 /*
00567  * Logarithmically-distributed random(ish) sequence.  Same as
00568  * sequence_geometric_next() but where each repetition of the sequence has
00569  * a random factor applied whose logarithm is uniformly distributed.  The
00570  * result is a random variable whose logarithm is uniformly distributed on
00571  * average, but in which neighbouring numbers are separated by a guaranteed
00572  * minimum ratio.
00573  */
00574 
00575 
00576 static double ran_flat_log_discrete(gsl_rng *rng, double a, double b, double ratio)
00577 {
00578         static double factor = 0.0;
00579         double x = sequence_geometric_next(a, b / ratio, ratio);
00580 
00581         if(x == a)
00582                 /* sequence has looped.  must happen first time through */
00583                 factor = ran_flat_log(rng, 1.0, ratio);
00584 
00585         return x * factor;
00586 }
00587 
00588 
00589 /* 
00590  * ============================================================================
00591  *
00592  *                          String Cusp Simulations
00593  *
00594  * ============================================================================
00595  */
00596 
00597 
00598 /*
00599  * \theta is the angle the line of sight makes with the cusp.  \theta^{2}
00600  * is distributed uniformly, and the high frequency cut-off of the
00601  * injection is \propto \theta^{-3}.  So we first solve for the limits of
00602  * \theta^{2} from the low- and high bounds of the frequency band of
00603  * interest, pick a uniformly-distributed number in that range, and convert
00604  * back to a high frequency cut-off.
00605  */
00606 
00607 
00608 static double random_string_cusp_fhigh(double flow, double fhigh, gsl_rng *rng)
00609 {
00610         const double thetasqmin = pow(fhigh, -2.0 / 3.0);
00611         const double thetasqmax = pow(flow, -2.0 / 3.0);
00612         const double thetasq = gsl_ran_flat(rng, thetasqmin, thetasqmax);
00613 
00614         return pow(thetasq, -3.0 / 2.0);
00615 }
00616 
00617 
00618 /*
00619  * Uniform sky location, and uniform polarization orientation.
00620  */
00621 
00622 
00623 static void random_location_and_polarization(double *ra, double *dec, double *psi, gsl_rng *rng)
00624 {
00625         *ra = gsl_ran_flat(rng, 0, LAL_TWOPI);
00626         *dec = asin(gsl_ran_flat(rng, -1, +1));
00627         *psi = gsl_ran_flat(rng, 0, LAL_TWOPI);
00628 }
00629 
00630 
00631 /*
00632  * Pick a random string cusp injection.
00633  */
00634 
00635 
00636 static SimBurst *random_string_cusp(double flow, double fhigh, double Alow, double Ahigh, gsl_rng *rng)
00637 {
00638         SimBurst *sim_burst = XLALCreateSimBurst();
00639 
00640         if(!sim_burst)
00641                 return NULL;
00642 
00643         strcpy(sim_burst->waveform, "StringCusp");
00644 
00645         /* high frequency cut-off and amplitude */
00646 
00647         sim_burst->frequency = random_string_cusp_fhigh(flow, fhigh, rng);
00648         sim_burst->amplitude = ran_flat_log(rng, Alow, Ahigh);
00649 
00650         /* sky location and wave frame orientation */
00651 
00652         random_location_and_polarization(&sim_burst->ra, &sim_burst->dec, &sim_burst->psi, rng);
00653 
00654         /* string cusp waveform generator makes a linearly polarized
00655          * waveform in the + polarization.  it ignores these parameters,
00656          * but just for consistency we populate them appropriately */
00657 
00658         sim_burst->pol_ellipse_e = 1.0;
00659         sim_burst->pol_ellipse_angle = 0.0;
00660 
00661         return sim_burst;
00662 }
00663 
00664 
00665 /* 
00666  * ============================================================================
00667  *
00668  *                            Directed Simulations
00669  *
00670  * ============================================================================
00671  */
00672 
00673 
00674 /*
00675  * Pick a random band- and time-limited white noise burst.
00676  */
00677 
00678 
00679 static SimBurst *random_directed_btlwnb(double ra, double dec, double minf, double maxf, double minband, double maxband, double mindur, double maxdur, double minEoverr2, double maxEoverr2, gsl_rng *rng)
00680 {
00681         REAL8TimeSeries *hplus, *hcross;
00682         SimBurst *sim_burst = XLALCreateSimBurst();
00683 
00684         if(!sim_burst)
00685                 return NULL;
00686 
00687         strcpy(sim_burst->waveform, "BTLWNB");
00688 
00689         /* sky location and wave frame orientation */
00690 
00691         sim_burst->ra = ra;
00692         sim_burst->dec = dec;
00693         sim_burst->psi = gsl_ran_flat(rng, 0, LAL_TWOPI);
00694 
00695         /* pick a waveform */
00696         /* FIXME:  this should be ULONG_MAX but metaio can't handle
00697          * unsigned ints yet */
00698 
00699         sim_burst->waveform_number = floor(gsl_ran_flat(rng, 0, LONG_MAX));
00700 
00701         /* centre frequency.  three steps between minf and maxf */
00702 
00703         sim_burst->frequency = ran_flat_log_discrete(rng, minf, maxf, pow(maxf / minf, 1.0 / 3.0));
00704 
00705         /* duration and bandwidth.  keep picking until a valid pair is
00706          * obtained (i.e. their product is >= 2 / \pi) */
00707 
00708         do {
00709                 sim_burst->duration = ran_flat_log(rng, mindur, maxdur);
00710                 sim_burst->bandwidth = ran_flat_log(rng, minband, maxband);
00711         } while(sim_burst->bandwidth * sim_burst->duration < LAL_2_PI);
00712 
00713         /* energy */
00714 
00715         sim_burst->egw_over_rsquared = ran_flat_log(rng, minEoverr2, maxEoverr2) * pow(sim_burst->frequency / 100.0, 4.0);
00716 
00717         /* populate the hrss column for convenience later */
00718         /* FIXME:  sample rate is hard-coded to 8192 Hz, which is what the
00719          * excess power pipeline's .ini file is configured for in CVS, but
00720          * because it can be easily changed this is not good */
00721 
00722         XLALGenerateSimBurst(&hplus, &hcross, sim_burst, 1.0 / 8192);
00723         sim_burst->hrss = XLALMeasureHrss(hplus, hcross);
00724         XLALDestroyREAL8TimeSeries(hplus);
00725         XLALDestroyREAL8TimeSeries(hcross);
00726 
00727         /* not sure if this makes sense.  these parameters are ignored by
00728          * the injection code, but post-processing tools sometimes wish to
00729          * know with what amplitude an injection should've been seen in an
00730          * instrument, and they use these to project the + and x amplitudes
00731          * onto the detector.  setting the eccentricity to 0 and the angle
00732          * to anything makes such tools believe the amplitude is equally
00733          * partitioned between the two polarizations which is true for
00734          * these injections. */
00735 
00736         sim_burst->pol_ellipse_e = 0.0;
00737         sim_burst->pol_ellipse_angle = 0.0;
00738 
00739         /* done */
00740 
00741         return sim_burst;
00742 }
00743 
00744 
00745 /* 
00746  * ============================================================================
00747  *
00748  *                           All-Sky sine-Gaussians
00749  *
00750  * ============================================================================
00751  */
00752 
00753 
00754 /*
00755  * the duration of a sine-Gaussian from its Q and centre frequency
00756  */
00757 
00758 
00759 static double duration_from_q_and_f(double Q, double f)
00760 {
00761         /* compute duration from Q and frequency */
00762         return Q / (sqrt(2.0) * LAL_PI * f);
00763 }
00764 
00765 
00766 /*
00767  * pick a sine-Gaussian
00768  */
00769 
00770 
00771 static SimBurst *random_all_sky_sineGaussian(double minf, double maxf, double q, double minhrsst, double maxhrsst, gsl_rng *rng)
00772 {
00773         SimBurst *sim_burst = XLALCreateSimBurst();
00774 
00775         if(!sim_burst)
00776                 return NULL;
00777 
00778         strcpy(sim_burst->waveform, "SineGaussian");
00779 
00780         /* sky location and wave frame orientation */
00781 
00782         random_location_and_polarization(&sim_burst->ra, &sim_burst->dec, &sim_burst->psi, rng);
00783 
00784         /* q and centre frequency.  three steps between minf and maxf */
00785 
00786         sim_burst->q = q;
00787         sim_burst->frequency = ran_flat_log_discrete(rng, minf, maxf, pow(maxf / minf, 1.0 / 3.0));
00788 
00789         /* hrss */
00790 
00791         sim_burst->hrss = ran_flat_log(rng, minhrsst, maxhrsst) / duration_from_q_and_f(sim_burst->q, sim_burst->frequency);
00792 
00793         /* hard-code for linearly polarized waveforms in the x
00794          * polarization.  induces LAL's sine-Gaussian generator to produce
00795          * linearly polarized sine-Gaussians (+ would be a cosine
00796          * Gaussian). */
00797 
00798         sim_burst->pol_ellipse_e = 1.0;
00799         sim_burst->pol_ellipse_angle = LAL_PI_2;
00800 
00801         /* done */
00802 
00803         return sim_burst;
00804 }
00805 
00806 
00807 /* 
00808  * ============================================================================
00809  *
00810  *                                   Output
00811  *
00812  * ============================================================================
00813  */
00814 
00815 
00816 static void write_xml(const ProcessTable *process_table_head, const ProcessParamsTable *process_params_table_head, const SearchSummaryTable *search_summary_head, const SimBurst *sim_burst, struct options options)
00817 {
00818         LALStatus status = blank_status;
00819         char fname[256];
00820         LIGOLwXMLStream xmlfp;
00821 
00822         memset(&xmlfp, 0, sizeof(xmlfp));
00823 
00824         if(options.user_tag)
00825                 snprintf(fname, sizeof(fname), "HL-INJECTIONS_%s-%d-%d.xml", options.user_tag, (int) (options.gps_start_time / LAL_INT8_C(1000000000)), (int) ((options.gps_end_time - options.gps_start_time) / LAL_INT8_C(1000000000)));
00826         else
00827                 snprintf(fname, sizeof(fname), "HL-INJECTIONS-%d-%d.xml", (int) (options.gps_start_time / LAL_INT8_C(1000000000)), (int) ((options.gps_end_time - options.gps_start_time) / LAL_INT8_C(1000000000)));
00828 
00829         LAL_CALL(LALOpenLIGOLwXMLFile(&status, &xmlfp, fname), &status);
00830 
00831         /* process table */
00832         if(XLALWriteLIGOLwXMLProcessTable(&xmlfp, process_table_head)) {
00833                 /* error occured.  ?? do anything else ?? */
00834                 exit(1);
00835         }
00836 
00837         /* process params table */
00838         if(XLALWriteLIGOLwXMLProcessParamsTable(&xmlfp, process_params_table_head)) {
00839                 /* error occured.  ?? do anything else ?? */
00840                 exit(1);
00841         }
00842 
00843         /* search summary table */
00844         if(XLALWriteLIGOLwXMLSearchSummaryTable(&xmlfp, search_summary_head)) {
00845                 /* error occured.  ?? do anything else ?? */
00846                 exit(1);
00847         }
00848 
00849         /* sim burst table */
00850         if(XLALWriteLIGOLwXMLSimBurstTable(&xmlfp, sim_burst)) {
00851                 /* error occured.  ?? do anything else ?? */
00852                 exit(1);
00853         }
00854 
00855         LAL_CALL(LALCloseLIGOLwXMLFile(&status, &xmlfp), &status);
00856 }
00857 
00858 
00859 /* 
00860  * ============================================================================
00861  *
00862  *                                Entry Point
00863  *
00864  * ============================================================================
00865  */
00866 
00867 
00868 int main(int argc, char *argv[])
00869 {
00870         struct options options;
00871         INT8 tinj;
00872         gsl_rng *rng;
00873         ProcessTable *process_table_head;
00874         ProcessParamsTable *process_params_table_head = NULL;
00875         SearchSummaryTable *search_summary_head;
00876         SimBurst *sim_burst_table_head = NULL;
00877         SimBurst **sim_burst = &sim_burst_table_head;
00878 
00879 
00880         /*
00881          * Initialize debug handler
00882          */
00883 
00884 
00885         lal_errhandler = LAL_ERR_EXIT;
00886         lalDebugLevel = LALINFO | LALWARNING | LALERROR | LALNMEMDBG | LALNMEMPAD | LALNMEMTRK;
00887 
00888 
00889         /*
00890          * Process table
00891          */
00892 
00893 
00894         process_table_head = XLALCreateProcessTableRow();
00895         if(XLALPopulateProcessTable(process_table_head, PROGRAM_NAME, CVS_REVISION, CVS_SOURCE, CVS_DATE, 0))
00896                 exit(1);
00897         XLALGPSTimeNow(&process_table_head->start_time);
00898         snprintf(process_table_head->ifos, sizeof(process_table_head->ifos), "H1,H2,L1");
00899 
00900 
00901         /*
00902          * Command line and process params table.
00903          */
00904 
00905 
00906         options = parse_command_line(&argc, &argv, process_table_head, &process_params_table_head);
00907         if(options.user_tag)
00908                 snprintf(process_table_head->comment, sizeof(process_table_head->comment), "%s", options.user_tag);
00909 
00910 
00911         /*
00912          * Search summary table
00913          */
00914 
00915 
00916         search_summary_head = XLALCreateSearchSummaryTableRow(process_table_head);
00917         if(options.user_tag)
00918                 snprintf(search_summary_head->comment, sizeof(search_summary_head->comment), "%s", options.user_tag);
00919         search_summary_head->nnodes = 1;
00920         search_summary_head->out_start_time = *XLALINT8NSToGPS(&search_summary_head->in_start_time, options.gps_start_time);
00921         search_summary_head->out_end_time = *XLALINT8NSToGPS(&search_summary_head->in_end_time, options.gps_end_time);
00922         snprintf(search_summary_head->ifos, sizeof(search_summary_head->ifos), process_table_head->ifos);
00923 
00924 
00925         /*
00926          * Initialize random number generator
00927          */
00928 
00929 
00930         rng = gsl_rng_alloc(gsl_rng_mt19937);
00931         if(options.seed)
00932                 gsl_rng_set(rng, options.seed);
00933 
00934 
00935         /*
00936          * Main loop
00937          */
00938 
00939 
00940         for(tinj = options.gps_start_time; tinj <= options.gps_end_time; tinj += options.time_step * 1e9) {
00941                 /*
00942                  * Progress bar.
00943                  */
00944 
00945                 XLALPrintInfo("%s: ", argv[0]);
00946                 XLALPrintProgressBar((tinj - options.gps_start_time) / (double) (options.gps_end_time - options.gps_start_time));
00947                 XLALPrintInfo(" complete\n");
00948 
00949                 /*
00950                  * Create an injection
00951                  */
00952 
00953                 switch(options.population) {
00954                 case POPULATION_TARGETED:
00955                         *sim_burst = random_directed_btlwnb(options.ra, options.dec, options.minf, options.maxf, options.minbandwidth, options.maxbandwidth, options.minduration, options.maxduration, options.minEoverr2, options.maxEoverr2, rng);
00956                         break;
00957 
00958                 case POPULATION_ALL_SKY_SINEGAUSSIAN:
00959                         *sim_burst = random_all_sky_sineGaussian(options.minf, options.maxf, options.q, options.minhrss, options.maxhrss, rng);
00960                         break;
00961 
00962                 case POPULATION_STRING_CUSP:
00963                         *sim_burst = random_string_cusp(options.minf, options.maxf, options.minA, options.maxA, rng);
00964                         break;
00965                 }
00966 
00967                 if(!*sim_burst) {
00968                         XLALPrintError("can't make injection");
00969                         exit(1);
00970                 }
00971 
00972                 /*
00973                  * Peak time at geocentre in GPS seconds
00974                  */
00975 
00976                 XLALINT8NSToGPS(&(*sim_burst)->time_geocent_gps, tinj);
00977 
00978                 /*
00979                  * Peak time at geocentre in GMST radians
00980                  */
00981 
00982                 (*sim_burst)->time_geocent_gmst = XLALGreenwichMeanSiderealTime(&(*sim_burst)->time_geocent_gps);
00983 
00984                 /*
00985                  * Move to next injection
00986                  */
00987 
00988                 sim_burst = &(*sim_burst)->next;
00989         }
00990 
00991         /* output */
00992 
00993         XLALGPSTimeNow(&process_table_head->end_time);
00994         search_summary_head->nevents = XLALSimBurstAssignIDs(sim_burst_table_head, process_table_head->process_id, 0);
00995         write_xml(process_table_head, process_params_table_head, search_summary_head, sim_burst_table_head, options);
00996 
00997         /* done */
00998 
00999         gsl_rng_free(rng);
01000         XLALDestroyProcessTable(process_table_head);
01001         XLALDestroyProcessParamsTable(process_params_table_head);
01002         XLALDestroySearchSummaryTable(search_summary_head);
01003         XLALDestroySimBurstTable(sim_burst_table_head);
01004         exit(0);
01005 }

Generated on Sun Jul 20 03:14:04 2008 for LAL by  doxygen 1.5.2