GENIEGenerator
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
gFNALExptEvGen.cxx
Go to the documentation of this file.
1 //______________________________________________________________________________
2 /*!
3 
4 \program gevgen_fnal
5 
6 \brief A GENIE event generation driver `customized' for the FNAL neutrino
7  experiments using flux drivers for file types used by those expts.
8  This program was adapted from gevgen_t2k used at T2K.
9 
10  This driver can use either the FNAL-supported neutrino flux ntuples
11  (generated by gNuMI, g4numi, g4numi_flugg), or "dk2nu" format,
12  or plain flux histograms for all neutrino species you are considering.
13  You can specify either a ROOT-based detailed detector geometry
14  description or a simple target mix. See below for further details.
15 
16  Users should note that the generic GENIE event generation driver
17  (gevgen) may still be a more appropriate tool to use for the simpler
18  event generation casesrequired for many 4-vector level / systematic
19  studies.
20  Please see the GENIE documentation (http://www.genie-mc.org) and
21  contact me <c.andreopoulos \at cern.ch> if in doubt.
22 
23  *** Synopsis :
24 
25  gevgen_fnal [-h]
26  [-r run#]
27  -f flux
28  -g geometry
29  [-t top_volume_name_at_geom || -t +Vol1-Vol2...]
30  [-m max_path_lengths_xml_file]
31  [-L length_units_at_geom]
32  [-D density_units_at_geom]
33  [-n n_of_events]
34  [-e exposure_in_POTs]
35  [-o output_event_file_prefix]
36  [-F fid_cut_string]
37  [-S nrays]
38  [-z zmin]
39  [-d debug flags]
40  [--seed random_number_seed]
41  --cross-sections xml_file
42 
43  // command line args handled by RunOpt:
44  [--event-generator-list list_name] // default "Default"
45  [--tune tune_name] // default "G18_02a_00_000"
46  [--xml-path path]
47  [--message-thresholds xml_file]
48  [--event-record-print-level level]
49  [--mc-job-status-refresh-rate rate]
50  [--cache-file root_file]
51  [--enable-bare-xsec-pre-calc]
52  [--disable-bare-xsec-pre-calc]
53  [--unphysical-event-mask mask]
54 
55  *** Options :
56 
57  [] Denotes an optional argument
58 
59  -h
60  Prints out the gevgen_fnal syntax and exits.
61  -r
62  Specifies the MC run number [default: 1000]
63  -g
64  Input 'geometry'.
65  This option can be used to specify any of:
66  1 > A ROOT file containing a ROOT/GEANT geometry description
67  [Examples]
68  - To use the master volume from the ROOT geometry stored
69  in the /some/path/nova-geom.root file, type:
70  '-g /some/path/nova-geom.root'
71  2 > A mix of target materials, each with its corresponding weight,
72  typed as a comma-separated list of nuclear PDG codes (in the
73  std PDG2006 convention: 10LZZZAAAI) with the weight fractions
74  in brackets, eg code1[fraction1],code2[fraction2],...
75  If that option is used (no detailed input geometry description)
76  then the interaction vertices are distributed in the detector
77  by the detector MC.
78  [Examples]
79  - To use a target mix of 95% O16 and 5% H type:
80  '-g 1000080160[0.95],1000010010[0.05]'
81  - To use a target which is 100% C12, type:
82  '-g 1000060120'
83  -t
84  Input 'top volume' for event generation -
85  can be used to force event generation in given sub-detector.
86  [default: the 'master volume' of the input geometry]
87 
88  You can also use the -t option to switch generation on/off at
89  multiple volumes as, for example, in:
90  `-t +Vol1-Vol2+Vol3-Vol4',
91  `-t "+Vol1 -Vol2 +Vol3 -Vol4"',
92  `-t -Vol2-Vol4+Vol1+Vol3',
93  `-t "-Vol2 -Vol4 +Vol1 +Vol3"'m
94  where:
95  "+Vol1" and "+Vol3" tells GENIE to `switch on' Vol1 and Vol3, while
96  "-Vol2" and "-Vol4" tells GENIE to `switch off' Vol2 and Vol4.
97  If the very first character is a '+', GENIE will neglect all volumes
98  except the ones explicitly turned on. Vice versa, if the very first
99  character is a `-', GENIE will keep all volumes except the ones
100  explicitly turned off (feature contributed by J.Holeczek).
101 
102  -m
103  An XML file (generated by gmxpl) with the max (density weighted)
104  path-lengths for each target material in the input ROOT geometry.
105  If no file is input, then the geometry will be scanned at MC job
106  initialization to determine those max path lengths.
107  Supplying this file can speed-up the MC job initialization.
108  -L
109  Input geometry length units, eg "m", "cm", "mm", ...
110  [default: mm]
111  -D
112  Input geometry density units, eg "g_cm3", "clhep_def_density_unit",...
113  [default: g_cm3]
114  -f
115  Input 'neutrino flux'.
116  This option can be used to specify any of:
117  1 > A g[4][numi|lbne][_flugg] beam simulation output file
118  and the detector location
119  The general sytax is:
120  -f /full/path/flux_file.root,detector,flavor1,flavor2...
121  [Notes]
122  - For more information on the flux ntuples, see the gNuMI doc.
123  - The original hbook ntuples need to be converted to a ROOT
124  format using the h2root ROOT utility.
125  - If flavors aren't specified then use default (12,-12,14,-14)
126  - See GNuMIFlux.xml or Dk2Nu.xml for all supported detector
127  locations
128  - The g[4][NuMI|LBNE][_flugg] flux ntuples are read via GENIE's
129  GNuMIFlux driver, and dk2nu file via an external
130  product w/ the GDk2NuFlux driver (if it can be found).
131  This customized GENIE event generation driver passes-through
132  the complete gNuMI input flux information (eg parent decay
133  kinematics / position etc) for each neutrino event it
134  generates (as an additional 'flux' branch is added on the
135  output event tree).
136  [Examples]
137  - To use the gNuMI flux ntuple flux.root at MINOS near
138  detector location, type:
139  '-f /path/flux.root,MINOS-NearDet'
140  1a> Similar to 1 above, but filename contains "dk2nu", then use
141  the GDk2NuFlux driver
142  1b> Similar to 1 above, but filename contains "gsimple", then
143  use GSimpleNtpFlux driver
144  2 > A set of histograms stored in a ROOT file.
145  The general syntax is:
146  -f /path/histogram_file.root,neutrino_code[histo_name],...
147  [Notes]
148  - The neutrino codes are the PDG ones.
149  - The 'neutrino_code[histogram_name]' part of the option can
150  be repeated multiple times (separated by commas), once for
151  each flux neutrino species you want to consider, eg
152  '-f somefile.root,12[nuehst],-12[nuebarhst],14[numuhst]'
153  - Code implicitly assumes the binning for multiple flavors
154  is the same for all histograms (no checks are made)
155  - Note that the relative normalization of the flux histograms
156  is taken into account and is reflected in the relative
157  frequency of flux neutrinos thrown by the flux driver
158  - Variable bin width flux histograms are modified to account
159  for incorrect sampling of TH1D::GetRandom() prior to ROOT
160  version 9.99.99
161  + notation such as "12[nuehst]WIDTH" forces bin width
162  adjustment no matter variable bin width histogram or not
163  + notation such as "12[nuehst]NOWIDTH" prevents bin width
164  adjustment
165  - By default this is a beam of zero radius
166  originating at (0,0,0) travelling in the direction (0,0,1)
167  To override this append a string of the **exact** form:
168  ";r=1.1,dir=(0.1,0.2,0.3),spot=(-1,-2,-3)"
169  use exactly this layout (in order) with numbers modified
170  - When using flux from histograms there is no branch with
171  neutrino parent information added in the output event
172  tree as your flux input contains no such information.
173  [Examples]
174  - To use the histogram 'h100' (representing the nu_mu flux)
175  and the histogram 'h300' (representing the nu_e flux) and
176  the histogram 'h301' (representing the nu_e_bar flux) from
177  the flux.root file in /path/ type:
178  '-f /path/flux.root,14[h100],12[h300],-12[h301]
179 
180  -e
181  Specifies how many POTs to generate.
182  -n
183  Specifies how many events to generate.
184 
185  -------
186  [Note on exposure / statistics]
187  Both -e and -n options can be used to set the exposure.
188  - If the input flux is a non-histogram driver then any of these
189  options can be used (one at a time).
190  - If the input flux is described with histograms then only the -n
191  option is available.
192  -------
193 
194  -F
195  Apply a fiducial cut (for now hard coded ... generalize)
196  Only used with ROOTGeomAnalyzer
197  if string starts with "-" then reverses sense (ie. anti-fiducial)
198  -S
199  Number of rays to use to scan geometry for max path length
200  Only used with ROOTGeomAnalyzer & { GNuMIFlux, GSimpleNtpFlux, GDk2NuFlux }
201  +N Use flux to scan geometry for max path length
202  -N Use N rays x N points on each face of a box
203  -z
204  Z from which to start flux ray in user world coordinates
205  Only use with ROOTGeomAnalyzer & { GNuMIFlux, GSimpleNtpFlux, GDk2NuFlux }
206  If left unset then flux originates on the flux window
207  [No longer attempts to determine z from geometry, generally
208  got this wrong]
209  -o
210  Sets the prefix of the output event file.
211  The output filename is built as:
212  [prefix].[run_number].[event_tree_format].[file_format]
213  The default output filename is:
214  gntp.[run_number].ghep.root
215  This cmd line arguments lets you override 'gntp'
216  --seed
217  Random number seed.
218  --cross-sections
219  Name (incl. full path) of an XML file with pre-computed
220  cross-section values used for constructing splines.
221 
222  --event-generator-list
223  List of event generators to load in event generation drivers.
224  [default: "Default"].
225  --tune
226  Specifies a GENIE comprehensive neutrino interaction model tune.
227  [default: "Default"].
228  --xml-path
229  A directory to load XML files from - overrides $GXMLPATH, and $GENIE/config
230  --message-thresholds
231  Allows users to customize the message stream thresholds.
232  The thresholds are specified using an XML file.
233  See $GENIE/config/Messenger.xml for the XML schema.
234  --event-record-print-level
235  Allows users to set the level of information shown when the event
236  record is printed in the screen. See GHepRecord::Print().
237  --mc-job-status-refresh-rate
238  Allows users to customize the refresh rate of the status file.
239  --cache-file
240  Allows users to specify a cache file so that the cache can be
241  re-used in subsequent MC jobs.
242  --unphysical-event-mask
243  Allows users to specify a 16-bit mask to allow certain types of
244  unphysical events to be written in the output file.
245  [default: all unphysical events are rejected]
246 
247 
248  *** Examples:
249 
250  (1) shell% gevgen_fnal
251  -r 1001
252  -f /data/mc_inputs/flux/flux_00001.root,MINOS-NearDet,12,-12
253  -g /data/mc_inputs/geom/minos.root
254  -L mm -D g_cm3
255  -e 5E+17
256  --cross-sections /data/xsec.xml
257 
258  Generate events (run number 1001) using the flux ntuple in
259  /data/mc_inputs/flux/v1/flux_00001.root
260  The job will load the MINOS near detector detector geometry
261  description from /data/mc_inputs/geom/minos.root and interpret it
262  using 'mm' as the length unit and 'g/cm^3' as the density unit.
263  The job will stop as it accumulates a sample corresponding to
264  5E+17 POT.
265  Pre-computed cross-section data are loaded from /data/xsec.xml
266 
267  (2) shell% gevgen_fnal
268  -r 1001
269  -f /data/t2k/flux/hst/flux.root,12[h100],-12[h101],14[h200]
270  -g 1000080160[0.95],1000010010[0.05]
271  -n 50000
272  --cross-sections /data/xsec.xml
273 
274  Please read the GENIE user manual for more information.
275 
276 \author Costas Andreopoulos <c.andreopoulos \at cern.ch>
277  University of Liverpool
278 
279  Robert Hatcher <rhatcher \at fnal.gov>
280  Fermi National Accelerator Laboratory
281 
282 \created August 20, 2008
283 
284 \cpright Copyright (c) 2003-2024, The GENIE Collaboration
285  For the full text of the license visit http://copyright.genie-mc.org
286 
287 */
288 //_________________________________________________________________________________________
289 
290 #include <cassert>
291 #include <cstdlib>
292 #include <csignal>
293 
294 #include <string>
295 #include <sstream>
296 #include <vector>
297 #include <map>
298 #include <algorithm> // for transform()
299 #include <fstream>
300 
301 #include <TSystem.h>
302 #include <TError.h> // for gErrorIgnoreLevel
303 #include <TTree.h>
304 #include <TFile.h>
305 #include <TH1D.h>
306 #include <TMath.h>
307 #include <TGeoVolume.h>
308 #include <TGeoShape.h>
309 
325 #include "Framework/Utils/AppInit.h"
326 #include "Framework/Utils/RunOpt.h"
330 
331 #ifdef __GENIE_FLUX_DRIVERS_ENABLED__
336 
337 //#include "Tools/Flux/GNuMIFlux.h"
338 //#include "Tools/Flux/GSimpleNtpFlux.h"
339 // #ifdef __DK2NU_FLUX_DRIVER_AVAILABLE__
340 // #include "dk2nu/tree/dk2nu.h"
341 // #include "dk2nu/tree/dkmeta.h"
342 // #include "dk2nu/tree/NuChoice.h"
343 // #include "dk2nu/genie/GDk2NuFlux.h"
344 // #endif
345 
346 #endif
347 
348 #ifdef __GENIE_GEOM_DRIVERS_ENABLED__
349 #include "Tools/Geometry/GeoUtils.h"
354 #endif
355 
356 using std::string;
357 using std::vector;
358 using std::map;
359 using std::ostringstream;
360 
361 using namespace genie;
362 
363 // Forward function declarations
364 //
365 void LoadExtraOptions (void);
366 void GetCommandLineArgs (int argc, char ** argv);
367 void PrintSyntax (void);
368 void CreateFidSelection (string fidcut, GeomAnalyzerI* geom_driver);
369 void CreateRockBoxSelection (string fidcut, GeomAnalyzerI* geom_driver);
370 void DetermineFluxDriver(string fopt);
371 void ParseFluxHst (string fopt);
372 void ParseFluxFileConfig(string fopt);
373 
374 // Default options (override them using the command line arguments):
375 //
376 string kDefOptGeomLUnits = "mm"; // default geometry length units
377 string kDefOptGeomDUnits = "g_cm3"; // default geometry density units
378 NtpMCFormat_t kDefOptNtpFormat = kNFGHEP; // default event tree format
379 string kDefOptEvFilePrefix = "gntp";
380 
381 // User-specified options:
382 //
383 Long_t gOptRunNu; // run number
384 bool gOptUsingRootGeom = false; // using root geom or target mix?
385 bool gOptUsingHistFlux = false; // using beam flux ntuples or flux from histograms?
386 string gOptFluxDriver = ""; // flux driver class to use
387 map<string,string> gOptFluxShortNames;
388 PDGCodeList gOptFluxPdg; // list of neutrino flavors to accept
389 map<int,double> gOptTgtMix; // target mix (tgt pdg -> wght frac) / if not using detailed root geom
390 
391 map<int,TH1D*> gOptFluxHst; // flux histos (nu pdg -> spectrum) / if not using beam sim ntuples
392 TVector3 gOptBeamDir = TVector3(0,0,1);
393 TVector3 gOptBeamSpot = TVector3(0,0,0);
394 Double_t gOptBeamRadius = -1; // meters
395 string gOptBeamRtDependence = "x";
396 
397 string gOptRootGeom; // input ROOT file with realistic detector geometry
398 string gOptRootGeomTopVol = ""; // input geometry top event generation volume
399 string gOptRootGeomMasterVol = ""; // (highest level of geometry)
400 double gOptGeomLUnits = 0; // input geometry length units
401 double gOptGeomDUnits = 0; // input geometry density units
402 string gOptExtMaxPlXml = ""; // max path lengths XML file for input geometry
403 bool gOptWriteMaxPlXml = false; // rather than read file, write the file
404  // triggered by leading '+' on given ofilename
405 string gOptFluxFile; // ROOT file with beam flux ntuple
406 string gOptDetectorLocation; // detector location (see GNuMIFlux.xml for supported locations))
407 int gOptNev; // number of events to generate
408 double gOptPOT; // exposure (in POT)
409 string gOptFidCut; // fiducial cut selection
410 int gOptNScan = 0; // # of geometry scan rays
411 double gOptZmin = -2.0e30; // starting z position [ if abs() < 1e30 ]
412 string gOptEvFilePrefix; // event file prefix
413 int gOptDebug = 0; // debug flags
414 long int gOptRanSeed; // random number seed
415 string gOptInpXSecFile; // cross-section splines
416 
417 bool gSigTERM = false; // was TERM signal sent?
418 
419 static void gsSIGTERMhandler(int /* s */)
420 {
421  gSigTERM = true;
422  std::cerr << "Caught SIGTERM" << std::endl;
423 }
424 
425 //____________________________________________________________________________
426 int main(int argc, char ** argv)
427 {
429  GetCommandLineArgs(argc,argv);
430 
431  if ( ! RunOpt::Instance()->Tune() ) {
432  LOG("gmkspl", pFATAL) << " No TuneId in RunOption";
433  exit(-1);
434  }
436 
437  // Initialization of random number generators, cross-section table,
438  // messenger thresholds, cache file
439  utils::app_init::MesgThresholds(RunOpt::Instance()->MesgThresholdFiles());
443 
444  // Set GHEP print level
445  int print_level = RunOpt::Instance()->EventRecordPrintLevel();
446  GHepRecord::SetPrintLevel(print_level);
447 
448  // *************************************************************************
449  // * Create / configure the geometry driver
450  // *************************************************************************
451  GeomAnalyzerI * geom_driver = 0;
452 
453  if(gOptUsingRootGeom) {
454  //
455  // *** Using a realistic root-based detector geometry description
456  //
457 
458  // creating & configuring a root geometry driver
461  TGeoVolume * topvol = rgeom->GetGeometry()->GetTopVolume();
462  if ( ! topvol ) {
463  LOG("gevgen_fnal", pFATAL) << "Null top ROOT geometry volume!";
464  exit(1);
465  }
466  // retrieve the master volume name
467  gOptRootGeomMasterVol = topvol->GetName();
468 
469  rgeom -> SetLengthUnits (gOptGeomLUnits);
470  rgeom -> SetDensityUnits (gOptGeomDUnits);
471  rgeom -> SetTopVolName (gOptRootGeomTopVol); // set user defined "topvol"
472 
473  // getting the bounding box dimensions along z so as to set the
474  // appropriate upstream generation surface for the NuMI flux driver
475 
476  // RWH 2010-07-16: do not try to automatically get zmin from geometry, rather
477  // by default let the flux start from the window. If the user wants to
478  // override this then they need to explicitly set a "zmin". Trying to use
479  // the geometry is fraught with problems in local vs. global coordinates and
480  // units where it can appear to work in some cases but it actually isn't really
481  // universally correct.
482  //was// TGeoShape * bounding_box = topvol->GetShape();
483  //was// bounding_box->GetAxisRange(3, zmin, zmax);
484  //was// zmin *= rgeom->LengthUnits();
485  //was// zmax *= rgeom->LengthUnits();
486 
487  // switch on/off volumes as requested
488  if ( (gOptRootGeomTopVol[0] == '+') || (gOptRootGeomTopVol[0] == '-') ) {
489  bool exhaust = (*gOptRootGeomTopVol.c_str() == '+');
491  }
492 
493  // casting to the GENIE geometry driver interface
494  geom_driver = dynamic_cast<GeomAnalyzerI *> (rgeom);
495 
496  // user specifid a fiducial volume cut ... parse that out
497  if ( gOptFidCut.find("rock") != std::string::npos )
499  else if ( gOptFidCut != "" ) CreateFidSelection(gOptFidCut,rgeom);
500 
501  }
502  else {
503  //
504  // *** Using a 'point' geometry with the specified target mix
505  // *** ( = a list of targets with their corresponding weight fraction)
506  //
507 
508  // creating & configuring a point geometry driver
511  // casting to the GENIE geometry driver interface
512  geom_driver = dynamic_cast<GeomAnalyzerI *> (pgeom);
513  }
514 
515  // *************************************************************************
516  // * Create / configure the flux driver
517  // *************************************************************************
518  GFluxI * flux_driver =
520  if ( ! flux_driver ) {
521  // couldn't get the requested flux driver ?
522  std::ostringstream s;
523  s << "Known FluxDrivers:\n";
524  const std::vector<std::string>& known =
526  std::vector<std::string>::const_iterator itr = known.begin();
527  for ( ; itr != known.end(); ++itr ) s << " " << (*itr) << "\n";
528  LOG("gevgen_fnal", pFATAL)
529  << "Failed to get any flux driver from GFluxDriverFactory\n"
530  << "when using \"" << gOptFluxDriver << "\"\n" << s.str();
531  exit(1);
532  }
533 
534  if ( ! gOptUsingHistFlux ) {
535  genie::flux::GFluxFileConfigI* flux_file_config =
536  dynamic_cast<genie::flux::GFluxFileConfigI*>(flux_driver);
537  if ( ! flux_file_config ) {
538  LOG("gevgen_fnal", pFATAL)
539  << "Failed to get GFluxFileConfigI driver from GFluxDriverFactory\n"
540  << "when using \"" << gOptFluxDriver << "\"";
541  exit(1);
542  }
543 
544  //
545  // *** Using the detailed ntuple neutrino flux description
546  //
548  flux_file_config->SetUpstreamZ(gOptZmin); // was "zmin" from bounding_box
549  flux_file_config->SetNumOfCycles(0);
550 
551  if ( gOptFluxPdg.size() > 0 ) {
552  // user specified list of neutrino PDGs
553  flux_file_config->SetFluxParticles(gOptFluxPdg);
554  std::ostringstream s;
555  PDGCodeList::const_iterator itr = gOptFluxPdg.begin();
556  for ( ; itr != gOptFluxPdg.end(); ++itr) s << (*itr) << " ";
557  LOG("gevgen_fnal", pNOTICE)
558  << "Limiting to nu PDGs: " << s.str();
559  }
560  }
561  else {
562  //
563  // *** Using fluxes from histograms (for all specified neutrino species)
564  //
565  genie::flux::GCylindTH1Flux * hist_flux_driver =
566  dynamic_cast<genie::flux::GCylindTH1Flux*>(flux_driver);
567  if ( ! hist_flux_driver ) {
568  LOG("gevgen_fnal", pFATAL)
569  << "Failed to get GCylinderTH1Flux driver from GFluxDriverFactory\n"
570  << "when using " << gOptFluxDriver;
571  exit(1);
572  }
573 
574  // creating & configuring a generic GCylindTH1Flux flux driver
575  hist_flux_driver->SetNuDirection (gOptBeamDir);
576  hist_flux_driver->SetBeamSpot (gOptBeamSpot);
577  hist_flux_driver->SetTransverseRadius (gOptBeamRadius);
578  //hist_flux_driver->SetRtDependence (gOptBeamRtDependence); // "x");
579 
580  map<int,TH1D*>::iterator it = gOptFluxHst.begin();
581  for( ; it != gOptFluxHst.end(); ++it) {
582  int pdg_code = it->first;
583  TH1D * spectrum = it->second;
584  hist_flux_driver->AddEnergySpectrum(pdg_code, spectrum);
585  // once the histogram has been added to the GCylindTH1Flux driver
586  // it is owned by the driver and it is up to the the driver
587  // to clean up (i.e. delete it).
588  // remove it from this map to avoid double deletion.
589  it->second = 0;
590  }
591  }
592 
593  // these might come in handy ... avoid repeated dynamic_cast<> calls
594  genie::flux::GFluxExposureI* fluxExposureI =
595  dynamic_cast<genie::flux::GFluxExposureI*>(flux_driver);
596  genie::flux::GFluxFileConfigI* fluxFileConfigI =
597  dynamic_cast<genie::flux::GFluxFileConfigI*>(flux_driver);
598 
599 
600  // *************************************************************************
601  // * Handle chicken/egg problem: geom analyzer vs. flux.
602  // * Need both at this point change geom scan defaults.
603  // *************************************************************************
605 
607  dynamic_cast<geometry::ROOTGeomAnalyzer *>(geom_driver);
608  if ( ! rgeom ) assert(0);
609 
610  rgeom -> SetDebugFlags(gOptDebug);
611 
612  // even if user doesn't specify gOptNScan configure to scan using flux
613  if ( gOptNScan >= 0 ) {
614  LOG("gevgen_fnal", pNOTICE)
615  << "Using ROOTGeomAnalyzer: geom scan using flux: nparticles=" << gOptNScan;
616  rgeom->SetScannerFlux(flux_driver);
617  if ( gOptNScan > 0 ) rgeom->SetScannerNParticles(gOptNScan);
618  } else {
619  int nabs = TMath::Abs(gOptNScan);
620  LOG("gevgen_fnal", pNOTICE)
621  << "Using ROOTGeomAnalyzer: geom scan using box: npoints=nrays=" << nabs;
622  rgeom->SetScannerNPoints(nabs);
623  rgeom->SetScannerNRays(nabs);
624  }
625  }
626 
627  // *************************************************************************
628  // * Create/configure the event generation driver
629  // *************************************************************************
630  GMCJDriver * mcj_driver = new GMCJDriver;
632  mcj_driver->UseFluxDriver(flux_driver);
633  mcj_driver->UseGeomAnalyzer(geom_driver);
634  if ( ( gOptExtMaxPlXml != "" ) && ! gOptWriteMaxPlXml ) {
635  mcj_driver->UseMaxPathLengths(gOptExtMaxPlXml);
636  }
637  mcj_driver->Configure();
638  mcj_driver->UseSplines();
639  mcj_driver->ForceSingleProbScale();
640 
641  if ( ( gOptExtMaxPlXml != "" ) && gOptWriteMaxPlXml ) {
643  dynamic_cast<geometry::ROOTGeomAnalyzer *>(geom_driver);
644  if ( rgeom ) {
645  const genie::PathLengthList& maxpath = rgeom->GetMaxPathLengths();
646  std::string maxplfile = gOptExtMaxPlXml;
647  maxpath.SaveAsXml(maxplfile);
648  // append extra info to file
649  std::ofstream mpfile(maxplfile.c_str(), std::ios_base::app);
650  mpfile
651  << std::endl
652  << "<!-- this file is only relevant for a setup compatible with:"
653  << std::endl
654  << "geom: " << gOptRootGeom << " top: \"" << gOptRootGeomTopVol << "\""
655  << std::endl
656  << "flux: " << gOptFluxFile
657  << std::endl
658  << "location: " << gOptDetectorLocation
659  << std::endl
660  << "fidcut: " << gOptFidCut
661  << std::endl
662  << "nscan: " << gOptNScan << " (0>= use flux, <0 use box |nscan| points/rays)"
663  << std::endl
664  << "zmin: " << gOptZmin << " (if |zmin| > 1e30, leave ray on flux window)"
665  << std::endl
666  << "-->"
667  << std::endl;
668  mpfile.close();
669  }
670  }
671 
672  // *************************************************************************
673  // * Prepare for writing the output event tree & status file
674  // *************************************************************************
675 
676  // Initialize an Ntuple Writer to save GHEP records into a TTree
679  ntpw.Initialize();
680 
681 
682  std::vector<TBranch*> extraBranches;
683  std::vector<std::string> branchNames;
684  std::vector<std::string> branchClassNames;
685  std::vector<void**> branchObjPointers;
686 
687  // Add custom branch(s) to the standard GENIE event tree so that
688  // info on the flux neutrino parent particle can be passed-through
689  if ( fluxFileConfigI ) {
690  fluxFileConfigI->GetBranchInfo(branchNames,branchClassNames,
691  branchObjPointers);
692  size_t nn = branchNames.size();
693  size_t nc = branchClassNames.size();
694  size_t np = branchObjPointers.size();
695  if ( nn != nc || nc != np ) {
696  LOG("gevgen_fnal", pERROR)
697  << "Inconsistent info back from \"" << gOptFluxDriver << "\" "
698  << "for branch info: " << nn << " " << nc << " " << np;
699  } else {
700  for (size_t ii = 0; ii < nn; ++ii) {
701  const char* bname = branchNames[ii].c_str();
702  const char* cname = branchClassNames[ii].c_str();
703  void**& optr = branchObjPointers[ii]; // note critical '&' !
704  if ( ! optr || ! *optr ) continue; // no pointer supplied, skip it
705  int split = 99; // 1
706  LOG("gevgen_fnal", pNOTICE)
707  << "Adding extra branch \"" << bname << "\" of type \""
708  << cname << "\" (" << optr << ") to output tree";
709  TBranch* bptr = ntpw.EventTree()->Branch(bname,cname,optr,32000,split);
710  extraBranches.push_back(bptr);
711 
712  if ( bptr ) {
713  // don't delete this !!! we're sharing
714  bptr->SetAutoDelete(false);
715  } else {
716  LOG("gevgen_fnal", pERROR)
717  << "FAILED to add extra branch \"" << bname << "\" of type \""
718  << cname << "\" to output tree";
719  }
720  } // loop over additions
721  } // same # of entries
722  } // of genie::flux::GFluxFileConfigI type
723 
724  // Create a MC job monitor for a periodically updated status file
725  GMCJMonitor mcjmonitor(gOptRunNu);
726  mcjmonitor.SetRefreshRate(RunOpt::Instance()->MCJobStatusRefreshRate());
727 
728  // *************************************************************************
729  // * Event generation loop
730  // *************************************************************************
731 
732  // define handler to allow signal to end job gracefully
733  signal(SIGTERM,gsSIGTERMhandler);
734 
735  int ievent = 0;
736  while ( ! gSigTERM )
737  {
738  LOG("gevgen_fnal", pINFO)
739  << " *** Generating event............ " << ievent;
740 
741  // In case the required statistics was expressed as 'number of events'
742  // then quit if that number has been generated
743  if ( ievent == gOptNev ) break;
744 
745  // In case the required statistics was expressed as 'number of POT'
746  // then exit the event loop if the requested POT has been generated.
747  if ( gOptPOT > 0 && fluxExposureI ) {
748  double fpot = fluxExposureI->GetTotalExposure(); // current POTs used
749  double psc = mcj_driver->GlobProbScale(); // interaction prob. scale
750  double pot = fpot / psc; // POTs for generated sample
751  if ( pot >= gOptPOT ) break;
752  }
753 
754  // Generate a single event using neutrinos coming from the specified flux
755  // and hitting the specified geometry or target mix
756  EventRecord * event = mcj_driver->GenerateEvent();
757 
758  // Check whether a null event was returned due to the flux driver reaching
759  // the end of the input flux ntuple - exit the event generation loop
760  if ( ! event && flux_driver->End() ) {
761  LOG("gevgen_fnal", pWARN)
762  << "** The flux driver read all the input flux entries: End()==true";
763  break;
764  }
765  if ( ! event ) {
766  LOG("gevgen_fnal", pERROR)
767  << "Got a null generated neutino event! Retrying ...";
768  continue;
769  }
770  if ( print_level >= 0 ) {
771  LOG("gevgen_fnal", pINFO)
772  << "Generated event: " << *event;
773  }
774  // A valid event was generated: flux info (parent decay/prod
775  // position/kinematics) for that simulated event should already
776  // be connected to the right output tree branch
777 
778  // Add event at the output ntuple, refresh the mc job monitor & clean-up
779  ntpw.AddEventRecord(ievent, event);
780  mcjmonitor.Update(ievent,event);
781  delete event;
782  ievent++;
783 
784  } //1
785 
786  // Copy metadata tree, if available
787  if ( fluxFileConfigI ) {
788  TTree* t1 = fluxFileConfigI->GetMetaDataTree();
789  if ( t1 ) {
790  TTree* t2 = (TTree*)t1->CloneTree();
791  t2->Write();
792  }
793  }
794 
795  LOG("gevgen_fnal", pINFO)
796  << "The GENIE MC job is done generating events - Cleaning up & exiting...";
797 
798  // *************************************************************************
799  // * Print job statistics &
800  // * calculate normalization factor for the generated sample
801  // *************************************************************************
803  // POT normalization will only be calculated if event generation was based
804  // on beam simulation ntuples (not just histograms) & a detailed detector
805  // geometry description.
806  // Get nunber of flux neutrinos read-in by flux driver, number of flux
807  // neutrinos actually thrown to the event generation driver and number
808  // of neutrino interactions actually generated
809  long int nflx = 0;
810  long int nflx_evg = mcj_driver-> NFluxNeutrinos();
811  double fpot = 0;
812  const char* exposureUnits = "(unknown units)";
813  if ( fluxExposureI ) {
814  fpot = fluxExposureI->GetTotalExposure(); // POTs used so far
815  nflx = fluxExposureI->NFluxNeutrinos();
816  //genie::flux::Exposure_t etype = fluxExposureI->GetExposureType();
817  //exposureUnits = genie::flux::GFluxExposureI::AsString(etype);
818  exposureUnits = fluxExposureI->GetExposureUnits();
819  }
820  if ( fluxFileConfigI ) {
821  fluxFileConfigI->PrintConfig();
822  }
823  double psc = mcj_driver->GlobProbScale(); // interaction prob. scale
824  if ( psc <= 0.0 ) {
825  LOG("gevgen_fnal", pFATAL) << "MCJobDriver GlobalProbScale was " << psc;
826  }
827  double pot = fpot / psc; // POT for generated sample
828  long int nev = ievent;
829 
830  LOG("gevgen_fnal", pNOTICE)
831  << "\n >> Interaction probability scaling factor: " << psc
832  << "\n >> using: " << gOptFluxDriver
833  << "\n >> N of flux v read-in by flux driver: " << nflx
834  << "\n >> N of flux v thrown to event gen driver: " << nflx_evg
835  << "\n >> N of generated v interactions: " << nev
836  << "\n ** Normalization for generated sample: " << pot
837  << " " << exposureUnits << " * detector";
838 
839  ntpw.EventTree()->SetWeight(pot); // store POT
840 
841  }
842 
843  // *************************************************************************
844  // * Save & clean-up
845  // *************************************************************************
846 
847  // Save the generated event tree & close the output file
848  ntpw.Save();
849 
850  // Clean-up
851  delete geom_driver;
852  delete flux_driver;
853  delete mcj_driver;
854  // this list should only be histograms that have (for some reason)
855  // not been handed over to the GCylindTH1Flux driver.
856  map<int,TH1D*>::iterator it = gOptFluxHst.begin();
857  for( ; it != gOptFluxHst.end(); ++it) {
858  TH1D * spectrum = it->second;
859  if(spectrum) delete spectrum;
860  }
861  gOptFluxHst.clear();
862 
863  LOG("gevgen_fnal", pNOTICE) << "Done!";
864 
865  return 0;
866 }
867 
868 //____________________________________________________________________________
869 void LoadExtraOptions(void)
870 {
871  /// potentially load extra libraries that might extend the list of
872  /// potential flux drivers, and how to map short names to classes ...
873  // we really should at this point read some external file to get
874  // an expandable list of libraries ... but for now just hard code it
875 
876  vector<string> extraLibs;
877 
878  ///***** this part should come from reading an external file
879  /// placeholder file $GENIE/config/FluxDriverExpansion.xml
880 
881  extraLibs.push_back("libdk2nuTree");
882  extraLibs.push_back("libdk2nuGenie");
883 
884  // what one might expect to find at the beginning of -f <arg>
885 
886  gOptFluxShortNames["histo"] = "genie::flux::GCylindTH1Flux";
887  gOptFluxShortNames["hist"] = "genie::flux::GCylindTH1Flux";
888 
889  gOptFluxShortNames["simple"] = "genie::flux::GSimpleNtpFlux";
890  gOptFluxShortNames["numi"] = "genie::flux::GNuMIFlux";
891  gOptFluxShortNames["dk2nu"] = "genie::flux::GDk2NuFlux";
892 
893  ///******* done with fake "read"
894 
895  // see if there are any libraries to load
896  // just let ROOT do it ... check error code on return
897  // but tweak ROOT's ERROR message output so we don't see huge messages
898  // for failures
899  // gErrorIgnoreLevel to kError (from TError.h)
900 
901  Int_t prevErrorLevel = gErrorIgnoreLevel;
902  gErrorIgnoreLevel = kFatal;
903  for (size_t ilib=0; ilib < extraLibs.size(); ++ilib) {
904  // library names should be like libXYZZY without extension [e.g. .so]
905  // but with the leading "lib"
906  int loadStatus = gSystem->Load(extraLibs[ilib].c_str());
907  const char* statWords = "failed to load";
908  if ( loadStatus == 0 ) { statWords = "successfully loaded"; }
909  else if ( loadStatus == 1 ) { statWords = "already had loaded"; }
910  else if ( loadStatus == -1 ) { statWords = "couldn't find library"; }
911  else if ( loadStatus == -2 ) { statWords = "mismatched version"; }
912 
913  LOG("gevgen_fnal",pNOTICE)
914  << statWords << " (" << loadStatus << ") " << extraLibs[ilib];
915  }
916  // restore the ROOT error message level
917  gErrorIgnoreLevel = prevErrorLevel;
918 
919 }
920 
921 //____________________________________________________________________________
922 void GetCommandLineArgs(int argc, char ** argv)
923 {
924  LOG("gevgen_fnal", pINFO) << "Parsing command line arguments";
925 
926  // Common run options. Set defaults and read.
929 
930  // Parse run options for this app
931 
932  CmdLnArgParser parser(argc,argv);
933 
934  // help?
935  bool help = parser.OptionExists('h');
936  if(help) {
937  PrintSyntax();
938  exit(0);
939  }
940 
941  // run number:
942  if ( parser.OptionExists('r') ) {
943  LOG("gevgen_fnal", pDEBUG) << "Reading MC run number";
944  gOptRunNu = parser.ArgAsLong('r');
945  } else {
946  LOG("gevgen_fnal", pDEBUG)
947  << "Unspecified run number - Using default";
948  gOptRunNu = 0;
949  } //-r
950 
951  //
952  // *** geometry
953  //
954 
955  string geom = "";
956  string lunits, dunits;
957  if( parser.OptionExists('g') ) {
958  LOG("gevgen_fnal", pDEBUG) << "Getting input geometry";
959  geom = parser.ArgAsString('g');
960 
961  // is it a ROOT file that contains a ROOT geometry?
962  bool accessible_geom_file =
963  ! (gSystem->AccessPathName(geom.c_str()));
964  if (accessible_geom_file) {
965  gOptRootGeom = geom;
966  gOptUsingRootGeom = true;
967  }
968  } else {
969  LOG("gevgen_fnal", pFATAL)
970  << "No geometry option specified - Exiting";
971  PrintSyntax();
972  exit(1);
973  } //-g
974 
975  if(gOptUsingRootGeom) {
976  // using a ROOT geometry - get requested geometry units
977 
978  // legth units:
979  if( parser.OptionExists('L') ) {
980  LOG("gevgen_fnal", pDEBUG)
981  << "Checking for input geometry length units";
982  lunits = parser.ArgAsString('L');
983  } else {
984  LOG("gevgen_fnal", pDEBUG) << "Using default geometry length units";
985  lunits = kDefOptGeomLUnits;
986  } // -L
987  // density units:
988  if( parser.OptionExists('D') ) {
989  LOG("gevgen_fnal", pDEBUG)
990  << "Checking for input geometry density units";
991  dunits = parser.ArgAsString('D');
992  } else {
993  LOG("gevgen_fnal", pDEBUG) << "Using default geometry density units";
994  dunits = kDefOptGeomDUnits;
995  } // -D
998 
999  // check whether an event generation volume name has been
1000  // specified -- default is the 'top volume'
1001  if( parser.OptionExists('t') ) {
1002  LOG("gevgen_fnal", pDEBUG) << "Checking for input volume name";
1003  gOptRootGeomTopVol = parser.ArgAsString('t');
1004  } else {
1005  LOG("gevgen_fnal", pDEBUG) << "Using the <master volume>";
1006  } // -t
1007 
1008  // check whether an XML file with the maximum (density weighted)
1009  // path lengths for each detector material is specified -
1010  // otherwise will compute the max path lengths at job init
1011  // if passed name starts with '+', then compute max at job init, but write out the result
1012  if ( parser.OptionExists('m') ) {
1013  LOG("gevgen_fnal", pDEBUG)
1014  << "Checking for maximum path lengths XML file";
1015  gOptExtMaxPlXml = parser.ArgAsString('m');
1016  gOptWriteMaxPlXml = false;
1017  if ( gOptExtMaxPlXml[0] == '+' ) {
1018  gOptExtMaxPlXml = gOptExtMaxPlXml.substr(1,std::string::npos);
1019  gOptWriteMaxPlXml = true;
1020  LOG("gevgen_fnal", pINFO)
1021  << "Will write maximum path lengths XML file: " << gOptExtMaxPlXml;
1022  }
1023  } else {
1024  LOG("gevgen_fnal", pDEBUG)
1025  << "Will compute the maximum path lengths at job init";
1026  gOptExtMaxPlXml = "";
1027  } // -m
1028 
1029  // fidcut:
1030  if( parser.OptionExists('F') ) {
1031  LOG("gevgen_fnal", pDEBUG) << "Using Fiducial cut?";
1032  gOptFidCut = parser.ArgAsString('F');
1033  } else {
1034  LOG("gevgen_fnal", pDEBUG) << "No fiducial volume cut";
1035  gOptFidCut = "";
1036  } //-F
1037 
1038  if(!gOptUsingHistFlux) {
1039  // how to scan the geometry (if relevant)
1040  if( parser.OptionExists('S') ) {
1041  LOG("gevgen_fnal", pDEBUG) << "Reading requested geom scan count";
1042  gOptNScan = parser.ArgAsInt('S');
1043  } else {
1044  LOG("gevgen_fnal", pDEBUG) << "No geom scan count was requested";
1045  gOptNScan = 0;
1046  } //-S
1047 
1048  // z for flux rays to start
1049  if( parser.OptionExists('z') ) {
1050  LOG("gevgen_fnal", pDEBUG) << "Reading requested zmin";
1051  gOptZmin = parser.ArgAsDouble('z');
1052  } else {
1053  LOG("gevgen_fnal", pDEBUG) << "No zmin was requested";
1054  gOptZmin = -2.0e30; // < -1.0e30 ==> leave it on flux window
1055  } //-z
1056 
1057  // debug flags
1058  if ( parser.OptionExists('d') ) {
1059  LOG("gevgen_fnal", pDEBUG) << "Reading debug flag value";
1060  gOptDebug = parser.ArgAsInt('d');
1061  } else {
1062  LOG("gevgen_fnal", pDEBUG) << "Unspecified debug flags - Using default";
1063  gOptDebug = 0;
1064  } //-d
1065 
1066  } // root geom && gnumi flux
1067 
1068  } // using root geom?
1069 
1070  else {
1071  // User has specified a target mix.
1072  // Decode the list of target pdf codes & their corresponding weight fraction
1073  // (specified as 'pdg_code_1[fraction_1],pdg_code_2[fraction_2],...')
1074  // See documentation on top section of this file.
1075  //
1076  gOptTgtMix.clear();
1077  vector<string> tgtmix = utils::str::Split(geom,",");
1078  if(tgtmix.size()==1) {
1079  int pdg = atoi(tgtmix[0].c_str());
1080  double wgt = 1.0;
1081  gOptTgtMix.insert(map<int, double>::value_type(pdg, wgt));
1082  } else {
1083  vector<string>::const_iterator tgtmix_iter = tgtmix.begin();
1084  for( ; tgtmix_iter != tgtmix.end(); ++tgtmix_iter) {
1085  string tgt_with_wgt = *tgtmix_iter;
1086  string::size_type open_bracket = tgt_with_wgt.find("[");
1087  string::size_type close_bracket = tgt_with_wgt.find("]");
1088  if (open_bracket ==string::npos ||
1089  close_bracket==string::npos)
1090  {
1091  LOG("gevgen_fnal", pFATAL)
1092  << "You made an error in specifying the target mix";
1093  PrintSyntax();
1094  exit(1);
1095  }
1096  string::size_type ibeg = 0;
1097  string::size_type iend = open_bracket;
1098  string::size_type jbeg = open_bracket+1;
1099  string::size_type jend = close_bracket;
1100  int pdg = atoi(tgt_with_wgt.substr(ibeg,iend-ibeg).c_str());
1101  double wgt = atof(tgt_with_wgt.substr(jbeg,jend-jbeg).c_str());
1102  LOG("gevgen_fnal", pDEBUG)
1103  << "Adding to target mix: pdg = " << pdg << ", wgt = " << wgt;
1104  gOptTgtMix.insert(map<int, double>::value_type(pdg, wgt));
1105 
1106  }// tgtmix_iter
1107  } // >1 materials in mix
1108  } // using tgt mix?
1109 
1110  //
1111  // *** flux
1112  //
1113  if ( parser.OptionExists('f') ) {
1114  LOG("gevgen_fnal", pDEBUG) << "Getting input flux";
1115  DetermineFluxDriver(parser.ArgAsString('f'));
1116  } else {
1117  LOG("gevgen_fnal", pFATAL) << "No flux info was specified - Exiting";
1118  PrintSyntax();
1119  exit(1);
1120  }
1121 
1122  // number of events to generate
1123  if( parser.OptionExists('n') ) {
1124  LOG("gevgen_fnal", pDEBUG)
1125  << "Reading limit on number of events to generate";
1126  gOptNev = parser.ArgAsInt('n');
1127  } else {
1128  LOG("gevgen_fnal", pDEBUG)
1129  << "Will keep on generating events till the flux driver stops";
1130  gOptNev = -1;
1131  } //-n
1132 
1133  // statistics to generate in terms of POT
1134  if( parser.OptionExists('e') ) {
1135  LOG("gevgen_fnal", pDEBUG) << "Reading requested exposure in POT";
1136  gOptPOT = parser.ArgAsDouble('e');
1137  } else {
1138  LOG("gevgen_fnal", pDEBUG) << "No POT exposure was requested";
1139  gOptPOT = -1;
1140  } //-e
1141 
1142  // event file prefix
1143  if( parser.OptionExists('o') ) {
1144  LOG("gevgen_fnal", pDEBUG) << "Reading the event filename prefix";
1145  gOptEvFilePrefix = parser.ArgAsString('o');
1146  } else {
1147  LOG("gevgen_fnal", pDEBUG)
1148  << "Will set the default event filename prefix";
1150  } //-o
1151 
1152 
1153  // random number seed
1154  if( parser.OptionExists("seed") ) {
1155  LOG("gevgen_fnal", pINFO) << "Reading random number seed";
1156  gOptRanSeed = parser.ArgAsLong("seed");
1157  } else {
1158  LOG("gevgen_fnal", pINFO) << "Unspecified random number seed - Using default";
1159  gOptRanSeed = -1;
1160  }
1161 
1162  // input cross-section file
1163  if( parser.OptionExists("cross-sections") ) {
1164  LOG("gevgen_fnal", pINFO) << "Reading cross-section file";
1165  gOptInpXSecFile = parser.ArgAsString("cross-sections");
1166  } else {
1167  LOG("gevgen_fnal", pINFO) << "Unspecified cross-section file";
1168  gOptInpXSecFile = "";
1169  }
1170 
1171 
1172  //
1173  // >>> perform 'sanity' checks on command line arguments
1174  //
1175 
1176  // Tthe 'exposure' may be set either as:
1177  // - a number of POTs
1178  // - a number of generated events
1179  // Only one of those options can be set.
1180  if(!gOptUsingHistFlux) {
1181  int nset=0;
1182  if(gOptPOT > 0) nset++;
1183  if(gOptNev > 0) nset++;
1184  if(nset==0) {
1185  LOG("gevgen_fnal", pFATAL)
1186  << "** To use a gNuMI flux ntuple you need to specify an exposure, "
1187  << "either via the -e or -n options";
1188  PrintSyntax();
1189  exit(1);
1190  }
1191  if(nset>1) {
1192  LOG("gevgen_fnal", pFATAL)
1193  << "You can not specify more than one of the -e or -n options";
1194  PrintSyntax();
1195  exit(1);
1196  }
1197  }
1198  // If we use a flux histograms (not flux ntuples) then -currently- the
1199  // only way to control exposure is via a number of events
1200  if(gOptUsingHistFlux) {
1201  if(gOptNev < 0) {
1202  LOG("gevgen_fnal", pFATAL)
1203  << "If you're using flux from histograms you need to specify the -n option";
1204  PrintSyntax();
1205  exit(1);
1206  }
1207  }
1208  // If we don't use a detailed ROOT detector geometry (but just a target mix)
1209  // then don't accept POT as a way to control job statistics (not enough info
1210  // is passed in the target mix to compute POT & the calculation can be easily
1211  // done offline)
1212  if(!gOptUsingRootGeom) {
1213  if(gOptPOT > 0) {
1214  LOG("gevgen_fnal", pFATAL)
1215  << "You may not use the -e option without a detector geometry description";
1216  exit(1);
1217  }
1218  }
1219 
1220  //
1221  // >>> print the command line options
1222  //
1223 
1224  PDGLibrary * pdglib = PDGLibrary::Instance();
1225 
1226  ostringstream gminfo;
1227  if (gOptUsingRootGeom) {
1228  gminfo << "Using ROOT geometry - file = " << gOptRootGeom
1229  << ", top volume = "
1230  << ((gOptRootGeomTopVol.size()==0) ? "<master volume>" : gOptRootGeomTopVol)
1231  << ", max{PL} file = "
1232  << ((gOptExtMaxPlXml.size()==0) ? "<none>" : gOptExtMaxPlXml)
1233  << ", length units = " << lunits
1234  << ", density units = " << dunits;
1235  } else {
1236  gminfo << "Using target mix: ";
1237  map<int,double>::const_iterator iter;
1238  for(iter = gOptTgtMix.begin(); iter != gOptTgtMix.end(); ++iter) {
1239  int pdg_code = iter->first;
1240  double wgt = iter->second;
1241  TParticlePDG * p = pdglib->Find(pdg_code);
1242  if(p) {
1243  string name = p->GetName();
1244  gminfo << "(" << name << ") -> " << 100*wgt << "% / ";
1245  }//p?
1246  }
1247  }
1248 
1249  ostringstream fluxinfo;
1250  if (gOptUsingHistFlux) {
1251  fluxinfo << "Using histograms: ";
1252  map<int,TH1D*>::const_iterator iter;
1253  for(iter = gOptFluxHst.begin(); iter != gOptFluxHst.end(); ++iter) {
1254  int pdg_code = iter->first;
1255  TH1D * spectrum = iter->second;
1256  TParticlePDG * p = pdglib->Find(pdg_code);
1257  if(p) {
1258  string name = p->GetName();
1259  fluxinfo << "(" << name << ") -> " << spectrum->GetName() << " / ";
1260  }//p?
1261  }
1262  } else {
1263  fluxinfo << "Using " << gOptFluxDriver << " flux driver- "
1264  << "file = " << gOptFluxFile
1265  << ", location = " << gOptDetectorLocation;
1266  }
1267 
1268  ostringstream exposure;
1269  if(gOptPOT > 0)
1270  exposure << "Number of POTs = " << gOptPOT;
1271  if(gOptNev > 0)
1272  exposure << "Number of events = " << gOptNev;
1273 
1274 
1275  LOG("gevgen_fnal", pNOTICE)
1276  << "\n\n"
1277  << utils::print::PrintFramedMesg("FNAL expt. event generation job configuration");
1278 
1279  LOG("gevgen_fnal", pNOTICE)
1280  << "\n - Run number: " << gOptRunNu
1281  << "\n - Random number seed: " << gOptRanSeed
1282  << "\n - Using cross-section file: " << gOptInpXSecFile
1283  << "\n - Flux @ " << fluxinfo.str()
1284  << "\n - Geometry @ " << gminfo.str()
1285  << "\n - Exposure @ " << exposure.str();
1286 
1287  LOG("gevgen_fnal", pNOTICE) << *RunOpt::Instance();
1288 }
1289 //____________________________________________________________________________
1290 void PrintSyntax(void)
1291 {
1292  LOG("gevgen_fnal", pFATAL)
1293  << "\n **Syntax**"
1294  << "\n gevgen_fnal [-h] [-r run#]"
1295  << "\n -f flux -g geometry"
1296  << "\n [-t top_volume_name_at_geom] [-m max_path_lengths_xml_file]"
1297  << "\n [-L length_units_at_geom] [-D density_units_at_geom]"
1298  << "\n [-n n_of_events] [-e exposure_in_POTs]"
1299  << "\n [-o output_event_file_prefix]"
1300  << "\n [-F fid_cut_string] [-S nrays_scan]"
1301  << "\n [-z zmin_start]"
1302  << "\n [--seed random_number_seed]"
1303  << "\n --cross-sections xml_file"
1305  << "\n"
1306  << " Please also read the detailed documentation at "
1307  << "$GENIE/src/Apps/gFNALExptEvGen.cxx"
1308  << "\n";
1309 }
1310 //____________________________________________________________________________
1311 void CreateFidSelection (string fidcut, GeomAnalyzerI* geom_driver)
1312 {
1313  ///
1314  /// User defined fiducial volume cut
1315  /// [0][M]<SHAPE>:val1,val2,...
1316  /// "0" means reverse the cut (i.e. exclude the volume)
1317  /// "M" means the coordinates are given in the ROOT geometry
1318  /// "master" system and need to be transformed to "top vol" system
1319  /// <SHAPE> can be any of "zcyl" "box" "zpoly" "sphere"
1320  /// [each takes different # of args]
1321  /// This must be followed by a ":" and a list of values separated by punctuation
1322  /// (allowed separators: commas , parentheses () braces {} or brackets [] )
1323  /// Value mapping:
1324  /// zcyl:x0,y0,radius,zmin,zmax - cylinder along z at (x0,y0) capped at z's
1325  /// xcyl:y0,z0,radius,xmin,ymax - cylinder along x
1326  /// ycyl:x0,z0,radius,ymin,ymax - cylinder along y
1327  /// gcyl:{x0,y0,z0}{dx,dy,dz},radius,{plane1}{plane2} -- generic cylinder w/ arbitrary orientation and caps
1328  /// {planeX} = 4 values to define plane and orientation
1329  /// box:xmin,ymin,zmin,xmax,ymax,zmax - box w/ upper & lower extremes
1330  /// zpoly:nfaces,x0,y0,r_in,phi,zmin,zmax - nfaces sided polygon in x-y plane
1331  // sphere:x0,y0,z0,radius - sphere of fixed radius at (x0,y0,z0)
1332  /// Examples:
1333  /// 1) 0mbox:0,0,0.25,1,1,8.75
1334  /// exclude (i.e. reverse) a box in master coordinates w/ corners (0,0,0.25) (1,1,8.75)
1335  /// 2) mzpoly:6,(2,-1),1.75,0,{0.25,8.75}
1336  /// six sided polygon in x-y plane, centered at x,y=(2,-1) w/ inscribed radius 1.75
1337  /// no rotation (so first face is in y-z plane +r from center, i.e. hex sits on point)
1338  /// limited to the z range of {0.25,8.75} in the master ROOT geom coordinates
1339  /// 3) zcyl:(3,4),5.5,-2,10
1340  /// a cylinder oriented parallel to the z axis in the "top vol" coordinates
1341  /// at x,y=(3,4) with radius 5.5 and z range of {-2,10}
1342  ///
1343  geometry::ROOTGeomAnalyzer * rgeom =
1344  dynamic_cast<geometry::ROOTGeomAnalyzer *>(geom_driver);
1345  if ( ! rgeom ) {
1346  LOG("gevgen_fnal", pWARN)
1347  << "Can not create GeomVolSelectorFiduction,"
1348  << " geometry driver is not ROOTGeomAnalyzer";
1349  return;
1350  }
1351 
1352  LOG("gevgen_fnal", pNOTICE) << "-F " << fidcut;
1353 
1356 
1357  fidsel->SetRemoveEntries(true); // drop segments that won't be considered
1358 
1359  // convert string to lowercase
1360  std::transform(fidcut.begin(),fidcut.end(),fidcut.begin(),::tolower);
1361 
1362  vector<string> strtok = genie::utils::str::Split(fidcut,":");
1363  if ( strtok.size() != 2 ) {
1364  LOG("gevgen_fnal", pWARN)
1365  << "Can not create GeomVolSelectorFiduction,"
1366  << " no \":\" separating type from values. nsplit=" << strtok.size();
1367  for ( unsigned int i=0; i < strtok.size(); ++i )
1368  LOG("gevgen_fnal",pNOTICE)
1369  << "strtok[" << i << "] = \"" << strtok[i] << "\"";
1370  return;
1371  }
1372 
1373  // parse out optional "x" and "m"
1374  string stype = strtok[0];
1375  bool reverse = ( stype.find("0") != string::npos );
1376  bool master = ( stype.find("m") != string::npos ); // action after values are set
1377 
1378  // parse out values
1379  vector<double> vals;
1380  vector<string> valstrs = genie::utils::str::Split(strtok[1]," ,;(){}[]");
1381  vector<string>::const_iterator iter = valstrs.begin();
1382  for ( ; iter != valstrs.end(); ++iter ) {
1383  const string& valstr1 = *iter;
1384  if ( valstr1 != "" ) vals.push_back(atof(valstr1.c_str()));
1385  }
1386  size_t nvals = vals.size();
1387 
1388  std::cout << "ivals = [";
1389  for (unsigned int i=0; i < nvals; ++i) {
1390  if (i>0) cout << ",";
1391  std::cout << vals[i];
1392  }
1393  std::cout << "]" << std::endl;
1394 
1395  // std::vector elements are required to be adjacent so we can treat address as ptr
1396 
1397  if ( stype.find("zcyl") != string::npos ) {
1398  // cylinder along z direction at (x0,y0) radius zmin zmax
1399  if ( nvals < 5 )
1400  LOG("gevgen_fnal", pFATAL) << "MakeZCylinder needs 5 values, not " << nvals
1401  << " fidcut=\"" << fidcut << "\"";
1402  fidsel->MakeZCylinder(vals[0],vals[1],vals[2],vals[3],vals[4]);
1403 
1404  } else if ( stype.find("xcyl") != string::npos ) {
1405  // cylinder along x direction at (y0,z0) radius xmin xmax
1406  if ( nvals < 5 )
1407  LOG("gevgen_fnal", pFATAL) << "MakeXCylinder needs 5 values, not " << nvals
1408  << " fidcut=\"" << fidcut << "\"";
1409  fidsel->MakeXCylinder(vals[0],vals[1],vals[2],vals[3],vals[4]);
1410 
1411  } else if ( stype.find("ycyl") != string::npos ) {
1412  // cylinder along y direction at (x0,z0) radius ymin ymax
1413  if ( nvals < 5 )
1414  LOG("gevgen_fnal", pFATAL) << "MakeYCylinder needs 5 values, not " << nvals
1415  << " fidcut=\"" << fidcut << "\"";
1416  fidsel->MakeYCylinder(vals[0],vals[1],vals[2],vals[3],vals[4]);
1417 
1418  } else if ( stype.find("gcyl") != string::npos ) {
1419  // cylinder along arbitrary direction at (x0,y0,z0) radius {plane1} {plane2}
1420  if ( nvals < 14 )
1421  LOG("gevgen_fnal", pFATAL) << "MakeYCylinder needs 14 values, not " << nvals
1422  << " fidcut=\"" << fidcut << "\"";
1423  Double_t base[3] = { vals[0], vals[1], vals[2] };
1424  Double_t axis[3] = { vals[3], vals[4], vals[5] };
1425  Double_t radius = vals[6];
1426  Double_t cap1[4] = { vals[ 7], vals[ 8], vals[ 9], vals[10] };
1427  Double_t cap2[4] = { vals[11], vals[12], vals[13], vals[14] };
1428 
1429  fidsel->MakeCylinder(base,axis,radius,cap1,cap2);
1430 
1431  } else if ( stype.find("box") != string::npos ) {
1432  // box (xmin,ymin,zmin) (xmax,ymax,zmax)
1433  if ( nvals < 6 )
1434  LOG("gevgen_fnal", pFATAL) << "MakeBox needs 6 values, not " << nvals
1435  << " fidcut=\"" << fidcut << "\"";
1436  double xyzmin[3] = { vals[0], vals[1], vals[2] };
1437  double xyzmax[3] = { vals[3], vals[4], vals[5] };
1438  fidsel->MakeBox(xyzmin,xyzmax);
1439 
1440  } else if ( stype.find("zpoly") != string::npos ) {
1441  // polygon along z direction nfaces at (x0,y0) radius phi zmin zmax
1442  if ( nvals < 7 )
1443  LOG("gevgen_fnal", pFATAL) << "MakeZPolygon needs 7 values, not " << nvals
1444  << " fidcut=\"" << fidcut << "\"";
1445  int nfaces = (int)vals[0];
1446  if ( nfaces < 3 )
1447  LOG("gevgen_fnal", pFATAL) << "MakeZPolygon needs nfaces>=3, not " << nfaces
1448  << " fidcut=\"" << fidcut << "\"";
1449  fidsel->MakeZPolygon(nfaces,vals[1],vals[2],vals[3],vals[4],vals[5],vals[6]);
1450 
1451  } else if ( stype.find("sphere") != string::npos ) {
1452  // sphere at (x0,y0,z0) radius
1453  if ( nvals < 4 )
1454  LOG("gevgen_fnal", pFATAL) << "MakeZSphere needs 4 values, not " << nvals
1455  << " fidcut=\"" << fidcut << "\"";
1456  fidsel->MakeSphere(vals[0],vals[1],vals[2],vals[3]);
1457 
1458  } else {
1459  LOG("gevgen_fnal", pFATAL)
1460  << "Can not create GeomVolSelectorFiduction for shape \"" << stype << "\"";
1461  }
1462 
1463  if ( master ) {
1464  fidsel->ConvertShapeMaster2Top(rgeom);
1465  LOG("gevgen_fnal", pNOTICE) << "Convert fiducial volume from master to topvol coords";
1466  }
1467  if ( reverse ) {
1468  fidsel->SetReverseFiducial(true);
1469  LOG("gevgen_fnal", pNOTICE) << "Reverse sense of fiducial volume cut";
1470  }
1471  rgeom->AdoptGeomVolSelector(fidsel);
1472 
1473 }
1474 //____________________________________________________________________________
1475 void CreateRockBoxSelection (string fidcut, GeomAnalyzerI* geom_driver)
1476 {
1477 
1478  if( fidcut.find_first_not_of(" \t\n") != 0) // trim any leading whitespace
1479  fidcut.erase( 0, fidcut.find_first_not_of(" \t\n") );
1480 
1481  // convert string to lowercase
1482  std::transform(fidcut.begin(),fidcut.end(),fidcut.begin(),::tolower);
1483 
1485  dynamic_cast<genie::geometry::ROOTGeomAnalyzer *>(geom_driver);
1486  if ( ! rgeom ) {
1487  LOG("gevgen_fnal", pWARN)
1488  << "Can not create GeomVolSelectorRockBox,"
1489  << " geometry driver is not ROOTGeomAnalyzer";
1490  return;
1491  }
1492 
1493  LOG("gevgen_fnal", pWARN) << "fiducial (rock) cut: " << fidcut;
1494 
1495  // for now, only fiducial no "rock box"
1498 
1499  vector<string> strtok = genie::utils::str::Split(fidcut,":");
1500  if ( strtok.size() != 2 ) {
1501  LOG("gevgen_fnal", pWARN)
1502  << "Can not create GeomVolSelectorRockBox,"
1503  << " no \":\" separating type from values. nsplit=" << strtok.size();
1504  for ( unsigned int i=0; i < strtok.size(); ++i )
1505  LOG("gevgen_fnal", pWARN)
1506  << "strtok[" << i << "] = \"" << strtok[i] << "\"";
1507  return;
1508  }
1509 
1510  string stype = strtok[0];
1511 
1512  // parse out values
1513  vector<double> vals;
1514  vector<string> valstrs = genie::utils::str::Split(strtok[1]," ,;(){}[]\t\n\r");
1515  vector<string>::const_iterator iter = valstrs.begin();
1516  for ( ; iter != valstrs.end(); ++iter ) {
1517  const string& valstr1 = *iter;
1518  if ( valstr1 != "" ) {
1519  double aval = atof(valstr1.c_str());
1520  LOG("gevgen_fnal", pWARN) << "rock value [" << vals.size() << "] "
1521  << aval;
1522  vals.push_back(aval);
1523  }
1524  }
1525  size_t nvals = vals.size();
1526 
1527  rocksel->SetRemoveEntries(true); // drop segments that won't be considered
1528 
1529  // assume coordinates are in the *master* (not "top volume") system
1530  // need to set fTopVolume to fWorldVolume
1531  //fTopVolume = fWorldVolume;
1532  //rgeom->SetTopVolName(fTopVolume.c_str());
1534  rgeom->SetTopVolName(gOptRootGeomMasterVol);
1535 
1536  if ( nvals < 6 ) {
1537  LOG("gevgen_fnal", pFATAL) << "rockbox needs at "
1538  << "least 6 values, found "
1539  << nvals << "in \""
1540  << strtok[1] << "\"";
1541  exit(1);
1542 
1543  }
1544  double xyzmin[3] = { vals[0], vals[1], vals[2] };
1545  double xyzmax[3] = { vals[3], vals[4], vals[5] };
1546 
1547  bool rockonly = true;
1548  double wallmin = 800.; // geometry in cm, ( 8 meter buffer)
1549  double dedx = 2.5 * 1.7e-3; // GeV/cm, rho=2.5, 1.7e-3 ~ rock like loss
1550  double fudge = 1.05;
1551 
1552  if ( nvals >= 7 ) rockonly = vals[6];
1553  if ( nvals >= 8 ) wallmin = vals[7];
1554  if ( nvals >= 9 ) dedx = vals[8];
1555  if ( nvals >= 10 ) fudge = vals[9];
1556 
1557  rocksel->SetRockBoxMinimal(xyzmin,xyzmax);
1558  rocksel->SetMinimumWall(wallmin);
1559  rocksel->SetDeDx(dedx/fudge);
1560 
1561  if ( nvals >= 11 ) rocksel->SetExpandFromInclusion((int)vals[10]);
1562 
1563  // if not rock-only then make a tiny exclusion bubble
1564  // call to MakeBox shouldn't be necessary
1565  // should be done by SetRockBoxMinimal but for some GENIE versions isn't
1566  if ( ! rockonly ) rocksel->MakeSphere(0,0,0,1.0e-10);
1567  else rocksel->MakeBox(xyzmin,xyzmax);
1568 
1569  rgeom->AdoptGeomVolSelector(rocksel);
1570 
1571 }
1572 
1573 //____________________________________________________________________________
1574 void DetermineFluxDriver(string fopt)
1575 {
1576  // based on the -f option string determine which flux driver to use
1577  // this may take some guessing
1578 
1579  // first look for strings that look like "<proto>:..."
1580  // or ....<proto>.root,....
1581  // where "<proto>" is a key the gOptFluxShortNames map
1582 
1583  map<string,string>::const_iterator mitr = gOptFluxShortNames.begin();
1584  map<string,string>::const_iterator mitr_end = gOptFluxShortNames.end();
1585  for ( ; mitr != mitr_end; ++mitr ) {
1586  string proto = mitr->first + string(":");
1587  string gproto = string("g") + proto;
1588  string protor = proto + ".root,";
1589  string full = mitr->second;
1590  if ( fopt.find(proto) == 0 ) {
1591  fopt.erase(0,proto.size());
1592  gOptFluxDriver = full;
1593  break;
1594  } else if ( fopt.find(gproto) == 0 ) {
1595  fopt.erase(0,gproto.size());
1596  gOptFluxDriver = full;
1597  break;
1598  } else if ( fopt.find(protor) != std::string::npos ) {
1599  gOptFluxDriver = full;
1600  break;
1601  }
1602  }
1603  // tested all cases where user might have specified explicitly
1604  // or been part of an extended file extension
1605  // this is where it gets messy
1606  if ( gOptFluxDriver == "" ) {
1607 
1608  // not specified ? guess from file name itself
1609  if ( fopt.find("gsimple") != std::string::npos ) {
1610  // put dk2nu after gsimple in case simple files are derived from dk2nu
1611  // then both are in name we should choose gsimple
1612  gOptFluxDriver = "genie::flux::GSimpleNtpFlux";
1613  } else if ( fopt.find("dk2nu") != std::string::npos ) {
1614  gOptFluxDriver = "genie::flux::GDk2NuFlux";
1615  } else {
1616  // does it look like the histogram format
1617  const char* hstrings[] = { ",12[", ",+12[", ",-12[",
1618  ",14[", ",+14[", ",-14[",
1619  ",16[", ",+16[", ",-16[" };
1620  size_t nh = sizeof(hstrings)/sizeof(const char*);
1621  for (size_t ih=0; ih<nh; ++ih) {
1622  if ( fopt.find(hstrings[ih]) != std::string::npos ) {
1623  // hey!
1624  gOptFluxDriver = "genie::flux::GCylindTH1Flux";
1625  break;
1626  }
1627  } // loop over possible histogram specifiers
1628  }
1629 
1630  // fall through default ... hope it works
1631  if ( gOptFluxDriver == "" ) {
1632  gOptFluxDriver = "genie::flux::GNuMIFlux";
1633  }
1634  }
1635 
1636  gOptUsingHistFlux = ( gOptFluxDriver == "genie::flux::GCylindTH1Flux" );
1637  if ( gOptUsingHistFlux ) ParseFluxHst(fopt);
1638  else ParseFluxFileConfig(fopt);
1639 }
1640 //____________________________________________________________________________
1641 void ParseFluxHst(string flux)
1642 {
1643  // Using flux from histograms
1644  // Extract the root file name & the list of histogram names & neutrino
1645  // species (specified as 'filename,histo1[species1],histo2[species2],...')
1646  // for variable width histograms, default to multiply in the width
1647  // histo1[species1]WIDTH = multiply in the width
1648  // histo1[species1]NOWIDTH = don't multiply in the width
1649  // possibly with configuration ";r=1.1,dir=(0.1,0.2,0.3),spot=(-1,-2,-3)"
1650  // appended
1651  // See documentation on top section of this file.
1652 
1653  vector<string> fluxtop = utils::str::Split(flux,";");
1654  string beamsettings;
1655  string histosettings;
1656  if ( fluxtop.size() == 1 ) {
1657  histosettings = fluxtop[0];
1658  LOG("gevgen_fnal", pNOTICE)
1659  << "ParseFluxHst: no settings for radius, direction, spot position"
1660  << " using defaults, r=" << gOptBeamRadius;
1661  } else if ( fluxtop.size() > 2 ) {
1662  LOG("gevgen_fnal", pFATAL)
1663  << "ParseFluxHst: too many ';' separated fields";
1664  PrintSyntax();
1665  exit(1);
1666  } else {
1667  // decide which is which depending on which has "=" in it
1668  string::size_type eqpos0 = fluxtop[0].find("=");
1669  string::size_type eqpos1 = fluxtop[1].find("=");
1670  if (eqpos0 == string::npos && eqpos1 != string::npos ) {
1671  histosettings = fluxtop[0];
1672  beamsettings = fluxtop[1];
1673  } else if (eqpos0 != string::npos && eqpos1 == string::npos ) {
1674  beamsettings = fluxtop[0];
1675  histosettings = fluxtop[1];
1676  } else {
1677  LOG("gevgen_fnal", pFATAL)
1678  << "ParseFluxHst: too many / not enough fields with '=' "
1679  << "\n" << fluxtop[0] << "\n" << fluxtop[1];
1680  PrintSyntax();
1681  exit(1);
1682  }
1683  // now parse the string we have
1684  double r=-1, dx=0, dy=0, dz=1, x=0, y=0, z=0;
1685  int nscan = sscanf(beamsettings.c_str(),
1686  "r=%lf,dir=(%lf,%lf,%lf),spot=(%lf,%lf,%lf)",
1687  &r,&dx,&dy,&dz,&x,&y,&z);
1688  cout << "nscan = " << nscan << endl;
1689  gOptBeamRadius = r;
1690  gOptBeamDir = TVector3(dx,dy,dz);
1691  gOptBeamSpot = TVector3(x,y,z);
1692 
1693  }
1694  LOG("gevgen_fnal", pNOTICE)
1695  << "ParseFluxHst: "
1696  << " using r=" << gOptBeamRadius
1697  << ",dir=("
1698  << gOptBeamDir.Px() << ","
1699  << gOptBeamDir.Py() << ","
1700  << gOptBeamDir.Pz() << ")"
1701  << ",spot=("
1702  << gOptBeamSpot.X() << ","
1703  << gOptBeamSpot.Y() << ","
1704  << gOptBeamSpot.Z() << ")";
1705 
1706 
1707  vector<string> fluxv = utils::str::Split(histosettings,",");
1708  if(fluxv.size()<2) {
1709  LOG("gevgen_fnal", pFATAL)
1710  << "You need to specify both a flux histogram ROOT file "
1711  << " _AND_ at least one histogram[pdg] mapping";
1712  PrintSyntax();
1713  exit(1);
1714  }
1715  gOptFluxFile = fluxv[0];
1716  bool accessible_flux_file = !(gSystem->AccessPathName(gOptFluxFile.c_str()));
1717  if (!accessible_flux_file) {
1718  LOG("gevgen_fnal", pFATAL)
1719  << "Can not access flux file: " << gOptFluxFile;
1720  PrintSyntax();
1721  exit(1);
1722  }
1723  // Extract energy spectra for all specified neutrino species
1724  TFile flux_file(gOptFluxFile.c_str(), "read");
1725  for(unsigned int inu=1; inu<fluxv.size(); inu++) {
1726  string nutype_and_histo = fluxv[inu];
1727  string::size_type open_bracket = nutype_and_histo.find("[");
1728  string::size_type close_bracket = nutype_and_histo.find("]");
1729  if (open_bracket ==string::npos ||
1730  close_bracket==string::npos)
1731  {
1732  LOG("gevgen_fnal", pFATAL)
1733  << "You made an error in specifying the flux histograms";
1734  PrintSyntax();
1735  exit(1);
1736  }
1737  string::size_type ibeg = 0;
1738  string::size_type iend = open_bracket;
1739  string::size_type jbeg = open_bracket+1;
1740  string::size_type jend = close_bracket;
1741  string nutype = nutype_and_histo.substr(ibeg,iend-ibeg);
1742  string histo = nutype_and_histo.substr(jbeg,jend-jbeg);
1743  string extra = nutype_and_histo.substr(jend+1,string::npos);
1744  std::transform(extra.begin(),extra.end(),extra.begin(),::toupper);
1745  LOG("gevgen_fnal", pNOTICE) // pDEBUG
1746  << " =======> nutype " << nutype << " histo " << histo << " extra " << extra;
1747  // access specified histogram from the input root file
1748  TH1D * ihst = (TH1D*) flux_file.Get(histo.c_str());
1749  if(!ihst) {
1750  LOG("gevgen_fnal", pFATAL)
1751  << "Can not find histogram: " << histo
1752  << " in flux file: " << gOptFluxFile;
1753  PrintSyntax();
1754  exit(1);
1755  }
1756 
1757  // Copy in the flux histogram from the root file
1758  // use Clone rather than assuming fix bin widths and rebooking
1759  TH1D* spectrum = (TH1D*)ihst->Clone();
1760  spectrum->SetNameTitle("spectrum","neutrino_flux");
1761  spectrum->SetDirectory(0);
1762  // get rid of original
1763  delete ihst;
1764 
1765  bool force_binwidth = false;
1766 #if ROOT_VERSION_CODE <= ROOT_VERSION(9,99,99)
1767  // GetRandom() sampling on variable bin width histograms does not
1768  // correctly account for bin widths for all versions of ROOT prior
1769  // to (currently forever). At some point this might change and
1770  // the necessity of this code snippet will go away
1771  TAxis* xaxis = spectrum->GetXaxis();
1772  if (xaxis->IsVariableBinSize()) force_binwidth = true;
1773 #endif
1774  if ( extra == "WIDTH" ) force_binwidth = true;
1775  if ( extra == "NOWIDTH" ) force_binwidth = false;
1776  if ( force_binwidth ) {
1777  LOG("gevgen_fnal", pNOTICE)
1778  << "multiplying by bin width for histogram " << histo
1779  << " as " << spectrum->GetName() << " for nutype " << nutype
1780  << " from " << gOptFluxFile;
1781  for(int ibin = 1; ibin <= spectrum->GetNbinsX(); ++ibin) {
1782  double data = spectrum->GetBinContent(ibin);
1783  double width = spectrum->GetBinWidth(ibin);
1784  spectrum->SetBinContent(ibin,data*width);
1785  }
1786  }
1787 
1788  // convert neutrino name -> pdg code
1789  int pdg = atoi(nutype.c_str());
1790  if(!pdg::IsNeutrino(pdg) && !pdg::IsAntiNeutrino(pdg)) {
1791  LOG("gevgen_fnal", pFATAL)
1792  << "Unknown neutrino type: " << nutype;
1793  PrintSyntax();
1794  exit(1);
1795  }
1796  // store flux neutrino code / energy spectrum
1797  LOG("gevgen_fnal", pDEBUG)
1798  << "Adding energy spectrum for flux neutrino: pdg = " << pdg;
1799  gOptFluxHst.insert(map<int, TH1D*>::value_type(pdg, spectrum));
1800  }//inu
1801 
1802  if(gOptFluxHst.size()<1) {
1803  LOG("gevgen_fnal", pFATAL)
1804  << "You have not specified any flux histogram!";
1805  PrintSyntax();
1806  exit(1);
1807  }
1808 
1809  flux_file.Close();
1810 }
1811 
1812 //____________________________________________________________________________
1813 void ParseFluxFileConfig(string flux)
1814 {
1815  // Using gnumi/gsimple/dk2nu beam flux ntuples
1816  // Extract beam flux (root) file name & detector location
1817  //
1818  vector<string> fluxv = utils::str::Split(flux,",");
1819  if(fluxv.size()<2) {
1820  LOG("gevgen_fnal", pFATAL)
1821  << "You need to specify both a flux ntuple ROOT file "
1822  << " _AND_ a detector location";
1823  PrintSyntax();
1824  exit(1);
1825  }
1826  gOptFluxFile = fluxv[0];
1827  gOptDetectorLocation = fluxv[1];
1828 
1829  for ( size_t j = 2; j < fluxv.size(); ++j ) {
1830  int ipdg = atoi(fluxv[j].c_str());
1831  gOptFluxPdg.push_back(ipdg);
1832  }
1833 
1834 }
1835 
1836 //____________________________________________________________________________
void MakeSphere(Double_t x0, Double_t y0, Double_t z0, Double_t radius)
void RandGen(long int seed)
Definition: AppInit.cxx:30
static void SetPrintLevel(int print_level)
Definition: GHepRecord.cxx:956
const char * GetExposureUnits() const
what units are returned by GetTotalExposure?
GENIE Interface for user-defined volume selector functors Trim path segments based on the intersectio...
void XSecTable(string inpfile, bool require_table)
Definition: AppInit.cxx:38
virtual GeomVolSelectorI * AdoptGeomVolSelector(GeomVolSelectorI *selector)
configure processing to perform path segment trimming
void CreateRockBoxSelection(string fidcut, GeomAnalyzerI *geom_driver)
PDGCodeList gOptFluxPdg
virtual void LoadBeamSimData(const std::vector< std::string > &filenames, const std::string &det_loc)=0
bool IsNeutrino(int pdgc)
Definition: PDGUtils.cxx:110
string gOptEvFilePrefix
Definition: gAtmoEvGen.cxx:310
#define pERROR
Definition: Messenger.h:59
int gOptDebug
void MakeZCylinder(Double_t x0, Double_t y0, Double_t radius, Double_t zmin, Double_t zmax)
string gOptRootGeomMasterVol
void SetEventGeneratorList(string listname)
Definition: GMCJDriver.cxx:66
virtual void PrintConfig()=0
print the current configuration
void DetermineFluxDriver(string fopt)
void ReadFromCommandLine(int argc, char **argv)
Definition: RunOpt.cxx:99
#define pFATAL
Definition: Messenger.h:56
double gOptGeomDUnits
Definition: gAtmoEvGen.cxx:303
string gOptDetectorLocation
void MakeCylinder(Double_t *base, Double_t *axis, Double_t radius, Double_t *cap1, Double_t *cap2)
virtual const PathLengthList & GetMaxPathLengths(void) const
static constexpr double s
Definition: Units.h:95
void MakeXCylinder(Double_t y0, Double_t z0, Double_t radius, Double_t xmin, Double_t xmax)
virtual void ConvertShapeMaster2Top(const ROOTGeomAnalyzer *rgeom)
void Update(int iev, const EventRecord *event)
Definition: GMCJMonitor.cxx:48
double gOptZmin
const std::vector< std::string > & AvailableFluxDrivers() const
int gOptNScan
static void gsSIGTERMhandler(int)
map< int, double > gOptTgtMix
Definition: gAtmoEvGen.cxx:299
string gOptFidCut
virtual void SetNumOfCycles(long int ncycle)
limit cycling through input files
virtual void SetScannerNParticles(int np)
virtual long int NFluxNeutrinos() const =0
of rays generated
int main(int argc, char **argv)
Definition: gAtmoEvGen.cxx:327
A list of PDG codes.
Definition: PDGCodeList.h:32
bool gOptWriteMaxPlXml
TVector3 gOptBeamDir
void CreateFidSelection(string fidcut, GeomAnalyzerI *geom_driver)
virtual void SetUpstreamZ(double z0)
void SetRefreshRate(int rate)
Definition: GMCJMonitor.cxx:43
string kDefOptEvFilePrefix
Definition: gAtmoEvGen.cxx:320
void UseFluxDriver(GFluxI *flux)
Definition: GMCJDriver.cxx:83
Object to be filled with the neutrino path-length, for all detector geometry materials, when starting from a position x and travelling along the direction of the neutrino 4-momentum.
string gOptBeamRtDependence
string lunits
Simple class to create &amp; update MC job status files and env. vars. This is used to be able to keep tr...
Definition: GMCJMonitor.h:31
A GENIE `MC Job Driver&#39;. 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 SaveAsXml(string filename) const
void ParseFluxFileConfig(string fopt)
const double e
TTree * EventTree(void)
Definition: NtpWriter.h:55
string gOptInpXSecFile
Definition: gAtmoEvGen.cxx:313
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
void MakeYCylinder(Double_t x0, Double_t z0, Double_t radius, Double_t ymin, Double_t ymax)
string gOptRootGeomTopVol
Definition: gAtmoEvGen.cxx:301
void ForceSingleProbScale(void)
Definition: GMCJDriver.cxx:172
virtual void SetScannerNRays(int nr)
A generic GENIE flux driver. Generates a &#39;cylindrical&#39; neutrino beam along the input direction...
static std::string RunOptSyntaxString(bool include_generator_specific)
Definition: RunOpt.cxx:157
string geom
bool IsAntiNeutrino(int pdgc)
Definition: PDGUtils.cxx:118
bool UseMaxPathLengths(string xml_filename)
Definition: GMCJDriver.cxx:99
double UnitFromString(string u)
Definition: UnitUtils.cxx:18
bool gSigTERM
GENIE interface for uniform flux exposure iterface.
int gOptNev
Definition: gAtmoEvGen.cxx:305
GENIE Interface for limiting vertex selection in the rock to a volume that depends (in part) on the n...
TVector3 gOptBeamSpot
void Configure(bool calc_prob_scales=true)
Definition: GMCJDriver.cxx:399
Double_t gOptBeamRadius
void MakeZPolygon(Int_t n, Double_t x0, Double_t y0, Double_t inradius, Double_t phi0deg, Double_t zmin, Double_t zmax)
A ROOT/GEANT4 geometry driver.
map< string, string > gOptFluxShortNames
void Save(void)
get the even tree
Definition: NtpWriter.cxx:225
#define pINFO
Definition: Messenger.h:62
void AddEventRecord(int ievent, const EventRecord *ev_rec)
save the event tree
Definition: NtpWriter.cxx:57
void SetTransverseRadius(double Rt)
void SetNuDirection(const TVector3 &direction)
string gOptRootGeom
Definition: gAtmoEvGen.cxx:300
void BuildTune()
build tune and inform XSecSplineList
Definition: RunOpt.cxx:92
#define pWARN
Definition: Messenger.h:60
EventRecord * GenerateEvent(void)
Definition: GMCJDriver.cxx:815
virtual void GetBranchInfo(std::vector< std::string > &branchNames, std::vector< std::string > &branchClassNames, std::vector< void ** > &branchObjPointers)
Long_t gOptRunNu
Definition: gAtmoEvGen.cxx:295
Generated Event Record. It is a GHepRecord object that can accept / be visited by EventRecordVisitorI...
Definition: EventRecord.h:37
void CustomizeFilenamePrefix(string prefix)
Definition: NtpWriter.cxx:133
void LoadExtraOptions(void)
void Initialize(void)
add event
Definition: NtpWriter.cxx:83
string kDefOptGeomLUnits
Definition: gAtmoEvGen.cxx:321
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:68
static RunOpt * Instance(void)
Definition: RunOpt.cxx:54
vector< string > Split(string input, string delim)
Definition: StringUtils.cxx:36
A utility class to facilitate creating the GENIE MC Ntuple from the output GENIE GHEP event records...
Definition: NtpWriter.h:39
Singleton class to load &amp; serve a TDatabasePDG.
Definition: PDGLibrary.h:35
virtual void SetScannerNPoints(int np)
set geometry driver&#39;s configuration options
virtual void SetScannerFlux(GFluxI *f)
static GFluxDriverFactory & Instance()
string gOptFluxFile
virtual bool End(void)=0
true if no more flux nu&#39;s can be thrown (eg reaching end of beam sim ntuples)
virtual double GetTotalExposure() const =0
NtpMCFormat_t kDefOptNtpFormat
Definition: gAtmoEvGen.cxx:319
string PrintFramedMesg(string mesg, unsigned int nl=1, const char f='*')
Definition: PrintUtils.cxx:164
A vector of EventGeneratorI objects.
virtual void SetFluxParticles(const PDGCodeList &particles)
specify list of flux neutrino species
void SetBeamSpot(const TVector3 &spot)
void MesgThresholds(string inpfile)
Definition: AppInit.cxx:99
void AddEnergySpectrum(int nu_pdgc, TH1D *spectrum)
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:86
Command line argument parser.
enum genie::ENtpMCFormat NtpMCFormat_t
#define pNOTICE
Definition: Messenger.h:61
void SetRockBoxMinimal(Double_t *xyzmin, Double_t *xyzmax)
void MakeBox(Double_t *xyzmin, Double_t *xyzmax)
void GetCommandLineArgs(int argc, char **argv)
Definition: gAtmoEvGen.cxx:563
double gOptPOT
Defines the GENIE Geometry Analyzer Interface.
Definition: GeomAnalyzerI.h:29
int EventRecordPrintLevel(void) const
Definition: RunOpt.h:51
string kDefOptGeomDUnits
Definition: gAtmoEvGen.cxx:322
void RecursiveExhaust(TGeoVolume *topvol, string volnames, bool exhaust)
Definition: GeoUtils.cxx:16
double gOptGeomLUnits
Definition: gAtmoEvGen.cxx:302
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils...
void UseGeomAnalyzer(GeomAnalyzerI *geom)
Definition: GMCJDriver.cxx:88
void UseSplines(bool useLogE=true)
Definition: GMCJDriver.cxx:93
The PointGeomAnalyzer class is the simplest implementation of the GeomAnalyserI interface and defines...
void PrintSyntax(void)
void ParseFluxHst(string fopt)
void CacheFile(string inpfile)
Definition: AppInit.cxx:117
void EnableBareXSecPreCalc(bool flag)
Definition: RunOpt.h:62
string gOptExtMaxPlXml
Definition: gAtmoEvGen.cxx:304
void push_back(int pdg_code)
Definition: PDGCodeList.cxx:58
bool gOptUsingRootGeom
Definition: gAtmoEvGen.cxx:298
map< int, TH1D * > gOptFluxHst
virtual TGeoManager * GetGeometry(void) const
genie::GFluxI * GetFluxDriver(const std::string &)
#define pDEBUG
Definition: Messenger.h:63
GENIE Interface for user-defined flux classes.
Definition: GFluxI.h:29
string gOptFluxDriver
long int gOptRanSeed
Definition: gAtmoEvGen.cxx:312
bool gOptUsingHistFlux