GENIEGenerator
Loading...
Searching...
No Matches
gT2KEvGen.cxx
Go to the documentation of this file.
1//________________________________________________________________________________________
2/*!
3
4\program gevgen_t2k
5
6\brief A GENIE event generation driver 'customized' for T2K.
7
8 This driver can handle the JPARC neutrino flux files generated by JNUBEAM and the
9 realistic detector geometries / target mix of T2K detectors.
10
11 T2K users should note that the generic GENIE event generation driver (gevgen) may
12 still be a more appropriate tool to use for the simpler event generation cases
13 required for many 4-vector level / systematic studies.
14 Please see the GENIE documentation (http://www.genie-mc.org) and contact me
15 <c.andreopoulos \at cern.ch> if in doubt.
16
17 *** Synopsis :
18
19 gevgen_t2k [-h]
20 [-r run#]
21 -f flux
22 -g geometry
23 [-p pot_normalization_of_flux_file]
24 [-t top_volume_name_at_geom || -t +Vol1-Vol2...]
25 [-P pre_gen_prob_file_name]
26 [-S] [output_name]
27 [-m max_path_lengths_xml_file]
28 [-L length_units_at_geom]
29 [-D density_units_at_geom]
30 [-n n_of_events]
31 [-c flux_cycles]
32 [-e, -E exposure_in_POTs]
33 [-o output_event_file_prefix]
34 [-R]
35 [--seed random_number_seed]
36 --cross-sections xml_file
37 [--tune genie_tune]
38 [--message-thresholds xml_file]
39 [--unphysical-event-mask mask]
40 [--event-record-print-level level]
41 [--mc-job-status-refresh-rate rate]
42 [--cache-file root_file]
43
44 *** Options :
45
46 [] Denotes an optional argument
47
48 -h
49 Prints out the gevgen_t2k syntax and exits.
50 -r
51 Specifies the MC run number [default: 1000].
52 -g
53 Input 'geometry'.
54 This option can be used to specify any of:
55 1 > A ROOT file containing a ROOT/GEANT geometry description
56 [Note]
57 - This is the standard option for generating events in the
58 nd280, 2km and INGRID detectors.
59 [Examples]
60 - To use the master volume from the ROOT geometry stored
61 in the nd280-geom.root file, type:
62 '-g /some/path/nd280-geom.root'
63 2 > A mix of target materials, each with its corresponding weight,
64 typed as a comma-separated list of nuclear PDG codes (in the
65 std PDG2006 convention: 10LZZZAAAI) with the weight fractions
66 in brackets, eg code1[fraction1],code2[fraction2],...
67 If that option is used (no detailed input geometry description)
68 then the interaction vertices are distributed in the detector
69 by the detector MC.
70 [Note]
71 - This is the standard option for generating events in the
72 SuperK detector.
73 [Examples]
74 - To use a target mix of 95% O16 and 5% H type:
75 '-g 1000080160[0.95],1000010010[0.05]'
76 - To use a target which is 100% C12, type:
77 '-g 1000060120'
78 -t
79 Input 'top volume' for event generation -
80 can be used to force event generation in given sub-detector.
81 [default: the 'master volume' of the input geometry]
82 You can also use the -t option to switch generation on/off at
83 multiple volumes as, for example, in:
84 `-t +Vol1-Vol2+Vol3-Vol4',
85 `-t "+Vol1 -Vol2 +Vol3 -Vol4"',
86 `-t -Vol2-Vol4+Vol1+Vol3',
87 `-t "-Vol2 -Vol4 +Vol1 +Vol3"'m
88 where:
89 "+Vol1" and "+Vol3" tells GENIE to `switch on' Vol1 and Vol3, while
90 "-Vol2" and "-Vol4" tells GENIE to `switch off' Vol2 and Vol4.
91 If the very first character is a '+', GENIE will neglect all volumes
92 except the ones explicitly turned on. Vice versa, if the very first
93 character is a `-', GENIE will keep all volumes except the ones
94 explicitly turned off (feature contributed by J.Holeczek).
95 -P
96 Use exact interaction probabilities (stored in the file specified
97 via this option) corresponding to the flux file input in this job.
98 For complex geometries this will dramatically speed up event generation!
99 If this option is chosen then no max path length file needs to be provided.
100 Instead of calculating the interaction probabilities on the fly
101 per job you can pre-generate them using a dedicated job (see the
102 -S option) and tell the event generator to use them via the -P option.
103 This is advised when processing flux files with more than ~100k entries
104 as the time to pre-calculate the interaction probabilities becomes
105 comparable to the event generation time. For smaller flux files there is
106 less book-keeping if you just calculate them per job and on the fly.
107 -S [output_name]
108 Pre-generate flux interaction probabilities and save to root
109 output file for use with future event generation jobs. With this
110 option no events are generated, it is just a preparatory step
111 before actual event generation. When using this option care
112 should be taken to run with same arguments used for the actual
113 event generation job (same ROOT geomtery, flux file, probe
114 particle types, etc...).
115 The output root file is specified as the input for event
116 generation using the [pre_gen_prob_file] optional value of the
117 -P option. The default output interaction probabilities file
118 name is constructed as: [FLUXFILENAME].[TOPVOL].flxprobs.root.
119 Specifying [output_name] will override this.
120 Introducing multiple functionality to the executable is not
121 desirable but is less error prone than duplicating a lot of the
122 functionality in a separate application.
123 -m
124 An XML file (generated by gmxpl) with the max (density weighted)
125 path-lengths for each target material in the input ROOT geometry.
126 If no file is input, then the geometry will be scanned at MC job
127 initialization to determine those max path lengths.
128 Supplying this file can speed-up the MC job initialization.
129 -L
130 Input geometry length units, eg 'm', 'cm', 'mm', ...
131 [default: 'mm']
132 Note that typically:
133 - nd280m uses: 'mm'
134 - ...
135 -D
136 Input geometry density units, eg 'g_cm3', 'clhep_def_density_unit',...
137 [default: 'g_cm3']
138 Note that typically:
139 - nd280m uses: 'clhep_def_density_unit'
140 - ...
141 -f
142 Input 'neutrino flux'.
143 This option can be used to specify any of:
144 1 > A JNUBEAM beam simulation output file and the detector location.
145 The general sytax is:
146 -f /full/path/flux_file.root,detector_location(,list_of_neutrino_codes)
147 [Notes]
148 - For more information on the flux ntuples, see (T2K internal):
149 http://jnusrv01.kek.jp/internal/t2k/nubeam/flux/
150 - The original HBOOK JNUBAM ntuples need to be converted to a ROOT
151 format using the h2root ROOT utility.
152 - The detector location can be any of 'sk' or the near detector
153 positions 'nd1',...,'nd6','nd13' simulated in JNUBEAM.
154 See the above JNUBEAM web page for more info.
155 - The neutrino codes are the PDG ones.
156 - The (comma separated) list of neutrino codes is optional.
157 It may be used for considering only certain neutrino flavours
158 (eg. nu_e only).
159 If no neutrino list is specified then GENIE will consider all
160 possible flavours (nu_e, nu_e_bar, nu_mu, nu_mu_bar).
161 See examples below.
162 - The JNUBEAM flux ntuples are read via GENIE's GJPARCNuFlux
163 driver. This customized GENIE event generation driver
164 passes-through the complete JPARC input flux information
165 (eg parent decay kinematics / position etc) for each neutrino
166 event it generates (an additional 'flux' branch is added at
167 the output event tree).
168 [Examples]
169 - To use the SuperK flux ntuple from the flux.root file,
170 type:
171 '-f /path/flux.root,sk'
172 - To do the same as above, but considering only nu_mu and
173 nu_mu_bar, type:
174 '-f /path/flux.root,sk,14,-14'
175 - To use the 2km flux ntuple [near detector position 'nd1'
176 in the JNUBEAM flux simulation] from the flux.root file,
177 type:
178 '-f /path/flux.root,nd1'
179 - To use the nd280 flux ntuple [near detector position 'nd5'
180 in the JNUBEAM flux simulation] from the flux.root file,
181 type:
182 '-f /path/flux.root,nd5'
183 - To do the same as above, but considering only nu_e
184 type:
185 '-f /path/flux.root,nd5,12'
186 2 > A list of JNUBEAM beam simulation output files and the detector location.
187 The general sytax is:
188 -f /full/path/flux_file_prefix@first_file_number@last_file_number,detector_location(,list_of_neutrino_codes)
189 [Notes]
190 - The ".root" is assumed.
191 - All the files in the series between flux_file_prefix.first_file_number.root
192 and flux_file_prefix.last_file_number.root should exist.
193 - It is important that you set the -p option correctly. See note below.
194 - Also see notes from option 1.
195 [Examples]
196 - To use the SuperK flux ntuples from the files: flux.0.root, flux.1.root, flux.2.root, flux.3.root
197 type:
198 '-f /path/flux.@0@3,sk'
199 - To use the nd280 flux ntuple [near detector position 'nd5'
200 in the JNUBEAM flux simulation] from the files in the series
201 flux.0.root file to flux.100.root, considering only nu_e, type:
202 '-f /path/flux.@0@100,nd5,12'
203 3 > A set of histograms stored in a ROOT file.
204 The general syntax is:
205 -f /path/histogram_file.root,neutrino_code[histo_name],...
206 [Notes]
207 - The neutrino codes are the PDG ones.
208 - The 'neutrino_code[histogram_name]' part of the option can be
209 repeated multiple times (separated by commas), once for each
210 flux neutrino species you want to consider, eg
211 '-f somefile.root,12[nuehst],-12[nuebarhst],14[numuhst]'
212 - Code implicitly assumes the binning for multiple flavors
213 is the same for all histograms (no checks are made)
214 - Note that the relative normalization of the flux histograms
215 is taken into account and is reflected in the relative
216 frequency of flux neutrinos thrown by the flux driver
217 - Variable bin width flux histograms are modified to account
218 for incorrect sampling of TH1D::GetRandom() prior to ROOT
219 version 9.99.99
220 + notation such as "12[nuehst]WIDTH" forces bin width
221 adjustment no matter variable bin width histogram or not
222 + notation such as "12[nuehst]NOWIDTH" prevents bin width
223 adjustment
224 - When using flux from histograms then there is no point in using
225 a 'detailed detector geometry description' as your flux input
226 contains no directional information for those flux neutrinos.
227 The neutrino direction is conventionally set to be +z {x=0,y=0}.
228 So, when using this option you must be using a simple 'target mix'
229 See the -g option for possible geometry settings.
230 If you want to use the detailed detector geometry description
231 then you should be feeding this driver with the JNUBEAM flux
232 simulation outputs.
233 - When using flux from histograms there is no branch with neutrino
234 parent information (code, decay mode, 4-momentum at prod & decay)
235 added in the output event tree as your flux input contains no
236 such information. If you want to be getting the neutrino parent
237 information written out as well then you should be feeding this
238 driver with the JNUBEAM flux simulation outputs.
239 - Note that the relative normalization of the flux histograms is
240 taken into account and is reflected in the relative frequency
241 of flux neutrinos thrown by the flux driver
242 [Examples]
243 - To use the histogram 'h100' (representing the nu_mu flux) and
244 the histogram 'h300' (representing the nu_e flux) and the
245 histogram 'h301' (representing the nu_e_bar flux) from the
246 flux.root file in /path/
247 type:
248 '-f /path/flux.root,14[h100],12[h300],-12[h301]
249 -p
250 POT normalization of the input flux file.
251 [default: The 'standard' JNUBEAM flux ntuple normalization of
252 1E+21 POT/detector for the near detectors and
253 1E+21 POT/cm2 for the far detector]
254 That will be used to interpret the flux weights & calculate the actual
255 POT normalization for the generated event sample.
256 [Note]
257 - If you are using the multiple JNUBEAM flux file entry method it is
258 very important that you set this. It should be set to the total POT
259 of all input flux files.
260 [Examples]
261 - If you have 10 standard JNUBEAM files use '-p 10E+21'
262 -c
263 Specifies how many times to cycle a JNUBEAM flux ntuple.
264 Due to the large rejection factor when generating unweighted events
265 in the full energy range (approximately ~500 flux neutrinos will be
266 rejected before getting an interaction in nd280), an option is
267 provided to recycle the flux ntuples for a number of times.
268 That option can be used to boost the generated statistics without
269 requiring enormous flux files.
270 See also 'Note on exposure / statistics' below.
271 -e
272 Specifies how many POTs to generate.
273 If that option is set, gevgen_t2k will work out how many times it has
274 to cycle through the input flux ntuple in order to accumulate the
275 requested statistics. The program will stop at the earliest complete
276 flux ntuple cycle after accumulating the required statistics, so the
277 actual statistics will 'slightly' overshoot that number.
278 See also 'Note on exposure / statistics' below.
279 -E
280 Specifies how many POTs to generate.
281 That option is similar to -e but the program will stop immediatelly
282 after the requested POT has been accumulated. That reduces the
283 generated POT overshoot of the requested POT, but the POT calculation
284 may not be as exact as with -e.
285 See also 'Note on exposure / statistics' below.
286 -n
287 Specifies how many events to generate.
288
289 --------------------------
290 [Note on setting the exposure / statistics]
291 All -c, -e (-E) and -n options can be used to set the exposure.
292 - If the input flux is a JNUBEAM ntuple then any of these options can
293 be used (only one at a time).
294 If no option is set, then the program will automatically set '-c 1'
295 - If the input flux is described with histograms then only the -n
296 option is available.
297 --------------------------
298
299 -o
300 Sets the prefix of the output event file.
301 The output filename is built as:
302 [prefix].[run_number].[event_tree_format].[file_format]
303 The default output filename is:
304 gntp.[run_number].ghep.root
305 This cmd line arguments lets you override 'gntp'
306 -R
307 Tell the flux driver to start looping over the flux ntuples with a
308 random offset. May be necessary to avoid biases introduced by always
309 starting at the same point when using very large flux input files.
310 --seed
311 Random number seed.
312 --cross-sections
313 Name (incl. full path) of an XML file with pre-computed
314 cross-section values used for constructing splines.
315 --tune
316 Specifies a GENIE comprehensive neutrino interaction model tune.
317 [default: "Default"].
318 --message-thresholds
319 Allows users to customize the message stream thresholds.
320 The thresholds are specified using an XML file.
321 See $GENIE/config/Messenger.xml for the XML schema.
322 --unphysical-event-mask
323 Allows users to specify a 16-bit mask to allow certain types of
324 unphysical events to be written in the output file.
325 [default: all unphysical events are rejected]
326 --event-record-print-level
327 Allows users to set the level of information shown when the event
328 record is printed in the screen. See GHepRecord::Print().
329 --mc-job-status-refresh-rate
330 Allows users to customize the refresh rate of the status file.
331 --cache-file
332 Allows users to specify a cache file so that the cache can be
333 re-used in subsequent MC jobs.
334
335 *** Examples:
336
337 (1) shell% gevgen_t2k
338 -r 1001
339 -f /data/t2k/flux/07a/jnubeam001.root,nd5
340 -g /data/t2k/geom/nd280.root
341 -L mm -D clhep_def_density_unit
342 -e 5E+17
343 --cross-sections /data/t2k/xsec/xsec.xml
344
345 Generate events (run number 1001) using the JNUBEAM flux ntuple in
346 /data/t2k/flux/07a/jnubeam001.root & picking up the flux entries for
347 the detector location nd5 (:nd280m). The job will load the nd280
348 detector geometry description from /data/t2k/geom/nd280.root and
349 use it thinking that the geometry length unit is 'mm' and the density
350 unit is 'clhep_def_density_unit' (g_cm3 / 0.62415185185E+19)
351 The job will stop on the first complete flux ntuple cycle after
352 generating 5E+17 POT. Pre-computed cross-sections for all relevant
353 initial states are loaded from /data/t2k/xsec/xsec.xml.
354
355 (2) shell% gevgen_t2k
356 -r 1001
357 -f /data/t2k/flux/07a/jnubeam001.root,nd5
358 -g /data/t2k/geom/nd280.root
359 -L mm -D clhep_def_density_unit
360 -c 100
361 --cross-sections /data/t2k/xsec/xsec.xml
362
363 As before, but now the job will stop after 100 flux ntuple cycles -
364 whatever POT & number of events that may correspond to.
365
366 (3) shell% gevgen_t2k
367 -r 1001
368 -f /data/t2k/flux/07a/jnubeam001.root,nd5
369 -g /data/t2k/geom/nd280.root
370 -L mm -D clhep_def_density_unit
371 -n 100000
372 --cross-sections /data/t2k/xsec/xsec.xml
373
374 As before, but now the job will stop after generating 100000 events -
375 whatever POT & number of flux ntuple cycles that may correspond to.
376
377 (4) shell% gevgen_t2k
378 -r 1001
379 -f /data/t2k/flux/07a/jnubeam001.root,nd5,-12,12
380 -g /data/t2k/geom/nd280.root
381 -L mm -D clhep_def_density_unit
382 -n 100000
383 --cross-sections /data/t2k/xsec/xsec.xml
384
385 As before, but now the job will consider flux nu_e and nu_e_bar only!
386
387 (5) shell% gevgen_t2k
388 -r 1001
389 -f /data/t2k/flux/07a/jnubeam001.root,sk
390 -g 1000080160[0.95],1000010010[0.05]
391 -n 50000
392 --cross-sections /data/t2k/xsec/xsec.xml
393
394 Generate events (run number 1001) using the JNUBEAM flux ntuple in
395 /data/t2k/flux/07a/jnubeam001.root & picking up the flux entries for
396 the SuperK detector location. This time, the job will not use any
397 detailed detector geometry description but just (95% O16 + 5% H)
398 target mix. The job will stop after generating 50000 events.
399
400 (6) shell% gevgen_t2k
401 -r 1001
402 -f /data/t2k/flux/hst/flux.root,12[h100],-12[h101],14[h200]
403 -g 1000080160[0.95],1000010010[0.05]
404 -n 50000
405 --cross-sections /data/t2k/xsec/xsec.xml
406
407 As before, but in this case the flux description is not based on a JNUBEAM
408 ntuple but a set of histograms at the /data/t2k/flux/hst/flux.root file:
409 The histogram named 'h100' will be used for the nu_e flux, 'h101' will
410 will be used for the nu_e_bar flux, and 'h200' for the nu_mu flux.
411
412
413 Please read the GENIE User Manual for more information.
414
415\author Costas Andreopoulos <c.andreopoulos \at cern.ch>
416 University of Liverpool
417
418\created February 05, 2008
419
420\cpright Copyright (c) 2003-2025, The GENIE Collaboration
421 For the full text of the license visit http://copyright.genie-mc.org
422
423*/
424//_________________________________________________________________________________________
425
426#include <cassert>
427#include <cstdlib>
428#include <string>
429#include <sstream>
430#include <vector>
431#include <map>
432
433#include <TSystem.h>
434#include <TTree.h>
435#include <TFile.h>
436#include <TH1D.h>
437#include <TMath.h>
438#include <TGeoVolume.h>
439#include <TGeoShape.h>
440#include <TList.h>
441#include <TObject.h>
442
464
465#ifdef __GENIE_FLUX_DRIVERS_ENABLED__
468#endif
469
470#ifdef __GENIE_GEOM_DRIVERS_ENABLED__
474#endif
475
476using std::string;
477using std::vector;
478using std::map;
479using std::ostringstream;
480
481using namespace genie;
482
483void GetCommandLineArgs (int argc, char ** argv);
484void PrintSyntax (void);
485
486// Default options (override them using the command line arguments):
487//
488string kDefOptGeomLUnits = "mm"; // default geometry length units
489string kDefOptGeomDUnits = "g_cm3"; // default geometry density units
490NtpMCFormat_t kDefOptNtpFormat = kNFGHEP; // default event tree format
491double kDefOptFluxNorm = 1E+21; // std JNUBEAM flux ntuple norm. (POT*detector [NDs] or POT*cm^2 [SK])
492string kDefOptEvFilePrefix = "gntp";
493
494// User-specified options:
495//
496Long_t gOptRunNu; // run number
497bool gOptUsingRootGeom = false; // using root geom or target mix?
498bool gOptUsingHistFlux = false; // using JNUBEAM flux ntuples or flux from histograms?
499map<int,double> gOptTgtMix; // target mix (tgt pdg -> wght frac) / if not using detailed root geom
500map<int,TH1D*> gOptFluxHst; // flux histos (nu pdg -> spectrum) / if not using beam sim ntuples
501string gOptRootGeom; // input ROOT file with realistic detector geometry
502string gOptRootGeomTopVol = ""; // input geometry top event generation volume
503double gOptGeomLUnits = 0; // input geometry length units
504double gOptGeomDUnits = 0; // input geometry density units
505string gOptExtMaxPlXml; // max path lengths XML file for input geometry
506string gOptFluxFile; // ROOT file with JNUBEAM flux ntuple
507string gOptDetectorLocation; // detector location ('sk','nd1','nd2',...)
508double gOptFluxNorm; // JNUBEAM flux ntuple normalization
509PDGCodeList gOptFluxNtpNuList(false); // JNUBEAM flux ntuple neutrinos to consider (can be used for considering certain flavours only)
510int gOptFluxNCycles; // number of flux ntuple cycles
511int gOptNev; // number of events to generate
512double gOptPOT; // exposure (in POT)
513bool gOptExitAtEndOfFullFluxCycles; // once POT >= requested_POT, stop at once or at the end of the flux cycle?
514string gOptEvFilePrefix; // event file prefix
515bool gOptUseFluxProbs = false; // use pre-calculated flux interaction probs instead of estimating them using the max paths
516bool gOptSaveFluxProbsFile = false; // special mode: no events generated, calculate and save flux interaction probs to root file
517string gOptFluxProbFileName; // filename for file containg flux probs
518string gOptSaveFluxProbsFileName; // output filename for pre-generated flux probabilities
519bool gOptRandomFluxOffset = false; // start looping over flux file from random start entry
520long int gOptRanSeed; // random number seed
521string gOptInpXSecFile; // cross-section splines
522
523//____________________________________________________________________________
524int main(int argc, char ** argv)
525{
526 // Parse command line arguments
527 GetCommandLineArgs(argc,argv);
528
529
530 if ( ! RunOpt::Instance()->Tune() ) {
531 LOG("gmkspl", pFATAL) << " No TuneId in RunOption";
532 exit(-1);
533 }
535
536 // Initialization of random number generators, cross-section table,
537 // messenger thresholds, cache file
538 utils::app_init::MesgThresholds(RunOpt::Instance()->MesgThresholdFiles());
542
543 // Set GHEP print level
544 GHepRecord::SetPrintLevel(RunOpt::Instance()->EventRecordPrintLevel());
545
546 // *************************************************************************
547 // * Create / configure the geometry driver
548 // *************************************************************************
549 GeomAnalyzerI * geom_driver = 0;
550 double zmin=0, zmax=0;
551
553 //
554 // *** Using a realistic root-based detector geometry description
555 //
556
557 // creating & configuring a root geometry driver
560 rgeom -> SetLengthUnits (gOptGeomLUnits);
561 rgeom -> SetDensityUnits (gOptGeomDUnits);
562 // getting the bounding box dimensions, before setting topvolume,
563 // along z so as to set the appropriate upstream generation surface
564 // for the JPARC flux driver
565 TGeoVolume * topvol = rgeom->GetGeometry()->GetTopVolume();
566 TGeoShape * bounding_box = topvol->GetShape();
567 bounding_box->GetAxisRange(3, zmin, zmax);
568 zmin *= rgeom->LengthUnits();
569 zmax *= rgeom->LengthUnits();
570 // now update to the requested topvolume for use in recursive exhaust method
571 rgeom -> SetTopVolName (gOptRootGeomTopVol);
572 topvol = rgeom->GetGeometry()->GetTopVolume();
573 if(!topvol) {
574 LOG("gevgen_t2k", pFATAL) << "Null top ROOT geometry volume!";
575 exit(1);
576 }
577 // switch on/off volumes as requested
578 if ( (gOptRootGeomTopVol[0] == '+') || (gOptRootGeomTopVol[0] == '-') ) {
579 bool exhaust = (*gOptRootGeomTopVol.c_str() == '+');
581 }
582
583 // casting to the GENIE geometry driver interface
584 geom_driver = dynamic_cast<GeomAnalyzerI *> (rgeom);
585 }
586 else {
587 //
588 // *** Using a 'point' geometry with the specified target mix
589 // *** ( = a list of targets with their corresponding weight fraction)
590 //
591
592 // creating & configuring a point geometry driver
595 // casting to the GENIE geometry driver interface
596 geom_driver = dynamic_cast<GeomAnalyzerI *> (pgeom);
597 }
598
599 // *************************************************************************
600 // * Create / configure the flux driver
601 // *************************************************************************
602 GFluxI * flux_driver = 0;
603
604 flux::GJPARCNuFlux * jparc_flux_driver = 0;
605 flux::GCylindTH1Flux * hst_flux_driver = 0;
606
607 if(!gOptUsingHistFlux) {
608 //
609 // *** Using the detailed JPARC neutrino flux desription by feeding-in
610 // *** the JNUBEAM flux simulation ntuples
611 //
612
613 // creating JPARC neutrino flux driver
614 jparc_flux_driver = new flux::GJPARCNuFlux;
615 // before loading the beam sim data set whether to use a random offset when looping
616 if(gOptRandomFluxOffset == false) jparc_flux_driver->DisableOffset();
617 // specify input JNUBEAM file & detector location
618 bool beam_sim_data_success = jparc_flux_driver->LoadBeamSimData(gOptFluxFile, gOptDetectorLocation);
619 if(!beam_sim_data_success) {
620 LOG("gevgen_t2k", pFATAL)
621 << "The flux driver has not been properly configured. Exiting";
622 exit(1);
623 }
624 // specify JNUBEAM normalization
625 jparc_flux_driver->SetFilePOT(gOptFluxNorm);
626 // specify upstream generation surface
627 jparc_flux_driver->SetUpstreamZ(zmin);
628 // specify which flavours to consider -
629 // if no neutrino list was specified then the MC job will consider all flavours
630 if( gOptFluxNtpNuList.size() > 0 ) {
631 jparc_flux_driver->SetFluxParticles(gOptFluxNtpNuList);
632 }
633
634 // casting to the GENIE flux driver interface
635 flux_driver = dynamic_cast<GFluxI *> (jparc_flux_driver);
636 }
637 else {
638 //
639 // *** Using fluxes from histograms (for all specified neutrino species)
640 //
641
642 // creating & configuring a generic GCylindTH1Flux flux driver
643 TVector3 bdir (0,0,1); // dir along +z
644 TVector3 bspot(0,0,0);
645 hst_flux_driver = new flux::GCylindTH1Flux;
646 hst_flux_driver->SetNuDirection (bdir);
647 hst_flux_driver->SetBeamSpot (bspot);
648 hst_flux_driver->SetTransverseRadius (-1);
649 map<int,TH1D*>::iterator it = gOptFluxHst.begin();
650 for( ; it != gOptFluxHst.end(); ++it) {
651 int pdg_code = it->first;
652 TH1D * spectrum = new TH1D(*(it->second));
653 hst_flux_driver->AddEnergySpectrum(pdg_code, spectrum);
654 }
655 // casting to the GENIE flux driver interface
656 flux_driver = dynamic_cast<GFluxI *> (hst_flux_driver);
657 }
658
659 // *************************************************************************
660 // * Create/configure the event generation driver
661 // *************************************************************************
662 GMCJDriver * mcj_driver = new GMCJDriver;
664 mcj_driver->UseFluxDriver(flux_driver);
665 mcj_driver->UseGeomAnalyzer(geom_driver);
667 // do not calculate probability scales if using pre-generated flux probs
668 bool calc_prob_scales = (gOptSaveFluxProbsFile || gOptUseFluxProbs) ? false : true;
669 mcj_driver->Configure(calc_prob_scales);
670 mcj_driver->UseSplines();
671 mcj_driver->ForceSingleProbScale();
672
673 // *************************************************************************
674 // * If specified use pre-calculated flux interaction probabilities instead
675 // * of preselecting based on max path lengths
676 // *************************************************************************
677
679
680 bool success = false;
681
682 // set flux probs output file name
684 // default output name is ${FLUFILENAME}.${TOPVOL}.flxprobs.root
685 string basename = gOptFluxFile.substr(gOptFluxFile.rfind("/")+1);
686 string name = basename.substr(0, basename.rfind("."));
687 if(gOptRootGeomTopVol.length()>0)
688 name += "."+gOptRootGeomTopVol+".flxprobs.root";
689 else
690 name += ".master.flxprobs.root";
691 // if specified override with cmd line option
693 // Tell the driver save pre-generated probabilities to an output file
694 mcj_driver->SaveFluxProbabilities(name);
695 }
696
697 // Either load pre-generated flux probabilities
698 if(gOptFluxProbFileName.size() > 0){
699 success = mcj_driver->LoadFluxProbabilities(gOptFluxProbFileName);
700 }
701 // Or pre-calculate them
702 else success = mcj_driver->PreCalcFluxProbabilities();
703
704 if(success){
705 LOG("gevgen_t2k", pNOTICE)
706 << "Successfully calculated/loaded flux interaction probabilities!";
707 // Print out a list of expected number of events per POT and per cycle
708 // based on the pre-generated flux interaction probabilities
709 map<int, double> sum_probs_map = mcj_driver->SumFluxIntProbs();
710 map<int, double>::const_iterator sum_probs_it = sum_probs_map.begin();
711 double ntot_per_pot = 0.0;
712 double ntot_per_cycle = 0.0;
713 double pscale = mcj_driver->GlobProbScale();
714 double pot_1cycle = jparc_flux_driver->POT_1cycle();
715 LOG("T2KProdInfo", pNOTICE) <<
716 "Expected event rates based on flux interaction probabilities:";
717 for(; sum_probs_it != sum_probs_map.end(); sum_probs_it++){
718 double sum_probs = sum_probs_it->second;
719 double nevts_per_cycle = sum_probs / pscale; // take into account rescale
720 double nevts_per_pot = sum_probs/pot_1cycle; // == (sum_probs*pscale)/(pot_1cycle*pscale)
721 ntot_per_pot += nevts_per_pot;
722 ntot_per_cycle += nevts_per_cycle;
723 LOG("T2KProdInfo", pNOTICE) <<
724 " PDG "<< sum_probs_it->first << ": " << nevts_per_cycle <<
725 " Events/Cycle, "<< nevts_per_pot << " Events/POT";
726 }
727 LOG("T2KProdInfo", pNOTICE) << " -----------";
728 LOG("T2KProdInfo", pNOTICE) <<
729 " All neutrino species: " << ntot_per_cycle <<
730 " Events/Cycle, "<< ntot_per_pot << " Events/POT";
731 LOG("T2KProdInfo", pNOTICE) <<
732 "N.B. This assumes input flux file corresponds to "<< pot_1cycle <<
733 "POT, ensure this is correct if using these numbers!";
734 }
735 else {
736 LOG("gevgen_t2k", pFATAL)
737 << "Failed to calculated/loaded flux interaction probabilities!";
738 return 1;
739 }
740
741 // Exit now if just pre-generating interaction probabilities
743 LOG("gevgen_t2k", pNOTICE)
744 << "Will not generate events - just pre-calculating flux interaction"
745 << "probabilities";
746 return 0;
747 }
748 } // Pre-calculated flux interaction probabilities
749
750 // *************************************************************************
751 // * Work out number of cycles for current exposure settings
752 // *************************************************************************
753
754 if(!gOptUsingHistFlux) {
755 // If a number of POT was requested, then work out how many flux ntuple
756 // cycles are required for accumulating those statistics
757 if(gOptPOT>0) {
758 double fpot_1c = jparc_flux_driver->POT_1cycle(); // flux POT / cycle
759 double psc = mcj_driver->GlobProbScale(); // interaction prob. scale
760 double pot_1c = fpot_1c / psc; // actual POT / cycle
761 int ncycles = (int) TMath::Max(1., TMath::Ceil(gOptPOT/pot_1c));
762
763 LOG("gevgen_t2k", pNOTICE)
764 << " *** POT/cycle: " << pot_1c;
765 LOG("gevgen_t2k", pNOTICE)
766 << " *** Requested POT will be accumulated in: "
767 << ncycles << " flux ntuple cycles";
768
769 jparc_flux_driver->SetNumOfCycles(ncycles);
770 }
771 // If a number of events was requested, then set the number of flux
772 // ntuple cycles to 'infinite'
773 else if(gOptNev>0) {
774 jparc_flux_driver->SetNumOfCycles(0);
775 }
776 // Just set the number of cycles to the requested value
777 else {
778 jparc_flux_driver->SetNumOfCycles(gOptFluxNCycles);
779 }
780 }
781
782 // *************************************************************************
783 // * Prepare for writing the output event tree & status file
784 // *************************************************************************
785
786 // Initialize an Ntuple Writer to save GHEP records into a TTree
789 ntpw.Initialize();
790
791 // Add a custom-branch at the standard GENIE event tree so that
792 // info on the flux neutrino parent particle can be passed-through
793 flux::GJPARCNuFluxPassThroughInfo * flux_info = 0;
794 if(!gOptUsingHistFlux) {
795 TBranch * flux = ntpw.EventTree()->Branch("flux",
796 "genie::flux::GJPARCNuFluxPassThroughInfo", &flux_info, 32000, 1);
797 assert(flux);
798 flux->SetAutoDelete(kFALSE);
799 }
800
801 // Create a MC job monitor for a periodically updated status file
802 GMCJMonitor mcjmonitor(gOptRunNu);
803 mcjmonitor.SetRefreshRate(RunOpt::Instance()->MCJobStatusRefreshRate());
804
805 // *************************************************************************
806 // * Event generation loop
807 // *************************************************************************
808
809 int ievent = 0;
810 while (1)
811 {
812 LOG("gevgen_t2k", pNOTICE)
813 << " *** Generating event............ " << ievent;
814
815 // In case the required statistics was expressed as 'number of events'
816 // then quit if that number has been generated
817 if(ievent == gOptNev) break;
818
819 // In case the required statistics was expressed as 'number of POT' and
820 // the user does not want to wait till the end of the flux cycle to exit
821 // the event loop, then quit if the requested POT has been generated.
822 // In this case the computed POT may not be as accurate as if the program
823 // was waiting for the current flux cycle to be completed.
825 double fpot = jparc_flux_driver->POT_curravg(); // current POT in flux file
826 double psc = mcj_driver->GlobProbScale(); // interaction prob. scale
827 double pot = fpot / psc; // POT for generated sample
828 if(pot >= gOptPOT) break;
829 }
830
831 // Generate a single event using neutrinos coming from the specified flux
832 // and hitting the specified geometry or target mix
833 EventRecord * event = mcj_driver->GenerateEvent();
834
835 // Check whether a null event was returned due to the flux driver reaching
836 // the end of the input flux ntuple - exit the event generation loop
837 if(!event && jparc_flux_driver->End()) {
838 LOG("gevgen_t2k", pNOTICE)
839 << "** The JPARC flux driver read all the input flux ntuple entries";
840 break;
841 }
842 if(!event) {
843 LOG("gevgen_t2k", pERROR)
844 << "Got a null generated neutino event! Retrying ...";
845 continue;
846 }
847 LOG("gevgen_t2k", pINFO)
848 << "Generated event: " << *event;
849
850 // A valid event was generated: extract flux info (parent decay/prod
851 // position/kinematics) for that simulated event so that it can be
852 // passed-through.
853 // Can only do so if I am generating events using the JNUBEAM flux
854 // ntuples, not simple histograms
855 if(!gOptUsingHistFlux) {
856 flux_info = new flux::GJPARCNuFluxPassThroughInfo(
857 jparc_flux_driver->PassThroughInfo());
858 LOG("gevgen_t2k", pINFO)
859 << "Pass-through flux info associated with generated event: "
860 << *flux_info;
861 }
862
863 // Add event at the output ntuple, refresh the mc job monitor & clean-up
864 ntpw.AddEventRecord(ievent, event);
865 mcjmonitor.Update(ievent,event);
866 delete event;
867 if(flux_info) delete flux_info;
868 ievent++;
869 } //1
870
871 LOG("gevgen_t2k", pNOTICE)
872 << "The GENIE MC job is done generaing events - Cleaning up & exiting...";
873
874 // *************************************************************************
875 // * Print job statistics &
876 // * calculate normalization factor for the generated sample
877 // *************************************************************************
879 {
880 // POT normalization will only be calculated if event generation was based
881 // on beam simulation outputs (not just histograms) & a detailed detector
882 // geometry description.
883 double fpot = jparc_flux_driver->POT_curravg(); // current POT in flux file
884 double psc = mcj_driver->GlobProbScale(); // interaction prob. scale
885 double pot = fpot / psc; // POT for generated sample
886 // Get nunber of flux neutrinos read-in by flux friver, number of flux
887 // neutrinos actually thrown to the event generation driver and number
888 // of neutrino interactions actually generated
889 long int nflx_evg = mcj_driver -> NFluxNeutrinos();
890 long int nflx = jparc_flux_driver -> NFluxNeutrinos();
891 long int nev = ievent;
892
893 LOG("gevgen_t2k", pNOTICE)
894 << "\n >> Actual JNUBEAM flux file normalization: " << fpot
895 << " POT * " << ((gOptDetectorLocation == "sk") ? "cm^2" : "det")
896 << "\n >> Interaction probability scaling factor: " << psc
897 << "\n >> N of flux v read-in by flux driver: " << nflx
898 << "\n >> N of flux v thrown to event gen driver: " << nflx_evg
899 << "\n >> N of generated v interactions: " << nev
900 << "\n ** Normalization for generated sample: " << pot
901 << " POT * " << ((gOptDetectorLocation == "sk") ? "cm^2" : "det");
902
903 ntpw.EventTree()->SetWeight(pot); // POT
904 }
905
906 // *************************************************************************
907 // * MC job meta-data
908 // *************************************************************************
909
911
912 metadata -> jnubeam_file = gOptFluxFile;
913 metadata -> detector_location = gOptDetectorLocation;
914 metadata -> geom_file = gOptRootGeom;
915 metadata -> geom_top_volume = gOptRootGeomTopVol;
916 metadata -> geom_length_units = gOptGeomLUnits;
917 metadata -> geom_density_units = gOptGeomDUnits;
918 metadata -> using_root_geom = gOptUsingRootGeom;
919 metadata -> target_mix = gOptTgtMix;
920 metadata -> using_hist_flux = gOptUsingHistFlux;
921 metadata -> flux_hists = gOptFluxHst;
922
923 ntpw.EventTree()->GetUserInfo()->Add(metadata);
924
925 // *************************************************************************
926 // * Save & clean-up
927 // *************************************************************************
928
929 // Save the generated event tree & close the output file
930 ntpw.Save();
931
932 // Clean-up
933 delete geom_driver;
934 delete flux_driver;
935 delete mcj_driver;
936 map<int,TH1D*>::iterator it = gOptFluxHst.begin();
937 for( ; it != gOptFluxHst.end(); ++it) {
938 TH1D * spectrum = it->second;
939 if(spectrum) delete spectrum;
940 }
941 gOptFluxHst.clear();
942
943 LOG("gevgen_t2k", pNOTICE) << "Done!";
944
945 return 0;
946}
947//____________________________________________________________________________
948void GetCommandLineArgs(int argc, char ** argv)
949{
950 LOG("gevgen_t2k", pINFO) << "Parsing command line arguments";
951
952 // Common run options. Set defaults and read.
955
956 // Parse run options for this app
957
958 CmdLnArgParser parser(argc,argv);
959
960 // help?
961 bool help = parser.OptionExists('h');
962 if(help) {
963 PrintSyntax();
964 exit(0);
965 }
966
967 // run number:
968 if( parser.OptionExists('r') ) {
969 LOG("gevgen_t2k", pDEBUG) << "Reading MC run number";
970 gOptRunNu = parser.ArgAsLong('r');
971 } else {
972 LOG("gevgen_t2k", pDEBUG) << "Unspecified run number - Using default";
973 gOptRunNu = 0;
974 } //-r
975
976 //
977 // *** geometry
978 //
979
980 string geom = "";
981 string lunits, dunits;
982 if( parser.OptionExists('g') ) {
983 LOG("gevgen_t2k", pDEBUG) << "Getting input geometry";
984 geom = parser.ArgAsString('g');
985
986 // is it a ROOT file that contains a ROOT geometry?
987 bool accessible_geom_file =
988 ! (gSystem->AccessPathName(geom.c_str()));
989 if (accessible_geom_file) {
991 gOptUsingRootGeom = true;
992 }
993 } else {
994 LOG("gevgen_t2k", pFATAL)
995 << "No geometry option specified - Exiting";
996 PrintSyntax();
997 exit(1);
998 } //-g
999
1000 if(gOptUsingRootGeom) {
1001 // using a ROOT geometry - get requested geometry units
1002
1003 // legth units:
1004 if( parser.OptionExists('L') ) {
1005 LOG("gevgen_t2k", pDEBUG)
1006 << "Checking for input geometry length units";
1007 lunits = parser.ArgAsString('L');
1008 } else {
1009 LOG("gevgen_t2k", pDEBUG) << "Using default geometry length units";
1011 } // -L
1012 // density units:
1013 if( parser.OptionExists('D') ) {
1014 LOG("gevgen_t2k", pDEBUG)
1015 << "Checking for input geometry density units";
1016 dunits = parser.ArgAsString('D');
1017 } else {
1018 LOG("gevgen_t2k", pDEBUG) << "Using default geometry density units";
1019 dunits = kDefOptGeomDUnits;
1020 } // -D
1023
1024 // check whether an event generation volume name has been
1025 // specified -- default is the 'top volume'
1026 if( parser.OptionExists('t') ) {
1027 LOG("gevgen_t2k", pDEBUG) << "Checking for input volume name";
1028 gOptRootGeomTopVol = parser.ArgAsString('t');
1029 } else {
1030 LOG("gevgen_t2k", pDEBUG) << "Using the <master volume>";
1031 } // -t
1032
1033 // check whether an XML file with the maximum (density weighted)
1034 // path lengths for each detector material is specified -
1035 // otherwise will compute the max path lengths at job init
1036 if( parser.OptionExists('m') ) {
1037 LOG("gevgen_t2k", pDEBUG)
1038 << "Checking for maximum path lengths XML file";
1039 gOptExtMaxPlXml = parser.ArgAsString('m');
1040 } else {
1041 LOG("gevgen_t2k", pDEBUG)
1042 << "Will compute the maximum path lengths at job init";
1043 gOptExtMaxPlXml = "";
1044 } // -m
1045 } // using root geom?
1046
1047 else {
1048 // User has specified a target mix.
1049 // Decode the list of target pdf codes & their corresponding weight fraction
1050 // (specified as 'pdg_code_1[fraction_1],pdg_code_2[fraction_2],...')
1051 // See documentation on top section of this file.
1052 //
1053 gOptTgtMix.clear();
1054 vector<string> tgtmix = utils::str::Split(geom,",");
1055 if(tgtmix.size()==1) {
1056 int pdg = atoi(tgtmix[0].c_str());
1057 double wgt = 1.0;
1058 gOptTgtMix.insert(map<int, double>::value_type(pdg, wgt));
1059 } else {
1060 vector<string>::const_iterator tgtmix_iter = tgtmix.begin();
1061 for( ; tgtmix_iter != tgtmix.end(); ++tgtmix_iter) {
1062 string tgt_with_wgt = *tgtmix_iter;
1063 string::size_type open_bracket = tgt_with_wgt.find("[");
1064 string::size_type close_bracket = tgt_with_wgt.find("]");
1065 if (open_bracket ==string::npos ||
1066 close_bracket==string::npos)
1067 {
1068 LOG("gevgen_t2k", pFATAL)
1069 << "You made an error in specifying the target mix";
1070 PrintSyntax();
1071 exit(1);
1072 }
1073 string::size_type ibeg = 0;
1074 string::size_type iend = open_bracket;
1075 string::size_type jbeg = open_bracket+1;
1076 string::size_type jend = close_bracket;
1077 int pdg = atoi(tgt_with_wgt.substr(ibeg,iend-ibeg).c_str());
1078 double wgt = atof(tgt_with_wgt.substr(jbeg,jend-jbeg).c_str());
1079 LOG("gevgen_t2k", pDEBUG)
1080 << "Adding to target mix: pdg = " << pdg << ", wgt = " << wgt;
1081 gOptTgtMix.insert(map<int, double>::value_type(pdg, wgt));
1082
1083 }// tgtmix_iter
1084 } // >1 materials in mix
1085 } // using tgt mix?
1086
1087 //
1088 // *** flux
1089 //
1090
1091 if( parser.OptionExists('f') ) {
1092 LOG("gevgen_t2k", pDEBUG) << "Getting input flux";
1093 string flux = parser.ArgAsString('f');
1094 gOptUsingHistFlux = (flux.find("[") != string::npos);
1095
1096 if(!gOptUsingHistFlux) {
1097 // Using JNUBEAM flux ntuples
1098 // Extract JNUBEAM flux (root) file name & detector location.
1099 // Also extract the list of JNUBEAM neutrinos to consider (if none
1100 // is specified here then I will consider all flavours)
1101 //
1102 vector<string> fluxv = utils::str::Split(flux,",");
1103 if(fluxv.size()<2) {
1104 LOG("gevgen_t2k", pFATAL)
1105 << "You need to specify both a flux ntuple ROOT file "
1106 << " _AND_ a detector location";
1107 PrintSyntax();
1108 exit(1);
1109 }
1110 gOptFluxFile = fluxv[0];
1111 gOptDetectorLocation = fluxv[1];
1112
1113 // Extract the list of neutrinos to consider (if any).
1114 //
1115 for(unsigned int inu = 2; inu < fluxv.size(); inu++)
1116 {
1117 gOptFluxNtpNuList.push_back( atoi(fluxv[inu].c_str()) );
1118 }
1119
1120 } else {
1121 // Using flux from histograms
1122 // Extract the root file name & the list of histogram names & neutrino
1123 // species (specified as 'filename,histo1[species1],histo2[species2],...')
1124 // for variable width histograms, default to multiply in the width
1125 // histo1[species1]WIDTH = multiply in the width
1126 // histo1[species1]NOWIDTH = don't multiply in the width
1127 // See documentation on top section of this file.
1128 //
1129 vector<string> fluxv = utils::str::Split(flux,",");
1130 if(fluxv.size()<2) {
1131 LOG("gevgen_t2k", pFATAL)
1132 << "You need to specify both a flux ntuple ROOT file "
1133 << " _AND_ a detector location";
1134 PrintSyntax();
1135 exit(1);
1136 }
1137 gOptFluxFile = fluxv[0];
1138 bool accessible_flux_file =
1139 !(gSystem->AccessPathName(gOptFluxFile.c_str()));
1140 if (!accessible_flux_file) {
1141 LOG("gevgen_t2k", pFATAL)
1142 << "Can not access flux file: " << gOptFluxFile;
1143 PrintSyntax();
1144 exit(1);
1145 }
1146 // Extract energy spectra for all specified neutrino species
1147 TFile flux_file(gOptFluxFile.c_str(), "read");
1148 for(unsigned int inu=1; inu<fluxv.size(); inu++) {
1149 string nutype_and_histo = fluxv[inu];
1150 string::size_type open_bracket = nutype_and_histo.find("[");
1151 string::size_type close_bracket = nutype_and_histo.find("]");
1152 if (open_bracket ==string::npos ||
1153 close_bracket==string::npos)
1154 {
1155 LOG("gevgen_t2k", pFATAL)
1156 << "You made an error in specifying the flux histograms";
1157 PrintSyntax();
1158 exit(1);
1159 }
1160 string::size_type ibeg = 0;
1161 string::size_type iend = open_bracket;
1162 string::size_type jbeg = open_bracket+1;
1163 string::size_type jend = close_bracket;
1164 string nutype = nutype_and_histo.substr(ibeg,iend-ibeg);
1165 string histo = nutype_and_histo.substr(jbeg,jend-jbeg);
1166 string extra = nutype_and_histo.substr(jend+1,string::npos);
1167 std::transform(extra.begin(),extra.end(),extra.begin(),::toupper);
1168 // access specified histogram from the input root file
1169 TH1D * ihst = (TH1D*) flux_file.Get(histo.c_str());
1170 if(!ihst) {
1171 LOG("gevgen_t2k", pFATAL)
1172 << "Can not find histogram: " << histo
1173 << " in flux file: " << gOptFluxFile;
1174 PrintSyntax();
1175 exit(1);
1176 }
1177 // create a local copy of the input histogram
1178 TString origname = ihst->GetName();
1179 TString tmpname; tmpname.Form("%s_", origname.Data());
1180 // Copy in the flux histogram from the root file
1181 // use Clone rather than assuming fix bin widths and rebooking
1182 TH1D* spectrum = (TH1D*)ihst->Clone();
1183 spectrum->SetNameTitle(tmpname.Data(),ihst->GetName());
1184 spectrum->SetDirectory(0);
1185
1186 // get rid of original
1187 delete ihst;
1188 // rename copy
1189 spectrum->SetName(origname.Data());
1190
1191 bool force_binwidth = false;
1192#if ROOT_VERSION_CODE <= ROOT_VERSION(9,99,99)
1193 // GetRandom() sampling on variable bin width histograms does not
1194 // correctly account for bin widths for all versions of ROOT prior
1195 // to (currently forever). At some point this might change and
1196 // the necessity of this code snippet will go away
1197 TAxis* xaxis = spectrum->GetXaxis();
1198 if (xaxis->IsVariableBinSize()) force_binwidth = true;
1199#endif
1200 if ( extra == "WIDTH" ) force_binwidth = true;
1201 if ( extra == "NOWIDTH" ) force_binwidth = false;
1202 if ( force_binwidth ) {
1203 LOG("gevgen_t2k", pNOTICE)
1204 << "multiplying by bin width for histogram " << histo
1205 << " as " << spectrum->GetName() << " for nutype " << nutype
1206 << " from " << gOptFluxFile;
1207 for(int ibin = 1; ibin <= spectrum->GetNbinsX(); ++ibin) {
1208 double data = spectrum->GetBinContent(ibin);
1209 double width = spectrum->GetBinWidth(ibin);
1210 spectrum->SetBinContent(ibin,data*width);
1211 }
1212 }
1213
1214 // convert neutrino name -> pdg code
1215 int pdg = atoi(nutype.c_str());
1217 LOG("gevgen_t2k", pFATAL)
1218 << "Unknown neutrino type: " << nutype;
1219 PrintSyntax();
1220 exit(1);
1221 }
1222 // store flux neutrino code / energy spectrum
1223 LOG("gevgen_t2k", pDEBUG)
1224 << "Adding energy spectrum for flux neutrino: pdg = " << pdg;
1225 gOptFluxHst.insert(map<int, TH1D*>::value_type(pdg, spectrum));
1226 }//inu
1227 if(gOptFluxHst.size()<1) {
1228 LOG("gevgen_t2k", pFATAL)
1229 << "You have not specified any flux histogram!";
1230 PrintSyntax();
1231 exit(1);
1232 }
1233 flux_file.Close();
1234 } // flux from histograms or from JNUBEAM ntuples?
1235
1236 } else {
1237 LOG("gevgen_t2k", pFATAL) << "No flux info was specified - Exiting";
1238 PrintSyntax();
1239 exit(1);
1240 }
1241
1242 // Use a random offset when looping over flux entries
1243 if(parser.OptionExists('R')) gOptRandomFluxOffset = true;
1244
1245 //
1246 // *** pre-calculated flux interaction probabilities
1247 //
1248
1249 // using pre-calculated flux interaction probabilities
1250 if( parser.OptionExists('P') ){
1251 gOptFluxProbFileName = parser.ArgAsString('P');
1252 if(gOptFluxProbFileName.length() > 0){
1253 gOptUseFluxProbs = true;
1254 bool accessible =
1255 !(gSystem->AccessPathName(gOptFluxProbFileName.c_str()));
1256 if(!accessible){
1257 LOG("gevgen_t2k", pFATAL)
1258 << "Can not access pre-calculated flux probabilities file: " << gOptFluxProbFileName;
1259 PrintSyntax();
1260 exit(1);
1261 }
1262 }
1263 else {
1264 LOG("gevgen_t2k", pFATAL)
1265 << "No flux interaction probabilites were specified - exiting";
1266 PrintSyntax();
1267 exit(1);
1268 }
1269 }
1270
1271 // pre-generating interaction probs and saving to output file
1272 if( parser.OptionExists('S') ){
1273 gOptSaveFluxProbsFile = true;
1275 }
1276
1277 // cannot save and run at the same time
1279 LOG("gevgen_t2k", pFATAL)
1280 << "Cannot specify both the -P and -S options at the same time!";
1281 exit(1);
1282 }
1283
1284 // only makes sense to be setting these options for a realistic flux
1286 LOG("gevgen_t2k", pFATAL)
1287 << "Using pre-calculated flux interaction probabilities only makes "
1288 << "sense when using JNUBEAM flux option!";
1289 exit(1);
1290 }
1291
1292 // flux file POT normalization
1293 // only relevant when using the JNUBEAM flux ntuples
1294 if( parser.OptionExists('p') ) {
1295 LOG("gevgen_t2k", pDEBUG) << "Reading flux file normalization";
1296 gOptFluxNorm = parser.ArgAsDouble('p');
1297 } else {
1298 LOG("gevgen_t2k", pDEBUG)
1299 << "Setting standard normalization for JNUBEAM flux ntuples";
1301 } //-p
1302
1303 // number of times to cycle through the JNUBEAM flux ntuple contents
1304 if( parser.OptionExists('c') ) {
1305 LOG("gevgen_t2k", pDEBUG) << "Reading number of flux ntuple cycles";
1306 gOptFluxNCycles = parser.ArgAsInt('c');
1307 } else {
1308 LOG("gevgen_t2k", pDEBUG)
1309 << "Setting standard number of cycles for JNUBEAM flux ntuples";
1310 gOptFluxNCycles = -1;
1311 } //-c
1312
1313 // limit on max number of events that can be generated
1314 if( parser.OptionExists('n') ) {
1315 LOG("gevgen_t2k", pDEBUG)
1316 << "Reading limit on number of events to generate";
1317 gOptNev = parser.ArgAsInt('n');
1318 } else {
1319 LOG("gevgen_t2k", pDEBUG)
1320 << "Will keep on generating events till the flux driver stops";
1321 gOptNev = -1;
1322 } //-n
1323
1324 // exposure (in POT)
1325 bool uppc_e = parser.OptionExists('E');
1326 char pot_args = 'e';
1327 bool pot_exit = true;
1328 if(uppc_e) {
1329 pot_args = 'E';
1330 pot_exit = false;
1331 }
1333 if( parser.OptionExists(pot_args) ) {
1334 LOG("gevgen_t2k", pDEBUG) << "Reading requested exposure in POT";
1335 gOptPOT = parser.ArgAsDouble(pot_args);
1336 } else {
1337 LOG("gevgen_t2k", pDEBUG) << "No POT exposure was requested";
1338 gOptPOT = -1;
1339 } //-e, -E
1340
1341 // event file prefix
1342 if( parser.OptionExists('o') ) {
1343 LOG("gevgen_t2k", pDEBUG) << "Reading the event filename prefix";
1344 gOptEvFilePrefix = parser.ArgAsString('o');
1345 } else {
1346 LOG("gevgen_t2k", pDEBUG)
1347 << "Will set the default event filename prefix";
1349 } //-o
1350
1351
1352 // random number seed
1353 if( parser.OptionExists("seed") ) {
1354 LOG("gevgen_t2k", pINFO) << "Reading random number seed";
1355 gOptRanSeed = parser.ArgAsLong("seed");
1356 } else {
1357 LOG("gevgen_t2k", pINFO) << "Unspecified random number seed - Using default";
1358 gOptRanSeed = -1;
1359 }
1360
1361 // input cross-section file
1362 if( parser.OptionExists("cross-sections") ) {
1363 LOG("gevgen_t2k", pINFO) << "Reading cross-section file";
1364 gOptInpXSecFile = parser.ArgAsString("cross-sections");
1365 } else {
1366 LOG("gevgen_t2k", pINFO) << "Unspecified cross-section file";
1367 gOptInpXSecFile = "";
1368 }
1369
1370 //
1371 // >>> perform 'sanity' checks on command line arguments
1372 //
1373
1374 // If we use a JNUBEAM flux ntuple, the 'exposure' may be set
1375 // either as:
1376 // - a number of POTs (whichever number of flux ntuple cycles that corresponds to)
1377 // - a number of generated events (whichever number of POTs that corresponds to)
1378 // - a number of flux ntuple cycles (whichever number of POTs that corresponds to)
1379 // Only one of those options can be set.
1380 if(!gOptUsingHistFlux) {
1381 int nset=0;
1382 if(gOptPOT > 0) nset++;
1383 if(gOptFluxNCycles > 0) nset++;
1384 if(gOptNev > 0) nset++;
1385 if(nset==0) {
1386 LOG("gevgen_t2k", pWARN)
1387 << "** To use a JNUBEAM flux ntuple you need to specify an exposure, "
1388 << "either via the -c, -e or -n options";
1389 LOG("gevgen_t2k", pWARN)
1390 << "** gevgen_t2k automatically sets the exposure via '-c 1'";
1391 gOptFluxNCycles = 1;
1392 }
1393 if(nset>1) {
1394 LOG("gevgen_t2k", pFATAL)
1395 << "You can not specify more than one of the -c, -e or -n options";
1396 PrintSyntax();
1397 exit(1);
1398 }
1399 }
1400 // If we use a flux histograms (not JNUBEAM flux ntuples) then -currently- the
1401 // only way to control exposure is via a number of events
1402 if(gOptUsingHistFlux) {
1403 if(gOptNev < 0) {
1404 LOG("gevgen_t2k", pFATAL)
1405 << "If you're using flux from histograms you need to specify the -n option";
1406 PrintSyntax();
1407 exit(1);
1408 }
1409 }
1410 // If we don't use a detailed ROOT detector geometry (but just a target mix) then
1411 // don't accept POT as a way to control job statistics (not enough info is passed
1412 // in the target mix to compute POT & the calculation can be easily done offline)
1413 if(!gOptUsingRootGeom) {
1414 if(gOptPOT > 0) {
1415 LOG("gevgen_t2k", pFATAL)
1416 << "You may not use the -e, -E options "
1417 << "without a detailed detector geometry description input";
1418 exit(1);
1419 }
1420 }
1421
1422 //
1423 // >>> print the command line options
1424 //
1425 PDGLibrary * pdglib = PDGLibrary::Instance();
1426
1427 ostringstream gminfo;
1428 if (gOptUsingRootGeom) {
1429 gminfo << "Using ROOT geometry - file: " << gOptRootGeom
1430 << ", top volume: "
1431 << ((gOptRootGeomTopVol.size()==0) ? "<master volume>" : gOptRootGeomTopVol)
1432 << ", max{PL} file: "
1433 << ((gOptExtMaxPlXml.size()==0) ? "<none>" : gOptExtMaxPlXml)
1434 << ", length units: " << lunits
1435 << ", density units: " << dunits;
1436 } else {
1437 gminfo << "Using target mix - ";
1438 map<int,double>::const_iterator iter;
1439 for(iter = gOptTgtMix.begin(); iter != gOptTgtMix.end(); ++iter) {
1440 int pdg_code = iter->first;
1441 double wgt = iter->second;
1442 TParticlePDG * p = pdglib->Find(pdg_code);
1443 if(p) {
1444 string name = p->GetName();
1445 gminfo << "(" << name << ") -> " << 100*wgt << "% / ";
1446 }//p?
1447 }
1448 }
1449
1450 ostringstream fluxinfo;
1451 if(gOptUsingHistFlux) {
1452 fluxinfo << "Using flux histograms - ";
1453 map<int,TH1D*>::const_iterator iter;
1454 for(iter = gOptFluxHst.begin(); iter != gOptFluxHst.end(); ++iter) {
1455 int pdg_code = iter->first;
1456 TH1D * spectrum = iter->second;
1457 TParticlePDG * p = pdglib->Find(pdg_code);
1458 if(p) {
1459 string name = p->GetName();
1460 fluxinfo << "(" << name << ") -> " << spectrum->GetName() << " / ";
1461 }//p?
1462 }
1463 } else {
1464 fluxinfo << "Using JNUBEAM flux ntuple - "
1465 << "file: " << gOptFluxFile
1466 << ", location: " << gOptDetectorLocation
1467 << ", pot norm: " << gOptFluxNorm;
1468
1469 if( gOptFluxNtpNuList.size() > 0 ) {
1470 fluxinfo << ", ** this job is considering only: ";
1471 for(unsigned int inu = 0; inu < gOptFluxNtpNuList.size(); inu++) {
1472 fluxinfo << gOptFluxNtpNuList[inu];
1473 if(inu < gOptFluxNtpNuList.size()-1) fluxinfo << ",";
1474 }
1475 }
1476
1477 }
1478
1479 ostringstream exposure;
1480 if(gOptPOT > 0)
1481 exposure << "Number of POTs = " << gOptPOT;
1482 if(gOptFluxNCycles > 0)
1483 exposure << "Number of flux cycles = " << gOptFluxNCycles;
1484 if(gOptNev > 0)
1485 exposure << "Number of events = " << gOptNev;
1486
1487 LOG("gevgen_t2k", pNOTICE)
1488 << "\n\n"
1489 << utils::print::PrintFramedMesg("T2K event generation job configuration");
1490
1491 LOG("gevgen_t2k", pNOTICE)
1492 << "\n - Run number: " << gOptRunNu
1493 << "\n - Random number seed: " << gOptRanSeed
1494 << "\n - Using cross-section file: " << gOptInpXSecFile
1495 << "\n - Flux @ " << fluxinfo.str()
1496 << "\n - Geometry @ " << gminfo.str()
1497 << "\n - Exposure @ " << exposure.str();
1498
1499 LOG("gevgen_t2k", pNOTICE) << *RunOpt::Instance();
1500}
1501//____________________________________________________________________________
1502void PrintSyntax(void)
1503{
1504 LOG("gevgen_t2k", pFATAL)
1505 << "\n **Syntax**"
1506 << "\n gevgen_t2k [-h] "
1507 << "\n [-r run#]"
1508 << "\n -f flux"
1509 << "\n -g geometry"
1510 << "\n [-p pot_normalization_of_flux_file]"
1511 << "\n [-t top_volume_name_at_geom]"
1512 << "\n [-P pre_gen_prob_file]"
1513 << "\n [-S] [output_name]"
1514 << "\n [-m max_path_lengths_xml_file]"
1515 << "\n [-L length_units_at_geom]"
1516 << "\n [-D density_units_at_geom]"
1517 << "\n [-n n_of_events]"
1518 << "\n [-c flux_cycles]"
1519 << "\n [-e, -E exposure_in_POTs]"
1520 << "\n [-o output_event_file_prefix]"
1521 << "\n [-R]"
1522 << "\n [--seed random_number_seed]"
1523 << "\n --cross-sections xml_file"
1524 << "\n [--event-generator-list list_name]"
1525 << "\n [--message-thresholds xml_file]"
1526 << "\n [--unphysical-event-mask mask]"
1527 << "\n [--event-record-print-level level]"
1528 << "\n [--mc-job-status-refresh-rate rate]"
1529 << "\n [--cache-file root_file]"
1530 << "\n"
1531 << " Please also read the detailed documentation at http://www.genie-mc.org"
1532 << " or look at the source code: $GENIE/src/Apps/gT2KEvGen.cxx"
1533 << "\n";
1534}
1535//____________________________________________________________________________
#define pNOTICE
Definition Messenger.h:61
#define pINFO
Definition Messenger.h:62
#define pERROR
Definition Messenger.h:59
#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
bool LoadFluxProbabilities(string filename)
double GlobProbScale(void) const
Definition GMCJDriver.h:76
bool PreCalcFluxProbabilities(void)
void ForceSingleProbScale(void)
void UseFluxDriver(GFluxI *flux)
map< int, double > SumFluxIntProbs(void) const
Definition GMCJDriver.h:78
void UseGeomAnalyzer(GeomAnalyzerI *geom)
void Configure(bool calc_prob_scales=true)
void SaveFluxProbabilities(string outfilename)
bool UseMaxPathLengths(string xml_filename)
void UseSplines(bool useLogE=true)
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)
TTree * EventTree(void)
Definition NtpWriter.h:55
A list of PDG codes.
Definition PDGCodeList.h:32
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
void EnableBareXSecPreCalc(bool flag)
Definition RunOpt.h:62
static RunOpt * Instance(void)
Definition RunOpt.cxx:54
A generic GENIE flux driver. Generates a 'cylindrical' neutrino beam along the input direction,...
void SetNuDirection(const TVector3 &direction)
void AddEnergySpectrum(int nu_pdgc, TH1D *spectrum)
void SetBeamSpot(const TVector3 &spot)
A GENIE flux driver encapsulating the JPARC neutrino flux. It reads-in the official JPARC neutrino fl...
void SetFilePOT(double pot)
flux file norm is in /N POT/det [ND] or /N POT/cm^2 [FD]. Specify N (typically 1E+21)
const GJPARCNuFluxPassThroughInfo & PassThroughInfo(void)
GJPARCNuFluxPassThroughInfo.
void SetNumOfCycles(int n)
set how many times to cycle through the ntuple (default: 1 / n=0 means 'infinite')
void SetUpstreamZ(double z0)
set flux neutrino initial z position (upstream of the detector)
bool End(void)
true if no more flux nu's can be thrown (eg reaching end of beam sim ntuples)
void DisableOffset(void)
switch off random offset, must be called before LoadBeamSimData to have any effect
void SetFluxParticles(const PDGCodeList &particles)
specify list of flux neutrino species
double POT_curravg(void)
current average POT
double POT_1cycle(void)
flux POT per cycle
bool LoadBeamSimData(string filename, string det_loc)
load a jnubeam root flux ntuple
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
Utility class to store MC job meta-data.
long int gOptRanSeed
double gOptGeomLUnits
string kDefOptEvFilePrefix
string kDefOptGeomLUnits
NtpMCFormat_t kDefOptNtpFormat
bool gOptUsingRootGeom
Long_t gOptRunNu
string kDefOptGeomDUnits
string gOptExtMaxPlXml
int gOptNev
map< int, double > gOptTgtMix
string gOptRootGeom
double gOptGeomDUnits
string gOptEvFilePrefix
string gOptInpXSecFile
string gOptRootGeomTopVol
string lunits
string geom
string gOptFluxFile
bool gOptUsingHistFlux
double gOptPOT
map< int, TH1D * > gOptFluxHst
string gOptDetectorLocation
PDGCodeList gOptFluxNtpNuList(false)
double kDefOptFluxNorm
double gOptFluxNorm
bool gOptSaveFluxProbsFile
string gOptFluxProbFileName
bool gOptRandomFluxOffset
bool gOptUseFluxProbs
void GetCommandLineArgs(int argc, char **argv)
int gOptFluxNCycles
void PrintSyntax(void)
bool gOptExitAtEndOfFullFluxCycles
string gOptSaveFluxProbsFileName
GENIE flux drivers.
Utilities for improving the code readability when using PDG codes.
bool IsNeutrino(int pdgc)
Definition PDGUtils.cxx:110
bool IsAntiNeutrino(int pdgc)
Definition PDGUtils.cxx:118
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)
double UnitFromString(string u)
Definition UnitUtils.cxx:18
THE MAIN GENIE PROJECT NAMESPACE
Definition AlgCmp.h:25
@ kNFGHEP
Definition NtpMCFormat.h:30
enum genie::ENtpMCFormat NtpMCFormat_t