GENIEGenerator
Loading...
Searching...
No Matches
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-2025, 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
35using std::string;
36using std::map;
37
38namespace genie {
39
40class EventRecord;
41class GFluxI;
42class GeomAnalyzerI;
43class GENIE;
44class GEVGPool;
45
47
48public :
49 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
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
86private:
87
88 // private methods:
89 void InitJob (void);
90 void InitEventGeneration (void);
91 void GetParticleLists (void);
92 void GetMaxPathLengthList (void);
93 void GetMaxFluxEnergy (void);
95 void BootstrapXSecSplines (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_
Generated Event Record. It is a GHepRecord object that can accept / be visited by EventRecordVisitorI...
Definition EventRecord.h:37
A pool of GEVGDriver objects with an initial state key.
Definition GEVGPool.h:37
GENIE Interface for user-defined flux classes.
Definition GFluxI.h:29
bool fUseLogE
[config] build splines = f(logE) (rather than f(E)) ?
Definition GMCJDriver.h:134
bool fPmaxLogBinning
[config] maximum interaction probability is computed in logarithmic energy bins
Definition GMCJDriver.h:124
PathLengthList fMaxPathLengths
[declared by the geom driver] maximum path length list
Definition GMCJDriver.h:116
bool fPreSelect
[config] set whether to pre-select events using max interaction paths
Definition GMCJDriver.h:138
bool LoadFluxProbabilities(string filename)
bool fForceInteraction
[config] force intearction?
Definition GMCJDriver.h:137
double fPmaxSafetyFactor
[config] safety factor to compute the maximum interaction probability
Definition GMCJDriver.h:126
double fEmax
[declared by the flux driver] maximum neutrino energy
Definition GMCJDriver.h:113
void SetXSecSplineNbins(int nbins)
double GlobProbScale(void) const
Definition GMCJDriver.h:76
bool fUseSplines
[config] compute all needed & not-loaded splines at init
Definition GMCJDriver.h:133
string fMaxPlXmlFilename
[config] input file with max density-weighted path lengths for all materials
Definition GMCJDriver.h:131
bool PreCalcFluxProbabilities(void)
void BootstrapXSecSplines(void)
void GenerateVertexPosition(void)
void InitEventGeneration(void)
void GenerateEventKinematics(void)
void ForceSingleProbScale(void)
int fBrFluxIndex
corresponding entry in flux input tree (set to address of branch:"FluxEntry")
Definition GMCJDriver.h:142
void UseFluxDriver(GFluxI *flux)
double fNFluxNeutrinos
[current] number of flux nuetrinos fired by the flux driver so far
Definition GMCJDriver.h:122
string fFluxIntFileName
whether to save pre-generated flux tree for use in later jobs
Definition GMCJDriver.h:146
const GeomAnalyzerI & GeomAnalyzer(void) const
Definition GMCJDriver.h:82
map< int, double > SumFluxIntProbs(void) const
Definition GMCJDriver.h:78
void UseGeomAnalyzer(GeomAnalyzerI *geom)
string fEventGenList
[config] list of event generators loaded by this driver (what used to be the $GEVGL setting)
Definition GMCJDriver.h:129
EventRecord * GenerateEvent1Try(void)
map< int, double > fCurCumulProbMap
[current] cummulative interaction probabilities
Definition GMCJDriver.h:121
bool GenerateFluxNeutrino(void)
double fGlobPmax
[computed at init] global interaction probability scale for given flux & geometry
Definition GMCJDriver.h:128
TBits * fUnphysEventMask
[config] controls whether unphysical events are returned (what used to be the $GUNPHYSMASK setting)
Definition GMCJDriver.h:130
double ComputeInteractionProbabilities(bool use_max_path_length)
int fSelTgtPdg
[current] selected target material PDG code
Definition GMCJDriver.h:120
TTree * fFluxIntTree
[computed-or-loaded] pre-computed flux interaction probabilities (expected tree name is "gFlxIntProbs...
Definition GMCJDriver.h:140
GeomAnalyzerI * GeomAnalyzerPtr(void) const
Definition GMCJDriver.h:84
void SetPmaxSafetyFactor(double sf)
void Configure(bool calc_prob_scales=true)
bool ComputePathLengths(void)
void GetMaxFluxEnergy(void)
GFluxI * fFluxDriver
[input] neutrino flux driver
Definition GMCJDriver.h:111
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
void GetParticleLists(void)
void ForceInteraction(void)
double fBrFluxIntProb
flux interaction probability (set to branch:"FluxIntProb")
Definition GMCJDriver.h:141
void SaveFluxProbabilities(string outfilename)
bool UseMaxPathLengths(string xml_filename)
void PopulateEventGenDriverPool(void)
void UseSplines(bool useLogE=true)
GeomAnalyzerI * fGeomAnalyzer
[input] detector geometry analyzer
Definition GMCJDriver.h:112
bool fGenerateUnweighted
[config] force single probability scale?
Definition GMCJDriver.h:136
int fBrFluxPDG
corresponding flux pdg code (set to address of branch: "FluxPDG")
Definition GMCJDriver.h:145
string fFluxIntTreeName
name for tree holding flux probabilities
Definition GMCJDriver.h:147
int SelectTargetMaterial(double R)
GFluxI * FluxDriverPtr(void) const
Definition GMCJDriver.h:83
TLorentzVector fCurVtx
[current] interaction vertex
Definition GMCJDriver.h:118
PathLengthList fCurPathLengths
[current] path length list for current flux neutrino
Definition GMCJDriver.h:117
map< int, TH1D * > fPmax
[computed at init] interaction probability scale /neutrino /energy for given geometry
Definition GMCJDriver.h:127
long int NFluxNeutrinos(void) const
Definition GMCJDriver.h:77
const GFluxI & FluxDriver(void) const
Definition GMCJDriver.h:81
double fBrFluxWeight
corresponding flux weight (set to address of branch: "FluxWeight")
Definition GMCJDriver.h:144
void ComputeEventProbability(void)
bool fUseExtMaxPl
[config] using external max path length estimate?
Definition GMCJDriver.h:132
void PreSelectEvents(bool preselect=true)
void BootstrapXSecSplineSummation(void)
EventRecord * fCurEvt
[current] generated event
Definition GMCJDriver.h:119
GEVGPool * fGPool
A pool of GEVGDrivers properly configured event generation drivers / one per init state.
Definition GMCJDriver.h:110
void SetUnphysEventMask(const TBits &mask)
void SetEventGeneratorList(string listname)
bool fKeepThrowingFluxNu
[config] keep firing flux neutrinos till one of them interacts
Definition GMCJDriver.h:135
PDGCodeList fNuList
[declared by the flux driver] list of neutrino codes
Definition GMCJDriver.h:114
void SetPmaxNbins(int nbins)
PDGCodeList fTgtList
[declared by the geom driver] list of target codes
Definition GMCJDriver.h:115
int fPmaxNbins
[config] number of bins in energy used in the maximum interaction probability
Definition GMCJDriver.h:125
void ComputeProbScales(void)
int fXSecSplineNbins
[config] number of bins in energy used in the xsec splines
Definition GMCJDriver.h:123
EventRecord * GenerateEvent(void)
double PreGenFluxInteractionProbability(void)
void GetMaxPathLengthList(void)
double fBrFluxEnu
corresponding flux P4 (set to address of branch:"FluxP4")
Definition GMCJDriver.h:143
void SetPmaxLogBinning(void)
TFile * fFluxIntProbFile
[input] pre-generated flux interaction probability file
Definition GMCJDriver.h:139
void KeepOnThrowingFluxNeutrinos(bool keep_on)
double InteractionProbability(double xsec, double pl, int A)
Defines the GENIE Geometry Analyzer Interface.
A list of PDG codes.
Definition PDGCodeList.h:32
Object to be filled with the neutrino path-length, for all detector geometry materials,...
string geom
GENIE flux drivers.
THE MAIN GENIE PROJECT NAMESPACE
Definition AlgCmp.h:25