18#ifndef _SIMPLE_NTP_FLUX_H_
19#define _SIMPLE_NTP_FLUX_H_
27#include <TLorentzVector.h>
28#include <TLorentzRotation.h>
66 void Print(
const Option_t* opt =
"")
const;
101 void Print(
const Option_t* opt =
"")
const;
140 void Print(
const Option_t* opt =
"")
const;
168 void Print(
const Option_t* opt =
"")
const;
221 void Clear (Option_t * opt);
272 const std::string& det_loc);
274 virtual void GetBranchInfo(std::vector<std::string>& branchNames,
275 std::vector<std::string>& branchClassNames,
276 std::vector<void**>& branchObjPointers);
293 void GetFluxWindow(TVector3& p1, TVector3& p2, TVector3& p3)
const;
304 void AddFile (TTree* fluxtree, TTree* metatree,
string fname);
GENIE Interface for user-defined flux classes.
GENIE interface for uniform flux exposure iterface.
PDGCodeList * fPdgCList
list of neutrino pdg-codes to generate
virtual void LoadBeamSimData(const std::vector< std::string > &filenames, const std::string &det_loc)=0
void Print(const Option_t *opt="") const
std::vector< Int_t > auxint
additional ints associated w/ entry
std::vector< Double_t > auxdbl
additional doubles associated w/ entry
friend ostream & operator<<(ostream &stream, const GSimpleNtpAux &info)
Double_t vtxx
x position in lab frame (meters)
Double_t E
energy in lab frame
Double_t dist
distance from hadron decay
Double_t px
x momentum in lab frame (GeV)
Double_t vtxt
time of ray start (seconds)
virtual ~GSimpleNtpEntry()
Double_t vtxz
z position in lab frame
Double_t py
y momentum in lab frame
Double_t vtxy
y position in lab frame
UInt_t metakey
key to meta data
friend ostream & operator<<(ostream &stream, const GSimpleNtpEntry &info)
Double_t pz
z momentum in lab frame
void Print(const Option_t *opt="") const
bool fIncludeVtxt
does fX4 include CurEntry.vtxt or 0
bool fEnd
end condition reached
const genie::flux::GSimpleNtpEntry * GetCurrentEntry(void)
GSimpleNtpEntry.
const TLorentzVector & Momentum(void)
returns the flux neutrino 4-momentum
void Clear(Option_t *opt)
reset state variables based on opt
void SetEntryReuse(long int nuse=1)
long int fNNeutrinos
number of flux neutrinos thrown so far
double GetDecayDist() const
dist (user units) from dk to current pos
GSimpleNtpNuMI * fCurNuMI
current "numi" branch extra info
const genie::flux::GSimpleNtpMeta * GetCurrentMeta(void)
GSimpleNtpMeta.
void GenerateWeighted(bool gen_weighted)
set whether to generate weighted or unweighted neutrinos
void SetGenWeighted(bool genwgt=false)
toggle whether GenerateNext() returns weight=1 flux (initial default false)
void SetIncludeVtxt(bool it=true)
bool GetIncludeVtxt()
should X4 include CurEntry.vtxt
virtual double GetTotalExposure() const
GFluxExposureI interface.
bool fGenWeighted
does GenerateNext() give weights?
double fMaxWeight
max flux neutrino weight in input file
double fEffPOTsPerNu
what a entry is worth ...
void GetFluxWindow(TVector3 &p1, TVector3 &p2, TVector3 &p3) const
3 points define a plane in beam coordinate
GSimpleNtpEntry * fCurEntryCopy
current entry
bool fAlreadyUnwgt
are input files already unweighted
void CalcEffPOTsPerNu(void)
bool OptionalAttachBranch(std::string bname)
virtual void GetBranchInfo(std::vector< std::string > &branchNames, std::vector< std::string > &branchClassNames, std::vector< void ** > &branchObjPointers)
void SetRequestedBranchList(string blist="entry,numi,aux")
TChain * fNuFluxTree
TTree // REF ONLY.
int PdgCode(void)
returns the flux neutrino pdg code
std::vector< std::string > GetFileList()
list of files currently part of chain
long int fNUse
how often to use same entry in a row
bool End(void)
true if no more flux nu's can be thrown (eg reaching end of beam sim ntuples)
GSimpleNtpAux * fCurAuxCopy
current "aux" branch extra info
double SumWeight(void) const
integrated weight for flux neutrinos looped so far
void AddFile(TTree *fluxtree, TTree *metatree, string fname)
TLorentzVector fX4
reconstituted position vector
double fMaxEv
maximum energy
virtual void LoadBeamSimData(const std::vector< string > &filenames, const std::string &det_loc)
double Weight(void)
returns the flux neutrino weight (if any)
long int NEntriesUsed(void) const
number of entries read from files
double MaxEnergy(void)
declare the max flux neutrino energy that can be generated (for init. purposes)
TLorentzVector fP4
reconstituted p4 vector
long int fNEntriesUsed
number of entries read from files
TChain * GetFluxTChain(void)
const genie::flux::GSimpleNtpNuMI * GetCurrentNuMI(void)
GSimpleNtpNuMI.
const genie::flux::GSimpleNtpAux * GetCurrentAux(void)
GSimpleNtpAux.
void ProcessMeta(void)
scan for max flux energy, weight
GSimpleNtpMeta * fCurMeta
current meta data
Int_t fIFileNumber
which file for the current entry
long int Index(void)
returns corresponding index for current flux neutrino (e.g. for a flux ntuple returns the current ent...
virtual long int NFluxNeutrinos() const
long int fIUse
current # of times an entry has been used
std::vector< string > fNuFluxFilePatterns
(potentially wildcarded) path(s)
Long64_t fIEntry
current flux ntuple entry
GSimpleNtpNuMI * fCurNuMICopy
current "numi" branch extra info
GSimpleNtpAux * fCurAux
current "aux" branch extra info
virtual TTree * GetMetaDataTree()
void PrintConfig()
print the current configuration
bool GenerateNext_weighted(void)
double UsedPOTs(void) const
bool fAllFilesMeta
do all files in chain have meta data
GSimpleNtpEntry * fCurEntry
current entry
void MoveToZ0(double z0)
move ray origin to user coord Z0
double fAccumPOTs
POTs used so far.
const TLorentzVector & Position(void)
returns the flux neutrino 4-position (note: expect SI rather than physical units)
void PrintCurrent(void)
print current entry from leaves
bool GenerateNext(void)
generate the next flux neutrino (return false in err)
double fWeight
current neutrino weight
int fNFiles
number of files in chain
double fSumWeight
sum of weights for nus thrown so far
TChain * fNuMetaTree
TTree // REF ONLY.
string fNuFluxBranchRequest
list of requested branches "entry,numi,au"
const PDGCodeList & FluxParticles(void)
declare list of flux neutrinos that can be generated (for init. purposes)
void SetMaxEnergy(double Ev)
specify maximum flx neutrino energy
Long64_t fNEntries
number of flux ntuple entries
Int_t ppmedium
tracking medium where parent was produced
Int_t ptype
parent type (PDG)
virtual ~GSimpleNtpNuMI()
Double_t vx
vertex position of hadron/muon decay
void Print(const Option_t *opt="") const
friend ostream & operator<<(ostream &stream, const GSimpleNtpNuMI &info)
Double_t pdpx
nu parent momentum at time of decay
Int_t tptype
parent particle type at target exit
Double_t pppx
nu parent momentum at production point
Double_t tpx
parent particle px at target exit
ostream & operator<<(ostream &stream, const genie::flux::GJPARCNuFluxPassThroughInfo &info)
THE MAIN GENIE PROJECT NAMESPACE