GENIEGenerator
Loading...
Searching...
No Matches
GSimpleNtpFlux.h
Go to the documentation of this file.
1//____________________________________________________________________________
2/*!
3
4\class genie::flux::GSimpleNtpFlux
5
6\brief A GENIE flux driver using a simple ntuple format
7
8\author Robert Hatcher <rhatcher \at fnal.gov>
9 Fermi National Accelerator Laboratory
10
11\created Jan 25, 2010
12
13\cpright Copyright (c) 2003-2025, The GENIE Collaboration
14 For the full text of the license visit http://copyright.genie-mc.org
15*/
16//____________________________________________________________________________
17
18#ifndef _SIMPLE_NTP_FLUX_H_
19#define _SIMPLE_NTP_FLUX_H_
20
21#include <string>
22#include <iostream>
23#include <vector>
24#include <set>
25
26#include <TVector3.h>
27#include <TLorentzVector.h>
28#include <TLorentzRotation.h>
29
34
35class TFile;
36class TChain;
37class TTree;
38class TBranch;
39
40using std::string;
41using std::ostream;
42
43namespace genie {
44namespace flux {
45
46class GSimpleNtpEntry;
47ostream & operator << (ostream & stream, const GSimpleNtpEntry & info);
48
49/// Small persistable C-struct -like classes that makes up the SimpleNtpFlux
50/// ntuple. This is only valid for a particular flux window (no reweighting,
51/// no coordinate transformation available).
52///
53/// Order elements from largest to smallest for ROOT alignment purposes
54
55/// GSimpleNtpEntry
56/// =========================
57/// This is the only required branch ("entry") of the "flux" tree
59 public:
61 /* allow default copy constructor ... for now nothing special
62 GSimpleNtpEntry(const GSimpleNtpEntry & info);
63 */
64 virtual ~GSimpleNtpEntry() { };
65 void Reset();
66 void Print(const Option_t* opt = "") const;
67 friend ostream & operator << (ostream & stream, const GSimpleNtpEntry & info);
68
69 Double_t wgt; ///< nu weight
70
71 Double_t vtxx; ///< x position in lab frame (meters)
72 Double_t vtxy; ///< y position in lab frame
73 Double_t vtxz; ///< z position in lab frame
74 Double_t vtxt; ///< time of ray start (seconds)
75 Double_t dist; ///< distance from hadron decay
76
77 Double_t px; ///< x momentum in lab frame (GeV)
78 Double_t py; ///< y momentum in lab frame
79 Double_t pz; ///< z momentum in lab frame
80 Double_t E; ///< energy in lab frame
81
82 Int_t pdg; ///< nu pdg-code
83 UInt_t metakey; ///< key to meta data
84
85 ClassDef(GSimpleNtpEntry,2)
86 };
87
88
89 class GSimpleNtpNuMI;
90 ostream & operator << (ostream & stream, const GSimpleNtpNuMI & info);
91
92/// GSimpleNtpNuMI
93/// =========================
94/// Additional elements for NuMI (allow SKZP reweighting and reference
95/// back to original GNuMI flux entries) as "numi" branch
97 public:
99 virtual ~GSimpleNtpNuMI() { };
100 void Reset();
101 void Print(const Option_t* opt = "") const;
102 friend ostream & operator << (ostream & stream, const GSimpleNtpNuMI & info);
103 Double_t tpx; ///< parent particle px at target exit
104 Double_t tpy;
105 Double_t tpz;
106 Double_t vx; ///< vertex position of hadron/muon decay
107 Double_t vy;
108 Double_t vz;
109 Double_t pdpx; ///< nu parent momentum at time of decay
110 Double_t pdpy;
111 Double_t pdpz;
112 Double_t pppx; ///< nu parent momentum at production point
113 Double_t pppy;
114 Double_t pppz;
115
116 Int_t ndecay; ///< decay mode
117 Int_t ptype; ///< parent type (PDG)
118 Int_t ppmedium; ///< tracking medium where parent was produced
119 Int_t tptype; ///< parent particle type at target exit
120
121 Int_t run; ///<
122 Int_t evtno; ///<
123 Int_t entryno; ///<
124
125 ClassDef(GSimpleNtpNuMI,3)
126 };
127
128
129 class GSimpleNtpAux;
130 ostream & operator << (ostream & stream, const GSimpleNtpAux & info);
131
132/// GSimpleNtpAux
133/// =========================
134/// Additional elements for expansion as "aux" branch
136 public:
138 virtual ~GSimpleNtpAux() { };
139 void Reset();
140 void Print(const Option_t* opt = "") const;
141 friend ostream & operator << (ostream & stream, const GSimpleNtpAux & info);
142
143 std::vector<Int_t> auxint; ///< additional ints associated w/ entry
144 std::vector<Double_t> auxdbl; ///< additional doubles associated w/ entry
145
146 ClassDef(GSimpleNtpAux,1)
147 };
148
149
150 class GSimpleNtpMeta;
151 ostream & operator << (ostream & stream, const GSimpleNtpMeta & info);
152
153/// GSimpleNtpMeta
154/// =========================
155/// A small persistable C-struct -like class that holds metadata
156/// about the the SimpleNtpFlux ntple.
157///
158 class GSimpleNtpMeta: public TObject {
159 public:
161 /* allow default copy constructor ... for now nothing special
162 GSimpleNtpMeta(const GSimpleNtpMeta & info);
163 */
164 virtual ~GSimpleNtpMeta();
165
166 void Reset();
167 void AddFlavor(Int_t nupdg);
168 void Print(const Option_t* opt = "") const;
169 friend ostream & operator << (ostream & stream, const GSimpleNtpMeta & info);
170
171 std::vector<Int_t> pdglist; ///< list of neutrino flavors
172
173 Double_t maxEnergy; ///< maximum energy
174 Double_t minWgt; ///< minimum weight
175 Double_t maxWgt; ///< maximum weight
176 Double_t protons; ///< represented number of protons-on-target
177
178 Double_t windowBase[3]; ///< x,y,z position of window base point
179 Double_t windowDir1[3]; ///< dx,dy,dz of window direction 1
180 Double_t windowDir2[3]; ///< dx,dy,dz of window direction 2
181
182 std::vector<std::string> auxintname; ///< tagname of aux ints associated w/ entry
183 std::vector<std::string> auxdblname; ///< tagname of aux doubles associated w/ entry
184 std::vector<std::string> infiles; ///< list of input files
185
186 Int_t seed; ///< random seed used in generation
187 UInt_t metakey; ///< index key to tie to individual entries
188
189 static UInt_t mxfileprint; ///< allow user to limit # of files to print
190
191 ClassDef(GSimpleNtpMeta,1)
192 };
193
194
195/// GSimpleNtpFlux:
196/// ==========
197/// An implementation of the GFluxI interface that provides NuMI flux
198///
200 : public genie::GFluxI
203{
204
205public :
208
209 // Methods implementing the GENIE GFluxI interface, required for integrating
210 // the NuMI neutrino flux simulations with the GENIE event generation drivers
211
212 const PDGCodeList & FluxParticles (void) { return *fPdgCList; }
213 double MaxEnergy (void) { return fMaxEv; }
214 bool GenerateNext (void);
215 int PdgCode (void) { return fCurEntry->pdg; }
216 double Weight (void) { return fWeight; }
217 const TLorentzVector & Momentum (void) { return fP4; }
218 const TLorentzVector & Position (void) { return fX4; }
219 bool End (void) { return fEnd; }
220 long int Index (void) { return fIEntry; }
221 void Clear (Option_t * opt);
222 void GenerateWeighted (bool gen_weighted);
223
224 // Methods specific to the NuMI flux driver,
225 // for configuration/initialization of the flux & event generation drivers
226 // and and for passing-through flux information (e.g. neutrino parent decay
227 // kinematics) not used by the generator but required by analyses/processing
228 // further downstream
229
230 //
231 // information about or actions on current entry
232 //
234 GetCurrentEntry(void) { return fCurEntry; } ///< GSimpleNtpEntry
236 GetCurrentNuMI(void) { return fCurNuMI; } ///< GSimpleNtpNuMI
238 GetCurrentAux(void) { return fCurAux; } ///< GSimpleNtpAux
240 GetCurrentMeta(void) { return fCurMeta; } ///< GSimpleNtpMeta
241
242 // allow access to main tree so we can call Branch() to retrieve extra stuff
243 TChain*
244 GetFluxTChain(void) { return fNuFluxTree; } ///<
245
246 double GetDecayDist() const; ///< dist (user units) from dk to current pos
247 void MoveToZ0(double z0); ///< move ray origin to user coord Z0
248
249 void SetIncludeVtxt(bool it=true) { fIncludeVtxt=it; }; ///< should X4 include CurEntry.vtxt
250 bool GetIncludeVtxt() { return fIncludeVtxt; }
251
252 //
253 // information about the current state
254 //
255 virtual double GetTotalExposure() const; ///< GFluxExposureI interface
256 virtual long int NFluxNeutrinos() const; ///< # of rays generated
257
258 double UsedPOTs(void) const; ///< # of protons-on-target used
259
260 long int NEntriesUsed(void) const { return fNEntriesUsed; } ///< number of entries read from files
261 double SumWeight(void) const { return fSumWeight; } ///< integrated weight for flux neutrinos looped so far
262
263 void PrintCurrent(void); ///< print current entry from leaves
264 void PrintConfig(); ///< print the current configuration
265
266 std::vector<std::string> GetFileList(); ///< list of files currently part of chain
267
268 //
269 // GFluxFileConfigI interface
270 //
271 virtual void LoadBeamSimData(const std::vector<string>& filenames,
272 const std::string& det_loc);
273 using GFluxFileConfigI::LoadBeamSimData; // inherit the rest
274 virtual void GetBranchInfo(std::vector<std::string>& branchNames,
275 std::vector<std::string>& branchClassNames,
276 std::vector<void**>& branchObjPointers);
277 virtual TTree* GetMetaDataTree();
278
279 //
280 // configuration of GSimpleNtpFlux
281 //
282
283 void SetRequestedBranchList(string blist="entry,numi,aux") { fNuFluxBranchRequest = blist; }
284
285 void SetMaxEnergy(double Ev); ///< specify maximum flx neutrino energy
286
287 void SetGenWeighted(bool genwgt=false) { fGenWeighted = genwgt; } ///< toggle whether GenerateNext() returns weight=1 flux (initial default false)
288
289 void SetEntryReuse(long int nuse=1); ///< # of times to use entry before moving to next
290
291 void ProcessMeta(void); ///< scan for max flux energy, weight
292
293 void GetFluxWindow(TVector3& p1, TVector3& p2, TVector3& p3) const; ///< 3 points define a plane in beam coordinate
294
295private:
296
297 // Private methods
298 //
299 bool GenerateNext_weighted (void);
300 void Initialize (void);
301 void SetDefaults (void);
302 void CleanUp (void);
303 void ResetCurrent (void);
304 void AddFile (TTree* fluxtree, TTree* metatree, string fname);
305 bool OptionalAttachBranch (std::string bname);
306 void CalcEffPOTsPerNu (void);
307 void ScanMeta (void);
308
309 // Private data members
310 //
311 double fMaxEv; ///< maximum energy
312 bool fEnd; ///< end condition reached
313
314 std::vector<string> fNuFluxFilePatterns; ///< (potentially wildcarded) path(s)
315 string fNuFluxBranchRequest; ///< list of requested branches "entry,numi,au"
316 TChain* fNuFluxTree; ///< TTree // REF ONLY
317 TChain* fNuMetaTree; ///< TTree // REF ONLY
318
319 int fNFiles; ///< number of files in chain
320 Long64_t fNEntries; ///< number of flux ntuple entries
321 Long64_t fIEntry; ///< current flux ntuple entry
322 Int_t fIFileNumber; ///< which file for the current entry
323
324 Double_t fFilePOTs; ///< # of protons-on-target represented by all files
325
326 double fWeight; ///< current neutrino weight
327 double fMaxWeight; ///< max flux neutrino weight in input file
328
329 long int fNUse; ///< how often to use same entry in a row
330 long int fIUse; ///< current # of times an entry has been used
331 double fSumWeight; ///< sum of weights for nus thrown so far
332 long int fNNeutrinos; ///< number of flux neutrinos thrown so far
333 long int fNEntriesUsed; ///< number of entries read from files
334 double fEffPOTsPerNu; ///< what a entry is worth ...
335 double fAccumPOTs; ///< POTs used so far
336
337 bool fGenWeighted; ///< does GenerateNext() give weights?
338 bool fAlreadyUnwgt; ///< are input files already unweighted
339 // i.e. are all entry "wgt" values = 1
340 bool fAllFilesMeta; ///< do all files in chain have meta data
341
342 GSimpleNtpEntry* fCurEntry; ///< current entry
343 GSimpleNtpNuMI* fCurNuMI; ///< current "numi" branch extra info
344 GSimpleNtpAux* fCurAux; ///< current "aux" branch extra info
345 TLorentzVector fP4; ///< reconstituted p4 vector
346 TLorentzVector fX4; ///< reconstituted position vector
347 GSimpleNtpMeta* fCurMeta; ///< current meta data
348
349 GSimpleNtpEntry* fCurEntryCopy; ///< current entry
350 GSimpleNtpNuMI* fCurNuMICopy; ///< current "numi" branch extra info
351 GSimpleNtpAux* fCurAuxCopy; ///< current "aux" branch extra info
352
353 bool fIncludeVtxt; ///< does fX4 include CurEntry.vtxt or 0
354
355};
356
357} // flux namespace
358} // genie namespace
359
360#endif // _SIMPLE_NTP_FLUX_H_
GENIE Interface for user-defined flux classes.
Definition GFluxI.h:29
A list of PDG codes.
Definition PDGCodeList.h:32
GENIE interface for uniform flux exposure iterface.
PDGCodeList * fPdgCList
list of neutrino pdg-codes to generate
virtual void LoadBeamSimData(const std::vector< std::string > &filenames, const std::string &det_loc)=0
void Print(const Option_t *opt="") const
std::vector< Int_t > auxint
additional ints associated w/ entry
std::vector< Double_t > auxdbl
additional doubles associated w/ entry
friend ostream & operator<<(ostream &stream, const GSimpleNtpAux &info)
Double_t vtxx
x position in lab frame (meters)
Double_t E
energy in lab frame
Double_t dist
distance from hadron decay
Double_t px
x momentum in lab frame (GeV)
Double_t vtxt
time of ray start (seconds)
Double_t vtxz
z position in lab frame
Double_t py
y momentum in lab frame
Double_t vtxy
y position in lab frame
UInt_t metakey
key to meta data
friend ostream & operator<<(ostream &stream, const GSimpleNtpEntry &info)
Double_t pz
z momentum in lab frame
void Print(const Option_t *opt="") const
bool fIncludeVtxt
does fX4 include CurEntry.vtxt or 0
bool fEnd
end condition reached
const genie::flux::GSimpleNtpEntry * GetCurrentEntry(void)
GSimpleNtpEntry.
const TLorentzVector & Momentum(void)
returns the flux neutrino 4-momentum
void Clear(Option_t *opt)
reset state variables based on opt
void SetEntryReuse(long int nuse=1)
long int fNNeutrinos
number of flux neutrinos thrown so far
double GetDecayDist() const
dist (user units) from dk to current pos
GSimpleNtpNuMI * fCurNuMI
current "numi" branch extra info
const genie::flux::GSimpleNtpMeta * GetCurrentMeta(void)
GSimpleNtpMeta.
void GenerateWeighted(bool gen_weighted)
set whether to generate weighted or unweighted neutrinos
void SetGenWeighted(bool genwgt=false)
toggle whether GenerateNext() returns weight=1 flux (initial default false)
void SetIncludeVtxt(bool it=true)
bool GetIncludeVtxt()
should X4 include CurEntry.vtxt
virtual double GetTotalExposure() const
GFluxExposureI interface.
bool fGenWeighted
does GenerateNext() give weights?
double fMaxWeight
max flux neutrino weight in input file
double fEffPOTsPerNu
what a entry is worth ...
void GetFluxWindow(TVector3 &p1, TVector3 &p2, TVector3 &p3) const
3 points define a plane in beam coordinate
GSimpleNtpEntry * fCurEntryCopy
current entry
bool fAlreadyUnwgt
are input files already unweighted
bool OptionalAttachBranch(std::string bname)
virtual void GetBranchInfo(std::vector< std::string > &branchNames, std::vector< std::string > &branchClassNames, std::vector< void ** > &branchObjPointers)
void SetRequestedBranchList(string blist="entry,numi,aux")
TChain * fNuFluxTree
TTree // REF ONLY.
int PdgCode(void)
returns the flux neutrino pdg code
std::vector< std::string > GetFileList()
list of files currently part of chain
long int fNUse
how often to use same entry in a row
bool End(void)
true if no more flux nu's can be thrown (eg reaching end of beam sim ntuples)
GSimpleNtpAux * fCurAuxCopy
current "aux" branch extra info
double SumWeight(void) const
integrated weight for flux neutrinos looped so far
void AddFile(TTree *fluxtree, TTree *metatree, string fname)
TLorentzVector fX4
reconstituted position vector
double fMaxEv
maximum energy
virtual void LoadBeamSimData(const std::vector< string > &filenames, const std::string &det_loc)
double Weight(void)
returns the flux neutrino weight (if any)
long int NEntriesUsed(void) const
number of entries read from files
double MaxEnergy(void)
declare the max flux neutrino energy that can be generated (for init. purposes)
TLorentzVector fP4
reconstituted p4 vector
long int fNEntriesUsed
number of entries read from files
const genie::flux::GSimpleNtpNuMI * GetCurrentNuMI(void)
GSimpleNtpNuMI.
const genie::flux::GSimpleNtpAux * GetCurrentAux(void)
GSimpleNtpAux.
void ProcessMeta(void)
scan for max flux energy, weight
GSimpleNtpMeta * fCurMeta
current meta data
Int_t fIFileNumber
which file for the current entry
long int Index(void)
returns corresponding index for current flux neutrino (e.g. for a flux ntuple returns the current ent...
virtual long int NFluxNeutrinos() const
long int fIUse
current # of times an entry has been used
std::vector< string > fNuFluxFilePatterns
(potentially wildcarded) path(s)
Long64_t fIEntry
current flux ntuple entry
GSimpleNtpNuMI * fCurNuMICopy
current "numi" branch extra info
GSimpleNtpAux * fCurAux
current "aux" branch extra info
virtual TTree * GetMetaDataTree()
void PrintConfig()
print the current configuration
bool fAllFilesMeta
do all files in chain have meta data
GSimpleNtpEntry * fCurEntry
current entry
void MoveToZ0(double z0)
move ray origin to user coord Z0
double fAccumPOTs
POTs used so far.
const TLorentzVector & Position(void)
returns the flux neutrino 4-position (note: expect SI rather than physical units)
void PrintCurrent(void)
print current entry from leaves
bool GenerateNext(void)
generate the next flux neutrino (return false in err)
double fWeight
current neutrino weight
int fNFiles
number of files in chain
double fSumWeight
sum of weights for nus thrown so far
TChain * fNuMetaTree
TTree // REF ONLY.
string fNuFluxBranchRequest
list of requested branches "entry,numi,au"
const PDGCodeList & FluxParticles(void)
declare list of flux neutrinos that can be generated (for init. purposes)
void SetMaxEnergy(double Ev)
specify maximum flx neutrino energy
Long64_t fNEntries
number of flux ntuple entries
Double_t windowBase[3]
x,y,z position of window base point
Double_t windowDir1[3]
dx,dy,dz of window direction 1
friend ostream & operator<<(ostream &stream, const GSimpleNtpMeta &info)
Double_t windowDir2[3]
dx,dy,dz of window direction 2
Double_t maxWgt
maximum weight
std::vector< std::string > infiles
list of input files
static UInt_t mxfileprint
allow user to limit # of files to print
Double_t minWgt
minimum weight
std::vector< Int_t > pdglist
list of neutrino flavors
UInt_t metakey
index key to tie to individual entries
Int_t seed
random seed used in generation
std::vector< std::string > auxdblname
tagname of aux doubles associated w/ entry
Double_t protons
represented number of protons-on-target
std::vector< std::string > auxintname
tagname of aux ints associated w/ entry
Double_t maxEnergy
maximum energy
void Print(const Option_t *opt="") const
Int_t ppmedium
tracking medium where parent was produced
Int_t ptype
parent type (PDG)
Double_t vx
vertex position of hadron/muon decay
void Print(const Option_t *opt="") const
friend ostream & operator<<(ostream &stream, const GSimpleNtpNuMI &info)
Double_t pdpx
nu parent momentum at time of decay
Int_t tptype
parent particle type at target exit
Double_t pppx
nu parent momentum at production point
Double_t tpx
parent particle px at target exit
GENIE flux drivers.
ostream & operator<<(ostream &stream, const genie::flux::GJPARCNuFluxPassThroughInfo &info)
THE MAIN GENIE PROJECT NAMESPACE
Definition AlgCmp.h:25