GENIEGenerator
Loading...
Searching...
No Matches
gAtmoEvGen.cxx
Go to the documentation of this file.
1//________________________________________________________________________________________
2/*!
3
4\program gevgen_atmo
5
6\brief A GENIE atmospheric neutrino event generation application.
7
8 *** Synopsis :
9
10 gevgen_atmo [-h]
11 [-r run#]
12 -f flux
13 -g geometry
14 [-R coordinate_rotation_matrix]
15 [-t geometry_top_volume_name]
16 [-m max_path_lengths_xml_file]
17 [-L geometry_length_units]
18 [-D geometry_density_units]
19 <-n n_of_events,
20 -e exposure_in_kton_x_yrs
21 -T exposure_in_seconds >
22 -E min_energy,max_energy
23 [-o output_event_file_prefix]
24 [--flux-ray-generation-surface-distance ]
25 [--flux-ray-generation-surface-radius ]
26 [--seed random_number_seed]
27 [--cross-sections xml_file]
28 [--event-generator-list list_name]
29 [--tune genie_tune]
30 [--message-thresholds xml_file]
31 [--unphysical-event-mask mask]
32 [--event-record-print-level level]
33 [--mc-job-status-refresh-rate rate]
34 [--cache-file root_file]
35
36 *** Options :
37
38 [] Denotes an optional argument.
39 <> Denotes a set of arguments out of which only one can be set.
40
41 -h
42 Prints out the syntax and exits
43 -r
44 Specifies the MC run number
45 [default: 100000000]
46 -f
47 Specifies the input flux files
48 The general syntax is: `-f simulation:/path/file.data[neutrino_code],...'
49 [Notes]
50 - The `simulation' string can be either `FLUKA', `BGLRS' or `HAKKM'.
51 See:
52 - $GENIE/src/Flux/GFLUKAAtmoFlux.h
53 - $GENIE/src/Flux/GBGLRSAtmoFlux.h
54 - $GENIE/src/Flux/GHAKKMAtmoFlux.h
55 - The neutrino codes are the PDG ones.
56 - The /path/file.data,neutrino_code part of the option can be
57 repeated multiple times (separated by commas), once for each
58 flux neutrino species you want to consider,
59 eg. '-f FLUKA:~/data/sdave_numu07.dat[14],~/data/sdave_nue07.dat[12]'
60 eg. '-f BGLRS:~/data/flux10_271003_z.kam_nue[12]'
61 eg. '-f HAKKM:~/data/kam-ally-20-12-solmax.d[12]'
62 -g
63 Input 'geometry'.
64 This option can be used to specify any of:
65 1 > A ROOT file containing a ROOT/GEANT geometry description
66 [Examples]
67 - To use the master volume from the ROOT geometry stored
68 in the nd280-geom.root file, type:
69 '-g /some/path/nd280-geom.root'
70 2 > A mix of target materials, each with its corresponding weight,
71 typed as a comma-separated list of nuclear pdg codes (in the
72 std PDG2006 convention: 10LZZZAAAI) with the weight fractions
73 in brackets, eg code1[fraction1],code2[fraction2],...
74 If that option is used (no detailed input geometry description)
75 then the interaction vertices are distributed in the detector
76 by the detector MC.
77 [Examples]
78 - To use a target mix of 89% O16 and 11% H, type:
79 '-g 1000080160[0.89],1000010010[0.11]'
80 - To use a target which is 100% C12, type:
81 '-g 1000060120'
82 -R
83 Input rotation matrix for transforming the flux neutrino coordinates
84 from the default Topocentric Horizontal (see GENIE manual) coordinate
85 system to the user-defined topocentric coordinate system.
86 The rotation is specified by the 3 Euler angles (phi, theta, psi).
87 The Euler angles are input as a comma separated list as:
88 `-R <convention>:phi,theta,psi',
89 where <convention> is either X (for X-convention), Y (for Y-convention),
90 X^-1 or Y^-1 (as previously, but using the inverse matrix).
91 By default, the X-convention (rotation about Z-axis, then about the
92 new X-axis, then about the Z-axis) is used.
93 Notes:
94 - (Extract from TRotation documentation)
95 "Euler angles usually define the rotation of the new coordinate
96 system with respect to the original system, however, the TRotation
97 class specifies the rotation of the object in the original system
98 (an active rotation). To recover the usual Euler rotations (ie.
99 rotate the system not the object), you must take the inverse of
100 the rotation."
101 Examples:
102 1. To set the Euler angles phi=3.14, theta=1.28, psi=1.0 using the
103 X-convention, type: `-R 3.14,1.28,1.0', or `-R X:3.14,1.28,1.0'
104 2. To set the Euler angles phi=3.14, theta=1.28, psi=1.0 using the
105 Y-convention, type: `-R Y:3.14,1.28,1.0'
106 3. To set the Euler angles phi=3.14, theta=1.28, psi=1.0 using the
107 Y-convention, and then use the inverse rotation matrix, type:
108 `-R Y^-1:3.14,1.28,1.0'
109 -L
110 Input geometry length units, eg 'm', 'cm', 'mm', ...
111 [default: 'mm']
112 -D
113 Input geometry density units, eg 'g_cm3', 'clhep_def_density_unit',...
114 [default: 'g_cm3']
115 -t
116 Input 'top volume' for event generation -
117 can be used to force event generation in given sub-detector
118 [default: the 'master volume' of the input geometry]
119 You can also use the -t option to switch generation on/off at
120 multiple volumes as, for example, in:
121 `-t +Vol1-Vol2+Vol3-Vol4',
122 `-t "+Vol1 -Vol2 +Vol3 -Vol4"',
123 `-t -Vol2-Vol4+Vol1+Vol3',
124 `-t "-Vol2 -Vol4 +Vol1 +Vol3"'
125 where:
126 "+Vol1" and "+Vol3" tells GENIE to `switch on' Vol1 and Vol3, while
127 "-Vol2" and "-Vol4" tells GENIE to `switch off' Vol2 and Vol4.
128 If the very first character is a '+', GENIE will neglect all volumes
129 except the ones explicitly turned on. Vice versa, if the very first
130 character is a `-', GENIE will keep all volumes except the ones
131 explicitly turned off (feature contributed by J.Holeczek).
132 -n
133 Specifies how many events to generate.
134 -e
135 Specifies requested exposure in terms of kton*yrs.
136 -T
137 Specifies requested exposure in terms of seconds.
138 -E
139 Specifies the neutrino energy in GeV.
140 Must be a comma-separated pair of numbers, eg `-E 0.3,70'
141 [default: 0.5,50]
142 --flux-ray-generation-surface-distance
143 --flux-ray-generation-surface-radius
144 See the User & Physics Manual for a graphical representation of the flux
145 ray generation surface: For a given zenith \theta and azimuthal angle \phi,
146 flux rays are produced within the area of a circle of radius Rt,
147 which is tangetial to a point P on a sphere of radius Rl, centred at the
148 detector. The point P has polar coordinates \theta and \phi.
149 The argument --flux-ray-generation-surface-distance sets Rl, while
150 the argument --flux-ray-generation-surface-distance sets Rt.
151 SI units are used.
152 -o
153 Sets the prefix of the output event file.
154 The output filename is built as:
155 [prefix].[run_number].[event_tree_format].[file_format]
156 The default output filename is:
157 gntp.[run_number].ghep.root
158 This cmd line arguments lets you override 'gntp'
159 --seed
160 Random number seed.
161 --cross-sections
162 Name (incl. full path) of an XML file with pre-computed
163 cross-section values used for constructing splines.
164 --tune
165 Specifies a GENIE comprehensive neutrino interaction model tune.
166 [default: "Default"].
167 --message-thresholds
168 Allows users to customize the message stream thresholds.
169 The thresholds are specified using an XML file.
170 See $GENIE/config/Messenger.xml for the XML schema.
171 Multiple files, delimited with a `:' can be specified.
172 --unphysical-event-mask
173 Allows users to specify a 16-bit mask to allow certain types of
174 unphysical events to be written in the output file.
175 [default: all unphysical events are rejected]
176 --event-record-print-level
177 Allows users to set the level of information shown when the event
178 record is printed in the screen. See GHepRecord::Print().
179 --mc-job-status-refresh-rate
180 Allows users to customize the refresh rate of the status file.
181 --cache-file
182 Allows users to specify a cache file so that the cache can be
183 re-used in subsequent MC jobs.
184
185 *** Examples:
186
187 (1) Generate 100k events (run number 999210) in the energy range 1-10 GeV
188 for nu_e and nu_mu only, using the sdave_numu07.dat FLUKA flux file for
189 nu_mu and the sdave_nue07.dat file for nu_e (files in /data/flx/).
190 Use the detector geometry in the /data/geo/SuperK.root file, where the
191 geometry length and density units are m and kgr/m^3. Generate events over
192 the entire geometry volume. Pre-computed cross-section data are loaded
193 from /data/xsec.xml.
194
195 % gevgen_atmo -r 999210 -n 100000 -E 1,10
196 -f FLUKA:/data/flx/sdave_numu07.dat[14],/data/flx/sdave_nue07.dat[12]
197 -g /data/geo/SuperK.root -L "m" -D "kg_m3"
198 --cross-sections /data/xsec.xml
199
200 (2) Like above but, instead of generating events in a realistic detector
201 geometry, use a simple target mix (88.79% O16 + 11.21% H, i.e. `water')
202
203 % gevgen_atmo -r 999210 -n 100000 -E 1,10
204 -f /data/flux/sdave_numu07.dat[14],/data/flux/sdave_nue07.dat[12]
205 -g 1000080160[0.8879],1000010010[0.1121]
206 --cross-sections /data/xsec.xml
207
208 ... to add more
209
210 Please read the GENIE User Manual for more information.
211
212\created August 20, 2010
213
214\author Costas Andreopoulos <c.andreopoulos \at cern.ch>
215 University of Liverpool
216
217 Torben Ferber <torben.ferber \at DESY.DE>
218 DESY
219
220 Hugh Gallagher <hugh.gallagher \at tufts.edu>
221 Tufts University
222
223 Tarak Thakore <tarak \at mailhost.tifr.res.in>
224 Tata Institute of Fundamental Research
225
226\cpright Copyright (c) 2003-2025, The GENIE Collaboration
227 For the full text of the license visit http://copyright.genie-mc.org
228
229*/
230//_________________________________________________________________________________________
231
232#include <cassert>
233#include <cstdlib>
234#include <cctype>
235#include <string>
236#include <vector>
237#include <sstream>
238#include <map>
239#include <iomanip>
240#include <cmath>
241
242#include <TRotation.h>
243#include <TMath.h>
244#include <TGeoShape.h>
245#include <TGeoBBox.h>
246
266
267#ifdef __GENIE_FLUX_DRIVERS_ENABLED__
271#endif
272
273#ifdef __GENIE_GEOM_DRIVERS_ENABLED__
277#endif
278
279using std::string;
280using std::vector;
281using std::map;
282using std::ostringstream;
283using std::setprecision;
284
285using namespace genie;
286using namespace genie::flux;
287
288void GetCommandLineArgs (int argc, char ** argv);
289void PrintSyntax (void);
290GAtmoFlux* GetFlux (void);
292
293// User-specified options:
294//
295Long_t gOptRunNu; // run number
296string gOptFluxSim; // flux simulation (FLUKA, BGLRS or HAKKM)
297map<int,string> gOptFluxFiles; // neutrino pdg code -> flux file map
298bool gOptUsingRootGeom = false; // using root geom or target mix?
299map<int,double> gOptTgtMix; // target mix (tgt pdg -> wght frac) / if not using detailed root geom
300string gOptRootGeom; // input ROOT file with realistic detector geometry
301string gOptRootGeomTopVol = ""; // input geometry top event generation volume
302double gOptGeomLUnits = 0; // input geometry length units
303double gOptGeomDUnits = 0; // input geometry density units
304string gOptExtMaxPlXml; // max path lengths XML file for input geometry
305int gOptNev = -1; // exposure - in terms of number of events
306double gOptKtonYrExposure = -1; // exposure - in terms of kton*yrs
307double gOptSecExposure = -1; // exposure - in terms of seconds
308double gOptEvMin; // minimum neutrino energy
309double gOptEvMax; // maximum neutrino energy
310string gOptEvFilePrefix; // event file prefix
311TRotation gOptRot; // coordinate rotation matrix: topocentric horizontal -> user-defined topocentric system
312long int gOptRanSeed; // random number seed
313string gOptInpXSecFile; // cross-section splines
314double gOptRL = -1; // distance of flux ray generation surface (m)
315double gOptRT = -1; // radius of flux ray generation surface (m)
316
317// Defaults:
318//
319NtpMCFormat_t kDefOptNtpFormat = kNFGHEP; // def event tree format
320string kDefOptEvFilePrefix = "gntp"; // def output prefix (override with -o)
321string kDefOptGeomLUnits = "mm"; // def geom length units (override with -L)
322string kDefOptGeomDUnits = "g_cm3"; // def geom density units (override with -D)
323double kDefOptEvMin = 0.5; // min neutrino energy (override with -E)
324double kDefOptEvMax = 50.0; // max neutrino energy (override with -E)
325
326//________________________________________________________________________________________
327int main(int argc, char** argv)
328{
329 GAtmoFlux* atmo_flux_driver;
330 double total_flux, expected_neutrinos;
331
332 // Parse command line arguments
333 GetCommandLineArgs(argc,argv);
334
335 if ( ! RunOpt::Instance()->Tune() ) {
336 LOG("gmkspl", pFATAL) << " No TuneId in RunOption";
337 exit(-1);
338 }
340
341 // Iinitialization of random number generators, cross-section table, messenger, cache etc...
342 utils::app_init::MesgThresholds(RunOpt::Instance()->MesgThresholdFiles());
346
347 // get geometry driver
348 GeomAnalyzerI * geom_driver = GetGeometry();
349
350 if (gOptRT < 0) {
351 gOptRT = 1000; // m
352 LOG("gevgen_atmo", pWARN) << "Warning! Flux surface transverse radius not specified so using default value of " << gOptRT << " meters!";
353 }
354
355 if (gOptRL < 0) {
356 gOptRL = 1000; // m
357 LOG("gevgen_atmo", pWARN) << "Warning! Flux surface longitudinal radius not specified so using default value of " << gOptRL << " meters!";
358 }
359
360 // get flux driver
361 atmo_flux_driver = GetFlux();
362
363 // Cast to GFluxI, the generic flux driver interface
364 GFluxI *flux_driver = dynamic_cast<GFluxI *>(atmo_flux_driver);
365
366 // create the GENIE monte carlo job driver
367 GMCJDriver* mcj_driver = new GMCJDriver;
369 mcj_driver->UseFluxDriver(flux_driver);
370 mcj_driver->UseGeomAnalyzer(geom_driver);
371 mcj_driver->Configure();
372 mcj_driver->UseSplines();
373 /* Note: For the method of calculating the total number of events using a
374 * livetime we *need* to force a single probability scale. So if you change
375 * the next line you need to make sure that the user didn't specify the -T
376 * option. */
377 mcj_driver->ForceSingleProbScale();
378
379 // initialize an ntuple writer
382 ntpw.Initialize();
383
384 // Create a MC job monitor for a periodically updated status file
385 GMCJMonitor mcjmonitor(gOptRunNu);
386 mcjmonitor.SetRefreshRate(RunOpt::Instance()->MCJobStatusRefreshRate());
387
388 // Set GHEP print level
389 GHepRecord::SetPrintLevel(RunOpt::Instance()->EventRecordPrintLevel());
390
391 total_flux = atmo_flux_driver->GetTotalFluxInEnergyRange();
392 LOG("gevgen_atmo", pNOTICE) << "Total atmospheric neutrino flux is " << setprecision(2) << total_flux << " neutrinos per m^2 per second.";
393 if (gOptSecExposure > 0) {
394 /* Calculate the expected value of the total number of neutrinos we need to
395 * throw. We do this by multiplying the total flux by the exposure time in
396 * seconds and the area of the flux surface. */
397 expected_neutrinos = total_flux*gOptSecExposure*atmo_flux_driver->GetFluxSurfaceArea();
398 LOG("gevgen_atmo", pNOTICE) << "Simulating an exposure of " << setprecision(0) << gOptSecExposure << " seconds which corresponds to a total of " << setprecision(0) << expected_neutrinos << " neutrinos.";
399 }
400
401 // event loop
402 for(int iev = 0; gOptNev > 0 ? iev < gOptNev : 1; iev++) {
403
404 // generate next event
405 EventRecord* event = mcj_driver->GenerateEvent();
406
407 // print-out
408 LOG("gevgen_atmo", pNOTICE) << "Generated event: " << *event;
409
410 // save the event, refresh the mc job monitor
411 ntpw.AddEventRecord(iev, event);
412 mcjmonitor.Update(iev,event);
413
414 // clean-up
415 delete event;
416
417 if (gOptSecExposure > 0 && mcj_driver->NFluxNeutrinos()/mcj_driver->GlobProbScale() > expected_neutrinos) {
418 break;
419 }
420 }
421
422 // save the event file
423 ntpw.Save();
424
425 // clean-up
426 delete geom_driver;
427 delete atmo_flux_driver;
428 delete mcj_driver;
429
430 return 0;
431}
432//________________________________________________________________________________________
434{
435 GeomAnalyzerI * geom_driver = 0;
436
437#ifdef __GENIE_GEOM_DRIVERS_ENABLED__
438
440 //
441 // *** Using a realistic root-based detector geometry description
442 //
443
444 // creating & configuring a root geometry driver
447 rgeom -> SetLengthUnits (gOptGeomLUnits);
448 rgeom -> SetDensityUnits (gOptGeomDUnits);
449 rgeom -> SetTopVolName (gOptRootGeomTopVol);
450 // getting the bounding box dimensions along z so as to set the
451 // appropriate upstream generation surface for the JPARC flux driver
452 TGeoVolume * topvol = rgeom->GetGeometry()->GetTopVolume();
453 if(!topvol) {
454 LOG("gevgen_atmo", pFATAL) << " ** Null top ROOT geometry volume!";
455 gAbortingInErr = true;
456 exit(1);
457 }
458
459 /* If flux generation surface isn't defined, get the bounding box for the
460 * geometry and set something appropriate. */
461 TGeoShape *bounding_box = topvol->GetShape();
462 TGeoBBox *box = (TGeoBBox *) bounding_box;
463 double dx = box->GetDX()*rgeom->LengthUnits();
464 double dy = box->GetDY()*rgeom->LengthUnits();
465 double dz = box->GetDZ()*rgeom->LengthUnits();
466
467 if (gOptRL < 0 && gOptRT < 0) {
468 gOptRL = TMath::Sqrt(dx*dx + dy*dy + dz*dz);
469 gOptRT = gOptRL;
470 LOG("gevgen_atmo", pNOTICE) << "Setting flux longitudinal and transverse radius to " << setprecision(2) << gOptRL << " meters based on bounding box of ROOT geometry.";
471 }
472
473 // switch on/off volumes as requested
474 if ( (gOptRootGeomTopVol[0] == '+') || (gOptRootGeomTopVol[0] == '-') ) {
475 bool exhaust = (*gOptRootGeomTopVol.c_str() == '+');
477 }
478
479 // casting to the GENIE geometry driver interface
480 geom_driver = dynamic_cast<GeomAnalyzerI *> (rgeom);
481 }
482 else {
483 //
484 // *** Using a 'point' geometry with the specified target mix
485 // *** ( = a list of targets with their corresponding weight fraction)
486 //
487
488 // creating & configuring a point geometry driver
491 // casting to the GENIE geometry driver interface
492 geom_driver = dynamic_cast<GeomAnalyzerI *> (pgeom);
493 }
494
495#else
496 LOG("gevgen_atmo", pFATAL) << "You need to enable the GENIE geometry drivers first!";
497 LOG("gevgen_atmo", pFATAL) << "Use --enable-geom-drivers at the configuration step.";
498 gAbortingInErr = true;
499 exit(1);
500#endif
501
502 return geom_driver;
503}
504//________________________________________________________________________________________
506{
507#ifdef __GENIE_FLUX_DRIVERS_ENABLED__
508
509 // Instantiate appropriate concrete flux driver
510 GAtmoFlux * atmo_flux_driver = 0;
511 if(gOptFluxSim == "FLUKA") {
512 GFLUKAAtmoFlux * fluka_flux = new GFLUKAAtmoFlux;
513 atmo_flux_driver = dynamic_cast<GAtmoFlux *>(fluka_flux);
514 } else
515 if(gOptFluxSim == "BGLRS") {
516 GBGLRSAtmoFlux * bartol_flux = new GBGLRSAtmoFlux;
517 atmo_flux_driver = dynamic_cast<GAtmoFlux *>(bartol_flux);
518 } else
519 if(gOptFluxSim == "HAKKM") {
520 GHAKKMAtmoFlux * honda_flux = new GHAKKMAtmoFlux;
521 atmo_flux_driver = dynamic_cast<GAtmoFlux *>(honda_flux);
522 } else {
523 LOG("gevgen_atmo", pFATAL) << "Unknown flux simulation: " << gOptFluxSim;
524 gAbortingInErr = true;
525 exit(1);
526 }
527 // Configure GAtmoFlux options (common to all concrete atmospheric flux drivers)
528 // set min/max energy:
529 atmo_flux_driver->ForceMinEnergy (gOptEvMin * units::GeV);
530 atmo_flux_driver->ForceMaxEnergy (gOptEvMax * units::GeV);
531 // set flux files:
532 map<int,string>::const_iterator file_iter = gOptFluxFiles.begin();
533 for( ; file_iter != gOptFluxFiles.end(); ++file_iter) {
534 int neutrino_code = file_iter->first;
535 string filename = file_iter->second;
536 atmo_flux_driver->AddFluxFile(neutrino_code, filename);
537 }
538
539 if (!atmo_flux_driver->LoadFluxData()) {
540 LOG("gevgen_atmo", pFATAL) << "Error loading flux data. Quitting...";
541 gAbortingInErr = true;
542 exit(1);
543 }
544
545 // configure flux generation surface:
546 atmo_flux_driver->SetRadii(gOptRL, gOptRT);
547 // set rotation for coordinate tranformation from the topocentric horizontal
548 // system to a user-defined coordinate system:
549 if(!gOptRot.IsIdentity()) {
550 atmo_flux_driver->SetUserCoordSystem(gOptRot);
551 }
552
553#else
554 LOG("gevgen_atmo", pFATAL) << "You need to enable the GENIE flux drivers first!";
555 LOG("gevgen_atmo", pFATAL) << "Use --enable-flux-drivers at the configuration step.";
556 gAbortingInErr = true;
557 exit(1);
558#endif
559
560 return atmo_flux_driver;
561}
562//________________________________________________________________________________________
563void GetCommandLineArgs(int argc, char ** argv)
564{
565// Get the command line arguments
566
568
569 LOG("gevgen_atmo", pNOTICE) << "Parsing command line arguments";
570
571 CmdLnArgParser parser(argc,argv);
572
573 // help?
574 bool help = parser.OptionExists('h');
575 if(help) {
576 PrintSyntax();
577 exit(0);
578 }
579
580 //
581 // *** run number
582 //
583 if( parser.OptionExists('r') ) {
584 LOG("gevgen_atmo", pDEBUG) << "Reading MC run number";
585 gOptRunNu = parser.ArgAsLong('r');
586 } else {
587 LOG("gevgen_atmo", pDEBUG) << "Unspecified run number - Using default";
588 gOptRunNu = 100000000;
589 } //-r
590
591 //
592 // *** exposure
593 //
594
595 // in number of events
596 bool have_required_statistics = false;
597 if( parser.OptionExists('n') ) {
598 LOG("gevgen_atmo", pDEBUG)
599 << "Reading number of events to generate";
600 gOptNev = parser.ArgAsInt('n');
601 have_required_statistics = true;
602 }//-n?
603 // or, in kton*yrs
604 if( parser.OptionExists('e') ) {
605 if(have_required_statistics) {
606 LOG("gevgen_atmo", pFATAL)
607 << "Can't request exposure both in terms of number of events and kton*yrs"
608 << "\nUse just one of the -n and -e options";
609 PrintSyntax();
610 gAbortingInErr = true;
611 exit(1);
612 }
613 LOG("gevgen_atmo", pDEBUG)
614 << "Reading requested exposure in kton*yrs";
615 gOptKtonYrExposure = parser.ArgAsDouble('e');
616 have_required_statistics = true;
617 }//-e?
618
619 if (parser.OptionExists('T')) {
620 if (have_required_statistics) {
621 LOG("gevgen_atmo", pFATAL)
622 << "Can't request exposure both in terms of number of events or kton*yrs and time"
623 << "\nUse just one of the -n, -e, or -T options";
624 PrintSyntax();
625 gAbortingInErr = true;
626 exit(1);
627 }
628 LOG("gevgen_atmo", pDEBUG)
629 << "Reading requested exposure in seconds";
630 gOptSecExposure = parser.ArgAsDouble('T');
631 have_required_statistics = true;
632 }
633
634 if(!have_required_statistics) {
635 LOG("gevgen_atmo", pFATAL)
636 << "You must request exposure either in terms of number of events and kton*yrs"
637 << "\nUse any of the -n, -e options";
638 PrintSyntax();
639 gAbortingInErr = true;
640 exit(1);
641 }
642
643 //
644 // *** event file prefix
645 //
646 if( parser.OptionExists('o') ) {
647 LOG("gevgen_atmo", pDEBUG) << "Reading the event filename prefix";
648 gOptEvFilePrefix = parser.ArgAsString('o');
649 } else {
650 LOG("gevgen_atmo", pDEBUG)
651 << "Will set the default event filename prefix";
653 } //-o
654
655 //
656 // *** neutrino energy range
657 //
658 if( parser.OptionExists('E') ) {
659 LOG("gevgen_atmo", pINFO) << "Reading neutrino energy range";
660 string nue = parser.ArgAsString('E');
661
662 // must be a comma separated set of values
663 if(nue.find(",") != string::npos) {
664 // split the comma separated list
665 vector<string> nurange = utils::str::Split(nue, ",");
666 assert(nurange.size() == 2);
667 double emin = atof(nurange[0].c_str());
668 double emax = atof(nurange[1].c_str());
669 assert(emax>emin && emin>=0.);
670 gOptEvMin = emin;
671 gOptEvMax = emax;
672 } else {
673 LOG("gevgen_atmo", pFATAL)
674 << "Invalid energy range. Use `-E emin,emax', eg `-E 0.5,100.";
675 PrintSyntax();
676 gAbortingInErr = true;
677 exit(1);
678 }
679 } else {
680 LOG("gevgen_atmo", pNOTICE)
681 << "No -e option. Using default energy range";
684 }
685
686 //
687 // *** flux files
688 //
689 // syntax:
690 // simulation:/path/file.data[neutrino_code],/path/file.data[neutrino_code],...
691 //
692 if( parser.OptionExists('f') ) {
693 LOG("gevgen_atmo", pDEBUG) << "Getting input flux files";
694 string flux = parser.ArgAsString('f');
695
696 // get flux simulation info (FLUKA,BGLRS) so as to instantiate the
697 // appropriate flux driver
698 string::size_type jsimend = flux.find_first_of(":",0);
699 if(jsimend==string::npos) {
700 LOG("gevgen_atmo", pFATAL)
701 << "You need to specify the flux file source";
702 PrintSyntax();
703 gAbortingInErr = true;
704 exit(1);
705 }
706 gOptFluxSim = flux.substr(0,jsimend);
707 for(string::size_type i=0; i<gOptFluxSim.size(); i++) {
708 gOptFluxSim[i] = toupper(gOptFluxSim[i]);
709 }
710 if((gOptFluxSim != "FLUKA") &&
711 (gOptFluxSim != "BGLRS") &&
712 (gOptFluxSim != "HAKKM")) {
713 LOG("gevgen_atmo", pFATAL)
714 << "The flux file source needs to be one of <FLUKA,BGLRS,HAKKM>";
715 PrintSyntax();
716 gAbortingInErr = true;
717 exit(1);
718 }
719 // now get the list of input files and the corresponding neutrino codes.
720 flux.erase(0,jsimend+1);
721 vector<string> fluxv = utils::str::Split(flux,",");
722 vector<string>::const_iterator fluxiter = fluxv.begin();
723 for( ; fluxiter != fluxv.end(); ++fluxiter) {
724 string filename_and_pdg = *fluxiter;
725 string::size_type open_bracket = filename_and_pdg.find("[");
726 string::size_type close_bracket = filename_and_pdg.find("]");
727 if (open_bracket ==string::npos ||
728 close_bracket==string::npos)
729 {
730 LOG("gevgen_atmo", pFATAL)
731 << "You made an error in specifying the flux info";
732 PrintSyntax();
733 gAbortingInErr = true;
734 exit(1);
735 }
736 string::size_type ibeg = 0;
737 string::size_type iend = open_bracket;
738 string::size_type jbeg = open_bracket+1;
739 string::size_type jend = close_bracket;
740 string flux_filename = filename_and_pdg.substr(ibeg,iend-ibeg);
741 string neutrino_pdg = filename_and_pdg.substr(jbeg,jend-jbeg);
742 gOptFluxFiles.insert(
743 map<int,string>::value_type(atoi(neutrino_pdg.c_str()), flux_filename));
744 }
745 if(gOptFluxFiles.size() == 0) {
746 LOG("gevgen_atmo", pFATAL)
747 << "You must specify at least one flux file!";
748 PrintSyntax();
749 gAbortingInErr = true;
750 exit(1);
751 }
752
753 } else {
754 LOG("gevgen_atmo", pFATAL) << "No flux info was specified! Use the -f option.";
755 PrintSyntax();
756 gAbortingInErr = true;
757 exit(1);
758 }
759
760 // *** options to fine tune the flux ray generation surface
761
762 if( parser.OptionExists("flux-ray-generation-surface-distance") ) {
763 LOG("gevgen_atmo", pINFO)
764 << "Reading distance of flux ray generation surface";
765 gOptRL = parser.ArgAsDouble("flux-ray-generation-surface-distance");
766 } else {
767 LOG("gevgen_atmo", pINFO)
768 << "Unspecified distance of flux ray generation surface - Using default";
769 }
770
771 if( parser.OptionExists("flux-ray-generation-surface-radius") ) {
772 LOG("gevgen_atmo", pINFO)
773 << "Reading radius of flux ray generation surface";
774 gOptRT = parser.ArgAsDouble("flux-ray-generation-surface-radius");
775 } else {
776 LOG("gevgen_atmo", pINFO)
777 << "Unspecified radius of flux ray generation surface - Using default";
778 }
779
780 //
781 // *** geometry
782 //
783 string geom = "";
784 string lunits, dunits;
785 if( parser.OptionExists('g') ) {
786 LOG("gevgen_atmo", pDEBUG) << "Getting input geometry";
787 geom = parser.ArgAsString('g');
788
789 // is it a ROOT file that contains a ROOT geometry?
790 bool accessible_geom_file =
792 if (accessible_geom_file) {
794 gOptUsingRootGeom = true;
795 }
796 } else {
797 LOG("gevgen_atmo", pFATAL)
798 << "No geometry option specified - Exiting";
799 PrintSyntax();
800 gAbortingInErr = true;
801 exit(1);
802 } //-g
803
805 // using a ROOT geometry - get requested geometry units
806
807 // legth units:
808 if( parser.OptionExists('L') ) {
809 LOG("gevgen_atmo", pDEBUG)
810 << "Checking for input geometry length units";
811 lunits = parser.ArgAsString('L');
812 } else {
813 LOG("gevgen_atmo", pDEBUG) << "Using default geometry length units";
815 } // -L
816 // density units:
817 if( parser.OptionExists('D') ) {
818 LOG("gevgen_atmo", pDEBUG)
819 << "Checking for input geometry density units";
820 dunits = parser.ArgAsString('D');
821 } else {
822 LOG("gevgen_atmo", pDEBUG) << "Using default geometry density units";
823 dunits = kDefOptGeomDUnits;
824 } // -D
827
828 // check whether an event generation volume name has been
829 // specified -- default is the 'top volume'
830 if( parser.OptionExists('t') ) {
831 LOG("gevgen_atmo", pDEBUG) << "Checking for input volume name";
832 gOptRootGeomTopVol = parser.ArgAsString('t');
833 } else {
834 LOG("gevgen_atmo", pDEBUG) << "Using the <master volume>";
835 } // -t
836
837 // check whether an XML file with the maximum (density weighted)
838 // path lengths for each detector material is specified -
839 // otherwise will compute the max path lengths at job init
840 if( parser.OptionExists('m') ) {
841 LOG("gevgen_atmo", pDEBUG)
842 << "Checking for maximum path lengths XML file";
843 gOptExtMaxPlXml = parser.ArgAsString('m');
844 } else {
845 LOG("gevgen_atmo", pDEBUG)
846 << "Will compute the maximum path lengths at job init";
847 gOptExtMaxPlXml = "";
848 } // -m
849 } // using root geom?
850
851 else {
852 // User has specified a target mix.
853 // Decode the list of target pdf codes & their corresponding weight fraction
854 // (specified as 'pdg_code_1[fraction_1],pdg_code_2[fraction_2],...')
855 // See documentation on top section of this file.
856 //
857 gOptTgtMix.clear();
858 vector<string> tgtmix = utils::str::Split(geom,",");
859 if(tgtmix.size()==1) {
860 int pdg = atoi(tgtmix[0].c_str());
861 double wgt = 1.0;
862 gOptTgtMix.insert(map<int, double>::value_type(pdg, wgt));
863 } else {
864 vector<string>::const_iterator tgtmix_iter = tgtmix.begin();
865 for( ; tgtmix_iter != tgtmix.end(); ++tgtmix_iter) {
866 string tgt_with_wgt = *tgtmix_iter;
867 string::size_type open_bracket = tgt_with_wgt.find("[");
868 string::size_type close_bracket = tgt_with_wgt.find("]");
869 if (open_bracket ==string::npos ||
870 close_bracket==string::npos)
871 {
872 LOG("gevgen_atmo", pFATAL)
873 << "You made an error in specifying the target mix";
874 PrintSyntax();
875 gAbortingInErr = true;
876 exit(1);
877 }
878 string::size_type ibeg = 0;
879 string::size_type iend = open_bracket;
880 string::size_type jbeg = open_bracket+1;
881 string::size_type jend = close_bracket;
882 int pdg = atoi(tgt_with_wgt.substr(ibeg,iend-ibeg).c_str());
883 double wgt = atof(tgt_with_wgt.substr(jbeg,jend-jbeg).c_str());
884 LOG("gevgen_atmo", pDEBUG)
885 << "Adding to target mix: pdg = " << pdg << ", wgt = " << wgt;
886 gOptTgtMix.insert(map<int, double>::value_type(pdg, wgt));
887
888 }// tgtmix_iter
889 } // >1 materials in mix
890 } // using tgt mix?
891
892 //
893 // Coordinate rotation matrix
894 //
895 gOptRot.SetToIdentity();
896 if( parser.OptionExists('R') ) {
897 string rotarg = parser.ArgAsString('R');
898 //get convention
899 string::size_type j = rotarg.find_first_of(":",0);
900 string convention = "";
901 if(j==string::npos) { convention = "X"; }
902 else { convention = rotarg.substr(0,j); }
903 //get angles phi,theta,psi
904 rotarg.erase(0,j+1);
905 vector<string> euler_angles = utils::str::Split(rotarg,",");
906 if(euler_angles.size() != 3) {
907 LOG("gevgen_atmo", pFATAL)
908 << "You didn't specify all 3 Euler angles using the -R option";
909 PrintSyntax();
910 gAbortingInErr = true;
911 exit(1);
912 }
913 double phi = atof(euler_angles[0].c_str());
914 double theta = atof(euler_angles[1].c_str());
915 double psi = atof(euler_angles[2].c_str());
916 //set Euler angles using appropriate convention
917 if(convention.find("X")!=string::npos ||
918 convention.find("x")!=string::npos)
919 {
920 LOG("gevgen_atmo", pNOTICE) << "Using X-convention for input Euler angles";
921 gOptRot.SetXEulerAngles(phi,theta,psi);
922 } else
923 if(convention.find("Y")!=string::npos ||
924 convention.find("y")!=string::npos)
925 {
926 LOG("gevgen_atmo", pNOTICE) << "Using Y-convention for input Euler angles";
927 gOptRot.SetYEulerAngles(phi,theta,psi);
928 } else {
929 LOG("gevgen_atmo", pFATAL)
930 << "Unknown Euler angle convention. Please use the X- or Y-convention";
931 PrintSyntax();
932 gAbortingInErr = true;
933 exit(1);
934 }
935 //invert?
936 if(convention.find("^-1")!=string::npos) {
937 LOG("gevgen_atmo", pNOTICE) << "Inverting rotation matrix";
938 gOptRot.Invert();
939 }
940 }
941
942 //
943 // *** random number seed
944 //
945 if( parser.OptionExists("seed") ) {
946 LOG("gevgen_atmo", pINFO) << "Reading random number seed";
947 gOptRanSeed = parser.ArgAsLong("seed");
948 } else {
949 LOG("gevgen_atmo", pINFO) << "Unspecified random number seed - Using default";
950 gOptRanSeed = -1;
951 }
952
953 //
954 // *** input cross-section file
955 //
956 if( parser.OptionExists("cross-sections") ) {
957 LOG("gevgen_atmo", pINFO) << "Reading cross-section file";
958 gOptInpXSecFile = parser.ArgAsString("cross-sections");
959 } else {
960 LOG("gevgen_atmo", pINFO) << "Unspecified cross-section file";
961 gOptInpXSecFile = "";
962 }
963
964 //
965 // print-out summary
966 //
967
968 PDGLibrary * pdglib = PDGLibrary::Instance();
969
970 ostringstream gminfo;
971 if (gOptUsingRootGeom) {
972 gminfo << "Using ROOT geometry - file: " << gOptRootGeom
973 << ", top volume: "
974 << ((gOptRootGeomTopVol.size()==0) ? "<master volume>" : gOptRootGeomTopVol)
975 << ", max{PL} file: "
976 << ((gOptExtMaxPlXml.size()==0) ? "<none>" : gOptExtMaxPlXml)
977 << ", length units: " << lunits
978 << ", density units: " << dunits;
979 } else {
980 gminfo << "Using target mix - ";
981 map<int,double>::const_iterator iter;
982 for(iter = gOptTgtMix.begin(); iter != gOptTgtMix.end(); ++iter) {
983 int pdg_code = iter->first;
984 double wgt = iter->second;
985 TParticlePDG * p = pdglib->Find(pdg_code);
986 if(p) {
987 string name = p->GetName();
988 gminfo << "(" << name << ") -> " << 100*wgt << "% / ";
989 }//p?
990 }
991 }
992
993 ostringstream fluxinfo;
994 fluxinfo << "Using " << gOptFluxSim << " flux files: ";
995 map<int,string>::const_iterator file_iter = gOptFluxFiles.begin();
996 for( ; file_iter != gOptFluxFiles.end(); ++file_iter) {
997 int neutrino_code = file_iter->first;
998 string filename = file_iter->second;
999 TParticlePDG * p = pdglib->Find(neutrino_code);
1000 if(p) {
1001 string name = p->GetName();
1002 fluxinfo << "(" << name << ") -> " << filename << " / ";
1003 }
1004 }
1005 fluxinfo << "Flux ray generation surface - Distance = "
1006 << gOptRL << " m, Radius = " << gOptRT << " m";
1007
1008 ostringstream expinfo;
1009 if(gOptNev > 0) { expinfo << gOptNev << " events"; }
1010 if(gOptKtonYrExposure > 0) { expinfo << gOptKtonYrExposure << " kton*yrs"; }
1011
1012 ostringstream rotation;
1013 rotation << "\t| " << gOptRot.XX() << " " << gOptRot.XY() << " " << gOptRot.XZ() << " |\n";
1014 rotation << "\t| " << gOptRot.YX() << " " << gOptRot.YY() << " " << gOptRot.YZ() << " |\n";
1015 rotation << "\t| " << gOptRot.ZX() << " " << gOptRot.ZY() << " " << gOptRot.ZZ() << " |\n";
1016
1017 LOG("gevgen_atmo", pNOTICE)
1018 << "\n\n"
1019 << utils::print::PrintFramedMesg("gevgen_atmo job configuration");
1020
1021 LOG("gevgen_atmo", pNOTICE)
1022 << "\n"
1023 << "\n @@ Run number: " << gOptRunNu
1024 << "\n @@ Random number seed: " << gOptRanSeed
1025 << "\n @@ Using cross-section file: " << gOptInpXSecFile
1026 << "\n @@ Geometry"
1027 << "\n\t" << gminfo.str()
1028 << "\n @@ Flux"
1029 << "\n\t" << fluxinfo.str()
1030 << "\n @@ Exposure"
1031 << "\n\t" << expinfo.str()
1032 << "\n @@ Cuts"
1033 << "\n\t Using energy range = (" << gOptEvMin << " GeV, " << gOptEvMax << " GeV)"
1034 << "\n @@ Coordinate transformation (Rotation THZ -> User-defined coordinate system)"
1035 << "\n" << rotation.str()
1036 << "\n\n";
1037
1038 //
1039 // final checks
1040 //
1041 if(gOptKtonYrExposure > 0) {
1042 LOG("gevgen_atmo", pFATAL)
1043 << "\n Option to set exposure in terms of kton*yrs not supported just yet!"
1044 << "\n Try the -n option instead";
1045 PrintSyntax();
1046 gAbortingInErr = true;
1047 exit(1);
1048 }
1049}
1050//________________________________________________________________________________________
1051void PrintSyntax(void)
1052{
1053 LOG("gevgen_atmo", pFATAL)
1054 << "\n **Syntax**"
1055 << "\n gevgen_atmo [-h]"
1056 << "\n [-r run#]"
1057 << "\n -f simulation:flux_file[neutrino_code],..."
1058 << "\n -g geometry"
1059 << "\n [-R coordinate_rotation_matrix]"
1060 << "\n [-t geometry_top_volume_name]"
1061 << "\n [-m max_path_lengths_xml_file]"
1062 << "\n [-L geometry_length_units]"
1063 << "\n [-D geometry_density_units]"
1064 << "\n <-n n_of_events,"
1065 << "\n -e exposure_in_kton_x_yrs"
1066 << "\n -T exposure_in_seconds>"
1067 << "\n -E min_energy,max_energy"
1068 << "\n [-o output_event_file_prefix]"
1069 << "\n [--flux-ray-generation-surface-distance]"
1070 << "\n [--flux-ray-generation-surface-radius]"
1071 << "\n [--seed random_number_seed]"
1072 << "\n --cross-sections xml_file"
1073 << "\n [--event-generator-list list_name]"
1074 << "\n [--message-thresholds xml_file]"
1075 << "\n [--unphysical-event-mask mask]"
1076 << "\n [--event-record-print-level level]"
1077 << "\n [--mc-job-status-refresh-rate rate]"
1078 << "\n [--cache-file root_file]"
1079 << "\n"
1080 << " Please also read the detailed documentation at http://www.genie-mc.org"
1081 << "\n";
1082}
1083//________________________________________________________________________________________
#define pNOTICE
Definition Messenger.h:61
#define pINFO
Definition Messenger.h:62
#define pFATAL
Definition Messenger.h:56
#define pDEBUG
Definition Messenger.h:63
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE,...
Definition Messenger.h:96
#define pWARN
Definition Messenger.h:60
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils.
int main()
Command line argument parser.
string ArgAsString(char opt)
double ArgAsDouble(char opt)
bool OptionExists(char opt)
was option set?
A vector of EventGeneratorI objects.
Generated Event Record. It is a GHepRecord object that can accept / be visited by EventRecordVisitorI...
Definition EventRecord.h:37
GENIE Interface for user-defined flux classes.
Definition GFluxI.h:29
static void SetPrintLevel(int print_level)
A GENIE ‘MC Job Driver’. Can be used for setting up complicated event generation cases involving deta...
Definition GMCJDriver.h:46
double GlobProbScale(void) const
Definition GMCJDriver.h:76
void ForceSingleProbScale(void)
void UseFluxDriver(GFluxI *flux)
void UseGeomAnalyzer(GeomAnalyzerI *geom)
void Configure(bool calc_prob_scales=true)
void UseSplines(bool useLogE=true)
long int NFluxNeutrinos(void) const
Definition GMCJDriver.h:77
void SetEventGeneratorList(string listname)
EventRecord * GenerateEvent(void)
Simple class to create & update MC job status files and env. vars. This is used to be able to keep tr...
Definition GMCJMonitor.h:31
void Update(int iev, const EventRecord *event)
void SetRefreshRate(int rate)
Defines the GENIE Geometry Analyzer Interface.
A utility class to facilitate creating the GENIE MC Ntuple from the output GENIE GHEP event records.
Definition NtpWriter.h:39
void Initialize(void)
add event
Definition NtpWriter.cxx:83
void Save(void)
get the even tree
void AddEventRecord(int ievent, const EventRecord *ev_rec)
save the event tree
Definition NtpWriter.cxx:57
void CustomizeFilenamePrefix(string prefix)
Singleton class to load & serve a TDatabasePDG.
Definition PDGLibrary.h:36
static PDGLibrary * Instance(void)
TParticlePDG * Find(int pdgc, bool must_exist=true)
void ReadFromCommandLine(int argc, char **argv)
Definition RunOpt.cxx:99
void BuildTune()
build tune and inform XSecSplineList
Definition RunOpt.cxx:92
static RunOpt * Instance(void)
Definition RunOpt.cxx:54
A base class for the FLUKA, BGLRS and ATMNC atmo. nu. flux drivers. The driver depends on data files ...
Definition GAtmoFlux.h:60
void AddFluxFile(int neutrino_pdg, string filename)
void ForceMaxEnergy(double emax)
double GetTotalFluxInEnergyRange(void)
double GetFluxSurfaceArea(void)
void ForceMinEnergy(double emin)
void SetRadii(double Rlongitudinal, double Rtransverse)
void SetUserCoordSystem(TRotation &rotation)
Rotation: Topocentric Horizontal -> User-defined Topocentric Coord System.
A flux driver for the Bartol Atmospheric Neutrino Flux.
A flux driver for the FLUKA 3-D Atmospheric Neutrino Flux.
A driver for the HAKKM 3-D atmospheric neutrino flux (commonly known as the ‘Honda flux’)
The PointGeomAnalyzer class is the simplest implementation of the GeomAnalyserI interface and defines...
A ROOT/GEANT4 geometry driver.
virtual TGeoManager * GetGeometry(void) const
virtual double LengthUnits(void) const
double kDefOptEvMax
GeomAnalyzerI * GetGeometry(void)
map< int, string > gOptFluxFiles
long int gOptRanSeed
double gOptGeomLUnits
double gOptEvMin
double gOptEvMax
double gOptKtonYrExposure
string kDefOptEvFilePrefix
string kDefOptGeomLUnits
NtpMCFormat_t kDefOptNtpFormat
bool gOptUsingRootGeom
double kDefOptEvMin
Long_t gOptRunNu
string kDefOptGeomDUnits
TRotation gOptRot
string gOptExtMaxPlXml
double gOptSecExposure
string gOptFluxSim
int gOptNev
map< int, double > gOptTgtMix
double gOptRT
string gOptRootGeom
void GetCommandLineArgs(int argc, char **argv)
void PrintSyntax(void)
GAtmoFlux * GetFlux(void)
double gOptGeomDUnits
double gOptRL
string gOptEvFilePrefix
string gOptInpXSecFile
string gOptRootGeomTopVol
string lunits
string geom
GENIE flux drivers.
Utilities for improving the code readability when using PDG codes.
static constexpr double GeV
Definition Units.h:28
void XSecTable(string inpfile, bool require_table)
Definition AppInit.cxx:38
void RandGen(long int seed)
Definition AppInit.cxx:30
void MesgThresholds(string inpfile)
Definition AppInit.cxx:99
void CacheFile(string inpfile)
Definition AppInit.cxx:117
void RecursiveExhaust(TGeoVolume *topvol, string volnames, bool exhaust)
Definition GeoUtils.cxx:16
string PrintFramedMesg(string mesg, unsigned int nl=1, const char f=' *')
vector< string > Split(string input, string delim)
bool FileExists(string filename)
double UnitFromString(string u)
Definition UnitUtils.cxx:18
THE MAIN GENIE PROJECT NAMESPACE
Definition AlgCmp.h:25
bool gAbortingInErr
Definition Messenger.cxx:34
@ kNFGHEP
Definition NtpMCFormat.h:30
enum genie::ENtpMCFormat NtpMCFormat_t