GENIEGenerator
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Private Member Functions | Private Attributes | List of all members
genie::flux::GSimpleNtpFlux Class Reference

A GENIE flux driver using a simple ntuple format. More...

#include <GSimpleNtpFlux.h>

Inheritance diagram for genie::flux::GSimpleNtpFlux:
Inheritance graph
[legend]
Collaboration diagram for genie::flux::GSimpleNtpFlux:
Collaboration graph
[legend]

Public Member Functions

 GSimpleNtpFlux ()
 
 ~GSimpleNtpFlux ()
 
const PDGCodeListFluxParticles (void)
 declare list of flux neutrinos that can be generated (for init. purposes) More...
 
double MaxEnergy (void)
 declare the max flux neutrino energy that can be generated (for init. purposes) More...
 
bool GenerateNext (void)
 generate the next flux neutrino (return false in err) More...
 
int PdgCode (void)
 returns the flux neutrino pdg code More...
 
double Weight (void)
 returns the flux neutrino weight (if any) More...
 
const TLorentzVector & Momentum (void)
 returns the flux neutrino 4-momentum More...
 
const TLorentzVector & Position (void)
 returns the flux neutrino 4-position (note: expect SI rather than physical units) More...
 
bool End (void)
 true if no more flux nu's can be thrown (eg reaching end of beam sim ntuples) More...
 
long int Index (void)
 returns corresponding index for current flux neutrino (e.g. for a flux ntuple returns the current entry number) More...
 
void Clear (Option_t *opt)
 reset state variables based on opt More...
 
void GenerateWeighted (bool gen_weighted)
 set whether to generate weighted or unweighted neutrinos More...
 
const
genie::flux::GSimpleNtpEntry
GetCurrentEntry (void)
 GSimpleNtpEntry. More...
 
const genie::flux::GSimpleNtpNuMIGetCurrentNuMI (void)
 GSimpleNtpNuMI. More...
 
const genie::flux::GSimpleNtpAuxGetCurrentAux (void)
 GSimpleNtpAux. More...
 
const genie::flux::GSimpleNtpMetaGetCurrentMeta (void)
 GSimpleNtpMeta. More...
 
TChain * GetFluxTChain (void)
 
double GetDecayDist () const
 dist (user units) from dk to current pos More...
 
void MoveToZ0 (double z0)
 move ray origin to user coord Z0 More...
 
void SetIncludeVtxt (bool it=true)
 
bool GetIncludeVtxt ()
 should X4 include CurEntry.vtxt More...
 
virtual double GetTotalExposure () const
 GFluxExposureI interface. More...
 
virtual long int NFluxNeutrinos () const
 

of rays generated

More...
 
double UsedPOTs (void) const
 

of protons-on-target used

More...
 
long int NEntriesUsed (void) const
 number of entries read from files More...
 
double SumWeight (void) const
 integrated weight for flux neutrinos looped so far More...
 
void PrintCurrent (void)
 print current entry from leaves More...
 
void PrintConfig ()
 print the current configuration More...
 
std::vector< std::string > GetFileList ()
 list of files currently part of chain More...
 
virtual void LoadBeamSimData (const std::vector< string > &filenames, const std::string &det_loc)
 
virtual void GetBranchInfo (std::vector< std::string > &branchNames, std::vector< std::string > &branchClassNames, std::vector< void ** > &branchObjPointers)
 
virtual TTree * GetMetaDataTree ()
 
void SetRequestedBranchList (string blist="entry,numi,aux")
 
void SetMaxEnergy (double Ev)
 specify maximum flx neutrino energy More...
 
void SetGenWeighted (bool genwgt=false)
 toggle whether GenerateNext() returns weight=1 flux (initial default false) More...
 
void SetEntryReuse (long int nuse=1)
 

of times to use entry before moving to next

More...
 
void ProcessMeta (void)
 scan for max flux energy, weight More...
 
void GetFluxWindow (TVector3 &p1, TVector3 &p2, TVector3 &p3) const
 3 points define a plane in beam coordinate More...
 
- Public Member Functions inherited from genie::GFluxI
virtual ~GFluxI ()
 
- Public Member Functions inherited from genie::flux::GFluxExposureI
 GFluxExposureI (genie::flux::Exposure_t etype)
 
virtual ~GFluxExposureI ()
 
const char * GetExposureUnits () const
 what units are returned by GetTotalExposure? More...
 
genie::flux::Exposure_t GetExposureType () const
 
- Public Member Functions inherited from genie::flux::GFluxFileConfigI
 GFluxFileConfigI ()
 
virtual ~GFluxFileConfigI ()
 
virtual void LoadBeamSimData (const std::vector< std::string > &filenames, const std::string &det_loc)=0
 
virtual void LoadBeamSimData (const std::set< std::string > &filenames, const std::string &det_loc)
 
virtual void LoadBeamSimData (const std::string &filename, const std::string &det_loc)
 
virtual void SetXMLFileBase (std::string xmlbasename="")
 
virtual std::string GetXMLFileBase () const
 
virtual void SetFluxParticles (const PDGCodeList &particles)
 specify list of flux neutrino species More...
 
virtual void SetUpstreamZ (double z0)
 
virtual void SetNumOfCycles (long int ncycle)
 limit cycling through input files More...
 

Private Member Functions

bool GenerateNext_weighted (void)
 
void Initialize (void)
 
void SetDefaults (void)
 
void CleanUp (void)
 
void ResetCurrent (void)
 
void AddFile (TTree *fluxtree, TTree *metatree, string fname)
 
bool OptionalAttachBranch (std::string bname)
 
void CalcEffPOTsPerNu (void)
 
void ScanMeta (void)
 

Private Attributes

double fMaxEv
 maximum energy More...
 
bool fEnd
 end condition reached More...
 
std::vector< string > fNuFluxFilePatterns
 (potentially wildcarded) path(s) More...
 
string fNuFluxBranchRequest
 list of requested branches "entry,numi,au" More...
 
TChain * fNuFluxTree
 TTree // REF ONLY. More...
 
TChain * fNuMetaTree
 TTree // REF ONLY. More...
 
int fNFiles
 number of files in chain More...
 
Long64_t fNEntries
 number of flux ntuple entries More...
 
Long64_t fIEntry
 current flux ntuple entry More...
 
Int_t fIFileNumber
 which file for the current entry More...
 
Double_t fFilePOTs
 

of protons-on-target represented by all files

More...
 
double fWeight
 current neutrino weight More...
 
double fMaxWeight
 max flux neutrino weight in input file More...
 
long int fNUse
 how often to use same entry in a row More...
 
long int fIUse
 current # of times an entry has been used More...
 
double fSumWeight
 sum of weights for nus thrown so far More...
 
long int fNNeutrinos
 number of flux neutrinos thrown so far More...
 
long int fNEntriesUsed
 number of entries read from files More...
 
double fEffPOTsPerNu
 what a entry is worth ... More...
 
double fAccumPOTs
 POTs used so far. More...
 
bool fGenWeighted
 does GenerateNext() give weights? More...
 
bool fAlreadyUnwgt
 are input files already unweighted More...
 
bool fAllFilesMeta
 do all files in chain have meta data More...
 
GSimpleNtpEntryfCurEntry
 current entry More...
 
GSimpleNtpNuMIfCurNuMI
 current "numi" branch extra info More...
 
GSimpleNtpAuxfCurAux
 current "aux" branch extra info More...
 
TLorentzVector fP4
 reconstituted p4 vector More...
 
TLorentzVector fX4
 reconstituted position vector More...
 
GSimpleNtpMetafCurMeta
 current meta data More...
 
GSimpleNtpEntryfCurEntryCopy
 current entry More...
 
GSimpleNtpNuMIfCurNuMICopy
 current "numi" branch extra info More...
 
GSimpleNtpAuxfCurAuxCopy
 current "aux" branch extra info More...
 
bool fIncludeVtxt
 does fX4 include CurEntry.vtxt or 0 More...
 

Additional Inherited Members

- Static Public Member Functions inherited from genie::flux::GFluxExposureI
static const char * AsString (genie::flux::Exposure_t etype)
 
static genie::flux::Exposure_t StringToEnum (const char *chars, int maxChar=0)
 
- Protected Member Functions inherited from genie::GFluxI
 GFluxI ()
 
- Protected Attributes inherited from genie::flux::GFluxFileConfigI
PDGCodeListfPdgCList
 list of neutrino pdg-codes to generate More...
 
PDGCodeListfPdgCListRej
 list of nu pdg-codes seen but rejected More...
 
std::string fXMLbasename
 XML file that might hold config param_sets. More...
 
long int fNCycles
 

times to cycle through the ntuple(s)

More...
 
long int fICycle
 
double fZ0
 

Detailed Description

A GENIE flux driver using a simple ntuple format.

GSimpleNtpFlux:

An implementation of the GFluxI interface that provides NuMI flux

Author
Robert Hatcher <rhatcher fnal.gov> Fermi National Accelerator Laboratory
Created:
Jan 25, 2010
License:
Copyright (c) 2003-2024, The GENIE Collaboration For the full text of the license visit http://copyright.genie-mc.org

Definition at line 199 of file GSimpleNtpFlux.h.

Constructor & Destructor Documentation

GSimpleNtpFlux::GSimpleNtpFlux ( )

Definition at line 67 of file GSimpleNtpFlux.cxx.

References Initialize().

67  :
69 {
70  this->Initialize();
71 }
GFluxExposureI(genie::flux::Exposure_t etype)
GSimpleNtpFlux::~GSimpleNtpFlux ( )

Definition at line 73 of file GSimpleNtpFlux.cxx.

References CleanUp().

74 {
75  this->CleanUp();
76 }

Member Function Documentation

void GSimpleNtpFlux::AddFile ( TTree *  fluxtree,
TTree *  metatree,
string  fname 
)
private

Definition at line 802 of file GSimpleNtpFlux.cxx.

References fAllFilesMeta, fNFiles, fNuFluxTree, fNuMetaTree, LOG, and pINFO.

Referenced by LoadBeamSimData().

803 {
804  // Add a file to the chain
805  if ( ! fluxtree ) return;
806 
807  ULong64_t nentries = fluxtree->GetEntries();
808 
809  int stat = fNuFluxTree->AddFile(fname.c_str());
810  if ( metatree ) fNuMetaTree->AddFile(fname.c_str());
811  else fAllFilesMeta = false;
812 
813  LOG("Flux",pINFO)
814  << "flux->AddFile() of " << nentries
815  << " " << ((metatree)?"[+meta]":"[no-meta]")
816  << " [status=" << stat << "]"
817  << " entries in file: " << fname;
818 
819  if ( stat != 1 ) return;
820 
821  fNFiles++;
822 
823 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
bool fAllFilesMeta
do all files in chain have meta data
#define pINFO
Definition: Messenger.h:62
int fNFiles
number of files in chain
TChain * fNuMetaTree
TTree // REF ONLY.
TChain * fNuFluxTree
TTree // REF ONLY.
void GSimpleNtpFlux::CalcEffPOTsPerNu ( void  )
private

Definition at line 362 of file GSimpleNtpFlux.cxx.

References fEffPOTsPerNu, fFilePOTs, fMaxWeight, fNEntries, and fNuFluxTree.

Referenced by LoadBeamSimData().

363 {
364  // do this if flux window changes or # of files changes
365 
366  if (!fNuFluxTree) return; // not yet fully configured
367 
368  fEffPOTsPerNu = ( (double)fFilePOTs / (double)fNEntries ) / fMaxWeight ;
369 }
Double_t fFilePOTs
of protons-on-target represented by all files
double fMaxWeight
max flux neutrino weight in input file
Long64_t fNEntries
number of flux ntuple entries
double fEffPOTsPerNu
what a entry is worth ...
TChain * fNuFluxTree
TTree // REF ONLY.
void GSimpleNtpFlux::CleanUp ( void  )
private

Definition at line 782 of file GSimpleNtpFlux.cxx.

References fCurAux, fCurEntry, fCurMeta, fCurNuMI, genie::flux::GFluxFileConfigI::fICycle, fIEntry, fIUse, genie::flux::GFluxFileConfigI::fNCycles, fNuFluxTree, fNuMetaTree, fNUse, genie::flux::GFluxFileConfigI::fPdgCList, genie::flux::GFluxFileConfigI::fPdgCListRej, LOG, pINFO, and pNOTICE.

Referenced by ~GSimpleNtpFlux().

783 {
784  LOG("Flux", pINFO) << "Cleaning up...";
785 
786  if (fPdgCList) delete fPdgCList;
787  if (fPdgCListRej) delete fPdgCListRej;
788  if (fCurEntry) delete fCurEntry;
789  if (fCurNuMI) delete fCurNuMI;
790  if (fCurAux) delete fCurAux;
791  if (fCurMeta) delete fCurMeta;
792 
793  if (fNuFluxTree) delete fNuFluxTree;
794  if (fNuMetaTree) delete fNuMetaTree;
795 
796  LOG("Flux", pNOTICE)
797  << " flux file cycles: " << fICycle << " of " << fNCycles
798  << ", entry " << fIEntry << " use: " << fIUse << " of " << fNUse;
799 }
GSimpleNtpAux * fCurAux
current &quot;aux&quot; branch extra info
GSimpleNtpMeta * fCurMeta
current meta data
long int fIUse
current # of times an entry has been used
PDGCodeList * fPdgCList
list of neutrino pdg-codes to generate
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pINFO
Definition: Messenger.h:62
long int fNCycles
times to cycle through the ntuple(s)
GSimpleNtpNuMI * fCurNuMI
current &quot;numi&quot; branch extra info
PDGCodeList * fPdgCListRej
list of nu pdg-codes seen but rejected
TChain * fNuMetaTree
TTree // REF ONLY.
GSimpleNtpEntry * fCurEntry
current entry
Long64_t fIEntry
current flux ntuple entry
long int fNUse
how often to use same entry in a row
#define pNOTICE
Definition: Messenger.h:61
TChain * fNuFluxTree
TTree // REF ONLY.
void GSimpleNtpFlux::Clear ( Option_t *  opt)
virtual

reset state variables based on opt

Implements genie::GFluxI.

Definition at line 685 of file GSimpleNtpFlux.cxx.

References fAccumPOTs, genie::flux::GFluxFileConfigI::fICycle, fNNeutrinos, fSumWeight, LOG, and pWARN.

686 {
687 // Clear the driver state
688 //
689  LOG("Flux", pWARN) << "GSimpleNtpFlux::Clear(" << opt << ") called";
690  // do it in all cases, but EVGDriver/GMCJDriver will pass "CycleHistory"
691 
692  fICycle = 0;
693 
694  fSumWeight = 0;
695  fNNeutrinos = 0;
696  fAccumPOTs = 0;
697 
698 }
long int fNNeutrinos
number of flux neutrinos thrown so far
double fSumWeight
sum of weights for nus thrown so far
#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
double fAccumPOTs
POTs used so far.
bool genie::flux::GSimpleNtpFlux::End ( void  )
inlinevirtual

true if no more flux nu's can be thrown (eg reaching end of beam sim ntuples)

Implements genie::GFluxI.

Definition at line 219 of file GSimpleNtpFlux.h.

References fEnd.

Referenced by GenerateNext().

219 { return fEnd; }
bool fEnd
end condition reached
const PDGCodeList& genie::flux::GSimpleNtpFlux::FluxParticles ( void  )
inlinevirtual

declare list of flux neutrinos that can be generated (for init. purposes)

Implements genie::GFluxI.

Definition at line 212 of file GSimpleNtpFlux.h.

References genie::flux::GFluxFileConfigI::fPdgCList.

212 { return *fPdgCList; }
PDGCodeList * fPdgCList
list of neutrino pdg-codes to generate
bool GSimpleNtpFlux::GenerateNext ( void  )
virtual

generate the next flux neutrino (return false in err)

Implements genie::GFluxI.

Definition at line 90 of file GSimpleNtpFlux.cxx.

References End(), fAlreadyUnwgt, fCurEntry, fGenWeighted, fMaxWeight, fP4, fWeight, fX4, GenerateNext_weighted(), GetCurrentEntry(), genie::RandomGen::Instance(), LOG, genie::utils::print::P4AsShortString(), genie::flux::GSimpleNtpEntry::pdg, pERROR, pNOTICE, genie::RandomGen::RndFlux(), Weight(), and genie::utils::print::X4AsString().

91 {
92 // Get next (unweighted) flux ntuple entry on the specified detector location
93 //
95  while ( true ) {
96  // Check for end of flux ntuple
97  bool end = this->End();
98  if ( end ) {
99  LOG("Flux", pNOTICE) << "GenerateNext signaled End() ";
100  return false;
101  }
102 
103  // Get next weighted flux ntuple entry
104  bool nextok = this->GenerateNext_weighted();
105  if ( fGenWeighted ) return nextok;
106  if ( ! nextok ) continue;
107  if ( fAlreadyUnwgt ) return true;
108 
109  /* RWH - debug purposes
110  if ( fNCycles == 0 ) {
111  LOG("Flux", pNOTICE)
112  << "Got flux entry: " << fIEntry
113  << " - Cycle: "<< fICycle << "/ infinite";
114  } else {
115  LOG("Flux", pNOTICE)
116  << "Got flux entry: "<< fIEntry
117  << " - Cycle: "<< fICycle << "/"<< fNCycles;
118  }
119  */
120 
121  // Get fractional weight & decide whether to accept curr flux neutrino
122  double f = this->Weight() / fMaxWeight;
123  //LOG("Flux", pNOTICE)
124  // << "Curr flux neutrino fractional weight = " << f;
125  if (f > 1.) {
126  fMaxWeight = this->Weight() * 1.01; // bump the weight
127  LOG("Flux", pERROR)
128  << "** Fractional weight = " << f
129  << " > 1 !! Bump fMaxWeight estimate to " << fMaxWeight
130  << GetCurrentEntry();
131  std::cout << std::flush;
132  }
133  double r = (f < 1.) ? rnd->RndFlux().Rndm() : 0;
134  bool accept = ( r < f );
135  if ( accept ) {
136 
137 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
138  LOG("Flux", pNOTICE)
139  << "Generated beam neutrino: "
140  << "\n pdg-code: " << fCurEntry->pdg
141  << "\n p4: " << utils::print::P4AsShortString(&(fP4))
142  << "\n x4: " << utils::print::X4AsString(&(fX4));
143 #endif
144 
145  fWeight = 1.;
146  return true;
147  }
148 
149  //LOG("Flux", pNOTICE)
150  // << "** Rejecting current flux neutrino based on the flux weight only";
151  }
152  return false;
153 }
string P4AsShortString(const TLorentzVector *p)
Definition: PrintUtils.cxx:45
#define pERROR
Definition: Messenger.h:59
static RandomGen * Instance()
Access instance.
Definition: RandomGen.cxx:71
bool End(void)
true if no more flux nu&#39;s can be thrown (eg reaching end of beam sim ntuples)
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition: RandomGen.h:29
const genie::flux::GSimpleNtpEntry * GetCurrentEntry(void)
GSimpleNtpEntry.
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
TLorentzVector fX4
reconstituted position vector
bool fGenWeighted
does GenerateNext() give weights?
double fMaxWeight
max flux neutrino weight in input file
double fWeight
current neutrino weight
GSimpleNtpEntry * fCurEntry
current entry
TLorentzVector fP4
reconstituted p4 vector
bool fAlreadyUnwgt
are input files already unweighted
#define pNOTICE
Definition: Messenger.h:61
TRandom3 & RndFlux(void) const
rnd number generator used by flux drivers
Definition: RandomGen.h:71
string X4AsString(const TLorentzVector *x)
Definition: PrintUtils.cxx:57
double Weight(void)
returns the flux neutrino weight (if any)
bool GSimpleNtpFlux::GenerateNext_weighted ( void  )
private

user might modify list via SetFluxParticles() in order to reject certain flavors, even if they're found in the file. So don't make a big fuss. Spit out a single message and then stop reporting that flavor as problematic.

Definition at line 155 of file GSimpleNtpFlux.cxx.

References genie::flux::GSimpleNtpEntry::E, genie::flux::GSimpleNtpNuMI::entryno, genie::flux::GSimpleNtpNuMI::evtno, genie::PDGCodeList::ExistsInPDGCodeList(), fAccumPOTs, fAllFilesMeta, fCurEntry, fCurMeta, fCurNuMI, fEffPOTsPerNu, fEnd, genie::flux::GFluxFileConfigI::fICycle, fIEntry, fIncludeVtxt, fIUse, fMaxEv, genie::flux::GFluxFileConfigI::fNCycles, fNEntries, fNEntriesUsed, fNNeutrinos, fNuFluxTree, fNuMetaTree, fP4, genie::flux::GFluxFileConfigI::fPdgCList, genie::flux::GFluxFileConfigI::fPdgCListRej, fSumWeight, fWeight, fX4, genie::flux::GFluxFileConfigI::fZ0, LOG, genie::flux::GSimpleNtpEntry::metakey, genie::flux::GSimpleNtpMeta::metakey, MoveToZ0(), genie::utils::print::P4AsShortString(), pDEBUG, genie::flux::GSimpleNtpEntry::pdg, pERROR, pFATAL, pINFO, genie::PDGCodeList::push_back(), pWARN, genie::flux::GSimpleNtpEntry::px, genie::flux::GSimpleNtpEntry::py, genie::flux::GSimpleNtpEntry::pz, ResetCurrent(), genie::flux::GSimpleNtpEntry::vtxt, genie::flux::GSimpleNtpEntry::vtxx, genie::flux::GSimpleNtpEntry::vtxy, genie::flux::GSimpleNtpEntry::vtxz, Weight(), genie::flux::GSimpleNtpEntry::wgt, and genie::utils::print::X4AsString().

Referenced by GenerateNext().

156 {
157 // Get next (weighted) flux ntuple entry
158 //
159 
160  // Check whether a flux ntuple has been loaded
161  if ( ! fNuFluxTree ) {
162  LOG("Flux", pFATAL)
163  << "The flux driver has not been properly configured";
164  // return false; // don't do this - creates an infinite loop!
165  exit(1);
166  }
167 
168  // Reuse an entry?
169  //std::cout << " ***** iuse " << fIUse << " nuse " << fNUse
170  // << " ientry " << fIEntry << " nentry " << fNEntries
171  // << " icycle " << fICycle << " ncycle " << fNCycles << std::endl;
172  if ( fIUse < fNUse && fIEntry >= 0 ) {
173  // Reuse this entry
174  fIUse++;
175  } else {
176  // Reset previously generated neutrino code / 4-p / 4-x
177  this->ResetCurrent();
178  // Move on, read next flux ntuple entry
179  ++fIEntry;
180  ++fNEntriesUsed; // count total # used
181  if ( fIEntry >= fNEntries ) {
182  // Ran out of entries @ the current cycle of this flux file
183  // Check whether more (or infinite) number of cycles is requested
184  if (fICycle < fNCycles || fNCycles == 0 ) {
185  fICycle++;
186  fIEntry=0;
187  } else {
188  LOG("Flux", pWARN)
189  << "No more entries in input flux neutrino ntuple, cycle "
190  << fICycle << " of " << fNCycles;
191  fEnd = true;
192  //assert(0);
193  return false;
194  }
195  }
196 
197  int nbytes = fNuFluxTree->GetEntry(fIEntry);
198  UInt_t metakey = fCurEntry->metakey;
199  if ( fAllFilesMeta && ( fCurMeta->metakey != metakey ) ) {
200  UInt_t oldkey = fCurMeta->metakey;
201 #ifdef USE_INDEX_FOR_META
202  int nbmeta = fNuMetaTree->GetEntryWithIndex(metakey);
203 #else
204  // unordered indices makes ROOT call Error() which might,
205  // if not DefaultErrorHandler, be fatal.
206  // so find the right one by a simple linear search.
207  // not a large burden since it only happens infrequently and
208  // the list is normally quite short.
209  int nmeta = fNuMetaTree->GetEntries();
210  int nbmeta = 0;
211  for (int imeta = 0; imeta < nmeta; ++imeta ) {
212  nbmeta = fNuMetaTree->GetEntry(imeta);
213  if ( fCurMeta->metakey == metakey ) break;
214  }
215  // next condition should never happen
216  if ( fCurMeta->metakey != metakey ) {
217  fCurMeta = 0; // didn't find it!?
218  LOG("Flux",pERROR) << "Failed to find right metakey=" << metakey
219  << " (was " << oldkey << ") out of " << nmeta
220  << " entries";
221  }
222 #endif
223  LOG("Flux",pDEBUG) << "Get meta " << metakey
224  << " (was " << oldkey << ") "
225  << fCurMeta->metakey
226  << " nb " << nbytes << " " << nbmeta;
227 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
228  LOG("Flux",pDEBUG) << "Get meta " << *fCurMeta;
229 #endif
230  }
231 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
232  Int_t ifile = fNuFluxTree->GetFileNumber();
233  LOG("Flux",pDEBUG)
234  << "got " << fNNeutrinos << " nu, using fIEntry " << fIEntry
235  << " ifile " << ifile << " nbytes " << nbytes
236  << *fCurEntry << *fCurMeta;
237 #endif
238 
239  fIUse = 1;
240 
241  // here we might want to do flavor oscillations or simple mappings
242  //RWH modify: fPdgC = fCurEntry->pdg;
243 
244  // Check neutrino pdg against declared list of neutrino species declared
245  // by the current instance of the NuMI neutrino flux driver.
246  // No undeclared neutrino species will be accepted at this point as GENIE
247  // has already been configured to handle the specified list.
248  // Make sure that the appropriate list of flux neutrino species was set at
249  // initialization via GSimpleNtpFlux::SetFluxParticles(const PDGCodeList &)
250 
251  // update the # POTs, sum of weights & number of neutrinos
252  // do this HERE (before rejecting flavors that users might be weeding out)
253  // in order to keep the POT accounting correct. This allows one to get
254  // the right normalization for generating only events from the intrinsic
255  // nu_e entries.
257  fSumWeight += this->Weight();
258  fNNeutrinos++;
259 
261  /// user might modify list via SetFluxParticles() in order to reject certain
262  /// flavors, even if they're found in the file. So don't make a big fuss.
263  /// Spit out a single message and then stop reporting that flavor as problematic.
264  int badpdg = fCurEntry->pdg;
265  if ( ! fPdgCListRej->ExistsInPDGCodeList(badpdg) ) {
266  fPdgCListRej->push_back(badpdg);
267  LOG("Flux", pWARN)
268  << "Encountered neutrino specie (" << badpdg
269  << ") that wasn't in SetFluxParticles() list, "
270  << "\nDeclared list of neutrino species: " << *fPdgCList;
271  }
272  return false;
273  }
274 
275  }
276 
277  // Update the curr neutrino p4/x4 lorentz vector
278  double Ev = fCurEntry->E;
279  fP4.SetPxPyPzE(fCurEntry->px,fCurEntry->py,fCurEntry->pz,Ev);
280  double t0 = ((fIncludeVtxt) ? fCurEntry->vtxt : 0 );
281  fX4.SetXYZT(fCurEntry->vtxx,fCurEntry->vtxy,fCurEntry->vtxz,t0);
282 
283  fWeight = fCurEntry->wgt;
284 
285  if (Ev > fMaxEv) {
286  LOG("Flux", pWARN)
287  << "Flux neutrino energy exceeds declared maximum neutrino energy"
288  << "\nEv = " << Ev << "(> Ev{max} = " << fMaxEv << ")";
289  }
290 
291  // if desired, move to user specified user coord z
292  if ( TMath::Abs(fZ0) < 1.0e30 ) this->MoveToZ0(fZ0);
293 
294 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
295  LOG("Flux", pINFO)
296  << "Generated neutrino: " << fIEntry
297 #ifdef NOT_YET
298  << " run " << fCurNuMI->evtno
299  << " evtno " << fCurNuMI->evtno
300  << " entryno " << fCurNuMI->entryno
301 #endif
302  << "\n pdg-code: " << fCurEntry->pdg << " wgt " << fCurEntry->wgt
303  << "\n p4: " << utils::print::P4AsShortString(&fP4)
304  << "\n x4: " << utils::print::X4AsString(&fX4);
305 #endif
306  if ( Ev > fMaxEv ) {
307  LOG("Flux", pFATAL)
308  << "Generated neutrino had E_nu = " << Ev << " > " << fMaxEv
309  << " maximum ";
310  assert(0);
311  }
312 
313  return true;
314 }
Double_t E
energy in lab frame
string P4AsShortString(const TLorentzVector *p)
Definition: PrintUtils.cxx:45
long int fNNeutrinos
number of flux neutrinos thrown so far
GSimpleNtpMeta * fCurMeta
current meta data
#define pERROR
Definition: Messenger.h:59
long int fIUse
current # of times an entry has been used
Double_t px
x momentum in lab frame (GeV)
double fSumWeight
sum of weights for nus thrown so far
Double_t vtxy
y position in lab frame
#define pFATAL
Definition: Messenger.h:56
bool ExistsInPDGCodeList(int pdg_code) const
void MoveToZ0(double z0)
move ray origin to user coord Z0
Double_t vtxt
time of ray start (seconds)
PDGCodeList * fPdgCList
list of neutrino pdg-codes to generate
bool fIncludeVtxt
does fX4 include CurEntry.vtxt or 0
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
TLorentzVector fX4
reconstituted position vector
bool fAllFilesMeta
do all files in chain have meta data
Double_t vtxz
z position in lab frame
#define pINFO
Definition: Messenger.h:62
long int fNEntriesUsed
number of entries read from files
long int fNCycles
times to cycle through the ntuple(s)
GSimpleNtpNuMI * fCurNuMI
current &quot;numi&quot; branch extra info
#define pWARN
Definition: Messenger.h:60
PDGCodeList * fPdgCListRej
list of nu pdg-codes seen but rejected
double fWeight
current neutrino weight
TChain * fNuMetaTree
TTree // REF ONLY.
GSimpleNtpEntry * fCurEntry
current entry
TLorentzVector fP4
reconstituted p4 vector
UInt_t metakey
index key to tie to individual entries
Long64_t fIEntry
current flux ntuple entry
Long64_t fNEntries
number of flux ntuple entries
Double_t vtxx
x position in lab frame (meters)
Double_t pz
z momentum in lab frame
double fEffPOTsPerNu
what a entry is worth ...
bool fEnd
end condition reached
string X4AsString(const TLorentzVector *x)
Definition: PrintUtils.cxx:57
double fMaxEv
maximum energy
UInt_t metakey
key to meta data
Double_t py
y momentum in lab frame
double Weight(void)
returns the flux neutrino weight (if any)
double fAccumPOTs
POTs used so far.
void push_back(int pdg_code)
Definition: PDGCodeList.cxx:58
TChain * fNuFluxTree
TTree // REF ONLY.
#define pDEBUG
Definition: Messenger.h:63
void GSimpleNtpFlux::GenerateWeighted ( bool  gen_weighted)
virtual

set whether to generate weighted or unweighted neutrinos

Implements genie::GFluxI.

Definition at line 700 of file GSimpleNtpFlux.cxx.

References fGenWeighted.

701 {
702 // Set whether to generate weighted rays
703 //
704  fGenWeighted = gen_weighted;
705 }
bool fGenWeighted
does GenerateNext() give weights?
void GSimpleNtpFlux::GetBranchInfo ( std::vector< std::string > &  branchNames,
std::vector< std::string > &  branchClassNames,
std::vector< void ** > &  branchObjPointers 
)
virtual

allow caller to copy current status / ntuple entry info in the output file by providing copies of internal info

Assumes that branch object pointers will not change which may require either a copy be made or, if using the class directly for reading the branch, one must force ROOT to not autodelete: myns::MyClassType* fCurrMyClass = new myns::MyClassType; myTree->SetBranchAddress("bname",&fCurMyClass); //? TBranch* b = myTree->GetBranch("bname"); //? b->SetAutoDelete(false);

ensure vectors are sized sufficiently (or use .push_back()) branchNames[i] = "bname" branchClassNames[i] = "myns::MyClassType" branchObjPointers[i] = (void**)

Reimplemented from genie::flux::GFluxFileConfigI.

Definition at line 564 of file GSimpleNtpFlux.cxx.

References fCurAux, fCurEntry, and fCurNuMI.

567 {
568  // allow flux driver to report back current status and/or ntuple entry
569  // info for possible recording in the output file by supplying
570  // the class name, and a pointer to the object that will be filled
571  // as well as a suggested name for the branch.
572 
573  if ( fCurEntry ) {
574  branchNames.push_back("simple");
575  branchClassNames.push_back("genie::flux::GSimpleNtpEntry");
576  branchObjPointers.push_back((void**)&fCurEntry);
577  }
578 
579  if ( fCurNuMI ) {
580  branchNames.push_back("numi");
581  branchClassNames.push_back("genie::flux::GSimpleNtpNuMI");
582  branchObjPointers.push_back((void**)&fCurNuMI);
583  }
584 
585  if ( fCurAux ) {
586  branchNames.push_back("aux");
587  branchClassNames.push_back("genie::flux::GSimpleNtpAux");
588  branchObjPointers.push_back((void**)&fCurAux);
589  }
590 }
GSimpleNtpAux * fCurAux
current &quot;aux&quot; branch extra info
GSimpleNtpNuMI * fCurNuMI
current &quot;numi&quot; branch extra info
GSimpleNtpEntry * fCurEntry
current entry
const genie::flux::GSimpleNtpAux* genie::flux::GSimpleNtpFlux::GetCurrentAux ( void  )
inline

GSimpleNtpAux.

Definition at line 238 of file GSimpleNtpFlux.h.

References fCurAux.

const genie::flux::GSimpleNtpEntry* genie::flux::GSimpleNtpFlux::GetCurrentEntry ( void  )
inline

GSimpleNtpEntry.

Definition at line 234 of file GSimpleNtpFlux.h.

References fCurEntry.

Referenced by GenerateNext().

const genie::flux::GSimpleNtpMeta* genie::flux::GSimpleNtpFlux::GetCurrentMeta ( void  )
inline

GSimpleNtpMeta.

Definition at line 240 of file GSimpleNtpFlux.h.

References fCurMeta.

const genie::flux::GSimpleNtpNuMI* genie::flux::GSimpleNtpFlux::GetCurrentNuMI ( void  )
inline

GSimpleNtpNuMI.

Definition at line 236 of file GSimpleNtpFlux.h.

References fCurNuMI.

double GSimpleNtpFlux::GetDecayDist ( ) const

dist (user units) from dk to current pos

Definition at line 316 of file GSimpleNtpFlux.cxx.

References genie::flux::GSimpleNtpEntry::dist, and fCurEntry.

Referenced by genie::flux::GFluxBlender::GenerateNext().

317 {
318  // return distance (user units) between dk point and start position
319  // these are in beam units
320  return fCurEntry->dist;
321 }
GSimpleNtpEntry * fCurEntry
current entry
Double_t dist
distance from hadron decay
std::vector< std::string > GSimpleNtpFlux::GetFileList ( )

list of files currently part of chain

Definition at line 1120 of file GSimpleNtpFlux.cxx.

References fNuFluxTree.

Referenced by PrintConfig().

1121 {
1122  std::vector<std::string> flist;
1123  TObjArray *fileElements=fNuFluxTree->GetListOfFiles();
1124  TIter next(fileElements);
1125  TChainElement *chEl=0;
1126  while (( chEl=(TChainElement*)next() )) {
1127  flist.push_back(chEl->GetTitle());
1128  }
1129  return flist;
1130 }
TChain * fNuFluxTree
TTree // REF ONLY.
TChain* genie::flux::GSimpleNtpFlux::GetFluxTChain ( void  )
inline

Definition at line 244 of file GSimpleNtpFlux.h.

References fNuFluxTree.

void GSimpleNtpFlux::GetFluxWindow ( TVector3 &  p1,
TVector3 &  p2,
TVector3 &  p3 
) const

3 points define a plane in beam coordinate

Definition at line 663 of file GSimpleNtpFlux.cxx.

References fCurMeta, genie::flux::GSimpleNtpMeta::windowBase, genie::flux::GSimpleNtpMeta::windowDir1, and genie::flux::GSimpleNtpMeta::windowDir2.

664 {
665  // return flux window points
666 
667  p0 = TVector3(fCurMeta->windowBase[0],
668  fCurMeta->windowBase[1],
669  fCurMeta->windowBase[2]);
670  p1 = p0 + TVector3(fCurMeta->windowDir1[0],
671  fCurMeta->windowDir1[1],
672  fCurMeta->windowDir1[2]);
673  p2 = p0 + TVector3(fCurMeta->windowDir2[0],
674  fCurMeta->windowDir2[1],
675  fCurMeta->windowDir2[2]);
676 
677 }
GSimpleNtpMeta * fCurMeta
current meta data
Double_t windowDir1[3]
dx,dy,dz of window direction 1
Double_t windowDir2[3]
dx,dy,dz of window direction 2
Double_t windowBase[3]
x,y,z position of window base point
bool genie::flux::GSimpleNtpFlux::GetIncludeVtxt ( )
inline

should X4 include CurEntry.vtxt

Definition at line 250 of file GSimpleNtpFlux.h.

References fIncludeVtxt.

250 { return fIncludeVtxt; }
bool fIncludeVtxt
does fX4 include CurEntry.vtxt or 0
TTree * GSimpleNtpFlux::GetMetaDataTree ( )
virtual

Reimplemented from genie::flux::GFluxFileConfigI.

Definition at line 591 of file GSimpleNtpFlux.cxx.

References fNuMetaTree.

591 { return fNuMetaTree; }
TChain * fNuMetaTree
TTree // REF ONLY.
double GSimpleNtpFlux::GetTotalExposure ( ) const
virtual

GFluxExposureI interface.

Implements genie::flux::GFluxExposureI.

Definition at line 78 of file GSimpleNtpFlux.cxx.

References UsedPOTs().

79 {
80  // complete the GFluxExposureI interface
81  return UsedPOTs();
82 }
double UsedPOTs(void) const
of protons-on-target used
long int genie::flux::GSimpleNtpFlux::Index ( void  )
inlinevirtual

returns corresponding index for current flux neutrino (e.g. for a flux ntuple returns the current entry number)

Implements genie::GFluxI.

Definition at line 220 of file GSimpleNtpFlux.h.

References fIEntry.

220 { return fIEntry; }
Long64_t fIEntry
current flux ntuple entry
void GSimpleNtpFlux::Initialize ( void  )
private

Definition at line 707 of file GSimpleNtpFlux.cxx.

References fAccumPOTs, fAllFilesMeta, fAlreadyUnwgt, fCurAux, fCurAuxCopy, fCurEntry, fCurEntryCopy, fCurMeta, fCurNuMI, fCurNuMICopy, fEffPOTsPerNu, fEnd, fFilePOTs, fGenWeighted, genie::flux::GFluxFileConfigI::fICycle, fIEntry, fIFileNumber, fIncludeVtxt, fIUse, fMaxEv, fMaxWeight, fNEntries, fNEntriesUsed, fNFiles, fNNeutrinos, fNuFluxBranchRequest, fNuFluxTree, fNuMetaTree, fNUse, fSumWeight, LOG, pINFO, ResetCurrent(), and SetDefaults().

708 {
709  LOG("Flux", pINFO) << "Initializing GSimpleNtpFlux driver";
710 
711  fMaxEv = 0;
712  fEnd = false;
713 
715  fCurNuMI = new GSimpleNtpNuMI;
716  fCurAux = new GSimpleNtpAux;
717  fCurMeta = new GSimpleNtpMeta;
718 
719  fCurEntryCopy = 0;
720  fCurNuMICopy = 0;
721  fCurAuxCopy = 0;
722 
723  fNuFluxTree = new TChain("flux");
724  fNuMetaTree = new TChain("meta");
725 
726  //fNuFluxFilePatterns = "";
727  fNuFluxBranchRequest = "entry,numi,aux"; // all branches
728 
729  fNFiles = 0;
730 
731  fNEntries = 0;
732  fIEntry = -1;
733  fIFileNumber = 0;
734  fICycle = 0;
735  fNUse = 1;
736  fIUse = 999999;
737 
738  fFilePOTs = 0;
739 
740  fMaxWeight = -1;
741 
742  fSumWeight = 0;
743  fNNeutrinos = 0;
744  fNEntriesUsed = 0;
745  fEffPOTsPerNu = 0;
746  fAccumPOTs = 0;
747 
748  fGenWeighted = false;
749  fAllFilesMeta = true;
750  fAlreadyUnwgt = false;
751 
752  fIncludeVtxt = true;
753 
754  this->SetDefaults();
755  this->ResetCurrent();
756 }
GSimpleNtpAux * fCurAux
current &quot;aux&quot; branch extra info
Int_t fIFileNumber
which file for the current entry
long int fNNeutrinos
number of flux neutrinos thrown so far
GSimpleNtpMeta * fCurMeta
current meta data
long int fIUse
current # of times an entry has been used
double fSumWeight
sum of weights for nus thrown so far
bool fIncludeVtxt
does fX4 include CurEntry.vtxt or 0
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
bool fAllFilesMeta
do all files in chain have meta data
Double_t fFilePOTs
of protons-on-target represented by all files
bool fGenWeighted
does GenerateNext() give weights?
GSimpleNtpNuMI * fCurNuMICopy
current &quot;numi&quot; branch extra info
#define pINFO
Definition: Messenger.h:62
long int fNEntriesUsed
number of entries read from files
GSimpleNtpNuMI * fCurNuMI
current &quot;numi&quot; branch extra info
int fNFiles
number of files in chain
double fMaxWeight
max flux neutrino weight in input file
TChain * fNuMetaTree
TTree // REF ONLY.
GSimpleNtpEntry * fCurEntry
current entry
bool fAlreadyUnwgt
are input files already unweighted
Long64_t fIEntry
current flux ntuple entry
Long64_t fNEntries
number of flux ntuple entries
GSimpleNtpAux * fCurAuxCopy
current &quot;aux&quot; branch extra info
long int fNUse
how often to use same entry in a row
GSimpleNtpEntry * fCurEntryCopy
current entry
double fEffPOTsPerNu
what a entry is worth ...
string fNuFluxBranchRequest
list of requested branches &quot;entry,numi,au&quot;
bool fEnd
end condition reached
double fMaxEv
maximum energy
double fAccumPOTs
POTs used so far.
TChain * fNuFluxTree
TTree // REF ONLY.
void GSimpleNtpFlux::LoadBeamSimData ( const std::vector< string > &  filenames,
const std::string &  det_loc 
)
virtual

Definition at line 385 of file GSimpleNtpFlux.cxx.

References AddFile(), CalcEffPOTsPerNu(), fAccumPOTs, fCurAux, fCurEntry, fCurNuMI, genie::flux::GFluxFileConfigI::fICycle, fIEntry, fIUse, fMaxWeight, fNEntries, fNEntriesUsed, fNNeutrinos, fNuFluxFilePatterns, fNuFluxTree, fSumWeight, genie::RandomGen::Instance(), LOG, OptionalAttachBranch(), pDEBUG, pERROR, pFATAL, pINFO, pNOTICE, ProcessMeta(), and genie::RandomGen::RndFlux().

387 {
388 // Loads a beam simulation root file into the GSimpleNtpFlux driver.
389 
390  fNuFluxFilePatterns = patterns;
391  std::vector<int> nfiles_from_pattern;
392 
393  // create a (sorted) set of file names
394  // this also helps ensure that the same file isn't listed multiple times
395  std::set<std::string> fnames;
396 
397  LOG("Flux", pINFO) << "LoadBeamSimData was passed " << patterns.size()
398  << " patterns";
399 
400  for (size_t ipatt = 0; ipatt < patterns.size(); ++ipatt ) {
401  string pattern = patterns[ipatt];
402  nfiles_from_pattern.push_back(0);
403  LOG("Flux", pINFO)
404  << "Loading flux tree from ROOT file pattern ["
405  << std::setw(3) << ipatt << "] \"" << pattern << "\"";
406 
407  // !WILDCARD only works for file name ... NOT directory
408  string dirname = gSystem->UnixPathName(gSystem->WorkingDirectory());
409  size_t slashpos = pattern.find_last_of("/");
410  size_t fbegin;
411  if ( slashpos != std::string::npos ) {
412  dirname = pattern.substr(0,slashpos);
413  LOG("Flux", pDEBUG) << "Look for flux using directory " << dirname;
414  fbegin = slashpos + 1;
415  } else { fbegin = 0; }
416 
417  void* dirp = gSystem->OpenDirectory(gSystem->ExpandPathName(dirname.c_str()));
418  if ( dirp ) {
419  std::string basename =
420  pattern.substr(fbegin,pattern.size()-fbegin);
421  TRegexp re(basename.c_str(),kTRUE);
422  const char* onefile;
423  while ( ( onefile = gSystem->GetDirEntry(dirp) ) ) {
424  TString afile = onefile;
425  if ( afile=="." || afile==".." ) continue;
426  if ( basename!=afile && afile.Index(re) == kNPOS ) continue;
427  std::string fullname = dirname + "/" + afile.Data();
428  fnames.insert(fullname);
429  nfiles_from_pattern[ipatt]++;
430  }
431  gSystem->FreeDirectory(dirp);
432  } // legal directory
433  } // loop over patterns
434 
435  size_t indx = 0;
436  std::set<string>::const_iterator sitr = fnames.begin();
437  for ( ; sitr != fnames.end(); ++sitr, ++indx) {
438  string filename = *sitr;
439  //std::cout << " [" << std::setw(3) << indx << "] \""
440  // << filename << "\"" << std::endl;
441  bool isok = true;
442  // this next test only works for local files, so we can't do that
443  // if we want to stream via xrootd
444  // ! (gSystem->AccessPathName(filename.c_str()));
445  if ( ! isok ) continue;
446  // open the file to see what it contains
447  LOG("Flux", pINFO) << "Load file " << filename;
448 
449  TFile* tf = TFile::Open(filename.c_str(),"READ");
450  TTree* etree = (TTree*)tf->Get("flux");
451  if ( etree ) {
452  TTree* mtree = (TTree*)tf->Get("meta");
453  // add the file to the chain
454  LOG("Flux", pDEBUG) << "AddFile " << filename
455  << " etree " << etree << " meta " << mtree;
456  this->AddFile(etree,mtree,filename);
457 
458  } // found a GSimpleNtpEntry "flux" tree
459  tf->Close();
460  delete tf;
461  } // loop over sorted file names
462 
463  // this will open all files and read headers!!
464  fNEntries = fNuFluxTree->GetEntries();
465 
466  if ( fNEntries == 0 ) {
467  LOG("Flux", pERROR)
468  << "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!";
469  LOG("Flux", pERROR)
470  << "Loaded flux tree contains " << fNEntries << " entries";
471  LOG("Flux", pERROR)
472  << "Was passed the file patterns: ";
473  for (size_t ipatt = 0; ipatt < patterns.size(); ++ipatt ) {
474  string pattern = patterns[ipatt];
475  LOG("Flux", pERROR)
476  << " [" << std::setw(3) << ipatt <<"] " << pattern;
477  }
478  LOG("Flux", pERROR)
479  << "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!";
480  } else {
481  LOG("Flux", pNOTICE)
482  << "Loaded flux tree contains " << fNEntries << " entries"
483  << " from " << fnames.size() << " unique files";
484  for (size_t ipatt = 0; ipatt < patterns.size(); ++ipatt ) {
485  string pattern = patterns[ipatt];
486  LOG("Flux", pINFO)
487  << " pattern: " << pattern << " yielded "
488  << nfiles_from_pattern[ipatt] << " files";
489  }
490  }
491 
492  int sba_status[3] = { -999, -999, -999 };
493  // "entry" branch isn't optional ... contains the neutrino info
494 #if ROOT_VERSION_CODE >= ROOT_VERSION(5,26,0)
495  sba_status[0] =
496 #endif
497  fNuFluxTree->SetBranchAddress("entry",&fCurEntry);
498 #if ROOT_VERSION_CODE >= ROOT_VERSION(5,26,0)
499  if ( sba_status[0] < 0 ) {
500  LOG("Flux", pFATAL)
501  << "flux chain has no \"entry\" branch " << sba_status[0];
502  assert(0);
503  }
504 #endif
505  //TBranch* bentry = fNuFluxTree->GetBranch("entry");
506  //bentry->SetAutoDelete(false);
507 
508  if ( OptionalAttachBranch("numi") ) {
509 #if ROOT_VERSION_CODE >= ROOT_VERSION(5,26,0)
510  sba_status[1] =
511 #endif
512  fNuFluxTree->SetBranchAddress("numi",&fCurNuMI);
513  //TBranch* bnumi = fNuFluxTree->GetBranch("numi");
514  //bnumi->SetAutoDelete(false);
515  } else { delete fCurNuMI; fCurNuMI = 0; }
516 
517  if ( OptionalAttachBranch("aux") ) {
518 #if ROOT_VERSION_CODE >= ROOT_VERSION(5,26,0)
519  sba_status[2] =
520 #endif
521  fNuFluxTree->SetBranchAddress("aux",&fCurAux);
522  //TBranch* baux = fNuFluxTree->GetBranch("aux");
523  //baux->SetAutoDelete(false);
524  } else { delete fCurAux; fCurAux = 0; }
525 
526  LOG("Flux", pDEBUG)
527  << " SetBranchAddress status: "
528  << " \"entry\"=" << sba_status[0]
529  << " \"numi\"=" << sba_status[1]
530  << " \"aux\"=" << sba_status[2];
531 
532  // attach requested branches
533 
534  if (fMaxWeight<=0) {
535  LOG("Flux", pDEBUG)
536  << "Run ProcessMeta() as part of LoadBeamSimData";
537  this->ProcessMeta();
538  }
539 
540  // current ntuple cycle # (flux ntuples may be recycled)
541  fICycle = 0;
542  // pick a starting entry index [0:fNEntries-1]
543  // pretend we just used up the the previous one
545  fIUse = 9999999;
546  fIEntry = rnd->RndFlux().Integer(fNEntries) - 1;
547  if ( config.find("no-offset-index") != string::npos ) {
548  LOG("Flux",pINFO) << "Config saw \"no-offset-index\"";
549  fIEntry = -1;
550  }
551  LOG("Flux",pINFO) << "Start with entry fIEntry=" << fIEntry;
552 
553  // don't count things we used to estimate max weight
554  fNEntriesUsed = 0;
555  fSumWeight = 0;
556  fNNeutrinos = 0;
557  fAccumPOTs = 0;
558 
559  LOG("Flux",pDEBUG) << "about to CalcEffPOTsPerNu";
560  this->CalcEffPOTsPerNu();
561 
562 }
GSimpleNtpAux * fCurAux
current &quot;aux&quot; branch extra info
std::vector< string > fNuFluxFilePatterns
(potentially wildcarded) path(s)
long int fNNeutrinos
number of flux neutrinos thrown so far
void AddFile(TTree *fluxtree, TTree *metatree, string fname)
#define pERROR
Definition: Messenger.h:59
long int fIUse
current # of times an entry has been used
static RandomGen * Instance()
Access instance.
Definition: RandomGen.cxx:71
double fSumWeight
sum of weights for nus thrown so far
#define pFATAL
Definition: Messenger.h:56
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition: RandomGen.h:29
void ProcessMeta(void)
scan for max flux energy, weight
bool OptionalAttachBranch(std::string bname)
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pINFO
Definition: Messenger.h:62
long int fNEntriesUsed
number of entries read from files
GSimpleNtpNuMI * fCurNuMI
current &quot;numi&quot; branch extra info
double fMaxWeight
max flux neutrino weight in input file
GSimpleNtpEntry * fCurEntry
current entry
Long64_t fIEntry
current flux ntuple entry
Long64_t fNEntries
number of flux ntuple entries
#define pNOTICE
Definition: Messenger.h:61
TRandom3 & RndFlux(void) const
rnd number generator used by flux drivers
Definition: RandomGen.h:71
double fAccumPOTs
POTs used so far.
TChain * fNuFluxTree
TTree // REF ONLY.
#define pDEBUG
Definition: Messenger.h:63
double genie::flux::GSimpleNtpFlux::MaxEnergy ( void  )
inlinevirtual

declare the max flux neutrino energy that can be generated (for init. purposes)

Implements genie::GFluxI.

Definition at line 213 of file GSimpleNtpFlux.h.

References fMaxEv.

213 { return fMaxEv; }
double fMaxEv
maximum energy
const TLorentzVector& genie::flux::GSimpleNtpFlux::Momentum ( void  )
inlinevirtual

returns the flux neutrino 4-momentum

Implements genie::GFluxI.

Definition at line 217 of file GSimpleNtpFlux.h.

References fP4.

217 { return fP4; }
TLorentzVector fP4
reconstituted p4 vector
void GSimpleNtpFlux::MoveToZ0 ( double  z0)

move ray origin to user coord Z0

Definition at line 323 of file GSimpleNtpFlux.cxx.

References e, fP4, fX4, genie::constants::kLightSpeed, LOG, genie::units::meter, pWARN, and genie::units::second.

Referenced by GenerateNext_weighted().

324 {
325  // move ray origin to specified user z0
326  // move beam coord entry correspondingly
327 
328  double pzusr = fP4.Pz();
329  if ( TMath::Abs(pzusr) < 1.0e-30 ) {
330  // neutrino is moving almost entirely in x-y plane
331  LOG("Flux", pWARN)
332  << "MoveToZ0(" << z0usr << ") not possible due to pz_usr (" << pzusr << ")";
333  return;
334  }
335 
336  // Find the 3-position shift needed to move to the requested z coordinate
337  double dz = z0usr - fX4.Z();
338  double scale = dz / pzusr;
339 
340  // LOG("Flux",pDEBUG)
341  // << "MoveToZ0: before x4=(" << fX4.X() << "," << fX4.Y() << "," << fX4.Z()
342  // << "," << fX4.T() << ") z0=" << z0usr << " pzusr=" << pzusr
343  // << " p4=(" << fP4.Px() << "," << fP4.Py() << "," << fP4.Pz() << "," << fP4.E() << ")";
344 
345  TVector3 dx3( scale*fP4.Vect() );
346 
347  // Find the corresponding time shift
348  double v = fP4.Beta() * constants::kLightSpeed
349  / ( units::meter / units::second );
350  double dt = fP4.P() * dz / ( pzusr * v );
351 
352  // Apply these shifts to update the neutrino 4-position
353  TLorentzVector dx4( dx3, dt );
354  fX4 += dx4;
355 
356  // LOG("Flux",pDEBUG)
357  // << "MoveToZ0: after x4=(" << fX4.X() << "," << fX4.Y() << "," << fX4.Z()
358  // << "," << fX4.T() << ")"
359 }
static constexpr double second
Definition: Units.h:37
const double e
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
TLorentzVector fX4
reconstituted position vector
#define pWARN
Definition: Messenger.h:60
TLorentzVector fP4
reconstituted p4 vector
static constexpr double meter
Definition: Units.h:35
long int genie::flux::GSimpleNtpFlux::NEntriesUsed ( void  ) const
inline

number of entries read from files

Definition at line 260 of file GSimpleNtpFlux.h.

References fNEntriesUsed.

long int GSimpleNtpFlux::NFluxNeutrinos ( void  ) const
virtual

of rays generated

< number of flux neutrinos looped so far

Implements genie::flux::GFluxExposureI.

Definition at line 84 of file GSimpleNtpFlux.cxx.

References fNNeutrinos.

85 {
86  ///< number of flux neutrinos looped so far
87  return fNNeutrinos;
88 }
long int fNNeutrinos
number of flux neutrinos thrown so far
bool GSimpleNtpFlux::OptionalAttachBranch ( std::string  bname)
private

Definition at line 827 of file GSimpleNtpFlux.cxx.

References fNuFluxBranchRequest, fNuFluxTree, LOG, and pINFO.

Referenced by LoadBeamSimData().

828 {
829 
830  if ( fNuFluxBranchRequest.find(name) == string::npos ) {
831  LOG("Flux", pINFO)
832  << "no request for \"" << name <<"\" branch ";
833  return false;
834  }
835 
836  if ( ( fNuFluxTree->GetBranch(name.c_str()) ) ) return true;
837 
838  LOG("Flux", pINFO)
839  << "no \"" << name << "\" branch in the \"flux\" tree";
840  return false;
841 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pINFO
Definition: Messenger.h:62
const char * name
string fNuFluxBranchRequest
list of requested branches &quot;entry,numi,au&quot;
TChain * fNuFluxTree
TTree // REF ONLY.
int genie::flux::GSimpleNtpFlux::PdgCode ( void  )
inlinevirtual

returns the flux neutrino pdg code

Implements genie::GFluxI.

Definition at line 215 of file GSimpleNtpFlux.h.

References fCurEntry, and genie::flux::GSimpleNtpEntry::pdg.

215 { return fCurEntry->pdg; }
GSimpleNtpEntry * fCurEntry
current entry
const TLorentzVector& genie::flux::GSimpleNtpFlux::Position ( void  )
inlinevirtual

returns the flux neutrino 4-position (note: expect SI rather than physical units)

Implements genie::GFluxI.

Definition at line 218 of file GSimpleNtpFlux.h.

References fX4.

218 { return fX4; }
TLorentzVector fX4
reconstituted position vector
void GSimpleNtpFlux::PrintConfig ( void  )
virtual

print the current configuration

Implements genie::flux::GFluxFileConfigI.

Definition at line 1076 of file GSimpleNtpFlux.cxx.

References fAccumPOTs, fAllFilesMeta, fAlreadyUnwgt, fEffPOTsPerNu, fFilePOTs, fGenWeighted, genie::flux::GFluxFileConfigI::fICycle, fIEntry, fIUse, fMaxEv, fMaxWeight, genie::flux::GFluxFileConfigI::fNCycles, fNEntries, fNEntriesUsed, fNFiles, fNNeutrinos, fNuFluxFilePatterns, fNUse, genie::flux::GFluxFileConfigI::fPdgCList, genie::flux::GFluxFileConfigI::fPdgCListRej, fSumWeight, genie::flux::GFluxFileConfigI::fZ0, GetFileList(), LOG, pNOTICE, and genie::units::s.

1077 {
1078 
1079  std::ostringstream s;
1080  PDGCodeList::const_iterator itr = fPdgCList->begin();
1081  for ( ; itr != fPdgCList->end(); ++itr) s << (*itr) << " ";
1082  s << "[rejected: ";
1083  itr = fPdgCListRej->begin();
1084  for ( ; itr != fPdgCListRej->end(); ++itr) s << (*itr) << " ";
1085  s << " ] ";
1086 
1087  std::ostringstream fpattout;
1088  for (size_t i = 0; i < fNuFluxFilePatterns.size(); ++i)
1089  fpattout << "\n [" << std::setw(3) << i << "] " << fNuFluxFilePatterns[i];
1090 
1091  std::ostringstream flistout;
1092  std::vector<std::string> flist = GetFileList();
1093  for (size_t i = 0; i < flist.size(); ++i)
1094  flistout << "\n [" << std::setw(3) << i << "] " << flist[i];
1095 
1096  LOG("Flux", pNOTICE)
1097  << "GSimpleNtpFlux Config:"
1098  << "\n Enu_max " << fMaxEv
1099  << "\n pdg-codes: " << s.str() << "\n "
1100  << "\"flux\" " << fNEntries << " entries, "
1101  << "\"meta\" " << fNFiles << " entries"
1102  << " (FilePOTs " << fFilePOTs << ") in files:"
1103  << flistout.str()
1104  << "\n from file patterns: "
1105  << fpattout.str()
1106  << "\n wgt max=" << fMaxWeight
1107  << "\n Z0 pushback " << fZ0
1108  << "\n used entry " << fIEntry << " " << fIUse << "/" << fNUse
1109  << " times, in " << fICycle << "/" << fNCycles << " cycles"
1110  << "\n SumWeight " << fSumWeight << " for " << fNNeutrinos << " neutrinos"
1111  << " with " << fNEntriesUsed << " entries read"
1112  << "\n EffPOTsPerNu " << fEffPOTsPerNu << " AccumPOTs " << fAccumPOTs
1113  << "\n GenWeighted \"" << (fGenWeighted?"true":"false") << "\""
1114  << " AlreadyUnwgt \"" << (fAlreadyUnwgt?"true":"false") << "\""
1115  << " AllFilesMeta \"" << (fAllFilesMeta?"true":"false") << "\"";
1116 
1117 }
std::vector< string > fNuFluxFilePatterns
(potentially wildcarded) path(s)
long int fNNeutrinos
number of flux neutrinos thrown so far
long int fIUse
current # of times an entry has been used
double fSumWeight
sum of weights for nus thrown so far
static constexpr double s
Definition: Units.h:95
PDGCodeList * fPdgCList
list of neutrino pdg-codes to generate
std::vector< std::string > GetFileList()
list of files currently part of chain
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
bool fAllFilesMeta
do all files in chain have meta data
Double_t fFilePOTs
of protons-on-target represented by all files
bool fGenWeighted
does GenerateNext() give weights?
long int fNEntriesUsed
number of entries read from files
long int fNCycles
times to cycle through the ntuple(s)
int fNFiles
number of files in chain
double fMaxWeight
max flux neutrino weight in input file
PDGCodeList * fPdgCListRej
list of nu pdg-codes seen but rejected
bool fAlreadyUnwgt
are input files already unweighted
Long64_t fIEntry
current flux ntuple entry
Long64_t fNEntries
number of flux ntuple entries
long int fNUse
how often to use same entry in a row
double fEffPOTsPerNu
what a entry is worth ...
#define pNOTICE
Definition: Messenger.h:61
double fMaxEv
maximum energy
double fAccumPOTs
POTs used so far.
void GSimpleNtpFlux::PrintCurrent ( void  )

print current entry from leaves

Definition at line 680 of file GSimpleNtpFlux.cxx.

References fCurEntry, LOG, and pNOTICE.

681 {
682  LOG("Flux", pNOTICE) << "CurrentEntry:" << *fCurEntry;
683 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
GSimpleNtpEntry * fCurEntry
current entry
#define pNOTICE
Definition: Messenger.h:61
void GSimpleNtpFlux::ProcessMeta ( void  )

scan for max flux energy, weight

Definition at line 594 of file GSimpleNtpFlux.cxx.

References fAllFilesMeta, fAlreadyUnwgt, fCurMeta, fFilePOTs, fIFileNumber, fMaxEv, fMaxWeight, fNFiles, fNuMetaTree, genie::flux::GFluxFileConfigI::fPdgCList, LOG, pDEBUG, pFATAL, pNOTICE, genie::PDGCodeList::push_back(), genie::flux::GSimpleNtpMeta::Reset(), and SetMaxEnergy().

Referenced by LoadBeamSimData().

595 {
596 
597  fAlreadyUnwgt = false;
598  fFilePOTs = 0;
599  double minwgt = +1.0e10;
600  double maxwgt = -1.0e10;
601  double maxenu = 0.0;
602 
603  // PDGLibrary* pdglib = PDGLibrary::Instance(); // get initialized now
604 
605  if ( fAllFilesMeta ) {
606  fNuMetaTree->SetBranchAddress("meta",&fCurMeta);
607 #ifdef USE_INDEX_FOR_META
608  int nindices = fNuMetaTree->BuildIndex("metakey"); // key used to tie entries to meta data
609  LOG("Flux", pDEBUG) << "ProcessMeta() BuildIndex nindices " << nindices;
610 #endif
611  int nmeta = fNuMetaTree->GetEntries();
612  for (int imeta = 0; imeta < nmeta; ++imeta ) {
613  fNuMetaTree->GetEntry(imeta);
614 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
615  LOG("Flux", pNOTICE) << "ProcessMeta() ifile " << imeta
616  << " (of " << fNFiles
617  << ") " << *fCurMeta;
618 #endif
619  minwgt = TMath::Min(minwgt,fCurMeta->minWgt);
620  maxwgt = TMath::Max(maxwgt,fCurMeta->maxWgt);
621  maxenu = TMath::Max(maxenu,fCurMeta->maxEnergy);
622  fFilePOTs += fCurMeta->protons;
623  for (size_t i = 0; i < fCurMeta->pdglist.size(); ++i)
624  fPdgCList->push_back(fCurMeta->pdglist[i]);
625  }
626  if ( minwgt == 1.0 && maxwgt == 1.0 ) fAlreadyUnwgt = true;
627  fMaxWeight = maxwgt;
628  this->SetMaxEnergy(maxenu);
629 
630  } else {
631  //
632  LOG("Flux", pFATAL) << "ProcessMeta() not all files have metadata";
633  // for now PUNT ... eventually could scan all the entries
634  }
635 
636 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
637  LOG("Flux", pNOTICE) << "ProcessMeta() Maximum flux weight = " << fMaxWeight
638  << ", energy = " << fMaxEv
639  << ", AlreadyUnwgt=" << (fAlreadyUnwgt?"true":"false");
640 #endif
641 
642  fCurMeta->Reset();
643  fIFileNumber = -999;
644 
645 }
Int_t fIFileNumber
which file for the current entry
GSimpleNtpMeta * fCurMeta
current meta data
#define pFATAL
Definition: Messenger.h:56
PDGCodeList * fPdgCList
list of neutrino pdg-codes to generate
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
bool fAllFilesMeta
do all files in chain have meta data
Double_t fFilePOTs
of protons-on-target represented by all files
void SetMaxEnergy(double Ev)
specify maximum flx neutrino energy
int fNFiles
number of files in chain
double fMaxWeight
max flux neutrino weight in input file
TChain * fNuMetaTree
TTree // REF ONLY.
bool fAlreadyUnwgt
are input files already unweighted
#define pNOTICE
Definition: Messenger.h:61
double fMaxEv
maximum energy
void push_back(int pdg_code)
Definition: PDGCodeList.cxx:58
#define pDEBUG
Definition: Messenger.h:63
void GSimpleNtpFlux::ResetCurrent ( void  )
private

Definition at line 771 of file GSimpleNtpFlux.cxx.

References fCurAux, fCurEntry, fCurNuMI, genie::flux::GSimpleNtpEntry::Reset(), genie::flux::GSimpleNtpNuMI::Reset(), and genie::flux::GSimpleNtpAux::Reset().

Referenced by GenerateNext_weighted(), and Initialize().

772 {
773 // reset running values of neutrino pdg-code, 4-position & 4-momentum
774 // and the input ntuple leaves
775 
776  if (fCurEntry) fCurEntry->Reset();
777  if (fCurNuMI) fCurNuMI->Reset();
778  if (fCurAux) fCurAux->Reset();
779  //allow caching//if (fCurMeta) fCurMeta->Reset();
780 }
GSimpleNtpAux * fCurAux
current &quot;aux&quot; branch extra info
GSimpleNtpNuMI * fCurNuMI
current &quot;numi&quot; branch extra info
GSimpleNtpEntry * fCurEntry
current entry
void genie::flux::GSimpleNtpFlux::ScanMeta ( void  )
private
void GSimpleNtpFlux::SetDefaults ( void  )
private

Definition at line 758 of file GSimpleNtpFlux.cxx.

References LOG, pINFO, SetEntryReuse(), genie::flux::GFluxFileConfigI::SetNumOfCycles(), and genie::flux::GFluxFileConfigI::SetUpstreamZ().

Referenced by Initialize().

759 {
760 
761 // - Set the default flux neutrino start z position at use flux window
762 // - Set number of cycles to 1
763 
764  LOG("Flux", pINFO) << "Setting default GSimpleNtpFlux driver options";
765 
766  this->SetUpstreamZ (-3.4e38); // way upstream ==> use flux window
767  this->SetNumOfCycles (0);
768  this->SetEntryReuse (1);
769 }
void SetEntryReuse(long int nuse=1)
of times to use entry before moving to next
virtual void SetNumOfCycles(long int ncycle)
limit cycling through input files
virtual void SetUpstreamZ(double z0)
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pINFO
Definition: Messenger.h:62
void GSimpleNtpFlux::SetEntryReuse ( long int  nuse = 1)

of times to use entry before moving to next

Definition at line 655 of file GSimpleNtpFlux.cxx.

References fNUse.

Referenced by SetDefaults().

656 {
657 // With nuse > 1 then the same entry in the file is used "nuse" times
658 // before moving on to the next entry in the ntuple
659 
660  fNUse = TMath::Max(1L, nuse);
661 }
long int fNUse
how often to use same entry in a row
void genie::flux::GSimpleNtpFlux::SetGenWeighted ( bool  genwgt = false)
inline

toggle whether GenerateNext() returns weight=1 flux (initial default false)

Definition at line 287 of file GSimpleNtpFlux.h.

References fGenWeighted.

void genie::flux::GSimpleNtpFlux::SetIncludeVtxt ( bool  it = true)
inline

Definition at line 249 of file GSimpleNtpFlux.h.

References fIncludeVtxt.

249 { fIncludeVtxt=it; }; ///< should X4 include CurEntry.vtxt
bool fIncludeVtxt
does fX4 include CurEntry.vtxt or 0
void GSimpleNtpFlux::SetMaxEnergy ( double  Ev)

specify maximum flx neutrino energy

Definition at line 647 of file GSimpleNtpFlux.cxx.

References fMaxEv, LOG, and pINFO.

Referenced by ProcessMeta().

648 {
649  fMaxEv = TMath::Max(0.,Ev);
650 
651  LOG("Flux", pINFO)
652  << "Declared maximum flux neutrino energy: " << fMaxEv;
653 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pINFO
Definition: Messenger.h:62
double fMaxEv
maximum energy
void genie::flux::GSimpleNtpFlux::SetRequestedBranchList ( string  blist = "entry,numi,aux")
inline

Definition at line 283 of file GSimpleNtpFlux.h.

References fNuFluxBranchRequest.

283 { fNuFluxBranchRequest = blist; }
string fNuFluxBranchRequest
list of requested branches &quot;entry,numi,au&quot;
double genie::flux::GSimpleNtpFlux::SumWeight ( void  ) const
inline

integrated weight for flux neutrinos looped so far

Definition at line 261 of file GSimpleNtpFlux.h.

References fSumWeight.

double GSimpleNtpFlux::UsedPOTs ( void  ) const

of protons-on-target used

Definition at line 372 of file GSimpleNtpFlux.cxx.

References fAccumPOTs, fNuFluxTree, LOG, and pWARN.

Referenced by GetTotalExposure().

373 {
374 // Compute current number of flux POTs
375 
376  if (!fNuFluxTree) {
377  LOG("Flux", pWARN)
378  << "The flux driver has not been properly configured";
379  return 0;
380  }
381  return fAccumPOTs;
382 }
#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
double fAccumPOTs
POTs used so far.
TChain * fNuFluxTree
TTree // REF ONLY.
double genie::flux::GSimpleNtpFlux::Weight ( void  )
inlinevirtual

returns the flux neutrino weight (if any)

Implements genie::GFluxI.

Definition at line 216 of file GSimpleNtpFlux.h.

References fWeight.

Referenced by GenerateNext(), and GenerateNext_weighted().

216 { return fWeight; }
double fWeight
current neutrino weight

Member Data Documentation

double genie::flux::GSimpleNtpFlux::fAccumPOTs
private

POTs used so far.

Definition at line 335 of file GSimpleNtpFlux.h.

Referenced by Clear(), GenerateNext_weighted(), Initialize(), LoadBeamSimData(), PrintConfig(), and UsedPOTs().

bool genie::flux::GSimpleNtpFlux::fAllFilesMeta
private

do all files in chain have meta data

Definition at line 340 of file GSimpleNtpFlux.h.

Referenced by AddFile(), GenerateNext_weighted(), Initialize(), PrintConfig(), and ProcessMeta().

bool genie::flux::GSimpleNtpFlux::fAlreadyUnwgt
private

are input files already unweighted

Definition at line 338 of file GSimpleNtpFlux.h.

Referenced by GenerateNext(), Initialize(), PrintConfig(), and ProcessMeta().

GSimpleNtpAux* genie::flux::GSimpleNtpFlux::fCurAux
private

current "aux" branch extra info

Definition at line 344 of file GSimpleNtpFlux.h.

Referenced by CleanUp(), GetBranchInfo(), GetCurrentAux(), Initialize(), LoadBeamSimData(), and ResetCurrent().

GSimpleNtpAux* genie::flux::GSimpleNtpFlux::fCurAuxCopy
private

current "aux" branch extra info

Definition at line 351 of file GSimpleNtpFlux.h.

Referenced by Initialize().

GSimpleNtpEntry* genie::flux::GSimpleNtpFlux::fCurEntry
private
GSimpleNtpEntry* genie::flux::GSimpleNtpFlux::fCurEntryCopy
private

current entry

Definition at line 349 of file GSimpleNtpFlux.h.

Referenced by Initialize().

GSimpleNtpMeta* genie::flux::GSimpleNtpFlux::fCurMeta
private

current meta data

Definition at line 347 of file GSimpleNtpFlux.h.

Referenced by CleanUp(), GenerateNext_weighted(), GetCurrentMeta(), GetFluxWindow(), Initialize(), and ProcessMeta().

GSimpleNtpNuMI* genie::flux::GSimpleNtpFlux::fCurNuMI
private

current "numi" branch extra info

Definition at line 343 of file GSimpleNtpFlux.h.

Referenced by CleanUp(), GenerateNext_weighted(), GetBranchInfo(), GetCurrentNuMI(), Initialize(), LoadBeamSimData(), and ResetCurrent().

GSimpleNtpNuMI* genie::flux::GSimpleNtpFlux::fCurNuMICopy
private

current "numi" branch extra info

Definition at line 350 of file GSimpleNtpFlux.h.

Referenced by Initialize().

double genie::flux::GSimpleNtpFlux::fEffPOTsPerNu
private

what a entry is worth ...

Definition at line 334 of file GSimpleNtpFlux.h.

Referenced by CalcEffPOTsPerNu(), GenerateNext_weighted(), Initialize(), and PrintConfig().

bool genie::flux::GSimpleNtpFlux::fEnd
private

end condition reached

Definition at line 312 of file GSimpleNtpFlux.h.

Referenced by End(), GenerateNext_weighted(), and Initialize().

Double_t genie::flux::GSimpleNtpFlux::fFilePOTs
private

of protons-on-target represented by all files

Definition at line 324 of file GSimpleNtpFlux.h.

Referenced by CalcEffPOTsPerNu(), Initialize(), PrintConfig(), and ProcessMeta().

bool genie::flux::GSimpleNtpFlux::fGenWeighted
private

does GenerateNext() give weights?

Definition at line 337 of file GSimpleNtpFlux.h.

Referenced by GenerateNext(), GenerateWeighted(), Initialize(), PrintConfig(), and SetGenWeighted().

Long64_t genie::flux::GSimpleNtpFlux::fIEntry
private

current flux ntuple entry

Definition at line 321 of file GSimpleNtpFlux.h.

Referenced by CleanUp(), GenerateNext_weighted(), Index(), Initialize(), LoadBeamSimData(), and PrintConfig().

Int_t genie::flux::GSimpleNtpFlux::fIFileNumber
private

which file for the current entry

Definition at line 322 of file GSimpleNtpFlux.h.

Referenced by Initialize(), and ProcessMeta().

bool genie::flux::GSimpleNtpFlux::fIncludeVtxt
private

does fX4 include CurEntry.vtxt or 0

Definition at line 353 of file GSimpleNtpFlux.h.

Referenced by GenerateNext_weighted(), GetIncludeVtxt(), Initialize(), and SetIncludeVtxt().

long int genie::flux::GSimpleNtpFlux::fIUse
private

current # of times an entry has been used

Definition at line 330 of file GSimpleNtpFlux.h.

Referenced by CleanUp(), GenerateNext_weighted(), Initialize(), LoadBeamSimData(), and PrintConfig().

double genie::flux::GSimpleNtpFlux::fMaxEv
private

maximum energy

Definition at line 311 of file GSimpleNtpFlux.h.

Referenced by GenerateNext_weighted(), Initialize(), MaxEnergy(), PrintConfig(), ProcessMeta(), and SetMaxEnergy().

double genie::flux::GSimpleNtpFlux::fMaxWeight
private

max flux neutrino weight in input file

Definition at line 327 of file GSimpleNtpFlux.h.

Referenced by CalcEffPOTsPerNu(), GenerateNext(), Initialize(), LoadBeamSimData(), PrintConfig(), and ProcessMeta().

Long64_t genie::flux::GSimpleNtpFlux::fNEntries
private

number of flux ntuple entries

Definition at line 320 of file GSimpleNtpFlux.h.

Referenced by CalcEffPOTsPerNu(), GenerateNext_weighted(), Initialize(), LoadBeamSimData(), and PrintConfig().

long int genie::flux::GSimpleNtpFlux::fNEntriesUsed
private

number of entries read from files

Definition at line 333 of file GSimpleNtpFlux.h.

Referenced by GenerateNext_weighted(), Initialize(), LoadBeamSimData(), NEntriesUsed(), and PrintConfig().

int genie::flux::GSimpleNtpFlux::fNFiles
private

number of files in chain

Definition at line 319 of file GSimpleNtpFlux.h.

Referenced by AddFile(), Initialize(), PrintConfig(), and ProcessMeta().

long int genie::flux::GSimpleNtpFlux::fNNeutrinos
private

number of flux neutrinos thrown so far

Definition at line 332 of file GSimpleNtpFlux.h.

Referenced by Clear(), GenerateNext_weighted(), Initialize(), LoadBeamSimData(), NFluxNeutrinos(), and PrintConfig().

string genie::flux::GSimpleNtpFlux::fNuFluxBranchRequest
private

list of requested branches "entry,numi,au"

Definition at line 315 of file GSimpleNtpFlux.h.

Referenced by Initialize(), OptionalAttachBranch(), and SetRequestedBranchList().

std::vector<string> genie::flux::GSimpleNtpFlux::fNuFluxFilePatterns
private

(potentially wildcarded) path(s)

Definition at line 314 of file GSimpleNtpFlux.h.

Referenced by LoadBeamSimData(), and PrintConfig().

TChain* genie::flux::GSimpleNtpFlux::fNuFluxTree
private
TChain* genie::flux::GSimpleNtpFlux::fNuMetaTree
private

TTree // REF ONLY.

Definition at line 317 of file GSimpleNtpFlux.h.

Referenced by AddFile(), CleanUp(), GenerateNext_weighted(), GetMetaDataTree(), Initialize(), and ProcessMeta().

long int genie::flux::GSimpleNtpFlux::fNUse
private

how often to use same entry in a row

Definition at line 329 of file GSimpleNtpFlux.h.

Referenced by CleanUp(), Initialize(), PrintConfig(), and SetEntryReuse().

TLorentzVector genie::flux::GSimpleNtpFlux::fP4
private

reconstituted p4 vector

Definition at line 345 of file GSimpleNtpFlux.h.

Referenced by GenerateNext(), GenerateNext_weighted(), Momentum(), and MoveToZ0().

double genie::flux::GSimpleNtpFlux::fSumWeight
private

sum of weights for nus thrown so far

Definition at line 331 of file GSimpleNtpFlux.h.

Referenced by Clear(), GenerateNext_weighted(), Initialize(), LoadBeamSimData(), PrintConfig(), and SumWeight().

double genie::flux::GSimpleNtpFlux::fWeight
private

current neutrino weight

Definition at line 326 of file GSimpleNtpFlux.h.

Referenced by GenerateNext(), GenerateNext_weighted(), and Weight().

TLorentzVector genie::flux::GSimpleNtpFlux::fX4
private

reconstituted position vector

Definition at line 346 of file GSimpleNtpFlux.h.

Referenced by GenerateNext(), GenerateNext_weighted(), MoveToZ0(), and Position().


The documentation for this class was generated from the following files: