GENIEGenerator
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GMCJDriver.h
Go to the documentation of this file.
1 //____________________________________________________________________________
2 /*!
3 
4 \class genie::GMCJDriver
5 
6 \brief A GENIE `MC Job Driver'. Can be used for setting up complicated event
7  generation cases involving detailed flux descriptions and detector
8  geometry descriptions.
9 
10 \author Costas Andreopoulos <c.andreopoulos \at cern.ch>
11  University of Liverpool
12 
13 \created May 25, 2005
14 
15 \cpright Copyright (c) 2003-2024, The GENIE Collaboration
16  For the full text of the license visit http://copyright.genie-mc.org
17 */
18 //____________________________________________________________________________
19 
20 #ifndef _GENIE_MC_JOB_DRIVER_H_
21 #define _GENIE_MC_JOB_DRIVER_H_
22 
23 #include <string>
24 #include <map>
25 
26 #include <TH1D.h>
27 #include <TLorentzVector.h>
28 #include <TFile.h>
29 #include <TTree.h>
30 #include <TBits.h>
31 
34 
35 using std::string;
36 using std::map;
37 
38 namespace genie {
39 
40 class EventRecord;
41 class GFluxI;
42 class GeomAnalyzerI;
43 class GENIE;
44 class GEVGPool;
45 
46 class GMCJDriver {
47 
48 public :
49  GMCJDriver();
50  ~GMCJDriver();
51 
52  // configure MC job
53  void SetEventGeneratorList (string listname);
54  void SetUnphysEventMask (const TBits & mask);
55  void UseFluxDriver (GFluxI * flux);
57  void UseSplines (bool useLogE = true);
58  bool UseMaxPathLengths (string xml_filename);
59  void KeepOnThrowingFluxNeutrinos (bool keep_on);
60  void SetXSecSplineNbins (int nbins);
61  void SetPmaxLogBinning (void);
62  void SetPmaxNbins (int nbins);
63  void SetPmaxSafetyFactor (double sf);
64  void ForceInteraction (void);
65  void ForceSingleProbScale (void);
66  void PreSelectEvents (bool preselect = true);
67  bool PreCalcFluxProbabilities (void);
68  bool LoadFluxProbabilities (string filename);
69  void SaveFluxProbabilities (string outfilename);
70  void Configure (bool calc_prob_scales = true);
71 
72  // generate single neutrino event for input flux & geometry
73  EventRecord * GenerateEvent (void);
74 
75  // info needed for computing the generated sample normalization
76  double GlobProbScale (void) const { return fGlobPmax; }
77  long int NFluxNeutrinos (void) const { return (long int) fNFluxNeutrinos; }
78  map<int, double> SumFluxIntProbs(void) const { return fSumFluxIntProbs; }
79 
80  // input flux and geometry drivers
81  const GFluxI & FluxDriver (void) const { return *fFluxDriver; }
82  const GeomAnalyzerI & GeomAnalyzer (void) const { return *fGeomAnalyzer; }
83  GFluxI * FluxDriverPtr (void) const { return fFluxDriver; }
84  GeomAnalyzerI * GeomAnalyzerPtr (void) const { return fGeomAnalyzer; }
85 
86 private:
87 
88  // private methods:
89  void InitJob (void);
90  void InitEventGeneration (void);
91  void GetParticleLists (void);
92  void GetMaxPathLengthList (void);
93  void GetMaxFluxEnergy (void);
94  void PopulateEventGenDriverPool (void);
95  void BootstrapXSecSplines (void);
96  void BootstrapXSecSplineSummation (void);
97  void ComputeProbScales (void);
99  bool GenerateFluxNeutrino (void);
100  bool ComputePathLengths (void);
101  double ComputeInteractionProbabilities (bool use_max_path_length);
102  int SelectTargetMaterial (double R);
103  void GenerateEventKinematics (void);
104  void GenerateVertexPosition (void);
105  void ComputeEventProbability (void);
106  double InteractionProbability (double xsec, double pl, int A);
108 
109  // private data members:
110  GEVGPool * fGPool; ///< A pool of GEVGDrivers properly configured event generation drivers / one per init state
111  GFluxI * fFluxDriver; ///< [input] neutrino flux driver
112  GeomAnalyzerI * fGeomAnalyzer; ///< [input] detector geometry analyzer
113  double fEmax; ///< [declared by the flux driver] maximum neutrino energy
114  PDGCodeList fNuList; ///< [declared by the flux driver] list of neutrino codes
115  PDGCodeList fTgtList; ///< [declared by the geom driver] list of target codes
116  PathLengthList fMaxPathLengths; ///< [declared by the geom driver] maximum path length list
117  PathLengthList fCurPathLengths; ///< [current] path length list for current flux neutrino
118  TLorentzVector fCurVtx; ///< [current] interaction vertex
119  EventRecord * fCurEvt; ///< [current] generated event
120  int fSelTgtPdg; ///< [current] selected target material PDG code
121  map<int,double> fCurCumulProbMap; ///< [current] cummulative interaction probabilities
122  double fNFluxNeutrinos; ///< [current] number of flux nuetrinos fired by the flux driver so far
123  int fXSecSplineNbins; ///< [config] number of bins in energy used in the xsec splines
124  bool fPmaxLogBinning; ///< [config] maximum interaction probability is computed in logarithmic energy bins
125  int fPmaxNbins; ///< [config] number of bins in energy used in the maximum interaction probability
126  double fPmaxSafetyFactor; ///< [config] safety factor to compute the maximum interaction probability
127  map<int,TH1D*> fPmax; ///< [computed at init] interaction probability scale /neutrino /energy for given geometry
128  double fGlobPmax; ///< [computed at init] global interaction probability scale for given flux & geometry
129  string fEventGenList; ///< [config] list of event generators loaded by this driver (what used to be the $GEVGL setting)
130  TBits * fUnphysEventMask; ///< [config] controls whether unphysical events are returned (what used to be the $GUNPHYSMASK setting)
131  string fMaxPlXmlFilename; ///< [config] input file with max density-weighted path lengths for all materials
132  bool fUseExtMaxPl; ///< [config] using external max path length estimate?
133  bool fUseSplines; ///< [config] compute all needed & not-loaded splines at init
134  bool fUseLogE; ///< [config] build splines = f(logE) (rather than f(E)) ?
135  bool fKeepThrowingFluxNu; ///< [config] keep firing flux neutrinos till one of them interacts
136  bool fGenerateUnweighted; ///< [config] force single probability scale?
137  bool fForceInteraction; ///< [config] force intearction?
138  bool fPreSelect; ///< [config] set whether to pre-select events using max interaction paths
139  TFile* fFluxIntProbFile; ///< [input] pre-generated flux interaction probability file
140  TTree* fFluxIntTree; ///< [computed-or-loaded] pre-computed flux interaction probabilities (expected tree name is "gFlxIntProbs")
141  double fBrFluxIntProb; ///< flux interaction probability (set to branch:"FluxIntProb")
142  int fBrFluxIndex; ///< corresponding entry in flux input tree (set to address of branch:"FluxEntry")
143  double fBrFluxEnu; ///< corresponding flux P4 (set to address of branch:"FluxP4")
144  double fBrFluxWeight; ///< corresponding flux weight (set to address of branch: "FluxWeight")
145  int fBrFluxPDG; ///< corresponding flux pdg code (set to address of branch: "FluxPDG")
146  string fFluxIntFileName; ///< whether to save pre-generated flux tree for use in later jobs
147  string fFluxIntTreeName; ///< name for tree holding flux probabilities
148  map<int, double> fSumFluxIntProbs; ///< map where the key is flux pdg code and the value is sum of fBrFluxWeight * fBrFluxIntProb for all these flux neutrinos
149 };
150 
151 } // genie namespace
152 #endif // _GENIE_MC_JOB_DRIVER_H_
TTree * fFluxIntTree
[computed-or-loaded] pre-computed flux interaction probabilities (expected tree name is &quot;gFlxIntProbs...
Definition: GMCJDriver.h:140
TBits * fUnphysEventMask
[config] controls whether unphysical events are returned (what used to be the $GUNPHYSMASK setting) ...
Definition: GMCJDriver.h:130
void SetUnphysEventMask(const TBits &mask)
Definition: GMCJDriver.cxx:74
void ForceInteraction(void)
Definition: GMCJDriver.cxx:162
GFluxI * FluxDriverPtr(void) const
Definition: GMCJDriver.h:83
void PopulateEventGenDriverPool(void)
Definition: GMCJDriver.cxx:564
bool PreCalcFluxProbabilities(void)
Definition: GMCJDriver.cxx:192
void InitJob(void)
Definition: GMCJDriver.cxx:448
map< int, double > SumFluxIntProbs(void) const
Definition: GMCJDriver.h:78
GEVGPool * fGPool
A pool of GEVGDrivers properly configured event generation drivers / one per init state...
Definition: GMCJDriver.h:110
GeomAnalyzerI * GeomAnalyzerPtr(void) const
Definition: GMCJDriver.h:84
void GetMaxPathLengthList(void)
Definition: GMCJDriver.cxx:536
void SetEventGeneratorList(string listname)
Definition: GMCJDriver.cxx:66
map< int, double > fCurCumulProbMap
[current] cummulative interaction probabilities
Definition: GMCJDriver.h:121
void InitEventGeneration(void)
Definition: GMCJDriver.cxx:807
double fBrFluxIntProb
flux interaction probability (set to branch:&quot;FluxIntProb&quot;)
Definition: GMCJDriver.h:141
PDGCodeList fTgtList
[declared by the geom driver] list of target codes
Definition: GMCJDriver.h:115
bool fPreSelect
[config] set whether to pre-select events using max interaction paths
Definition: GMCJDriver.h:138
bool GenerateFluxNeutrino(void)
bool ComputePathLengths(void)
GFluxI * fFluxDriver
[input] neutrino flux driver
Definition: GMCJDriver.h:111
const GFluxI & FluxDriver(void) const
Definition: GMCJDriver.h:81
bool fForceInteraction
[config] force intearction?
Definition: GMCJDriver.h:137
void ComputeEventProbability(void)
double fBrFluxEnu
corresponding flux P4 (set to address of branch:&quot;FluxP4&quot;)
Definition: GMCJDriver.h:143
bool fUseExtMaxPl
[config] using external max path length estimate?
Definition: GMCJDriver.h:132
A list of PDG codes.
Definition: PDGCodeList.h:32
double fBrFluxWeight
corresponding flux weight (set to address of branch: &quot;FluxWeight&quot;)
Definition: GMCJDriver.h:144
void UseFluxDriver(GFluxI *flux)
Definition: GMCJDriver.cxx:83
string fFluxIntFileName
whether to save pre-generated flux tree for use in later jobs
Definition: GMCJDriver.h:146
PathLengthList fCurPathLengths
[current] path length list for current flux neutrino
Definition: GMCJDriver.h:117
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.
void KeepOnThrowingFluxNeutrinos(bool keep_on)
Definition: GMCJDriver.cxx:120
void GetMaxFluxEnergy(void)
Definition: GMCJDriver.cxx:554
bool fGenerateUnweighted
[config] force single probability scale?
Definition: GMCJDriver.h:136
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
map< int, TH1D * > fPmax
[computed at init] interaction probability scale /neutrino /energy for given geometry ...
Definition: GMCJDriver.h:127
EventRecord * fCurEvt
[current] generated event
Definition: GMCJDriver.h:119
TFile * fFluxIntProbFile
[input] pre-generated flux interaction probability file
Definition: GMCJDriver.h:139
bool fKeepThrowingFluxNu
[config] keep firing flux neutrinos till one of them interacts
Definition: GMCJDriver.h:135
int fSelTgtPdg
[current] selected target material PDG code
Definition: GMCJDriver.h:120
static constexpr double A
Definition: Units.h:74
void SetPmaxSafetyFactor(double sf)
Definition: GMCJDriver.cxx:154
void ForceSingleProbScale(void)
Definition: GMCJDriver.cxx:172
double fNFluxNeutrinos
[current] number of flux nuetrinos fired by the flux driver so far
Definition: GMCJDriver.h:122
double PreGenFluxInteractionProbability(void)
int fBrFluxPDG
corresponding flux pdg code (set to address of branch: &quot;FluxPDG&quot;)
Definition: GMCJDriver.h:145
bool fUseLogE
[config] build splines = f(logE) (rather than f(E)) ?
Definition: GMCJDriver.h:134
bool fUseSplines
[config] compute all needed &amp; not-loaded splines at init
Definition: GMCJDriver.h:133
string geom
bool UseMaxPathLengths(string xml_filename)
Definition: GMCJDriver.cxx:99
void Configure(bool calc_prob_scales=true)
Definition: GMCJDriver.cxx:399
bool fPmaxLogBinning
[config] maximum interaction probability is computed in logarithmic energy bins
Definition: GMCJDriver.h:124
int fBrFluxIndex
corresponding entry in flux input tree (set to address of branch:&quot;FluxEntry&quot;)
Definition: GMCJDriver.h:142
string fFluxIntTreeName
name for tree holding flux probabilities
Definition: GMCJDriver.h:147
void SetPmaxLogBinning(void)
Definition: GMCJDriver.cxx:137
string fEventGenList
[config] list of event generators loaded by this driver (what used to be the $GEVGL setting) ...
Definition: GMCJDriver.h:129
const GeomAnalyzerI & GeomAnalyzer(void) const
Definition: GMCJDriver.h:82
void BootstrapXSecSplines(void)
Definition: GMCJDriver.cxx:601
double ComputeInteractionProbabilities(bool use_max_path_length)
double fGlobPmax
[computed at init] global interaction probability scale for given flux &amp; geometry ...
Definition: GMCJDriver.h:128
TLorentzVector fCurVtx
[current] interaction vertex
Definition: GMCJDriver.h:118
EventRecord * GenerateEvent(void)
Definition: GMCJDriver.cxx:815
long int NFluxNeutrinos(void) const
Definition: GMCJDriver.h:77
string fMaxPlXmlFilename
[config] input file with max density-weighted path lengths for all materials
Definition: GMCJDriver.h:131
Generated Event Record. It is a GHepRecord object that can accept / be visited by EventRecordVisitorI...
Definition: EventRecord.h:37
int fXSecSplineNbins
[config] number of bins in energy used in the xsec splines
Definition: GMCJDriver.h:123
double InteractionProbability(double xsec, double pl, int A)
void SetPmaxNbins(int nbins)
Definition: GMCJDriver.cxx:145
void PreSelectEvents(bool preselect=true)
Definition: GMCJDriver.cxx:184
void GenerateEventKinematics(void)
map< int, double > fSumFluxIntProbs
map where the key is flux pdg code and the value is sum of fBrFluxWeight * fBrFluxIntProb for all the...
Definition: GMCJDriver.h:148
int SelectTargetMaterial(double R)
void GenerateVertexPosition(void)
void SetXSecSplineNbins(int nbins)
Definition: GMCJDriver.cxx:128
GeomAnalyzerI * fGeomAnalyzer
[input] detector geometry analyzer
Definition: GMCJDriver.h:112
void ComputeProbScales(void)
Definition: GMCJDriver.cxx:670
bool LoadFluxProbabilities(string filename)
Definition: GMCJDriver.cxx:343
EventRecord * GenerateEvent1Try(void)
Definition: GMCJDriver.cxx:844
void BootstrapXSecSplineSummation(void)
Definition: GMCJDriver.cxx:628
double fPmaxSafetyFactor
[config] safety factor to compute the maximum interaction probability
Definition: GMCJDriver.h:126
Defines the GENIE Geometry Analyzer Interface.
Definition: GeomAnalyzerI.h:29
PathLengthList fMaxPathLengths
[declared by the geom driver] maximum path length list
Definition: GMCJDriver.h:116
int fPmaxNbins
[config] number of bins in energy used in the maximum interaction probability
Definition: GMCJDriver.h:125
PDGCodeList fNuList
[declared by the flux driver] list of neutrino codes
Definition: GMCJDriver.h:114
void UseGeomAnalyzer(GeomAnalyzerI *geom)
Definition: GMCJDriver.cxx:88
void UseSplines(bool useLogE=true)
Definition: GMCJDriver.cxx:93
void SaveFluxProbabilities(string outfilename)
Definition: GMCJDriver.cxx:390
void GetParticleLists(void)
Definition: GMCJDriver.cxx:519
GENIE Interface for user-defined flux classes.
Definition: GFluxI.h:29
A pool of GEVGDriver objects with an initial state key.
Definition: GEVGPool.h:37
double fEmax
[declared by the flux driver] maximum neutrino energy
Definition: GMCJDriver.h:113