GENIEGenerator
Loading...
Searching...
No Matches
GNuMIFlux.h
Go to the documentation of this file.
1//____________________________________________________________________________
2/*!
3
4\class genie::flux::GNuMIFlux
5
6\brief A GENIE flux driver encapsulating the NuMI neutrino flux.
7 It reads-in the official GNUMI neutrino flux ntuples.
8 Supports both geant3 and geant4 formats.
9
10\ref http://www.hep.utexas.edu/~zarko/wwwgnumi/v19/
11
12\author Robert Hatcher <rhatcher \at fnal.gov>
13 Fermi National Accelerator Laboratory
14
15\created Jun 27, 2008
16
17\cpright Copyright (c) 2003-2025, The GENIE Collaboration
18 For the full text of the license visit http://copyright.genie-mc.org
19
20*/
21//____________________________________________________________________________
22
23#ifndef _GNUMI_NEUTRINO_FLUX_H_
24#define _GNUMI_NEUTRINO_FLUX_H_
25
26#include <string>
27#include <iostream>
28#include <vector>
29#include <set>
30
31#include <TVector3.h>
32#include <TLorentzVector.h>
33#include <TLorentzRotation.h>
34
39
40class TFile;
41class TChain;
42class TTree;
43class TBranch;
44
45// MakeClass created classes for handling NuMI flux files
46class g3numi;
47class g4numi;
48class flugg;
49
50using std::string;
51using std::ostream;
52
53namespace genie {
54namespace flux {
55
57ostream & operator << (ostream & stream, const GNuMIFluxPassThroughInfo & info);
58
59/// GNuMIFluxPassThroughInfo:
60/// =========================
61/// A small persistable C-struct -like class that mirrors (some of) the
62/// structure of the gnumi ntuples. This can then be stored as an extra
63/// branch of the output event tree -alongside with the generated event
64/// branch- for use further upstream in the analysis chain - e.g. beam
65/// reweighting etc.
66/// To do future x-y reweighting users must retain the info found in:
67// Ntype Vx Vy Vz
68// pdPx pdPy pdPz
69// ppdxdz ppdydz pppz ppenergy ptype
70// muparpx muparpy muparpz mupare Necm
71// Nimpwt
72///
73class GNuMIFluxPassThroughInfo: public TObject {
74public:
76 /* allow default copy constructor ... for now nothing special
77 GNuMIFluxPassThroughInfo(const GNuMIFluxPassThroughInfo & info);
78 */
80
81 void MakeCopy(const g3numi*); ///< pull in from g3 ntuple
82 void MakeCopy(const g4numi*); ///< pull in from g4 ntuple
83 void MakeCopy(const flugg*); ///< pull in from flugg ntuple
84
85 void ResetCopy(); // reset portion copied from ntuple
86 void ResetCurrent(); // reset generated xy positioned info
87 void ConvertPartCodes();
88 void Print(const Option_t* opt = "") const;
89
90 int CalcEnuWgt(const TLorentzVector& xyz, double& enu, double& wgt_xy) const;
91
92 friend ostream & operator << (ostream & stream, const GNuMIFluxPassThroughInfo & info);
93
94 int pcodes; // 0=original GEANT particle codes, 1=converted to PDG
95 int units; // 0=original GEANT cm, 1=meters
96
97 // Values for GNuMIFlux chosen x-y-z position, not from flux ntuple
98 int fgPdgC; ///< generated nu pdg-code
99 double fgXYWgt; ///< generated nu x-y weight
100 /// not the same as GNuMIFlux::Weight()
101 /// which include importance wgt and deweighting
102 TLorentzVector fgP4; ///< generated nu 4-momentum beam coord
103 TLorentzVector fgX4; ///< generated nu 4-position beam coord
104 TLorentzVector fgP4User; ///< generated nu 4-momentum user coord
105 TLorentzVector fgX4User; ///< generated nu 4-position user coord
106
107 // Values copied from gnumi ntuples (generally maintained variable names)
108 // see http://www.hep.utexas.edu/~zarko/wwwgnumi/v19/[/v19/output_gnumi.html]
109
110 Int_t run;
111 Int_t evtno;
112 Double_t ndxdz;
113 Double_t ndydz;
114 Double_t npz;
115 Double_t nenergy;
116 Double_t ndxdznea;
117 Double_t ndydznea;
118 Double_t nenergyn;
119 Double_t nwtnear;
120 Double_t ndxdzfar;
121 Double_t ndydzfar;
122 Double_t nenergyf;
123 Double_t nwtfar;
124 Int_t norig;
125 Int_t ndecay;
126 Int_t ntype;
127 Double_t vx;
128 Double_t vy;
129 Double_t vz;
130 Double_t pdpx;
131 Double_t pdpy;
132 Double_t pdpz;
133 Double_t ppdxdz;
134 Double_t ppdydz;
135 Double_t pppz;
136 Double_t ppenergy;
137 Int_t ppmedium;
138 Int_t ptype; // converted to PDG
139 Double_t ppvx;
140 Double_t ppvy;
141 Double_t ppvz;
142 Double_t muparpx;
143 Double_t muparpy;
144 Double_t muparpz;
145 Double_t mupare;
146 Double_t necm;
147 Double_t nimpwt;
148 Double_t xpoint;
149 Double_t ypoint;
150 Double_t zpoint;
151 Double_t tvx;
152 Double_t tvy;
153 Double_t tvz;
154 Double_t tpx;
155 Double_t tpy;
156 Double_t tpz;
157 Int_t tptype; // converted to PDG
158 Int_t tgen;
159 Int_t tgptype; // converted to PDG
160 Double_t tgppx;
161 Double_t tgppy;
162 Double_t tgppz;
163 Double_t tprivx;
164 Double_t tprivy;
165 Double_t tprivz;
166 Double_t beamx;
167 Double_t beamy;
168 Double_t beamz;
169 Double_t beampx;
170 Double_t beampy;
171 Double_t beampz;
172
173#ifndef SKIP_MINERVA_MODS
174 //=========================================
175 // The following was inserted by MINERvA
176 //=========================================
177 int getProcessID(TString sval);
178 int getVolID(TString sval);
179
180 static const unsigned int MAX_N_TRAJ = 10; ///< Maximum number of trajectories to store
181
183 Bool_t overflow;
187
203
207 //END of minerva additions
208#endif
209
210ClassDef(GNuMIFluxPassThroughInfo,5)
211};
212
213/// GNuMIFlux:
214/// ==========
215/// An implementation of the GFluxI interface that provides NuMI flux
216///
218 : public genie::GFluxI
221{
222
223public :
224 GNuMIFlux();
225 ~GNuMIFlux();
226
227 // Methods implementing the GENIE GFluxI interface, required for integrating
228 // the NuMI neutrino flux simulations with the GENIE event generation drivers
229
230 const PDGCodeList & FluxParticles (void) { return *fPdgCList; }
231 double MaxEnergy (void) { return fMaxEv; }
232 bool GenerateNext (void);
233 int PdgCode (void) { return fCurEntry->fgPdgC; }
234 double Weight (void) { return fWeight; }
235 const TLorentzVector & Momentum (void) { return fCurEntry->fgP4User; }
236 const TLorentzVector & Position (void) { return fCurEntry->fgX4User; }
237 bool End (void) { return fEnd; }
238 long int Index (void) { return fIEntry; }
239 void Clear (Option_t * opt);
240 void GenerateWeighted (bool gen_weighted);
241
242 // Methods specific to the NuMI flux driver,
243 // for configuration/initialization of the flux & event generation drivers
244 // and and for passing-through flux information (e.g. neutrino parent decay
245 // kinematics) not used by the generator but required by analyses/processing
246 // further downstream
247
248 //
249 // information about or actions on current entry
250 //
252 PassThroughInfo(void) { return *fCurEntry; } ///< GNuMIFluxPassThroughInfo
253 Long64_t GetEntryNumber() { return fIEntry; } ///< index in chain
254
255 double GetDecayDist() const; ///< dist (user units) from dk to current pos
256 void MoveToZ0(double z0); ///< move ray origin to user coord Z0
257
258 //
259 // information about the current state
260 //
261 virtual double GetTotalExposure() const; // GFluxExposureI interface
262 virtual long int NFluxNeutrinos() const; ///< # of rays generated
263
264 double POT_curr(void); ///< current average POT (RWH?)
265 double UsedPOTs(void) const; ///< # of protons-on-target used
266
267 double SumWeight(void) const { return fSumWeight; } ///< integrated weight for flux neutrinos looped so far
268
269 void PrintCurrent(void); ///< print current entry from leaves
270 void PrintConfig(); ///< print the current configuration
271
272 std::vector<std::string> GetFileList(); ///< list of files currently part of chain
273
274 //
275 // GFluxFileConfigI interface
276 //
277 virtual void LoadBeamSimData(const std::vector<std::string>& filenames,
278 const std::string& det_loc);
279 using GFluxFileConfigI::LoadBeamSimData; // inherit the rest
280 virtual void GetBranchInfo(std::vector<std::string>& branchNames,
281 std::vector<std::string>& branchClassNames,
282 std::vector<void**>& branchObjPointers);
283 virtual TTree* GetMetaDataTree();
284
285 //
286 // configuration of GNuMIFlux
287 //
288
289 bool LoadConfig(string cfg); ///< load a named configuration
290 void SetMaxEnergy(double Ev); ///< specify maximum flx neutrino energy
291
292 void SetGenWeighted(bool genwgt=false) { fGenWeighted = genwgt; } ///< toggle whether GenerateNext() returns weight=1 flux (initial default false)
293
294 void SetEntryReuse(long int nuse=1); ///< # of times to use entry before moving to next
295
296 void SetTreeName(string name); ///< set input tree name (default: "h10")
297 void ScanForMaxWeight(void); ///< scan for max flux weight (before generating unweighted flux neutrinos)
298 void SetMaxWgtScan(double fudge = 1.05, long int nentries = 2500000) ///< configuration when estimating max weight
299 { fMaxWgtFudge = fudge; fMaxWgtEntries = nentries; }
300 void SetMaxEFudge(double fudge = 1.05) ///< extra fudge factor in estimating maximum energy
301 { fMaxEFudge = fudge; }
302 void SetApplyWindowTiltWeight(bool apply = true) ///< apply wgt due to tilt of flux window relative to beam
303 { fApplyTiltWeight = apply; }
304
305
306 // GNuMIFlux uses "cm" as the length unit consistently internally (this is
307 // the length units used by both the g3 and g4 ntuples). User interactions
308 // setting the beam-to-detector coordinate transform, flux window, and the
309 // returned position might need to be in other units. Use:
310 // double scale = genie::utils::units::UnitFromString("cm");
311 // ( #include "Utils/UnitUtils.h for declaration )
312 // to get the correct scale factor to pass in. This should get set
313 // FIRST before setting detector position/rotation
314
315 void SetLengthUnits(double user_units); ///< Set units assumed by user
316 double LengthUnits(void) const; ///< Return user units
317
318 // set relative orientation of user coords vs. beam system, i.e.
319 // x3_user = ( beamrot * x3_beam ) + x0beam_user
320 // p3_user = beamrot * p3_beam
321
322 ///< beam (0,0,0) relative to user frame, beam direction in user frame
323 void SetBeamRotation(TRotation beamrot);
324 void SetBeamCenter(TVector3 beam0);
325 TRotation GetBeamRotation() const; ///< rotation to apply from beam->user
326 TVector3 GetBeamCenter() const; ///< beam origin in user frame
327
328 // configure a flux window (or point) where E_nu and weight are evaluated
329
330 typedef enum EStdFluxWindow {
331 kMinosNearDet, // appropriate for Near Detector
332 kMinosFarDet, // appropriate for Far Detector
333 kMinosNearRock, // appropriate for Near rock generation
334 kMinosFarRock, // appropriate for Far rock generation
335 kMinosNearCenter, // point location mimic near value in ntuple
336 kMinosFarCenter // point location mimic far value in ntuple
338
339 // set both flux window in user coord and coordinate transform
340 // for some standard conditions
341 bool SetFluxWindow(StdFluxWindow_t stdwindow, double padding=0); ///< return false if unhandled
342
343 // rwh: potential upgrade: allow flux window set/get in beam coords
344 // as optional flag to *etFluxWindow
345 void SetFluxWindow(TVector3 p1, TVector3 p2, TVector3 p3); ///< 3 points define a plane (by default in user coordinates)
346 void GetFluxWindow(TVector3& p1, TVector3& p2, TVector3& p3) const; ///< 3 points define a plane in beam coordinate
347
348 /// force weights at MINOS detector "center" as found in ntuple
349 void UseFluxAtNearDetCenter(void);
350 void UseFluxAtFarDetCenter(void);
351
352 //
353 // Actual coordinate transformations b=beam, u=user (e.g. detector)
354 //
355 void Beam2UserPos(const TLorentzVector& beamxyz,
356 TLorentzVector& usrxyz ) const;
357 void Beam2UserDir(const TLorentzVector& beamdir,
358 TLorentzVector& usrdir ) const;
359 void Beam2UserP4 (const TLorentzVector& beamp4,
360 TLorentzVector& usrp4 ) const;
361 void User2BeamPos(const TLorentzVector& usrxyz,
362 TLorentzVector& beamxyz ) const;
363 void User2BeamDir(const TLorentzVector& usrdir,
364 TLorentzVector& beamdir ) const;
365 void User2BeamP4 (const TLorentzVector& usrp4,
366 TLorentzVector& beamp4 ) const;
367
368 TVector3 FluxWindowNormal() { return fWindowNormal; }
369
370private:
371
372 // Private methods
373 //
374 bool GenerateNext_weighted (void);
375 void Initialize (void);
376 void SetDefaults (void);
377 void CleanUp (void);
378 void ResetCurrent (void);
379 void AddFile (TTree* tree, string fname);
380 void CalcEffPOTsPerNu (void);
381
382 // Private data members
383 //
384 double fMaxEv; ///< maximum energy
385 bool fEnd; ///< end condition reached
386
387 std::vector<string> fNuFluxFilePatterns; ///< (potentially wildcarded) path(s)
388 string fNuFluxTreeName; ///< Tree name "h10" (g3) or "nudata" (g4)
389 TChain* fNuFluxTree; ///< TTree in g3numi or g4numi // REF ONLY!
390 string fNuFluxGen; ///< "g3numi" "g4numi" or "flugg"
391 g3numi* fG3NuMI; ///< g3numi ntuple
392 g4numi* fG4NuMI; ///< g4numi ntuple
393 flugg* fFlugg; ///< flugg ntuple
394 int fNFiles; ///< number of files in chain
395 Long64_t fNEntries; ///< number of flux ntuple entries
396 Long64_t fIEntry; ///< current flux ntuple entry
397 Long64_t fNuTot; ///< cummulative # of entries (=fNEntries)
398 Long64_t fFilePOTs; ///< # of protons-on-target represented by all files
399
400 double fWeight; ///< current neutrino weight, =1 if generating unweighted entries
401 double fMaxWeight; ///< max flux neutrino weight in input file
402 double fMaxWgtFudge; ///< fudge factor for estimating max wgt
403 long int fMaxWgtEntries; ///< # of entries in estimating max wgt
404 double fMaxEFudge; ///< fudge factor for estmating max enu (0=> use fixed 120GeV)
405
406 long int fNUse; ///< how often to use same entry in a row
407 long int fIUse; ///< current # of times an entry has been used
408 double fSumWeight; ///< sum of weights for nus thrown so far
409 long int fNNeutrinos; ///< number of flux neutrinos thrown so far
410 double fEffPOTsPerNu; ///< what a entry is worth ...
411 double fAccumPOTs; ///< POTs used so far
412
413 bool fGenWeighted; ///< does GenerateNext() give weights?
414 bool fApplyTiltWeight; ///< wgt due to window normal not || beam
415 bool fDetLocIsSet; ///< is a flux location (near/far) set?
416 int fUseFluxAtDetCenter; ///< use flux at near (-1) or far (+1) det center from ntuple?
417
418 double fLengthUnits; ///< units for coord in user exchanges
419 double fLengthScaleB2U; ///< scale factor beam (cm) --> user
420 double fLengthScaleU2B; ///< scale factor beam user --> (cm)
421
422 TLorentzVector fBeamZero; ///< beam origin in user coords
423 TLorentzRotation fBeamRot; ///< rotation applied beam --> user coord
424 TLorentzRotation fBeamRotInv;
425
426 TVector3 fFluxWindowPtUser[3]; ///< user points of flux window
427 TLorentzVector fFluxWindowBase; ///< base point for flux window - beam coord
428 TLorentzVector fFluxWindowDir1; ///< extent for flux window (direction 1)
429 TLorentzVector fFluxWindowDir2; ///< extent for flux window (direction 2)
432 TVector3 fWindowNormal; ///< normal direction for flux window
433
434 TLorentzVector fgX4dkvtx; ///< decay 4-position beam coord
435
436 GNuMIFluxPassThroughInfo* fCurEntry; ///< copy of current ntuple entry info (owned structure)
437
438};
439
440//#define GNUMI_TEST_XY_WGT
441#ifdef GNUMI_TEST_XY_WGT
442class xypartials {
443 // intermediate partial info from xy reweighting for comparison w/ f77 version
444 friend ostream & operator << (ostream & stream, const xypartials & info);
445public:
446 xypartials() { ; }
447 void ReadStream(std::ifstream& myfile);
448 int Compare(const xypartials& other) const;
449 void Print(const Option_t* opt = "") const;
450 static xypartials& GetStaticInstance(); // copy used by CalcEnuWgt()
451 // actual data
452 double xdet, ydet, zdet;
453 double parent_mass, parentp, parent_energy;
454 double gamma, beta_mag, enuzr, rad;
455 double costh_pardet, theta_pardet, emrat, eneu;
456 double sangdet, wgt;
457 double betanu[3], p_nu[3], partial1, p_dcm_nu[4];
458 double muparent_px, muparent_py, muparent_pz;
459 double gammamp, betamp[3], partial2, p_pcm_mp[3], p_pcm;
460 double costhmu, wgt_ratio;
461 int ptype, ntype;
462
463};
464#endif
465
466} // flux namespace
467} // genie namespace
468
469#endif // _GNUMI_NEUTRINO_FLUX_H_
Definition flugg.h:15
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
TLorentzVector fgP4User
generated nu 4-momentum user coord
Definition GNuMIFlux.h:104
TLorentzVector fgX4
generated nu 4-position beam coord
Definition GNuMIFlux.h:103
friend ostream & operator<<(ostream &stream, const GNuMIFluxPassThroughInfo &info)
int fgPdgC
generated nu pdg-code
Definition GNuMIFlux.h:98
void MakeCopy(const g3numi *)
pull in from g3 ntuple
TLorentzVector fgP4
generated nu 4-momentum beam coord
Definition GNuMIFlux.h:102
TLorentzVector fgX4User
generated nu 4-position user coord
Definition GNuMIFlux.h:105
int CalcEnuWgt(const TLorentzVector &xyz, double &enu, double &wgt_xy) const
static const unsigned int MAX_N_TRAJ
Maximum number of trajectories to store.
Definition GNuMIFlux.h:180
void Print(const Option_t *opt="") const
TLorentzVector fFluxWindowDir2
extent for flux window (direction 2)
Definition GNuMIFlux.h:429
int PdgCode(void)
returns the flux neutrino pdg code
Definition GNuMIFlux.h:233
void UseFluxAtFarDetCenter(void)
double fLengthScaleB2U
scale factor beam (cm) --> user
Definition GNuMIFlux.h:419
TVector3 fWindowNormal
normal direction for flux window
Definition GNuMIFlux.h:432
string fNuFluxGen
"g3numi" "g4numi" or "flugg"
Definition GNuMIFlux.h:390
g4numi * fG4NuMI
g4numi ntuple
Definition GNuMIFlux.h:392
GNuMIFluxPassThroughInfo * fCurEntry
copy of current ntuple entry info (owned structure)
Definition GNuMIFlux.h:436
bool fEnd
end condition reached
Definition GNuMIFlux.h:385
TRotation GetBeamRotation() const
rotation to apply from beam->user
void SetApplyWindowTiltWeight(bool apply=true)
Definition GNuMIFlux.h:302
double GetDecayDist() const
dist (user units) from dk to current pos
const TLorentzVector & Momentum(void)
returns the flux neutrino 4-momentum
Definition GNuMIFlux.h:235
double POT_curr(void)
current average POT (RWH?)
TVector3 fFluxWindowPtUser[3]
user points of flux window
Definition GNuMIFlux.h:426
TChain * fNuFluxTree
TTree in g3numi or g4numi // REF ONLY!
Definition GNuMIFlux.h:389
long int Index(void)
returns corresponding index for current flux neutrino (e.g. for a flux ntuple returns the current ent...
Definition GNuMIFlux.h:238
TLorentzVector fBeamZero
beam origin in user coords
Definition GNuMIFlux.h:422
void SetBeamCenter(TVector3 beam0)
Long64_t fIEntry
current flux ntuple entry
Definition GNuMIFlux.h:396
virtual long int NFluxNeutrinos() const
long int fNUse
how often to use same entry in a row
Definition GNuMIFlux.h:406
void CalcEffPOTsPerNu(void)
void SetMaxWgtScan(double fudge=1.05, long int nentries=2500000)
Definition GNuMIFlux.h:298
bool GenerateNext(void)
generate the next flux neutrino (return false in err)
void Beam2UserPos(const TLorentzVector &beamxyz, TLorentzVector &usrxyz) const
Long64_t fNuTot
cummulative # of entries (=fNEntries)
Definition GNuMIFlux.h:397
TLorentzRotation fBeamRot
rotation applied beam --> user coord
Definition GNuMIFlux.h:423
bool fDetLocIsSet
is a flux location (near/far) set?
Definition GNuMIFlux.h:415
void ScanForMaxWeight(void)
scan for max flux weight (before generating unweighted flux neutrinos)
void SetMaxEnergy(double Ev)
specify maximum flx neutrino energy
bool SetFluxWindow(StdFluxWindow_t stdwindow, double padding=0)
return false if unhandled
void SetLengthUnits(double user_units)
Set units assumed by user.
void AddFile(TTree *tree, string fname)
std::vector< string > fNuFluxFilePatterns
(potentially wildcarded) path(s)
Definition GNuMIFlux.h:387
void User2BeamP4(const TLorentzVector &usrp4, TLorentzVector &beamp4) const
virtual TTree * GetMetaDataTree()
enum genie::flux::GNuMIFlux::EStdFluxWindow StdFluxWindow_t
double fLengthUnits
units for coord in user exchanges
Definition GNuMIFlux.h:418
void SetTreeName(string name)
set input tree name (default: "h10")
double fMaxWgtFudge
fudge factor for estimating max wgt
Definition GNuMIFlux.h:402
int fNFiles
number of files in chain
Definition GNuMIFlux.h:394
double fWeight
current neutrino weight, =1 if generating unweighted entries
Definition GNuMIFlux.h:400
TLorentzVector fgX4dkvtx
decay 4-position beam coord
Definition GNuMIFlux.h:434
bool GenerateNext_weighted(void)
double fSumWeight
sum of weights for nus thrown so far
Definition GNuMIFlux.h:408
double SumWeight(void) const
integrated weight for flux neutrinos looped so far
Definition GNuMIFlux.h:267
const TLorentzVector & Position(void)
returns the flux neutrino 4-position (note: expect SI rather than physical units)
Definition GNuMIFlux.h:236
void MoveToZ0(double z0)
move ray origin to user coord Z0
double Weight(void)
returns the flux neutrino weight (if any)
Definition GNuMIFlux.h:234
void Beam2UserDir(const TLorentzVector &beamdir, TLorentzVector &usrdir) const
TLorentzRotation fBeamRotInv
Definition GNuMIFlux.h:424
Long64_t GetEntryNumber()
index in chain
Definition GNuMIFlux.h:253
virtual void GetBranchInfo(std::vector< std::string > &branchNames, std::vector< std::string > &branchClassNames, std::vector< void ** > &branchObjPointers)
double fMaxWeight
max flux neutrino weight in input file
Definition GNuMIFlux.h:401
TVector3 FluxWindowNormal()
Definition GNuMIFlux.h:368
double MaxEnergy(void)
declare the max flux neutrino energy that can be generated (for init. purposes)
Definition GNuMIFlux.h:231
long int fIUse
current # of times an entry has been used
Definition GNuMIFlux.h:407
void GetFluxWindow(TVector3 &p1, TVector3 &p2, TVector3 &p3) const
3 points define a plane in beam coordinate
double UsedPOTs(void) const
void Beam2UserP4(const TLorentzVector &beamp4, TLorentzVector &usrp4) const
void PrintCurrent(void)
print current entry from leaves
int fUseFluxAtDetCenter
use flux at near (-1) or far (+1) det center from ntuple?
Definition GNuMIFlux.h:416
void SetEntryReuse(long int nuse=1)
bool fApplyTiltWeight
wgt due to window normal not || beam
Definition GNuMIFlux.h:414
bool End(void)
true if no more flux nu's can be thrown (eg reaching end of beam sim ntuples)
Definition GNuMIFlux.h:237
long int fNNeutrinos
number of flux neutrinos thrown so far
Definition GNuMIFlux.h:409
void User2BeamPos(const TLorentzVector &usrxyz, TLorentzVector &beamxyz) const
void PrintConfig()
print the current configuration
void UseFluxAtNearDetCenter(void)
force weights at MINOS detector "center" as found in ntuple
void SetGenWeighted(bool genwgt=false)
toggle whether GenerateNext() returns weight=1 flux (initial default false)
Definition GNuMIFlux.h:292
double fEffPOTsPerNu
what a entry is worth ...
Definition GNuMIFlux.h:410
double fAccumPOTs
POTs used so far.
Definition GNuMIFlux.h:411
string fNuFluxTreeName
Tree name "h10" (g3) or "nudata" (g4)
Definition GNuMIFlux.h:388
g3numi * fG3NuMI
g3numi ntuple
Definition GNuMIFlux.h:391
double fMaxEv
maximum energy
Definition GNuMIFlux.h:384
const GNuMIFluxPassThroughInfo & PassThroughInfo(void)
GNuMIFluxPassThroughInfo.
Definition GNuMIFlux.h:252
double fLengthScaleU2B
scale factor beam user --> (cm)
Definition GNuMIFlux.h:420
void GenerateWeighted(bool gen_weighted)
set whether to generate weighted or unweighted neutrinos
TLorentzVector fFluxWindowBase
base point for flux window - beam coord
Definition GNuMIFlux.h:427
bool LoadConfig(string cfg)
load a named configuration
const PDGCodeList & FluxParticles(void)
declare list of flux neutrinos that can be generated (for init. purposes)
Definition GNuMIFlux.h:230
void SetBeamRotation(TRotation beamrot)
< beam (0,0,0) relative to user frame, beam direction in user frame
double fMaxEFudge
fudge factor for estmating max enu (0=> use fixed 120GeV)
Definition GNuMIFlux.h:404
bool fGenWeighted
does GenerateNext() give weights?
Definition GNuMIFlux.h:413
void Clear(Option_t *opt)
reset state variables based on opt
virtual void LoadBeamSimData(const std::vector< std::string > &filenames, const std::string &det_loc)
std::vector< std::string > GetFileList()
list of files currently part of chain
void SetMaxEFudge(double fudge=1.05)
Definition GNuMIFlux.h:300
virtual double GetTotalExposure() const
TLorentzVector fFluxWindowDir1
extent for flux window (direction 1)
Definition GNuMIFlux.h:428
flugg * fFlugg
flugg ntuple
Definition GNuMIFlux.h:393
Long64_t fNEntries
number of flux ntuple entries
Definition GNuMIFlux.h:395
void User2BeamDir(const TLorentzVector &usrdir, TLorentzVector &beamdir) const
TVector3 GetBeamCenter() const
beam origin in user frame
double LengthUnits(void) const
Return user units.
ostream & operator<<(ostream &stream, const TClonesArray *particle_list)
GENIE flux drivers.
ostream & operator<<(ostream &stream, const genie::flux::GJPARCNuFluxPassThroughInfo &info)
THE MAIN GENIE PROJECT NAMESPACE
Definition AlgCmp.h:25