22 #ifndef _GJPARC_NEUTRINO_FLUX_H_
23 #define _GJPARC_NEUTRINO_FLUX_H_
28 #include <TLorentzVector.h>
46 class GJPARCNuFluxPassThroughInfo;
48 ostream &
operator << (ostream & stream,
const GJPARCNuFluxPassThroughInfo & info);
68 long int Index (
void);
69 void Clear (Option_t * opt);
229 #endif // _GJPARC_NEUTRINO_FLUX_H_
TTree * fNuFluxTree
input flux ntuple
double POT_curravg(void)
current average POT
int fgPdgC
running generated nu pdg-code
double fSumWeightTot1c
total sum of weights for neutrinos to be thrown / cycle
bool LoadBeamSimData(string filename, string det_loc)
load a jnubeam root flux ntuple
bool End(void)
true if no more flux nu's can be thrown (eg reaching end of beam sim ntuples)
bool fLoadedNeutrino
set to true when GenerateNext_weighted has been called successfully
void SetUpstreamZ(double z0)
set flux neutrino initial z position (upstream of the detector)
int DLocName2Id(string name)
string fDetLoc
input detector location ('sk','nd1','nd2',...)
TChain * fNuFluxSumChain
input summary ntuple for flux file. This tree is only present for later flux versions > 10a ...
virtual ~GJPARCNuFluxPassThroughInfo()
double POT_1cycle(void)
flux POT per cycle
bool GenerateNext_weighted(void)
bool fNuFluxUsingTree
are we using a TTree or a TChain to view the input flux file?
const TLorentzVector & Position(void)
returns the flux neutrino 4-position (note: expect SI rather than physical units) ...
int PdgCode(void)
returns the flux neutrino pdg code
double fFilePOT
file POT normalization, typically 1E+21
TChain * fNuFluxChain
input flux ntuple
double Weight(void)
returns the flux neutrino weight (if any)
void GenerateWeighted(bool gen_weighted=true)
set whether to generate weighted or unweighted neutrinos
long int fOffset
start looping at entry fOffset
long int fNNeutrinos
number of flux neutrinos thrown so far
void SetFluxParticles(const PDGCodeList &particles)
specify list of flux neutrino species
GJPARCNuFluxPassThroughInfo()
TTree * fNuFluxSumTree
input summary ntuple for flux file. This tree is only present for later flux versions > 10a ...
A GENIE flux driver encapsulating the JPARC neutrino flux. It reads-in the official JPARC neutrino fl...
void SetMaxEnergy(double Ev)
specify maximum flx neutrino energy
ostream & operator<<(ostream &stream, const genie::flux::GJPARCNuFluxPassThroughInfo &info)
long int fNNeutrinosTot1c
total number of flux neutrinos to be thrown / cycle
long int NFluxNeutrinos(void) const
number of flux neutrinos looped so far
TLorentzVector fgP4
running generated nu 4-momentum
double MaxEnergy(void)
declare the max flux neutrino energy that can be generated (for init. purposes)
int fDetLocId
input detector location id (fDetLoc -> jnubeam idfd)
double SumWeight(void) const
intergated weight for flux neutrinos looped so far
void SetFilePOT(double pot)
flux file norm is in /N POT/det [ND] or /N POT/cm^2 [FD]. Specify N (typically 1E+21) ...
int fNCycles
how many times to cycle through the flux ntuple
int fNDetLocIdFound
per cycle keep track of the number of fDetLocId are found - if this is zero will exit job ...
void DisableOffset(void)
switch off random offset, must be called before LoadBeamSimData to have any effect ...
long int fIEntry
current flux ntuple entry
bool fIsNDLoc
input location is a 'near' detector location?
GJPARCNuFluxPassThroughInfo * fPassThroughInfo
bool fGenerateWeighted
generate weighted/deweighted flux neutrinos (default is false)
TFile * fNuFluxFile
input flux file
double fMaxWeight
max flux neutrino weight in input file for the specified detector location
void Clear(Option_t *opt)
reset state variables based on opt
friend ostream & operator<<(ostream &stream, const GJPARCNuFluxPassThroughInfo &info)
long int Index(void)
returns corresponding index for current flux neutrino (e.g. for a flux ntuple returns the current ent...
bool fIsFDLoc
input location is a 'far' detector location?
double fZ0
configurable starting z position for each flux neutrino (in detector coord system) ...
PDGCodeList * fPdgCList
list of neutrino pdg-codes
bool GenerateNext(void)
generate the next flux neutrino (return false in err)
long int fNEntries
number of flux ntuple entries
void SetNumOfCycles(int n)
set how many times to cycle through the ntuple (default: 1 / n=0 means 'infinite') ...
TLorentzVector fgX4
running generated nu 4-position
bool fUseRandomOffset
whether set random starting point when looping over flux ntuples
void RandomOffset(void)
choose a random offset as starting entry in flux ntuple
const TLorentzVector & Momentum(void)
returns the flux neutrino 4-momentum
double fNorm
current flux ntuple normalisation
double fSumWeight
sum of weights for neutrinos thrown so far
long int fEntriesThisCycle
keep track of number of entries used so far for this cycle
const GJPARCNuFluxPassThroughInfo & PassThroughInfo(void)
GJPARCNuFluxPassThroughInfo.
double fMaxEv
maximum energy
const PDGCodeList & FluxParticles(void)
declare list of flux neutrinos that can be generated (for init. purposes)
GENIE Interface for user-defined flux classes.