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::GMCJDriver Class Reference

A GENIE `MC Job Driver'. Can be used for setting up complicated event generation cases involving detailed flux descriptions and detector geometry descriptions. More...

#include <GMCJDriver.h>

Collaboration diagram for genie::GMCJDriver:
Collaboration graph
[legend]

Public Member Functions

 GMCJDriver ()
 
 ~GMCJDriver ()
 
void SetEventGeneratorList (string listname)
 
void SetUnphysEventMask (const TBits &mask)
 
void UseFluxDriver (GFluxI *flux)
 
void UseGeomAnalyzer (GeomAnalyzerI *geom)
 
void UseSplines (bool useLogE=true)
 
bool UseMaxPathLengths (string xml_filename)
 
void KeepOnThrowingFluxNeutrinos (bool keep_on)
 
void SetXSecSplineNbins (int nbins)
 
void SetPmaxLogBinning (void)
 
void SetPmaxNbins (int nbins)
 
void SetPmaxSafetyFactor (double sf)
 
void ForceInteraction (void)
 
void ForceSingleProbScale (void)
 
void PreSelectEvents (bool preselect=true)
 
bool PreCalcFluxProbabilities (void)
 
bool LoadFluxProbabilities (string filename)
 
void SaveFluxProbabilities (string outfilename)
 
void Configure (bool calc_prob_scales=true)
 
EventRecordGenerateEvent (void)
 
double GlobProbScale (void) const
 
long int NFluxNeutrinos (void) const
 
map< int, double > SumFluxIntProbs (void) const
 
const GFluxIFluxDriver (void) const
 
const GeomAnalyzerIGeomAnalyzer (void) const
 
GFluxIFluxDriverPtr (void) const
 
GeomAnalyzerIGeomAnalyzerPtr (void) const
 

Private Member Functions

void InitJob (void)
 
void InitEventGeneration (void)
 
void GetParticleLists (void)
 
void GetMaxPathLengthList (void)
 
void GetMaxFluxEnergy (void)
 
void PopulateEventGenDriverPool (void)
 
void BootstrapXSecSplines (void)
 
void BootstrapXSecSplineSummation (void)
 
void ComputeProbScales (void)
 
EventRecordGenerateEvent1Try (void)
 
bool GenerateFluxNeutrino (void)
 
bool ComputePathLengths (void)
 
double ComputeInteractionProbabilities (bool use_max_path_length)
 
int SelectTargetMaterial (double R)
 
void GenerateEventKinematics (void)
 
void GenerateVertexPosition (void)
 
void ComputeEventProbability (void)
 
double InteractionProbability (double xsec, double pl, int A)
 
double PreGenFluxInteractionProbability (void)
 

Private Attributes

GEVGPoolfGPool
 A pool of GEVGDrivers properly configured event generation drivers / one per init state. More...
 
GFluxIfFluxDriver
 [input] neutrino flux driver More...
 
GeomAnalyzerIfGeomAnalyzer
 [input] detector geometry analyzer More...
 
double fEmax
 [declared by the flux driver] maximum neutrino energy More...
 
PDGCodeList fNuList
 [declared by the flux driver] list of neutrino codes More...
 
PDGCodeList fTgtList
 [declared by the geom driver] list of target codes More...
 
PathLengthList fMaxPathLengths
 [declared by the geom driver] maximum path length list More...
 
PathLengthList fCurPathLengths
 [current] path length list for current flux neutrino More...
 
TLorentzVector fCurVtx
 [current] interaction vertex More...
 
EventRecordfCurEvt
 [current] generated event More...
 
int fSelTgtPdg
 [current] selected target material PDG code More...
 
map< int, double > fCurCumulProbMap
 [current] cummulative interaction probabilities More...
 
double fNFluxNeutrinos
 [current] number of flux nuetrinos fired by the flux driver so far More...
 
int fXSecSplineNbins
 [config] number of bins in energy used in the xsec splines More...
 
bool fPmaxLogBinning
 [config] maximum interaction probability is computed in logarithmic energy bins More...
 
int fPmaxNbins
 [config] number of bins in energy used in the maximum interaction probability More...
 
double fPmaxSafetyFactor
 [config] safety factor to compute the maximum interaction probability More...
 
map< int, TH1D * > fPmax
 [computed at init] interaction probability scale /neutrino /energy for given geometry More...
 
double fGlobPmax
 [computed at init] global interaction probability scale for given flux & geometry More...
 
string fEventGenList
 [config] list of event generators loaded by this driver (what used to be the $GEVGL setting) More...
 
TBits * fUnphysEventMask
 [config] controls whether unphysical events are returned (what used to be the $GUNPHYSMASK setting) More...
 
string fMaxPlXmlFilename
 [config] input file with max density-weighted path lengths for all materials More...
 
bool fUseExtMaxPl
 [config] using external max path length estimate? More...
 
bool fUseSplines
 [config] compute all needed & not-loaded splines at init More...
 
bool fUseLogE
 [config] build splines = f(logE) (rather than f(E)) ? More...
 
bool fKeepThrowingFluxNu
 [config] keep firing flux neutrinos till one of them interacts More...
 
bool fGenerateUnweighted
 [config] force single probability scale? More...
 
bool fForceInteraction
 [config] force intearction? More...
 
bool fPreSelect
 [config] set whether to pre-select events using max interaction paths More...
 
TFile * fFluxIntProbFile
 [input] pre-generated flux interaction probability file More...
 
TTree * fFluxIntTree
 [computed-or-loaded] pre-computed flux interaction probabilities (expected tree name is "gFlxIntProbs") More...
 
double fBrFluxIntProb
 flux interaction probability (set to branch:"FluxIntProb") More...
 
int fBrFluxIndex
 corresponding entry in flux input tree (set to address of branch:"FluxEntry") More...
 
double fBrFluxEnu
 corresponding flux P4 (set to address of branch:"FluxP4") More...
 
double fBrFluxWeight
 corresponding flux weight (set to address of branch: "FluxWeight") More...
 
int fBrFluxPDG
 corresponding flux pdg code (set to address of branch: "FluxPDG") More...
 
string fFluxIntFileName
 whether to save pre-generated flux tree for use in later jobs More...
 
string fFluxIntTreeName
 name for tree holding flux probabilities More...
 
map< int, double > fSumFluxIntProbs
 map where the key is flux pdg code and the value is sum of fBrFluxWeight * fBrFluxIntProb for all these flux neutrinos More...
 

Detailed Description

A GENIE `MC Job Driver'. Can be used for setting up complicated event generation cases involving detailed flux descriptions and detector geometry descriptions.

Author
Costas Andreopoulos <c.andreopoulos cern.ch> University of Liverpool
Created:
May 25, 2005
License:
Copyright (c) 2003-2024, The GENIE Collaboration For the full text of the license visit http://copyright.genie-mc.org

Definition at line 46 of file GMCJDriver.h.

Constructor & Destructor Documentation

GMCJDriver::GMCJDriver ( )

Definition at line 43 of file GMCJDriver.cxx.

44 {
45  this->InitJob();
46 }
void InitJob(void)
Definition: GMCJDriver.cxx:448
GMCJDriver::~GMCJDriver ( )

Definition at line 48 of file GMCJDriver.cxx.

49 {
51  if (fGPool) delete fGPool;
52 
53  map<int,TH1D*>::iterator pmax_iter = fPmax.begin();
54  for( ; pmax_iter != fPmax.end(); ++pmax_iter) {
55  TH1D * pmax = pmax_iter->second;
56  if(pmax) {
57  delete pmax; pmax = 0;
58  }
59  }
60  fPmax.clear();
61 
62  if(fFluxIntTree) delete fFluxIntTree;
64 }
TTree * fFluxIntTree
[computed-or-loaded] pre-computed flux interaction probabilities (expected tree name is &quot;gFlxIntProbs...
Definition: GMCJDriver.h:140
TBits * fUnphysEventMask
[config] controls whether unphysical events are returned (what used to be the $GUNPHYSMASK setting) ...
Definition: GMCJDriver.h:130
GEVGPool * fGPool
A pool of GEVGDrivers properly configured event generation drivers / one per init state...
Definition: GMCJDriver.h:110
map< int, TH1D * > fPmax
[computed at init] interaction probability scale /neutrino /energy for given geometry ...
Definition: GMCJDriver.h:127
TFile * fFluxIntProbFile
[input] pre-generated flux interaction probability file
Definition: GMCJDriver.h:139

Member Function Documentation

void GMCJDriver::BootstrapXSecSplines ( void  )
private

Definition at line 601 of file GMCJDriver.cxx.

References genie::InitialState::AsString(), genie::GEVGDriver::CreateSplines(), LOG, pINFO, and pNOTICE.

602 {
603 // Bootstrap cross section spline generation by the event generation drivers
604 // that handle each initial state.
605 
606  if(!fUseSplines) return;
607 
608  LOG("GMCJDriver", pNOTICE)
609  << "Asking event generation drivers to compute all needed xsec splines";
610 
611  PDGCodeList::const_iterator nuiter;
612  PDGCodeList::const_iterator tgtiter;
613  for(nuiter = fNuList.begin(); nuiter != fNuList.end(); ++nuiter){
614  for(tgtiter = fTgtList.begin(); tgtiter != fTgtList.end(); ++tgtiter) {
615  int target_pdgc = *tgtiter;
616  int neutrino_pdgc = *nuiter;
617  InitialState init_state(target_pdgc, neutrino_pdgc);
618  LOG("GMCJDriver", pINFO)
619  << "Computing all splines needed for init-state: "
620  << init_state.AsString();
621  GEVGDriver * evgdriver = fGPool->FindDriver(init_state);
622  evgdriver->CreateSplines(-1,-1,fUseLogE);
623  } // targets
624  } // neutrinos
625  LOG("GMCJDriver", pINFO) << "Finished creating cross section splines\n";
626 }
GEVGPool * fGPool
A pool of GEVGDrivers properly configured event generation drivers / one per init state...
Definition: GMCJDriver.h:110
PDGCodeList fTgtList
[declared by the geom driver] list of target codes
Definition: GMCJDriver.h:115
GEVGDriver * FindDriver(const InitialState &init) const
Definition: GEVGPool.cxx:44
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
bool fUseLogE
[config] build splines = f(logE) (rather than f(E)) ?
Definition: GMCJDriver.h:134
bool fUseSplines
[config] compute all needed &amp; not-loaded splines at init
Definition: GMCJDriver.h:133
GENIE Event Generation Driver. A minimalist user interface object for generating neutrino interaction...
Definition: GEVGDriver.h:54
#define pINFO
Definition: Messenger.h:62
void CreateSplines(int nknots=-1, double emax=-1, bool inLogE=true)
Definition: GEVGDriver.cxx:577
#define pNOTICE
Definition: Messenger.h:61
PDGCodeList fNuList
[declared by the flux driver] list of neutrino codes
Definition: GMCJDriver.h:114
Initial State information.
Definition: InitialState.h:48
void GMCJDriver::BootstrapXSecSplineSummation ( void  )
private

Definition at line 628 of file GMCJDriver.cxx.

References genie::GEVGDriver::CreateXSecSumSpline(), LOG, genie::Range1D_t::max, genie::Range1D_t::min, pFATAL, pNOTICE, and genie::GEVGDriver::ValidEnergyRange().

629 {
630 // Sum-up the cross section splines for all the interaction that can be
631 // simulated for each initial state
632 
633  LOG("GMCJDriver", pNOTICE)
634  << "Summing-up splines to get total cross section for each init state";
635 
636  GEVGPool::iterator diter;
637  for(diter = fGPool->begin(); diter != fGPool->end(); ++diter) {
638  string init_state = diter->first;
639  GEVGDriver * evgdriver = diter->second;
640  assert(evgdriver);
641  LOG("GMCJDriver", pNOTICE)
642  << "**** Summing xsec splines for init-state = " << init_state;
643 
644  Range1D_t rE = evgdriver->ValidEnergyRange();
645  if (fEmax>rE.max || fEmax<rE.min)
646  LOG("GMCJDriver",pFATAL)
647  << " rE (validEnergyRange) [" << rE.min << "," << rE.max << "] "
648  << " fEmax " << fEmax;
649  assert(fEmax<rE.max && fEmax>rE.min);
650 
651  // decide the energy range for the sum spline - extend the spline a little
652  // bit above the maximum beam energy (but below the maximum valid energy)
653  // to avoid the evaluation of the cubic spline around the viscinity of
654  // knots with zero y values (although the GENIE Spline object handles it)
655  double min = rE.min;
656  double max = TMath::Min(rE.max, fEmax);
657 
658  // Because of edge issue (see GENIE docdb 297) these lines are commented out
659  // double dE = fEmax/10.;
660  // double max = (fEmax+dE < rE.max) ? fEmax+dE : rE.max;
661 
662  // in the logaritmic binning is important to have a narrow binning to
663  // describe better the glashow resonance peak.
664  evgdriver->CreateXSecSumSpline(fXSecSplineNbins,min,max,true);
665  }
666  LOG("GMCJDriver", pNOTICE)
667  << "Finished summing all interaction xsec splines per initial state";
668 }
void CreateXSecSumSpline(int nk, double Emin, double Emax, bool inlogE=true)
Definition: GEVGDriver.cxx:440
GEVGPool * fGPool
A pool of GEVGDrivers properly configured event generation drivers / one per init state...
Definition: GMCJDriver.h:110
A simple [min,max] interval for doubles.
Definition: Range1.h:42
#define pFATAL
Definition: Messenger.h:56
Range1D_t ValidEnergyRange(void) const
Definition: GEVGDriver.cxx:674
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
GENIE Event Generation Driver. A minimalist user interface object for generating neutrino interaction...
Definition: GEVGDriver.h:54
int fXSecSplineNbins
[config] number of bins in energy used in the xsec splines
Definition: GMCJDriver.h:123
double max
Definition: Range1.h:53
double min
Definition: Range1.h:52
#define pNOTICE
Definition: Messenger.h:61
double fEmax
[declared by the flux driver] maximum neutrino energy
Definition: GMCJDriver.h:113
void GMCJDriver::ComputeEventProbability ( void  )
private

Definition at line 1267 of file GMCJDriver.cxx.

References genie::units::A, genie::pdg::IonPdgCodeToA(), genie::GHepParticle::P4(), and genie::GHepParticle::Pdg().

1268 {
1269 // Compute event probability for the given flux neutrino & detector geometry
1270 
1271  // get interaction cross section
1272  double xsec = fCurEvt->XSec();
1273 
1274  // get path length in detector along v direction for specified target material
1275  PathLengthList::const_iterator pliter = fCurPathLengths.find(fSelTgtPdg);
1276  double path_length = pliter->second;
1277 
1278  // get target material mass number
1280 
1281  // calculate interaction probability
1282  double P = this->InteractionProbability(xsec, path_length, A);
1283 
1284  //
1285  // get weight for selected event
1286  //
1287 
1288  GHepParticle * nu = fCurEvt->Probe();
1289  int nu_pdg = nu->Pdg();
1290  double Ev = nu->P4()->Energy();
1291 
1292  double weight = 1.0;
1293  if(!fGenerateUnweighted) {
1294  map<int,TH1D*>::const_iterator pmax_iter = fPmax.find(nu_pdg);
1295  assert(pmax_iter != fPmax.end());
1296  TH1D * pmax_hst = pmax_iter->second;
1297  assert(pmax_hst);
1298  double pmax = pmax_hst->GetBinContent(pmax_hst->FindBin(Ev));
1299  assert(pmax>0);
1300  weight = pmax/fGlobPmax;
1301  }
1302 
1303  // set probability & update weight
1304  fCurEvt->SetProbability(P);
1305  fCurEvt->SetWeight(weight * fCurEvt->Weight());
1306 }
virtual void SetWeight(double wght)
Definition: GHepRecord.h:130
const TLorentzVector * P4(void) const
Definition: GHepParticle.h:78
virtual void SetProbability(double prob)
Definition: GHepRecord.h:131
int IonPdgCodeToA(int pdgc)
Definition: PDGUtils.cxx:63
virtual double Weight(void) const
Definition: GHepRecord.h:124
virtual GHepParticle * Probe(void) const
Definition: GHepRecord.cxx:284
int Pdg(void) const
Definition: GHepParticle.h:63
PathLengthList fCurPathLengths
[current] path length list for current flux neutrino
Definition: GMCJDriver.h:117
bool fGenerateUnweighted
[config] force single probability scale?
Definition: GMCJDriver.h:136
map< int, TH1D * > fPmax
[computed at init] interaction probability scale /neutrino /energy for given geometry ...
Definition: GMCJDriver.h:127
EventRecord * fCurEvt
[current] generated event
Definition: GMCJDriver.h:119
int fSelTgtPdg
[current] selected target material PDG code
Definition: GMCJDriver.h:120
static constexpr double A
Definition: Units.h:74
double fGlobPmax
[computed at init] global interaction probability scale for given flux &amp; geometry ...
Definition: GMCJDriver.h:128
double InteractionProbability(double xsec, double pl, int A)
virtual double XSec(void) const
Definition: GHepRecord.h:126
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
double GMCJDriver::ComputeInteractionProbabilities ( bool  use_max_path_length)
private

Definition at line 1110 of file GMCJDriver.cxx.

References genie::units::A, genie::InitialState::AsString(), genie::units::cm2, genie::Spline::Evaluate(), genie::pdg::IonPdgCodeToA(), LOG, pDEBUG, pFATAL, pNOTICE, and genie::GEVGDriver::XSecSumSpline().

1111 {
1112  LOG("GMCJDriver", pNOTICE)
1113  << "Computing relative interaction probabilities for each material";
1114 
1115  // current flux neutrino code & 4-p
1116  int nupdg = fFluxDriver->PdgCode();
1117  const TLorentzVector & nup4 = fFluxDriver->Momentum();
1118 
1119  fCurCumulProbMap.clear();
1120 
1121  const PathLengthList & path_length_list =
1122  (use_max_path_length) ? fMaxPathLengths : fCurPathLengths;
1123 
1124  double probsum=0;
1125  PathLengthList::const_iterator pliter;
1126 
1127  for(pliter = path_length_list.begin();
1128  pliter != path_length_list.end(); ++pliter) {
1129  int mpdg = pliter->first; // material PDG code
1130  double pl = pliter->second; // density x path-length
1131  int A = pdg::IonPdgCodeToA(mpdg);
1132  double xsec = 0.; // sum of xsecs for all modelled processes for given init state
1133  double prob = 0.; // interaction probability
1134  double probn = 0.; // normalized interaction probability
1135 
1136  // find the GEVGDriver object that is handling the current init state
1137  InitialState init_state(mpdg, nupdg);
1138  GEVGDriver * evgdriver = fGPool->FindDriver(init_state);
1139  if(!evgdriver) {
1140  LOG("GMCJDriver", pFATAL)
1141  << "\n * The MC Job driver isn't properly configured!"
1142  << "\n * No event generation driver could be found for init state: "
1143  << init_state.AsString();
1144  exit(1);
1145  }
1146  // compute the interaction xsec and probability (if path-length>0)
1147  if(pl>0.) {
1148  const Spline * totxsecspl = evgdriver->XSecSumSpline();
1149  if(!totxsecspl) {
1150  LOG("GMCJDriver", pFATAL)
1151  << "\n * The MC Job driver isn't properly configured!"
1152  << "\n * Couldn't retrieve total cross section spline for init state: "
1153  << init_state.AsString();
1154  exit(1);
1155  } else {
1156  xsec = totxsecspl->Evaluate( nup4.Energy() );
1157  }
1158  prob = this->InteractionProbability(xsec,pl,A);
1159  LOG("GMCJDriver", pDEBUG)
1160  << " (xsec, pl, A)=(" << xsec << "," << pl << "," << A << ")";
1161 
1162  // scale the interaction probability to the maximum one so as not
1163  // to have to throw few billions of flux neutrinos before getting
1164  // an interaction...
1165  double pmax = 0;
1166  if(fForceInteraction) pmax = 1.;
1167  else if(fGenerateUnweighted) pmax = fGlobPmax;
1168  else {
1169  map<int,TH1D*>::const_iterator pmax_iter = fPmax.find(nupdg);
1170  assert(pmax_iter != fPmax.end());
1171  TH1D * pmax_hst = pmax_iter->second;
1172  assert(pmax_hst);
1173  int ie = pmax_hst->FindBin(nup4.Energy());
1174  pmax = pmax_hst->GetBinContent(ie);
1175  }
1176  assert(pmax>0);
1177  LOG("GMCJDriver", pDEBUG)
1178  << "Pmax=" << pmax;
1179  probn = prob/pmax;
1180  }
1181 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
1182  LOG("GMCJDriver", pNOTICE)
1183  << "tgt: " << mpdg << " -> TotXSec = "
1184  << xsec/units::cm2 << " cm^2, Norm.Prob = " << 100*probn << "%";
1185 #endif
1186 
1187  probsum += probn;
1188  fCurCumulProbMap.insert(map<int,double>::value_type(mpdg,probsum));
1189  }
1190  return probsum;
1191 }
GEVGPool * fGPool
A pool of GEVGDrivers properly configured event generation drivers / one per init state...
Definition: GMCJDriver.h:110
map< int, double > fCurCumulProbMap
[current] cummulative interaction probabilities
Definition: GMCJDriver.h:121
A numeric analysis tool class for interpolating 1-D functions.
Definition: Spline.h:58
#define pFATAL
Definition: Messenger.h:56
int IonPdgCodeToA(int pdgc)
Definition: PDGUtils.cxx:63
GFluxI * fFluxDriver
[input] neutrino flux driver
Definition: GMCJDriver.h:111
const Spline * XSecSumSpline(void) const
Definition: GEVGDriver.h:87
GEVGDriver * FindDriver(const InitialState &init) const
Definition: GEVGPool.cxx:44
double Evaluate(double x) const
Definition: Spline.cxx:363
bool fForceInteraction
[config] force intearction?
Definition: GMCJDriver.h:137
PathLengthList fCurPathLengths
[current] path length list for current flux neutrino
Definition: GMCJDriver.h:117
Object to be filled with the neutrino path-length, for all detector geometry materials, when starting from a position x and travelling along the direction of the neutrino 4-momentum.
bool fGenerateUnweighted
[config] force single probability scale?
Definition: GMCJDriver.h:136
map< int, TH1D * > fPmax
[computed at init] interaction probability scale /neutrino /energy for given geometry ...
Definition: GMCJDriver.h:127
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
static constexpr double A
Definition: Units.h:74
static constexpr double cm2
Definition: Units.h:69
GENIE Event Generation Driver. A minimalist user interface object for generating neutrino interaction...
Definition: GEVGDriver.h:54
double fGlobPmax
[computed at init] global interaction probability scale for given flux &amp; geometry ...
Definition: GMCJDriver.h:128
double InteractionProbability(double xsec, double pl, int A)
virtual int PdgCode(void)=0
returns the flux neutrino pdg code
virtual const TLorentzVector & Momentum(void)=0
returns the flux neutrino 4-momentum
#define pNOTICE
Definition: Messenger.h:61
PathLengthList fMaxPathLengths
[declared by the geom driver] maximum path length list
Definition: GMCJDriver.h:116
Initial State information.
Definition: InitialState.h:48
#define pDEBUG
Definition: Messenger.h:63
bool GMCJDriver::ComputePathLengths ( void  )
private

Definition at line 1081 of file GMCJDriver.cxx.

References LOG, pFATAL, and pNOTICE.

1082 {
1083 // Ask the geometry driver to compute (pathLength x density x weight frac.)
1084 // for all detector materials for the neutrino generated by the flux driver
1085 // and make sure that things look ok...
1086 
1087  fCurPathLengths.clear();
1088 
1089  const TLorentzVector & nup4 = fFluxDriver -> Momentum ();
1090  const TLorentzVector & nux4 = fFluxDriver -> Position ();
1091 
1093 
1094  LOG("GMCJDriver", pNOTICE) << fCurPathLengths;
1095 
1096  if(fCurPathLengths.size() == 0) {
1097  LOG("GMCJDriver", pFATAL)
1098  << "\n *** Geometry driver error ***"
1099  << "\n Got an empty PathLengthList - No material found in geometry?";
1100  return false;
1101  }
1102 
1103  if(fCurPathLengths.AreAllZero()) {
1104  LOG("GMCJDriver", pNOTICE)
1105  << "current flux v doesn't cross any geometry material...";
1106  }
1107  return true;
1108 }
bool AreAllZero(void) const
#define pFATAL
Definition: Messenger.h:56
GFluxI * fFluxDriver
[input] neutrino flux driver
Definition: GMCJDriver.h:111
PathLengthList fCurPathLengths
[current] path length list for current flux neutrino
Definition: GMCJDriver.h:117
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
GeomAnalyzerI * fGeomAnalyzer
[input] detector geometry analyzer
Definition: GMCJDriver.h:112
virtual const PathLengthList & ComputePathLengths(const TLorentzVector &x, const TLorentzVector &p)=0
#define pNOTICE
Definition: Messenger.h:61
void GMCJDriver::ComputeProbScales ( void  )
private

Definition at line 670 of file GMCJDriver.cxx.

References genie::units::A, genie::InitialState::AsString(), genie::Spline::Evaluate(), genie::pdg::IonPdgCodeToA(), LOG, pDEBUG, pINFO, pNOTICE, pWARN, and genie::GEVGDriver::XSecSumSpline().

671 {
672 // Computing interaction probability scales.
673 // To minimize the numbers or trials before choosing a neutrino+target init
674 // state for generating an event (note: each trial means selecting a flux
675 // neutrino, navigating it through the detector to compute path lengths,
676 // computing interaction probabilities for each material and so on...)
677 // a set of probability scales can be used for different neutrino species
678 // and at different energy bins.
679 // A global probability scale is also being constructed for keeping the correct
680 // proportions between differect flux neutrino species or flux neutrinos of
681 // different energies.
682 
683  LOG("GMCJDriver", pNOTICE)
684  << "Computing the max. interaction probability (probability scale)";
685 
686  // clean up global probability scale and maximum probabilties per neutrino
687  // type & energy bin
688  {
689  fGlobPmax = 0;
690  map<int,TH1D*>::iterator pmax_iter = fPmax.begin();
691  for( ; pmax_iter != fPmax.end(); ++pmax_iter) {
692  TH1D * pmax = pmax_iter->second;
693  if(pmax) {
694  delete pmax; pmax = 0;
695  }
696  }
697  fPmax.clear();
698  }
699 
700  double * ebins;
701  if (fPmaxLogBinning) {
702  double emin = 0.1;
703  double de = (TMath::Log10(fEmax) - TMath::Log10(emin)) / fPmaxNbins;
704  ebins = new double[fPmaxNbins+1];
705  for(int i=0; i<=fPmaxNbins; i++) ebins[i] = TMath::Power(10., TMath::Log10(emin) + i * de);
706  }
707  else {
708  // for maximum interaction probability vs E /for given geometry/ I will
709  // be using 300 bins up to the maximum energy for the input flux
710  double de = fEmax/fPmaxNbins;//djk june 5, 2013
711  double emin = 0.0;
712  double emax = fEmax + de;
713  fPmaxNbins++;
714  ebins = new double[fPmaxNbins+1];
715  for(int i=0; i<=fPmaxNbins; i++) ebins[i] = emin + i * (emax-emin)/fPmaxNbins;
716  }
717 
718  PDGCodeList::const_iterator nuiter;
719  PDGCodeList::const_iterator tgtiter;
720 
721  // loop over all neutrino types generated by the flux driver
722  for(nuiter = fNuList.begin(); nuiter != fNuList.end(); ++nuiter) {
723  int neutrino_pdgc = *nuiter;
724  TH1D * pmax_hst = new TH1D("pmax_hst",
725  "max interaction probability vs E | geom",fPmaxNbins,ebins);
726  pmax_hst->SetDirectory(0);
727 
728  // loop over energy bins
729  for(int ie = 1; ie <= pmax_hst->GetNbinsX(); ie++) {
730  double EvLow = pmax_hst->GetBinCenter(ie) - 0.5*pmax_hst->GetBinWidth(ie);
731  double EvHigh = pmax_hst->GetBinCenter(ie) + 0.5*pmax_hst->GetBinWidth(ie);
732  //double Ev = pmax_hst->GetBinCenter(ie);
733 
734  // loop over targets in input geometry, form initial state and compute
735  // the sum of maximum interaction probabilities at the current energy bin
736  //
737  for(tgtiter = fTgtList.begin(); tgtiter != fTgtList.end(); ++tgtiter) {
738  int target_pdgc = *tgtiter;
739 
740  InitialState init_state(target_pdgc, neutrino_pdgc);
741 
742  LOG("GMCJDriver", pDEBUG)
743  << "Computing Pmax for init-state: " << init_state.AsString() << " E from " << EvLow << "-" << EvHigh;
744 
745  // get the appropriate driver
746  GEVGDriver * evgdriver = fGPool->FindDriver(init_state);
747 
748  // get xsec sum over all modelled processes for given neutrino+target)
749  double sxsecLow = evgdriver->XSecSumSpline()->Evaluate(EvLow);
750  double sxsecHigh = evgdriver->XSecSumSpline()->Evaluate(EvHigh);
751 
752  // get max{path-length x density}
753  double plmax = fMaxPathLengths.PathLength(target_pdgc);
754 
755  // compute/store the max interaction probabiity (for given energy)
756  int A = pdg::IonPdgCodeToA(target_pdgc);
757  double pmaxLow = this->InteractionProbability(sxsecLow, plmax, A);
758  double pmaxHigh = this->InteractionProbability(sxsecHigh, plmax, A);
759 
760  double pmax = pmaxHigh;
761  if ( pmaxLow > pmaxHigh){
762  pmax = pmaxLow;
763  LOG("GMCJDriver", pWARN)
764  << "Lower energy neutrinos have a higher probability of interacting than those at higher energy."
765  << " pmaxLow(E=" << EvLow << ")=" << pmaxLow << " and " << " pmaxHigh(E=" << EvHigh << ")=" << pmaxHigh;
766  }
767 
768  pmax_hst->SetBinContent(ie, pmax_hst->GetBinContent(ie) + pmax);
769 
770  LOG("GMCJDriver", pDEBUG)
771  << "Pmax[" << init_state.AsString() << ", Ev from " << EvLow << "-" << EvHigh << "] = " << pmax;
772  } // targets
773 
774  pmax_hst->SetBinContent(ie, fPmaxSafetyFactor * pmax_hst->GetBinContent(ie));
775 
776  LOG("GMCJDriver", pINFO)
777  << "Pmax[nu=" << neutrino_pdgc << ", Ev from " << EvLow << "-" << EvHigh << "] = "
778  << pmax_hst->GetBinContent(ie);
779  } // E
780 
781  fPmax.insert(map<int,TH1D*>::value_type(neutrino_pdgc,pmax_hst));
782  } // nu
783 
784  delete ebins;
785 
786  // Compute global probability scale
787  // Sum Probabilities {
788  // all neutrinos, all targets, @ max path length, @ max energy}
789  //
790  {
791  for(nuiter = fNuList.begin(); nuiter != fNuList.end(); ++nuiter) {
792  int neutrino_pdgc = *nuiter;
793  map<int,TH1D*>::const_iterator pmax_iter = fPmax.find(neutrino_pdgc);
794  assert(pmax_iter != fPmax.end());
795  TH1D * pmax_hst = pmax_iter->second;
796  assert(pmax_hst);
797 // double pmax = pmax_hst->GetBinContent(pmax_hst->FindBin(fEmax));
798  double pmax = pmax_hst->GetMaximum();
799  assert(pmax>0);
800 // fGlobPmax += pmax;
801  fGlobPmax = TMath::Max(pmax, fGlobPmax); // ?;
802  }
803  LOG("GMCJDriver", pNOTICE) << "*** Probability scale = " << fGlobPmax;
804  }
805 }
double PathLength(int pdgc) const
GEVGPool * fGPool
A pool of GEVGDrivers properly configured event generation drivers / one per init state...
Definition: GMCJDriver.h:110
PDGCodeList fTgtList
[declared by the geom driver] list of target codes
Definition: GMCJDriver.h:115
int IonPdgCodeToA(int pdgc)
Definition: PDGUtils.cxx:63
const Spline * XSecSumSpline(void) const
Definition: GEVGDriver.h:87
GEVGDriver * FindDriver(const InitialState &init) const
Definition: GEVGPool.cxx:44
double Evaluate(double x) const
Definition: Spline.cxx:363
map< int, TH1D * > fPmax
[computed at init] interaction probability scale /neutrino /energy for given geometry ...
Definition: GMCJDriver.h:127
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
static constexpr double A
Definition: Units.h:74
bool fPmaxLogBinning
[config] maximum interaction probability is computed in logarithmic energy bins
Definition: GMCJDriver.h:124
GENIE Event Generation Driver. A minimalist user interface object for generating neutrino interaction...
Definition: GEVGDriver.h:54
#define pINFO
Definition: Messenger.h:62
double fGlobPmax
[computed at init] global interaction probability scale for given flux &amp; geometry ...
Definition: GMCJDriver.h:128
#define pWARN
Definition: Messenger.h:60
double InteractionProbability(double xsec, double pl, int A)
#define pNOTICE
Definition: Messenger.h:61
double fPmaxSafetyFactor
[config] safety factor to compute the maximum interaction probability
Definition: GMCJDriver.h:126
PathLengthList fMaxPathLengths
[declared by the geom driver] maximum path length list
Definition: GMCJDriver.h:116
int fPmaxNbins
[config] number of bins in energy used in the maximum interaction probability
Definition: GMCJDriver.h:125
PDGCodeList fNuList
[declared by the flux driver] list of neutrino codes
Definition: GMCJDriver.h:114
Initial State information.
Definition: InitialState.h:48
#define pDEBUG
Definition: Messenger.h:63
double fEmax
[declared by the flux driver] maximum neutrino energy
Definition: GMCJDriver.h:113
void GMCJDriver::Configure ( bool  calc_prob_scales = true)

Definition at line 399 of file GMCJDriver.cxx.

References LOG, pNOTICE, and genie::utils::print::PrintFramedMesg().

Referenced by main().

400 {
401  LOG("GMCJDriver", pNOTICE)
402  << utils::print::PrintFramedMesg("Configuring GMCJDriver");
403 
404  // Get the list of neutrino types from the input flux driver and the list
405  // of target materials from the input geometry driver
406  this->GetParticleLists();
407 
408  // Ask the input GFluxI for the max. neutrino energy (to compute Pmax)
409  this->GetMaxFluxEnergy();
410 
411  // Create all possible initial states and for each one initialize,
412  // configure & store an GEVGDriver event generation driver object.
413  // Once an 'initial state' has been selected from the input flux / geom,
414  // the responsibility for generating the neutrino interaction will be
415  // delegated to one of these drivers.
417 
418  // If the user wants to use cross section splines in order to speed things
419  // up, then coordinate spline creation from all GEVGDriver objects pushed
420  // into GEVGPool. This will create all xsec splines needed for all (enabled)
421  // processes that can be simulated involving the particles in the input flux
422  // and geometry.
423  // Spline creation will be skipped for every spline that has been pre-loaded
424  // into the the XSecSplineList.
425  // Once more it is noted that computing cross section splines is a huge
426  // overhead. The user is encouraged to generate them in advance and load
427  // them into the XSecSplineList
428  this->BootstrapXSecSplines();
429 
430  // Create cross section splines describing the total interaction xsec
431  // for a given initial state (Create them by summing all xsec splines
432  // for each possible initial state)
434 
435  if(calc_prob_scales && !fForceInteraction){
436  // Ask the input geometry driver to compute the max. path length for each
437  // material in the list of target materials (or load a precomputed list)
438  this->GetMaxPathLengthList();
439 
440  // Compute the max. interaction probability to scale all interaction
441  // probabilities to be computed by this driver
442  this->ComputeProbScales();
443  }
444  if (fForceInteraction) fGlobPmax = 1.;
445  LOG("GMCJDriver", pNOTICE) << "Finished configuring GMCJDriver\n\n";
446 }
void PopulateEventGenDriverPool(void)
Definition: GMCJDriver.cxx:564
void GetMaxPathLengthList(void)
Definition: GMCJDriver.cxx:536
bool fForceInteraction
[config] force intearction?
Definition: GMCJDriver.h:137
void GetMaxFluxEnergy(void)
Definition: GMCJDriver.cxx:554
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
void BootstrapXSecSplines(void)
Definition: GMCJDriver.cxx:601
double fGlobPmax
[computed at init] global interaction probability scale for given flux &amp; geometry ...
Definition: GMCJDriver.h:128
void ComputeProbScales(void)
Definition: GMCJDriver.cxx:670
string PrintFramedMesg(string mesg, unsigned int nl=1, const char f='*')
Definition: PrintUtils.cxx:164
void BootstrapXSecSplineSummation(void)
Definition: GMCJDriver.cxx:628
#define pNOTICE
Definition: Messenger.h:61
void GetParticleLists(void)
Definition: GMCJDriver.cxx:519
const GFluxI& genie::GMCJDriver::FluxDriver ( void  ) const
inline

Definition at line 81 of file GMCJDriver.h.

References fFluxDriver.

81 { return *fFluxDriver; }
GFluxI * fFluxDriver
[input] neutrino flux driver
Definition: GMCJDriver.h:111
GFluxI* genie::GMCJDriver::FluxDriverPtr ( void  ) const
inline

Definition at line 83 of file GMCJDriver.h.

References fFluxDriver.

83 { return fFluxDriver; }
GFluxI * fFluxDriver
[input] neutrino flux driver
Definition: GMCJDriver.h:111
void GMCJDriver::ForceInteraction ( void  )

Definition at line 162 of file GMCJDriver.cxx.

References LOG, and pNOTICE.

163 {
164 // Force interaction in detector volume. That generates unweighted events.
165 //
166  fForceInteraction = true;
167 
168  LOG("GMCJDriver", pNOTICE)
169  << "GMCJDriver will generate weighted events forcing the interaction. ";
170 }
bool fForceInteraction
[config] force intearction?
Definition: GMCJDriver.h:137
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pNOTICE
Definition: Messenger.h:61
void GMCJDriver::ForceSingleProbScale ( void  )

Definition at line 172 of file GMCJDriver.cxx.

References LOG, and pNOTICE.

Referenced by main().

173 {
174 // Use a single probability scale. That generates unweighted events.
175 // (Note that generating unweighted event kinematics is a different thing)
176 //
177  fGenerateUnweighted = true;
178 
179  LOG("GMCJDriver", pNOTICE)
180  << "GMCJDriver will generate un-weighted events. "
181  << "Note: That does not force unweighted event kinematics!";
182 }
bool fGenerateUnweighted
[config] force single probability scale?
Definition: GMCJDriver.h:136
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pNOTICE
Definition: Messenger.h:61
EventRecord * GMCJDriver::GenerateEvent ( void  )

Definition at line 815 of file GMCJDriver.cxx.

References LOG, pINFO, and pNOTICE.

Referenced by main().

816 {
817  LOG("GMCJDriver", pNOTICE) << "Generating next event...";
818 
819  this->InitEventGeneration();
820 
821  while(1) {
822  bool flux_end = fFluxDriver->End();
823  if(flux_end) {
824  LOG("GMCJDriver", pNOTICE)
825  << "No more neutrinos can be thrown by the flux driver";
826  return 0;
827  }
828 
829  EventRecord * event = this->GenerateEvent1Try();
830  if(event) return event;
831 
832  if(fKeepThrowingFluxNu) {
833  LOG("GMCJDriver", pNOTICE)
834  << "Flux neutrino didn't interact - Trying the next one...";
835  continue;
836  }
837  break;
838  } // (w(1)
839 
840  LOG("GMCJDriver", pINFO) << "Returning NULL event!";
841  return 0;
842 }
void InitEventGeneration(void)
Definition: GMCJDriver.cxx:807
GFluxI * fFluxDriver
[input] neutrino flux driver
Definition: GMCJDriver.h:111
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
bool fKeepThrowingFluxNu
[config] keep firing flux neutrinos till one of them interacts
Definition: GMCJDriver.h:135
#define pINFO
Definition: Messenger.h:62
Generated Event Record. It is a GHepRecord object that can accept / be visited by EventRecordVisitorI...
Definition: EventRecord.h:37
EventRecord * GenerateEvent1Try(void)
Definition: GMCJDriver.cxx:844
virtual bool End(void)=0
true if no more flux nu&#39;s can be thrown (eg reaching end of beam sim ntuples)
#define pNOTICE
Definition: Messenger.h:61
EventRecord * GMCJDriver::GenerateEvent1Try ( void  )
private

Definition at line 844 of file GMCJDriver.cxx.

References genie::gAbortingInErr, genie::RandomGen::Instance(), genie::controls::kASmallNum, LOG, genie::utils::print::P4AsString(), pDEBUG, genie::utils::res::PdgCode(), pERROR, pFATAL, pNOTICE, pWARN, genie::RandomGen::RndEvg(), and genie::utils::print::X4AsString().

845 {
846 // attempt generating a neutrino interaction by firing a single flux neutrino
847 //
848  RandomGen * rnd = RandomGen::Instance();
849 
850  double Pno=0, Psum=0;
851  double R = rnd->RndEvg().Rndm();
852  LOG("GMCJDriver", pDEBUG) << "Rndm [0,1] = " << R;
853 
854  // Generate a neutrino using the input GFluxI & get current pdgc/p4/x4
855  bool flux_ok = this->GenerateFluxNeutrino();
856  if(!flux_ok) {
857  LOG("GMCJDriver", pERROR)
858  << "** Rejecting current flux neutrino (flux driver err)";
859  return 0;
860  }
861 
862  if (fForceInteraction) {
863  bool pl_ok = this->ComputePathLengths();
864  if(!pl_ok) {
865  LOG("GMCJDriver", pFATAL) << "** Cannot calculate path lenths!";
866  exit(1);
867  }
869  LOG("GMCJDriver", pNOTICE)
870  << "** Rejecting current flux neutrino (misses generation volume)";
871  return 0;
872  }
873  Psum = this->ComputeInteractionProbabilities(false);
874  LOG("GMCJDriver", pNOTICE)
875  << "The interaction probability is: " << Psum;
876  R *= Psum;
877  }
878  else {
879  // Compute the interaction probabilities assuming max. path lengths
880  // and decide whether the neutrino would interact --
881  // Many flux neutrinos should be rejected here, drastically reducing
882  // the number of neutrinos that I need to propagate through the
883  // actual detector geometry (this is skipped when using
884  // pre-calculated flux interaction probabilities)
885  if(fPreSelect) {
886  LOG("GMCJDriver", pNOTICE)
887  << "Computing interaction probabilities for max. path lengths";
888 
889  Psum = this->ComputeInteractionProbabilities(true /* <- max PL*/);
890  Pno = 1-Psum;
891  LOG("GMCJDriver", pNOTICE)
892  << "The no-interaction probability (max. path lengths) is: "
893  << 100*Pno << " %";
894  if(Pno<0.) {
895  LOG("GMCJDriver", pFATAL)
896  << "Negative no-interaction probability! (P = " << 100*Pno << " %)"
897  << " Particle E=" << fFluxDriver->Momentum().E() << " type=" << fFluxDriver->PdgCode() << "Psum=" << Psum;
898  gAbortingInErr=true;
899  exit(1);
900  }
901  if(R>=1-Pno) {
902  LOG("GMCJDriver", pNOTICE)
903  << "** Rejecting current flux neutrino";
904  return 0;
905  }
906  } // preselect
907 
908  bool pl_ok = false;
909 
910 
911  // If possible use pre-generated flux neutrino interaction probabilities
912  if(fFluxIntTree){
913  Psum = this->PreGenFluxInteractionProbability();
914  }
915  // Else compute them in the usual manner
916  else {
917  // Compute (pathLength x density x weight fraction) for all materials
918  // in the input geometry, for the neutrino generated by the flux driver
919  pl_ok = this->ComputePathLengths();
920  if(!pl_ok) {
921  LOG("GMCJDriver", pERROR)
922  << "** Rejecting current flux neutrino (err computing path-lengths)";
923  return 0;
924  }
926  LOG("GMCJDriver", pNOTICE)
927  << "** Rejecting current flux neutrino (misses generation volume)";
928  return 0;
929  }
930  Psum = this->ComputeInteractionProbabilities(false /* <- actual PL */);
931  }
932 
933 
934  if(TMath::Abs(Psum) < controls::kASmallNum){
935  LOG("GMCJDriver", pNOTICE)
936  << "** Rejecting current flux neutrino (has null interaction probability)";
937  return 0;
938  }
939 
940  // Now decide whether the current neutrino interacts
941  Pno = 1-Psum;
942  LOG("GMCJDriver", pNOTICE)
943  << "The actual 'no interaction' probability is: " << 100*Pno << " %";
944  if(Pno<0.) {
945  LOG("GMCJDriver", pFATAL)
946  << "Negative no interactin probability! (P = " << 100*Pno << " %)";
947 
948  // print info about what caused the problem
949  int nupdg = fFluxDriver -> PdgCode ();
950  const TLorentzVector & nup4 = fFluxDriver -> Momentum ();
951  const TLorentzVector & nux4 = fFluxDriver -> Position ();
952 
953  LOG("GMCJDriver", pWARN)
954  << "\n [-] Problematic neutrino: "
955  << "\n |----o PDG-code : " << nupdg
956  << "\n |----o 4-momentum : " << utils::print::P4AsString(&nup4)
957  << "\n |----o 4-position : " << utils::print::X4AsString(&nux4)
958  << "\n Emax : " << fEmax;
959 
960  LOG("GMCJDriver", pWARN)
961  << "\n Problematic path lengths:" << fCurPathLengths;
962 
963  LOG("GMCJDriver", pWARN)
964  << "\n Maximum path lengths:" << fMaxPathLengths;
965 
966  exit(1);
967  }
968  if(R>=1-Pno) {
969  LOG("GMCJDriver", pNOTICE)
970  << "** Rejecting current flux neutrino";
971  return 0;
972  }
973 
974  //
975  // The flux neutrino interacts!
976  //
977 
978  // Calculate path lengths for first time and check potential mismatch if
979  // used pre-generated flux interaction probabilities
980  if(fFluxIntTree){
981  pl_ok = this->ComputePathLengths();
982  if(!pl_ok) {
983  LOG("GMCJDriver", pFATAL) << "** Cannot calculate path lenths!";
984  exit(1);
985  }
986  double Psum_curr = this->ComputeInteractionProbabilities(false /* <- actual PL */);
987  bool mismatch = TMath::Abs(Psum-Psum_curr) > controls::kASmallNum;
988  if(mismatch){
989  LOG("GMCJDriver", pFATAL) <<
990  "** Mismatch between pre-calculated and current interaction "<<
991  "probabilities!";
992  exit(1);
993  }
994  }
995  }
996 
997  // Select a target material
998  fSelTgtPdg = this->SelectTargetMaterial(R);
999  if(fSelTgtPdg==0) {
1000  LOG("GMCJDriver", pERROR)
1001  << "** Rejecting current flux neutrino (failed to select tgt!)";
1002  return 0;
1003  }
1004 
1005  // Ask the GEVGDriver object to select and generate an interaction and
1006  // its kinematics for the selected initial state & neutrino 4-momentum
1007  this->GenerateEventKinematics();
1008  if(!fCurEvt) {
1009  LOG("GMCJDriver", pWARN)
1010  << "** Couldn't generate kinematics for selected interaction";
1011  return 0;
1012  }
1013 
1014  // Generate an 'interaction position' in the selected material (in the
1015  // detector coord system), along the direction of nup4 & set it
1016  this->GenerateVertexPosition();
1017 
1018  // Set the event probability (probability for this event to happen given
1019  // the detector setup & the selected flux neutrino)
1020  // Note for users:
1021  // The above probability is stored at GHepRecord::Probability()
1022  // For normalization purposes make sure that you take into account the
1023  // GHepRecord::Weight() -if event generation is weighted-, and
1024  // GFluxI::Weight() -if beam simulation is weighted-.
1025  if(fForceInteraction) {
1026  double weight = 1. - TMath::Exp(-Psum);
1027  fCurEvt->SetProbability(Psum);
1028  fCurEvt->SetWeight(weight * fCurEvt->Weight());
1029  }
1030  else {
1031  this->ComputeEventProbability();
1032  }
1033 
1034  return fCurEvt;
1035 }
TTree * fFluxIntTree
[computed-or-loaded] pre-computed flux interaction probabilities (expected tree name is &quot;gFlxIntProbs...
Definition: GMCJDriver.h:140
virtual void SetWeight(double wght)
Definition: GHepRecord.h:130
bool AreAllZero(void) const
#define pERROR
Definition: Messenger.h:59
static RandomGen * Instance()
Access instance.
Definition: RandomGen.cxx:71
virtual void SetProbability(double prob)
Definition: GHepRecord.h:131
#define pFATAL
Definition: Messenger.h:56
bool fPreSelect
[config] set whether to pre-select events using max interaction paths
Definition: GMCJDriver.h:138
bool GenerateFluxNeutrino(void)
bool ComputePathLengths(void)
GFluxI * fFluxDriver
[input] neutrino flux driver
Definition: GMCJDriver.h:111
string P4AsString(const TLorentzVector *p)
Definition: PrintUtils.cxx:27
virtual double Weight(void) const
Definition: GHepRecord.h:124
bool fForceInteraction
[config] force intearction?
Definition: GMCJDriver.h:137
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition: RandomGen.h:29
void ComputeEventProbability(void)
PathLengthList fCurPathLengths
[current] path length list for current flux neutrino
Definition: GMCJDriver.h:117
EventRecord * fCurEvt
[current] generated event
Definition: GMCJDriver.h:119
int PdgCode(Resonance_t res, int Q)
(resonance id, charge) -&gt; PDG code
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
int fSelTgtPdg
[current] selected target material PDG code
Definition: GMCJDriver.h:120
double PreGenFluxInteractionProbability(void)
static const double kASmallNum
Definition: Controls.h:40
TRandom3 & RndEvg(void) const
rnd number generator used by the event generation drivers
Definition: RandomGen.h:74
double ComputeInteractionProbabilities(bool use_max_path_length)
#define pWARN
Definition: Messenger.h:60
void GenerateEventKinematics(void)
int SelectTargetMaterial(double R)
void GenerateVertexPosition(void)
virtual int PdgCode(void)=0
returns the flux neutrino pdg code
virtual const TLorentzVector & Momentum(void)=0
returns the flux neutrino 4-momentum
#define pNOTICE
Definition: Messenger.h:61
string X4AsString(const TLorentzVector *x)
Definition: PrintUtils.cxx:57
PathLengthList fMaxPathLengths
[declared by the geom driver] maximum path length list
Definition: GMCJDriver.h:116
bool gAbortingInErr
Definition: Messenger.cxx:34
#define pDEBUG
Definition: Messenger.h:63
double fEmax
[declared by the flux driver] maximum neutrino energy
Definition: GMCJDriver.h:113
void GMCJDriver::GenerateEventKinematics ( void  )
private

Definition at line 1215 of file GMCJDriver.cxx.

References genie::InitialState::AsString(), genie::GEVGDriver::GenerateEvent(), LOG, pFATAL, pNOTICE, and genie::GEVGDriver::SetUnphysEventMask().

1216 {
1217  int nupdg = fFluxDriver->PdgCode();
1218  const TLorentzVector & nup4 = fFluxDriver->Momentum();
1219 
1220  // Find the GEVGDriver object that generates interactions for the
1221  // given initial state (neutrino + target)
1222  InitialState init_state(fSelTgtPdg, nupdg);
1223  GEVGDriver * evgdriver = fGPool->FindDriver(init_state);
1224  if(!evgdriver) {
1225  LOG("GMCJDriver", pFATAL)
1226  << "No GEVGDriver object for init state: " << init_state.AsString();
1227  exit(1);
1228  }
1229 
1230  // propagate current unphysical event mask
1231  evgdriver->SetUnphysEventMask(*fUnphysEventMask);
1232 
1233  // Ask the GEVGDriver object to select and generate an interaction for
1234  // the selected initial state & neutrino 4-momentum
1235  LOG("GMCJDriver", pNOTICE)
1236  << "Asking the selected GEVGDriver object to generate an event";
1237  fCurEvt = evgdriver->GenerateEvent(nup4);
1238 }
TBits * fUnphysEventMask
[config] controls whether unphysical events are returned (what used to be the $GUNPHYSMASK setting) ...
Definition: GMCJDriver.h:130
GEVGPool * fGPool
A pool of GEVGDrivers properly configured event generation drivers / one per init state...
Definition: GMCJDriver.h:110
#define pFATAL
Definition: Messenger.h:56
GFluxI * fFluxDriver
[input] neutrino flux driver
Definition: GMCJDriver.h:111
GEVGDriver * FindDriver(const InitialState &init) const
Definition: GEVGPool.cxx:44
EventRecord * fCurEvt
[current] generated event
Definition: GMCJDriver.h:119
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
int fSelTgtPdg
[current] selected target material PDG code
Definition: GMCJDriver.h:120
GENIE Event Generation Driver. A minimalist user interface object for generating neutrino interaction...
Definition: GEVGDriver.h:54
virtual int PdgCode(void)=0
returns the flux neutrino pdg code
virtual const TLorentzVector & Momentum(void)=0
returns the flux neutrino 4-momentum
#define pNOTICE
Definition: Messenger.h:61
void SetUnphysEventMask(const TBits &mask)
Definition: GEVGDriver.cxx:219
Initial State information.
Definition: InitialState.h:48
EventRecord * GenerateEvent(const TLorentzVector &nu4p)
Definition: GEVGDriver.cxx:228
bool GMCJDriver::GenerateFluxNeutrino ( void  )
private

Definition at line 1037 of file GMCJDriver.cxx.

References LOG, genie::utils::print::P4AsString(), genie::utils::res::PdgCode(), pERROR, pFATAL, pNOTICE, and genie::utils::print::X4AsString().

1038 {
1039 // Ask the neutrino flux driver to generate a flux neutrino and make sure
1040 // that things look ok...
1041 //
1042  LOG("GMCJDriver", pNOTICE) << "Generating a flux neutrino";
1043 
1044  bool ok = fFluxDriver->GenerateNext();
1045  if(!ok) {
1046  LOG("GMCJDriver", pERROR)
1047  << "*** The flux driver couldn't generate a flux neutrino!!";
1048  return false;
1049  }
1050 
1051  fNFluxNeutrinos++;
1052  int nupdg = fFluxDriver -> PdgCode ();
1053  const TLorentzVector & nup4 = fFluxDriver -> Momentum ();
1054  const TLorentzVector & nux4 = fFluxDriver -> Position ();
1055 
1056  LOG("GMCJDriver", pNOTICE)
1057  << "\n [-] Generated flux neutrino: "
1058  << "\n |----o PDG-code : " << nupdg
1059  << "\n |----o 4-momentum : " << utils::print::P4AsString(&nup4)
1060  << "\n |----o 4-position : " << utils::print::X4AsString(&nux4);
1061 
1062  if(nup4.Energy() > fEmax) {
1063  LOG("GMCJDriver", pFATAL)
1064  << "\n *** Flux driver error ***"
1065  << "\n Generated flux v with E = " << nup4.Energy() << " GeV"
1066  << "\n Max v energy (declared by flux driver) = " << fEmax << " GeV"
1067  << "\n My interaction probability scaling is invalidated!!";
1068  return false;
1069  }
1070  if(!fNuList.ExistsInPDGCodeList(nupdg)) {
1071  LOG("GMCJDriver", pFATAL)
1072  << "\n *** Flux driver error ***"
1073  << "\n Generated flux v with pdg = " << nupdg
1074  << "\n It does not belong to the declared list of flux neutrinos"
1075  << "\n I was not configured to handle this!!";
1076  return false;
1077  }
1078  return true;
1079 }
#define pERROR
Definition: Messenger.h:59
#define pFATAL
Definition: Messenger.h:56
GFluxI * fFluxDriver
[input] neutrino flux driver
Definition: GMCJDriver.h:111
string P4AsString(const TLorentzVector *p)
Definition: PrintUtils.cxx:27
bool ExistsInPDGCodeList(int pdg_code) const
int PdgCode(Resonance_t res, int Q)
(resonance id, charge) -&gt; PDG code
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
double fNFluxNeutrinos
[current] number of flux nuetrinos fired by the flux driver so far
Definition: GMCJDriver.h:122
virtual bool GenerateNext(void)=0
generate the next flux neutrino (return false in err)
#define pNOTICE
Definition: Messenger.h:61
string X4AsString(const TLorentzVector *x)
Definition: PrintUtils.cxx:57
PDGCodeList fNuList
[declared by the flux driver] list of neutrino codes
Definition: GMCJDriver.h:114
double fEmax
[declared by the flux driver] maximum neutrino energy
Definition: GMCJDriver.h:113
void GMCJDriver::GenerateVertexPosition ( void  )
private

Definition at line 1240 of file GMCJDriver.cxx.

References genie::constants::kLightSpeed, LOG, genie::units::meter, pNOTICE, and genie::units::second.

1241 {
1242  // Generate an 'interaction position' in the selected material, along
1243  // the direction of nup4
1244  LOG("GMCJDriver", pNOTICE)
1245  << "Asking the geometry analyzer to generate a vertex";
1246 
1247  const TLorentzVector & p4 = fFluxDriver->Momentum ();
1248  const TLorentzVector & x4 = fFluxDriver->Position ();
1249 
1250  const TVector3 & vtx = fGeomAnalyzer->GenerateVertex(x4, p4, fSelTgtPdg);
1251 
1252  TVector3 origin(x4.X(), x4.Y(), x4.Z());
1253  origin-=vtx; // computes vector dr = origin - vtx
1254 
1255  double dL = origin.Mag();
1256  double v = p4.Beta() * kLightSpeed /(units::meter/units::second);
1257  double dt = dL/v;
1258 
1259  LOG("GMCJDriver", pNOTICE)
1260  << "|vtx - origin|: dL = " << dL << " m, dt = " << dt << " sec";
1261 
1262  fCurVtx.SetXYZT(vtx.x(), vtx.y(), vtx.z(), x4.T() + dt);
1263 
1265 }
virtual const TLorentzVector & Position(void)=0
returns the flux neutrino 4-position (note: expect SI rather than physical units) ...
GFluxI * fFluxDriver
[input] neutrino flux driver
Definition: GMCJDriver.h:111
virtual const TVector3 & GenerateVertex(const TLorentzVector &x, const TLorentzVector &p, int tgtpdg)=0
static constexpr double second
Definition: Units.h:37
EventRecord * fCurEvt
[current] generated event
Definition: GMCJDriver.h:119
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
int fSelTgtPdg
[current] selected target material PDG code
Definition: GMCJDriver.h:120
TLorentzVector fCurVtx
[current] interaction vertex
Definition: GMCJDriver.h:118
static constexpr double meter
Definition: Units.h:35
virtual void SetVertex(double x, double y, double z, double t)
Definition: GHepRecord.cxx:827
GeomAnalyzerI * fGeomAnalyzer
[input] detector geometry analyzer
Definition: GMCJDriver.h:112
virtual const TLorentzVector & Momentum(void)=0
returns the flux neutrino 4-momentum
#define pNOTICE
Definition: Messenger.h:61
const GeomAnalyzerI& genie::GMCJDriver::GeomAnalyzer ( void  ) const
inline

Definition at line 82 of file GMCJDriver.h.

References fGeomAnalyzer.

82 { return *fGeomAnalyzer; }
GeomAnalyzerI * fGeomAnalyzer
[input] detector geometry analyzer
Definition: GMCJDriver.h:112
GeomAnalyzerI* genie::GMCJDriver::GeomAnalyzerPtr ( void  ) const
inline

Definition at line 84 of file GMCJDriver.h.

References fGeomAnalyzer.

84 { return fGeomAnalyzer; }
GeomAnalyzerI * fGeomAnalyzer
[input] detector geometry analyzer
Definition: GMCJDriver.h:112
void GMCJDriver::GetMaxFluxEnergy ( void  )
private

Definition at line 554 of file GMCJDriver.cxx.

References LOG, and pNOTICE.

555 {
556  LOG("GMCJDriver", pNOTICE)
557  << "Querying the flux driver for the maximum energy of flux neutrinos";
559 
560  LOG("GMCJDriver", pNOTICE)
561  << "Maximum flux neutrino energy = " << fEmax << " GeV";
562 }
virtual double MaxEnergy(void)=0
declare the max flux neutrino energy that can be generated (for init. purposes)
GFluxI * fFluxDriver
[input] neutrino flux driver
Definition: GMCJDriver.h:111
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pNOTICE
Definition: Messenger.h:61
double fEmax
[declared by the flux driver] maximum neutrino energy
Definition: GMCJDriver.h:113
void GMCJDriver::GetMaxPathLengthList ( void  )
private

Definition at line 536 of file GMCJDriver.cxx.

References LOG, and pNOTICE.

537 {
538  if(fUseExtMaxPl) {
539  LOG("GMCJDriver", pNOTICE)
540  << "Loading external max path-length list for input geometry from "
543 
544  } else {
545  LOG("GMCJDriver", pNOTICE)
546  << "Querying the geometry driver to compute the max path-length list";
548  }
549  // Print maximum path lengths & neutrino energy
550  LOG("GMCJDriver", pNOTICE)
551  << "Maximum path length list: " << fMaxPathLengths;
552 }
bool fUseExtMaxPl
[config] using external max path length estimate?
Definition: GMCJDriver.h:132
virtual const PathLengthList & ComputeMaxPathLengths(void)=0
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
XmlParserStatus_t LoadFromXml(string filename)
string fMaxPlXmlFilename
[config] input file with max density-weighted path lengths for all materials
Definition: GMCJDriver.h:131
GeomAnalyzerI * fGeomAnalyzer
[input] detector geometry analyzer
Definition: GMCJDriver.h:112
#define pNOTICE
Definition: Messenger.h:61
PathLengthList fMaxPathLengths
[declared by the geom driver] maximum path length list
Definition: GMCJDriver.h:116
void GMCJDriver::GetParticleLists ( void  )
private

Definition at line 519 of file GMCJDriver.cxx.

References LOG, and pNOTICE.

520 {
521  // Get the list of flux neutrinos from the flux driver
522  LOG("GMCJDriver", pNOTICE)
523  << "Asking the flux driver for its list of neutrinos";
525 
526  LOG("GMCJDriver", pNOTICE) << "Flux particles: " << fNuList;
527 
528  // Get the list of target materials from the geometry driver
529  LOG("GMCJDriver", pNOTICE)
530  << "Asking the geometry driver for its list of targets";
532 
533  LOG("GMCJDriver", pNOTICE) << "Target materials: " << fTgtList;
534 }
PDGCodeList fTgtList
[declared by the geom driver] list of target codes
Definition: GMCJDriver.h:115
GFluxI * fFluxDriver
[input] neutrino flux driver
Definition: GMCJDriver.h:111
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
virtual const PDGCodeList & FluxParticles(void)=0
declare list of flux neutrinos that can be generated (for init. purposes)
GeomAnalyzerI * fGeomAnalyzer
[input] detector geometry analyzer
Definition: GMCJDriver.h:112
#define pNOTICE
Definition: Messenger.h:61
PDGCodeList fNuList
[declared by the flux driver] list of neutrino codes
Definition: GMCJDriver.h:114
virtual const PDGCodeList & ListOfTargetNuclei(void)=0
double genie::GMCJDriver::GlobProbScale ( void  ) const
inline

Definition at line 76 of file GMCJDriver.h.

References fGlobPmax.

Referenced by main().

76 { return fGlobPmax; }
double fGlobPmax
[computed at init] global interaction probability scale for given flux &amp; geometry ...
Definition: GMCJDriver.h:128
void GMCJDriver::InitEventGeneration ( void  )
private

Definition at line 807 of file GMCJDriver.cxx.

808 {
809  fCurPathLengths.clear();
810  fCurEvt = 0;
811  fSelTgtPdg = 0;
812  fCurVtx.SetXYZT(0.,0.,0.,0.);
813 }
PathLengthList fCurPathLengths
[current] path length list for current flux neutrino
Definition: GMCJDriver.h:117
EventRecord * fCurEvt
[current] generated event
Definition: GMCJDriver.h:119
int fSelTgtPdg
[current] selected target material PDG code
Definition: GMCJDriver.h:120
TLorentzVector fCurVtx
[current] interaction vertex
Definition: GMCJDriver.h:118
void GMCJDriver::InitJob ( void  )
private

Definition at line 448 of file GMCJDriver.cxx.

References genie::AlgConfigPool::Instance(), genie::Messenger::Instance(), and genie::GHepFlags::NFlags().

449 {
450  fEventGenList = "Default"; // <-- set of event generators to be loaded by this driver
451 
452  fUnphysEventMask = new TBits(GHepFlags::NFlags()); //<-- unphysical event mask
453  //fUnphysEventMask->ResetAllBits(true);
454  for(unsigned int i = 0; i < GHepFlags::NFlags(); i++) {
455  fUnphysEventMask->SetBitNumber(i, true);
456  }
457 
458  fFluxDriver = 0; // <-- flux driver
459  fGeomAnalyzer = 0; // <-- geometry driver
460  fGPool = 0; // <-- pool of GEVGDriver event generation drivers
461  fEmax = 0; // <-- maximum neutrino energy
462  fMaxPlXmlFilename = ""; // <-- XML file with external path lengths
463  fUseExtMaxPl = false;
464  fUseSplines = false;
465  fNFluxNeutrinos = 0; // <-- number of flux neutrinos thrown so far
466 
467  fXSecSplineNbins = 100; // <-- number of energy bins used in the xsec splines
468  fPmaxLogBinning = false; // <-- maximum interaction probability is computed in logarithmic energy bins
469  fPmaxNbins = 300; // <-- number of energy bins used in the maximum interaction probability
470  fPmaxSafetyFactor = 1.2; // <-- safety factor to compute maximum interaction probability per neutrino & per energy bin
471  fGlobPmax = 0; // <-- maximum interaction probability (global prob scale)
472  fPmax.clear(); // <-- maximum interaction probability per neutrino & per energy bin
473 
474  fForceInteraction = false; // <-- default opt to not force the interaction
475  fGenerateUnweighted = false; // <-- default opt to generate weighted events
476  fPreSelect = true; // <-- default to use pre-selection based on maximum path lengths
477 
478  fSelTgtPdg = 0;
479  fCurEvt = 0;
480  fCurVtx.SetXYZT(0.,0.,0.,0.);
481 
482  fFluxIntProbFile = 0;
483  fFluxIntTreeName = "gFlxIntProb";
484  fFluxIntFileName = "";
485  fFluxIntTree = 0;
486  fBrFluxIntProb = -1.;
487  fBrFluxIndex = -1;
488  fBrFluxEnu = -1.;
489  fBrFluxWeight = -1.;
490  fBrFluxPDG = 0;
491  fSumFluxIntProbs.clear();
492 
493  // Throw as many flux neutrinos as necessary till one has interacted
494  // so that GenerateEvent() never returns NULL (except when in error)
495  this->KeepOnThrowingFluxNeutrinos(true);
496 
497  // Allow the selected GEVGDriver to go into recursive mode and regenerate
498  // an interaction that turns out to be unphysical.
499  //TBits unphysmask(GHepFlags::NFlags());
500  //unphysmask.ResetAllBits(false);
501  //this->FilterUnphysical(unphysmask);
502 
503  // Force early initialization of singleton objects that are typically
504  // would be initialized at their first use later on.
505  // This is purely cosmetic and I do it to send the banner and some prolific
506  // initialization printout at the top.
507  assert( Messenger::Instance() );
508  assert( AlgConfigPool::Instance() );
509 
510  // Clear the target and neutrino lists
511  fNuList.clear();
512  fTgtList.clear();
513 
514  // Clear the maximum path length list
515  fMaxPathLengths.clear();
516  fCurPathLengths.clear();
517 }
TTree * fFluxIntTree
[computed-or-loaded] pre-computed flux interaction probabilities (expected tree name is &quot;gFlxIntProbs...
Definition: GMCJDriver.h:140
TBits * fUnphysEventMask
[config] controls whether unphysical events are returned (what used to be the $GUNPHYSMASK setting) ...
Definition: GMCJDriver.h:130
GEVGPool * fGPool
A pool of GEVGDrivers properly configured event generation drivers / one per init state...
Definition: GMCJDriver.h:110
double fBrFluxIntProb
flux interaction probability (set to branch:&quot;FluxIntProb&quot;)
Definition: GMCJDriver.h:141
PDGCodeList fTgtList
[declared by the geom driver] list of target codes
Definition: GMCJDriver.h:115
bool fPreSelect
[config] set whether to pre-select events using max interaction paths
Definition: GMCJDriver.h:138
GFluxI * fFluxDriver
[input] neutrino flux driver
Definition: GMCJDriver.h:111
bool fForceInteraction
[config] force intearction?
Definition: GMCJDriver.h:137
double fBrFluxEnu
corresponding flux P4 (set to address of branch:&quot;FluxP4&quot;)
Definition: GMCJDriver.h:143
bool fUseExtMaxPl
[config] using external max path length estimate?
Definition: GMCJDriver.h:132
double fBrFluxWeight
corresponding flux weight (set to address of branch: &quot;FluxWeight&quot;)
Definition: GMCJDriver.h:144
string fFluxIntFileName
whether to save pre-generated flux tree for use in later jobs
Definition: GMCJDriver.h:146
PathLengthList fCurPathLengths
[current] path length list for current flux neutrino
Definition: GMCJDriver.h:117
void KeepOnThrowingFluxNeutrinos(bool keep_on)
Definition: GMCJDriver.cxx:120
bool fGenerateUnweighted
[config] force single probability scale?
Definition: GMCJDriver.h:136
map< int, TH1D * > fPmax
[computed at init] interaction probability scale /neutrino /energy for given geometry ...
Definition: GMCJDriver.h:127
EventRecord * fCurEvt
[current] generated event
Definition: GMCJDriver.h:119
static unsigned int NFlags(void)
Definition: GHepFlags.h:76
TFile * fFluxIntProbFile
[input] pre-generated flux interaction probability file
Definition: GMCJDriver.h:139
int fSelTgtPdg
[current] selected target material PDG code
Definition: GMCJDriver.h:120
double fNFluxNeutrinos
[current] number of flux nuetrinos fired by the flux driver so far
Definition: GMCJDriver.h:122
int fBrFluxPDG
corresponding flux pdg code (set to address of branch: &quot;FluxPDG&quot;)
Definition: GMCJDriver.h:145
bool fUseSplines
[config] compute all needed &amp; not-loaded splines at init
Definition: GMCJDriver.h:133
static Messenger * Instance(void)
Definition: Messenger.cxx:49
bool fPmaxLogBinning
[config] maximum interaction probability is computed in logarithmic energy bins
Definition: GMCJDriver.h:124
int fBrFluxIndex
corresponding entry in flux input tree (set to address of branch:&quot;FluxEntry&quot;)
Definition: GMCJDriver.h:142
string fFluxIntTreeName
name for tree holding flux probabilities
Definition: GMCJDriver.h:147
string fEventGenList
[config] list of event generators loaded by this driver (what used to be the $GEVGL setting) ...
Definition: GMCJDriver.h:129
double fGlobPmax
[computed at init] global interaction probability scale for given flux &amp; geometry ...
Definition: GMCJDriver.h:128
TLorentzVector fCurVtx
[current] interaction vertex
Definition: GMCJDriver.h:118
string fMaxPlXmlFilename
[config] input file with max density-weighted path lengths for all materials
Definition: GMCJDriver.h:131
int fXSecSplineNbins
[config] number of bins in energy used in the xsec splines
Definition: GMCJDriver.h:123
map< int, double > fSumFluxIntProbs
map where the key is flux pdg code and the value is sum of fBrFluxWeight * fBrFluxIntProb for all the...
Definition: GMCJDriver.h:148
GeomAnalyzerI * fGeomAnalyzer
[input] detector geometry analyzer
Definition: GMCJDriver.h:112
double fPmaxSafetyFactor
[config] safety factor to compute the maximum interaction probability
Definition: GMCJDriver.h:126
PathLengthList fMaxPathLengths
[declared by the geom driver] maximum path length list
Definition: GMCJDriver.h:116
int fPmaxNbins
[config] number of bins in energy used in the maximum interaction probability
Definition: GMCJDriver.h:125
PDGCodeList fNuList
[declared by the flux driver] list of neutrino codes
Definition: GMCJDriver.h:114
static AlgConfigPool * Instance()
double fEmax
[declared by the flux driver] maximum neutrino energy
Definition: GMCJDriver.h:113
double GMCJDriver::InteractionProbability ( double  xsec,
double  pl,
int  A 
)
private

Definition at line 1308 of file GMCJDriver.cxx.

References genie::units::cm2, genie::units::gram, genie::units::kilogram, genie::constants::kNA, and genie::units::m2.

1309 {
1310 // P = Na (Avogadro number, atoms/mole) *
1311 // 1/A (1/mass number, mole/gr) *
1312 // xsec (total interaction cross section, cm^2) *
1313 // pL (density-weighted path-length, gr/cm^2)
1314 //
1315  xsec = xsec / units::cm2;
1317 
1318  return kNA*(xsec*pL)/A;
1319 }
static constexpr double gram
Definition: Units.h:140
static constexpr double kilogram
Definition: Units.h:36
static constexpr double m2
Definition: Units.h:72
static constexpr double A
Definition: Units.h:74
static constexpr double cm2
Definition: Units.h:69
void GMCJDriver::KeepOnThrowingFluxNeutrinos ( bool  keep_on)

Definition at line 120 of file GMCJDriver.cxx.

References genie::utils::print::BoolAsYNString(), LOG, and pNOTICE.

121 {
122  LOG("GMCJDriver", pNOTICE)
123  << "Keep on throwing flux neutrinos till one interacts? : "
124  << utils::print::BoolAsYNString(keep_on);
125  fKeepThrowingFluxNu = keep_on;
126 }
string BoolAsYNString(bool b)
Definition: PrintUtils.cxx:108
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
bool fKeepThrowingFluxNu
[config] keep firing flux neutrinos till one of them interacts
Definition: GMCJDriver.h:135
#define pNOTICE
Definition: Messenger.h:61
bool GMCJDriver::LoadFluxProbabilities ( string  filename)

Definition at line 343 of file GMCJDriver.cxx.

References LOG, pERROR, pNOTICE, and pWARN.

Referenced by main().

344 {
345 // Load a pre-generated set of flux interaction probabilities from an external
346 // file. This is recommended when using large flux files (>100k entries) as
347 // for these the time to calculate the interaction probabilities can exceed
348 // ~20 minutes. After loading the input tree we call PreCalcFluxProbabilities
349 // to check that has successfully loaded
350 //
351  if(fFluxIntProbFile){
352  LOG("GMCJDriver", pWARN)
353  << "Can't load flux interaction prob file as one is already loaded";
354  return false;
355  }
356 
357  fFluxIntProbFile = new TFile(filename.c_str(), "OPEN");
358 
359  if(fFluxIntProbFile){
360  fFluxIntTree = dynamic_cast<TTree*>(fFluxIntProbFile->Get(fFluxIntTreeName.c_str()));
361  if(fFluxIntTree){
362  bool set_addresses =
363  fFluxIntTree->SetBranchAddress("FluxIntProb", &fBrFluxIntProb) >= 0 &&
364  fFluxIntTree->SetBranchAddress("FluxIndex", &fBrFluxIndex) >= 0 &&
365  fFluxIntTree->SetBranchAddress("FluxPDG", &fBrFluxPDG) >= 0 &&
366  fFluxIntTree->SetBranchAddress("FluxWeight", &fBrFluxWeight) >= 0 &&
367  fFluxIntTree->SetBranchAddress("FluxEnu", &fBrFluxEnu) >= 0;
368  if(set_addresses){
369  // Finally check that can use them
370  if(this->PreCalcFluxProbabilities()) {
371  LOG("GMCJDriver", pNOTICE)
372  << "Successfully loaded pre-generated flux interaction probabilities";
373  return true;
374  }
375  }
376  // If cannot load then delete tree
377  LOG("GMCJDriver", pERROR) <<
378  "Cannot find expected branches in input flux probability tree!";
379  delete fFluxIntTree; fFluxIntTree = 0;
380  }
381  else LOG("GMCJDriver", pERROR)
382  << "Cannot find tree: "<< fFluxIntTreeName.c_str();
383  }
384 
385  LOG("GMCJDriver", pWARN)
386  << "Unable to load flux interaction probabilities file";
387  return false;
388 }
TTree * fFluxIntTree
[computed-or-loaded] pre-computed flux interaction probabilities (expected tree name is &quot;gFlxIntProbs...
Definition: GMCJDriver.h:140
bool PreCalcFluxProbabilities(void)
Definition: GMCJDriver.cxx:192
#define pERROR
Definition: Messenger.h:59
double fBrFluxIntProb
flux interaction probability (set to branch:&quot;FluxIntProb&quot;)
Definition: GMCJDriver.h:141
double fBrFluxEnu
corresponding flux P4 (set to address of branch:&quot;FluxP4&quot;)
Definition: GMCJDriver.h:143
double fBrFluxWeight
corresponding flux weight (set to address of branch: &quot;FluxWeight&quot;)
Definition: GMCJDriver.h:144
TFile * fFluxIntProbFile
[input] pre-generated flux interaction probability file
Definition: GMCJDriver.h:139
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
int fBrFluxPDG
corresponding flux pdg code (set to address of branch: &quot;FluxPDG&quot;)
Definition: GMCJDriver.h:145
int fBrFluxIndex
corresponding entry in flux input tree (set to address of branch:&quot;FluxEntry&quot;)
Definition: GMCJDriver.h:142
string fFluxIntTreeName
name for tree holding flux probabilities
Definition: GMCJDriver.h:147
#define pWARN
Definition: Messenger.h:60
#define pNOTICE
Definition: Messenger.h:61
long int genie::GMCJDriver::NFluxNeutrinos ( void  ) const
inline

Definition at line 77 of file GMCJDriver.h.

References fNFluxNeutrinos.

Referenced by main().

77 { return (long int) fNFluxNeutrinos; }
double fNFluxNeutrinos
[current] number of flux nuetrinos fired by the flux driver so far
Definition: GMCJDriver.h:122
void GMCJDriver::PopulateEventGenDriverPool ( void  )
private

Definition at line 564 of file GMCJDriver.cxx.

References genie::InitialState::AsString(), genie::GEVGDriver::Configure(), LOG, pDEBUG, pNOTICE, genie::GEVGDriver::SetEventGeneratorList(), and genie::GEVGDriver::UseSplines().

565 {
566  LOG("GMCJDriver", pDEBUG)
567  << "Creating GEVGPool & adding a GEVGDriver object per init-state";
568 
569  if (fGPool) delete fGPool;
570  fGPool = new GEVGPool;
571 
572  PDGCodeList::const_iterator nuiter;
573  PDGCodeList::const_iterator tgtiter;
574 
575  for(nuiter = fNuList.begin(); nuiter != fNuList.end(); ++nuiter) {
576  for(tgtiter = fTgtList.begin(); tgtiter != fTgtList.end(); ++tgtiter) {
577 
578  int target_pdgc = *tgtiter;
579  int neutrino_pdgc = *nuiter;
580 
581  InitialState init_state(target_pdgc, neutrino_pdgc);
582 
583  LOG("GMCJDriver", pNOTICE)
584  << "\n\n ---- Creating a GEVGDriver object configured for init-state: "
585  << init_state.AsString() << " ----\n\n";
586 
587  GEVGDriver * evgdriver = new GEVGDriver;
588  evgdriver->SetEventGeneratorList(fEventGenList); // specify list of generators
589  evgdriver->Configure(init_state);
590  evgdriver->UseSplines(); // check if all splines needed are loaded
591 
592  LOG("GMCJDriver", pDEBUG) << "Adding new GEVGDriver object to GEVGPool";
593  fGPool->insert( GEVGPool::value_type(init_state.AsString(), evgdriver) );
594  } // targets
595  } // neutrinos
596 
597  LOG("GMCJDriver", pNOTICE)
598  << "All necessary GEVGDriver object were pushed into GEVGPool\n";
599 }
GEVGPool * fGPool
A pool of GEVGDrivers properly configured event generation drivers / one per init state...
Definition: GMCJDriver.h:110
PDGCodeList fTgtList
[declared by the geom driver] list of target codes
Definition: GMCJDriver.h:115
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
void SetEventGeneratorList(string listname)
Definition: GEVGDriver.cxx:348
GENIE Event Generation Driver. A minimalist user interface object for generating neutrino interaction...
Definition: GEVGDriver.h:54
string fEventGenList
[config] list of event generators loaded by this driver (what used to be the $GEVGL setting) ...
Definition: GMCJDriver.h:129
void Configure(int nu_pdgc, int Z, int A)
Definition: GEVGDriver.cxx:137
void UseSplines(void)
Definition: GEVGDriver.cxx:508
#define pNOTICE
Definition: Messenger.h:61
PDGCodeList fNuList
[declared by the flux driver] list of neutrino codes
Definition: GMCJDriver.h:114
Initial State information.
Definition: InitialState.h:48
#define pDEBUG
Definition: Messenger.h:63
A pool of GEVGDriver objects with an initial state key.
Definition: GEVGPool.h:37
bool GMCJDriver::PreCalcFluxProbabilities ( void  )

Definition at line 192 of file GMCJDriver.cxx.

References genie::controls::kASmallNum, LOG, pFATAL, pNOTICE, and pWARN.

Referenced by main().

193 {
194 // Loop over complete set of flux entries satisfying input config options
195 // (such as neutrino type) and save the interaction probability in a tree
196 // relating flux index (entry number in input flux tree) to interaction
197 // probability. If a pre-generated flux interaction probability tree has
198 // already been loaded then just returns true. Also save tree to a TFile
199 // for use in later jobs if flag is set
200 //
201  bool success = true;
202 
203  bool save_to_file = fFluxIntProbFile == 0 && fFluxIntFileName.size()>0;
204 
205  // Clear map storing sum(fBrFluxWeight*fBrFluxIntProb) for each neutrino pdg
206  fSumFluxIntProbs.clear();
207 
208  // check if already loaded flux interaction probs using LoadFluxProbTree
209  if(fFluxIntTree){
210  LOG("GMCJDriver", pNOTICE) <<
211  "Skipping pre-generation of flux interaction probabilities - "<<
212  "using pre-generated file";
213  success = true;
214  }
215  // otherwise create them on the fly now
216  else {
217 
218  if(save_to_file){
219  fFluxIntProbFile = new TFile(fFluxIntFileName.c_str(), "CREATE");
220  if(fFluxIntProbFile->IsZombie()){
221  LOG("GMCJDriver", pFATAL) << "Cannot overwrite an existing file. Exiting!";
222  exit(1);
223  }
224  }
225 
226  // Create the tree to store flux probs
227  fFluxIntTree = new TTree(fFluxIntTreeName.c_str(),
228  "Tree storing pre-calculated flux interaction probs");
229  fFluxIntTree->Branch("FluxIndex", &fBrFluxIndex, "FluxIndex/I");
230  fFluxIntTree->Branch("FluxIntProb", &fBrFluxIntProb, "FluxIntProb/D");
231  fFluxIntTree->Branch("FluxEnu", &fBrFluxEnu, "FluxEnu/D");
232  fFluxIntTree->Branch("FluxWeight", &fBrFluxWeight, "FluxWeight/D");
233  fFluxIntTree->Branch("FluxPDG", &fBrFluxPDG, "FluxPDG/I");
234  // Associate to file otherwise get std::bad_alloc when writing large trees
235  if(save_to_file) fFluxIntTree->SetDirectory(fFluxIntProbFile);
236 
238 
239  fGlobPmax = 1.0; // Force ComputeInteractionProbabilities to return absolute value
240 
241  // Loop over flux entries and calculate interaction probabilities
242  TStopwatch stopwatch;
243  stopwatch.Start();
244  long int first_index = -1;
245  bool first_loop = true;
246  // loop until at end of flux ntuple
247  while(fFluxDriver->End() == false){
248 
249  // get the next flux neutrino
250  bool gotnext = fFluxDriver->GenerateNext();
251  if(!gotnext){
252  LOG("GMCJDriver", pWARN) << "*** Couldn't generate next flux ray! ";
253  continue;
254  }
255 
256  // stop if completed a full cycle (this check is necessary as fluxdriver
257  // may be set to loop over more than one cycle before reaching end)
258  bool already_been_here = first_loop ? false : first_index == fFluxDriver->Index();
259  if(already_been_here) break;
260 
261  // compute the path lengths for current flux neutrino
262  if(this->ComputePathLengths() == false){ success = false; break;}
263 
264  // compute and store the interaction probability
265  double psum = this->ComputeInteractionProbabilities(false /*Based on actual PLs*/);
266  assert(psum+controls::kASmallNum > 0.);
267  fBrFluxIntProb = psum;
272  fFluxIntTree->Fill();
273 
274  // store the first index so know when have cycled exactly once
275  if(first_loop){
276  first_index = fFluxDriver->Index();
277  first_loop = false;
278  }
279  } // flux loop
280  stopwatch.Stop();
281  LOG("GMCJDriver", pNOTICE)
282  << "Finished pre-calculating flux interaction probabilities. "
283  << "Total CPU time to process "<< fFluxIntTree->GetEntries()
284  << " entries: "<< stopwatch.CpuTime();
285 
286  // reset the flux driver so can be used at next stage. N.B. This
287  // should also reset flux driver to throw de-weighted flux neutrinos
288  fFluxDriver->Clear("CycleHistory");
289  }
290 
291  // If successfully calculated/loaded interaction probabilities then set global
292  // probability scale and, if requested, save tree to output file
293  if(success){
294  fGlobPmax = 0.0;
295  double safety_factor = 1.01;
296  for(int i = 0; i< fFluxIntTree->GetEntries(); i++){
297  fFluxIntTree->GetEntry(i);
298  // Check have non-negative probabilities
299  assert(fBrFluxIntProb+controls::kASmallNum > 0.0);
300  assert(fBrFluxWeight+controls::kASmallNum > 0.0);
301  // Update the global maximum
302  fGlobPmax = TMath::Max(fGlobPmax, fBrFluxIntProb*safety_factor);
303  // Update the sum of fBrFluxIntProb*fBrFluxWeight for different species
304  if(fSumFluxIntProbs.find(fBrFluxPDG) == fSumFluxIntProbs.end()){
306  }
308  }
309  LOG("GMCJDriver", pNOTICE) <<
310  "Updated global probability scale to fGlobPmax = "<< fGlobPmax;
311 
312  if(save_to_file){
313  LOG("GMCJDriver", pNOTICE) <<
314  "Saving pre-generated interaction probabilities to file: "<<
315  fFluxIntProbFile->GetName();
316  fFluxIntProbFile->cd();
317  fFluxIntTree->Write();
318  }
319 
320  // Also build index for use later
321  if(fFluxIntTree->BuildIndex("FluxIndex") != fFluxIntTree->GetEntries()){
322  LOG("GMCJDriver", pFATAL) <<
323  "Cannot build index using branch \"FluxIndex\" for flux prob tree!";
324  exit(1);
325  }
326 
327  // Now that have pre-generated flux probabilities need to trun off event
328  // preselection as this is only advantages when using max path lengths
329  this->PreSelectEvents(false);
330 
331  LOG("GMCJDriver", pNOTICE) << "Successfully generated/loaded pre-calculate flux interaction probabilities";
332  }
333  // Otherwise clean up
334  else if(fFluxIntTree){
335  delete fFluxIntTree;
336  fFluxIntTree = 0;
337  }
338 
339  // Return whether have successfully pre-calculated flux interaction probabilities
340  return success;
341 }
TTree * fFluxIntTree
[computed-or-loaded] pre-computed flux interaction probabilities (expected tree name is &quot;gFlxIntProbs...
Definition: GMCJDriver.h:140
double fBrFluxIntProb
flux interaction probability (set to branch:&quot;FluxIntProb&quot;)
Definition: GMCJDriver.h:141
#define pFATAL
Definition: Messenger.h:56
bool ComputePathLengths(void)
GFluxI * fFluxDriver
[input] neutrino flux driver
Definition: GMCJDriver.h:111
virtual void GenerateWeighted(bool gen_weighted)=0
set whether to generate weighted or unweighted neutrinos
double fBrFluxEnu
corresponding flux P4 (set to address of branch:&quot;FluxP4&quot;)
Definition: GMCJDriver.h:143
double fBrFluxWeight
corresponding flux weight (set to address of branch: &quot;FluxWeight&quot;)
Definition: GMCJDriver.h:144
string fFluxIntFileName
whether to save pre-generated flux tree for use in later jobs
Definition: GMCJDriver.h:146
virtual long int Index(void)=0
returns corresponding index for current flux neutrino (e.g. for a flux ntuple returns the current ent...
TFile * fFluxIntProbFile
[input] pre-generated flux interaction probability file
Definition: GMCJDriver.h:139
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
int fBrFluxPDG
corresponding flux pdg code (set to address of branch: &quot;FluxPDG&quot;)
Definition: GMCJDriver.h:145
int fBrFluxIndex
corresponding entry in flux input tree (set to address of branch:&quot;FluxEntry&quot;)
Definition: GMCJDriver.h:142
string fFluxIntTreeName
name for tree holding flux probabilities
Definition: GMCJDriver.h:147
static const double kASmallNum
Definition: Controls.h:40
double ComputeInteractionProbabilities(bool use_max_path_length)
double fGlobPmax
[computed at init] global interaction probability scale for given flux &amp; geometry ...
Definition: GMCJDriver.h:128
#define pWARN
Definition: Messenger.h:60
virtual void Clear(Option_t *opt)=0
reset state variables based on opt
virtual bool GenerateNext(void)=0
generate the next flux neutrino (return false in err)
void PreSelectEvents(bool preselect=true)
Definition: GMCJDriver.cxx:184
map< int, double > fSumFluxIntProbs
map where the key is flux pdg code and the value is sum of fBrFluxWeight * fBrFluxIntProb for all the...
Definition: GMCJDriver.h:148
virtual int PdgCode(void)=0
returns the flux neutrino pdg code
virtual const TLorentzVector & Momentum(void)=0
returns the flux neutrino 4-momentum
virtual bool End(void)=0
true if no more flux nu&#39;s can be thrown (eg reaching end of beam sim ntuples)
virtual double Weight(void)=0
returns the flux neutrino weight (if any)
#define pNOTICE
Definition: Messenger.h:61
double GMCJDriver::PreGenFluxInteractionProbability ( void  )
private

Definition at line 1321 of file GMCJDriver.cxx.

References genie::controls::kASmallNum, LOG, pERROR, and pFATAL.

1322 {
1323 // Return the pre-computed interaction probability for the current flux
1324 // neutrino index (entry number in flux file). Exit if not possible as
1325 // using meaningless interaction probability leads to incorrect physics
1326 //
1327  if(!fFluxIntTree){
1328  LOG("GMCJDriver", pERROR) <<
1329  "Cannot get pre-computed flux interaction probability as no tree!";
1330  exit(1);
1331  }
1332 
1333  assert(fFluxDriver->Index() >= 0); // Check trying to find meaningfull index
1334 
1335  // Check if can find relevant entry and no mismatch in energies -->
1336  // using correct pre-gen interaction prob file
1337  bool found_entry = fFluxIntTree->GetEntryWithIndex(fFluxDriver->Index()) > 0;
1338  bool enu_match = false;
1339  if(found_entry){
1340  double rel_err = fBrFluxEnu-fFluxDriver->Momentum().E();
1341  if(fBrFluxEnu > controls::kASmallNum) rel_err /= fBrFluxEnu;
1342  enu_match = TMath::Abs(rel_err)<controls::kASmallNum;
1343  if(enu_match == false){
1344  LOG("GMCJDriver", pERROR) <<
1345  "Mismatch between: Enu_curr = "<< fFluxDriver->Momentum().E() <<
1346  ", Enu_pre_gen = "<< fBrFluxEnu;
1347  }
1348  }
1349  else {
1350  LOG("GMCJDriver", pERROR) << "Cannot find flux entry in interaction prob tree!";
1351  }
1352 
1353  // Exit if not successful
1354  bool success = found_entry && enu_match;
1355  if(!success){
1356  LOG("GMCJDriver", pFATAL) <<
1357  "Cannot find pre-generated interaction probability! Check you "<<
1358  "are using the correct pre-generated interaction prob file " <<
1359  "generated using current flux input file with same input " <<
1360  "config (same geom TopVol, neutrino species list)";
1361  exit(1);
1362  }
1363  assert(fGlobPmax+controls::kASmallNum>0.0);
1364  return fBrFluxIntProb/fGlobPmax;
1365 }
TTree * fFluxIntTree
[computed-or-loaded] pre-computed flux interaction probabilities (expected tree name is &quot;gFlxIntProbs...
Definition: GMCJDriver.h:140
#define pERROR
Definition: Messenger.h:59
double fBrFluxIntProb
flux interaction probability (set to branch:&quot;FluxIntProb&quot;)
Definition: GMCJDriver.h:141
#define pFATAL
Definition: Messenger.h:56
GFluxI * fFluxDriver
[input] neutrino flux driver
Definition: GMCJDriver.h:111
double fBrFluxEnu
corresponding flux P4 (set to address of branch:&quot;FluxP4&quot;)
Definition: GMCJDriver.h:143
virtual long int Index(void)=0
returns corresponding index for current flux neutrino (e.g. for a flux ntuple returns the current ent...
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
static const double kASmallNum
Definition: Controls.h:40
double fGlobPmax
[computed at init] global interaction probability scale for given flux &amp; geometry ...
Definition: GMCJDriver.h:128
virtual const TLorentzVector & Momentum(void)=0
returns the flux neutrino 4-momentum
void GMCJDriver::PreSelectEvents ( bool  preselect = true)

Definition at line 184 of file GMCJDriver.cxx.

185 {
186 // Set whether to pre-select events based on a max-path lengths file. This
187 // should be turned off if using pre-generated interaction probabilities
188 // calculated from a given flux file.
189  fPreSelect = preselect;
190 }
bool fPreSelect
[config] set whether to pre-select events using max interaction paths
Definition: GMCJDriver.h:138
void GMCJDriver::SaveFluxProbabilities ( string  outfilename)

Definition at line 390 of file GMCJDriver.cxx.

Referenced by main().

391 {
392 // Configue the flux driver to save the calculated flux interaction
393 // probabilities to the specified output file name for use in later jobs. See
394 // the LoadFluxProbTree method for how they are fed into a later job.
395 //
396  fFluxIntFileName = outfilename;
397 }
string fFluxIntFileName
whether to save pre-generated flux tree for use in later jobs
Definition: GMCJDriver.h:146
int GMCJDriver::SelectTargetMaterial ( double  R)
private

Definition at line 1193 of file GMCJDriver.cxx.

References LOG, pERROR, and pNOTICE.

1194 {
1195 // Pick a target material using the pre-computed interaction probabilities
1196 // for a flux neutrino that has already been determined that interacts
1197 
1198  LOG("GMCJDriver", pNOTICE) << "Selecting target material";
1199  int tgtpdg = 0;
1200  map<int,double>::const_iterator probiter = fCurCumulProbMap.begin();
1201  for( ; probiter != fCurCumulProbMap.end(); ++probiter) {
1202  double prob = probiter->second;
1203  if(R<prob) {
1204  tgtpdg = probiter->first;
1205  LOG("GMCJDriver", pNOTICE)
1206  << "Selected target material = " << tgtpdg;
1207  return tgtpdg;
1208  }
1209  }
1210  LOG("GMCJDriver", pERROR)
1211  << "Could not select target material for an interacting neutrino";
1212  return 0;
1213 }
#define pERROR
Definition: Messenger.h:59
map< int, double > fCurCumulProbMap
[current] cummulative interaction probabilities
Definition: GMCJDriver.h:121
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pNOTICE
Definition: Messenger.h:61
void GMCJDriver::SetEventGeneratorList ( string  listname)

Definition at line 66 of file GMCJDriver.cxx.

References LOG, and pNOTICE.

Referenced by main().

67 {
68  LOG("GMCJDriver", pNOTICE)
69  << "Setting event generator list: " << listname;
70 
71  fEventGenList = listname;
72 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
string fEventGenList
[config] list of event generators loaded by this driver (what used to be the $GEVGL setting) ...
Definition: GMCJDriver.h:129
#define pNOTICE
Definition: Messenger.h:61
void GMCJDriver::SetPmaxLogBinning ( void  )

Definition at line 137 of file GMCJDriver.cxx.

References LOG, and pNOTICE.

138 {
139  fPmaxLogBinning = true;
140 
141  LOG("GMCJDriver", pNOTICE)
142  << "Pmax will be store in histogram with logarithmic energy bins";
143 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
bool fPmaxLogBinning
[config] maximum interaction probability is computed in logarithmic energy bins
Definition: GMCJDriver.h:124
#define pNOTICE
Definition: Messenger.h:61
void GMCJDriver::SetPmaxNbins ( int  nbins)

Definition at line 145 of file GMCJDriver.cxx.

References LOG, and pNOTICE.

146 {
147  fPmaxNbins = nbins;
148 
149  LOG("GMCJDriver", pNOTICE)
150  << "Number of bins in energy used for maximum int. probability: "
151  << fPmaxNbins;
152 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pNOTICE
Definition: Messenger.h:61
int fPmaxNbins
[config] number of bins in energy used in the maximum interaction probability
Definition: GMCJDriver.h:125
void GMCJDriver::SetPmaxSafetyFactor ( double  sf)

Definition at line 154 of file GMCJDriver.cxx.

References LOG, and pNOTICE.

155 {
156  fPmaxSafetyFactor = sf;
157 
158  LOG("GMCJDriver", pNOTICE)
159  << "Pmax safety factor = " << fPmaxSafetyFactor;
160 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pNOTICE
Definition: Messenger.h:61
double fPmaxSafetyFactor
[config] safety factor to compute the maximum interaction probability
Definition: GMCJDriver.h:126
void GMCJDriver::SetUnphysEventMask ( const TBits &  mask)

Definition at line 74 of file GMCJDriver.cxx.

References LOG, genie::GHepFlags::NFlags(), and pNOTICE.

75 {
76  *fUnphysEventMask = mask;
77 
78  LOG("GMCJDriver", pNOTICE)
79  << "Setting unphysical event mask (bits: " << GHepFlags::NFlags() - 1
80  << " -> 0) : " << *fUnphysEventMask;
81 }
TBits * fUnphysEventMask
[config] controls whether unphysical events are returned (what used to be the $GUNPHYSMASK setting) ...
Definition: GMCJDriver.h:130
static unsigned int NFlags(void)
Definition: GHepFlags.h:76
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pNOTICE
Definition: Messenger.h:61
void GMCJDriver::SetXSecSplineNbins ( int  nbins)

Definition at line 128 of file GMCJDriver.cxx.

References LOG, and pNOTICE.

129 {
130  fXSecSplineNbins = nbins;
131 
132  LOG("GMCJDriver", pNOTICE)
133  << "Number of bins in energy used for the xsec splines: "
134  << fXSecSplineNbins;
135 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
int fXSecSplineNbins
[config] number of bins in energy used in the xsec splines
Definition: GMCJDriver.h:123
#define pNOTICE
Definition: Messenger.h:61
map<int, double> genie::GMCJDriver::SumFluxIntProbs ( void  ) const
inline

Definition at line 78 of file GMCJDriver.h.

References fSumFluxIntProbs.

Referenced by main().

78 { return fSumFluxIntProbs; }
map< int, double > fSumFluxIntProbs
map where the key is flux pdg code and the value is sum of fBrFluxWeight * fBrFluxIntProb for all the...
Definition: GMCJDriver.h:148
void GMCJDriver::UseFluxDriver ( GFluxI flux)

Definition at line 83 of file GMCJDriver.cxx.

Referenced by main().

84 {
85  fFluxDriver = flux_driver;
86 }
GFluxI * fFluxDriver
[input] neutrino flux driver
Definition: GMCJDriver.h:111
void GMCJDriver::UseGeomAnalyzer ( GeomAnalyzerI geom)

Definition at line 88 of file GMCJDriver.cxx.

Referenced by main().

89 {
90  fGeomAnalyzer = geom_analyzer;
91 }
GeomAnalyzerI * fGeomAnalyzer
[input] detector geometry analyzer
Definition: GMCJDriver.h:112
bool GMCJDriver::UseMaxPathLengths ( string  xml_filename)

Definition at line 99 of file GMCJDriver.cxx.

References LOG, and pWARN.

Referenced by main().

100 {
101 // If you supply the maximum path lengths for all materials in your geometry
102 // (eg for ROOT/GEANT geometries they can be computed running GENIE's gmxpl
103 // application, see $GENIE/src/stdapp/gMaxPathLengths.cxx ) you can speed up
104 // the driver init phase by quite a bit (especially for complex geometries).
105 
106  fMaxPlXmlFilename = xml_filename;
107 
108  bool is_accessible = !(gSystem->AccessPathName(fMaxPlXmlFilename.c_str()));
109 
110  if ( is_accessible ) fUseExtMaxPl = true;
111  else {
112  fUseExtMaxPl = false;
113  LOG("GMCJDriver", pWARN)
114  << "UseMaxPathLengths could not find file: \"" << xml_filename << "\"";
115  }
116  return fUseExtMaxPl;
117 
118 }
bool fUseExtMaxPl
[config] using external max path length estimate?
Definition: GMCJDriver.h:132
#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
string fMaxPlXmlFilename
[config] input file with max density-weighted path lengths for all materials
Definition: GMCJDriver.h:131
void GMCJDriver::UseSplines ( bool  useLogE = true)

Definition at line 93 of file GMCJDriver.cxx.

Referenced by main().

94 {
95  fUseSplines = true;
96  fUseLogE = useLogE;
97 }
bool fUseLogE
[config] build splines = f(logE) (rather than f(E)) ?
Definition: GMCJDriver.h:134
bool fUseSplines
[config] compute all needed &amp; not-loaded splines at init
Definition: GMCJDriver.h:133

Member Data Documentation

double genie::GMCJDriver::fBrFluxEnu
private

corresponding flux P4 (set to address of branch:"FluxP4")

Definition at line 143 of file GMCJDriver.h.

int genie::GMCJDriver::fBrFluxIndex
private

corresponding entry in flux input tree (set to address of branch:"FluxEntry")

Definition at line 142 of file GMCJDriver.h.

double genie::GMCJDriver::fBrFluxIntProb
private

flux interaction probability (set to branch:"FluxIntProb")

Definition at line 141 of file GMCJDriver.h.

int genie::GMCJDriver::fBrFluxPDG
private

corresponding flux pdg code (set to address of branch: "FluxPDG")

Definition at line 145 of file GMCJDriver.h.

double genie::GMCJDriver::fBrFluxWeight
private

corresponding flux weight (set to address of branch: "FluxWeight")

Definition at line 144 of file GMCJDriver.h.

map<int,double> genie::GMCJDriver::fCurCumulProbMap
private

[current] cummulative interaction probabilities

Definition at line 121 of file GMCJDriver.h.

EventRecord* genie::GMCJDriver::fCurEvt
private

[current] generated event

Definition at line 119 of file GMCJDriver.h.

PathLengthList genie::GMCJDriver::fCurPathLengths
private

[current] path length list for current flux neutrino

Definition at line 117 of file GMCJDriver.h.

TLorentzVector genie::GMCJDriver::fCurVtx
private

[current] interaction vertex

Definition at line 118 of file GMCJDriver.h.

double genie::GMCJDriver::fEmax
private

[declared by the flux driver] maximum neutrino energy

Definition at line 113 of file GMCJDriver.h.

string genie::GMCJDriver::fEventGenList
private

[config] list of event generators loaded by this driver (what used to be the $GEVGL setting)

Definition at line 129 of file GMCJDriver.h.

GFluxI* genie::GMCJDriver::fFluxDriver
private

[input] neutrino flux driver

Definition at line 111 of file GMCJDriver.h.

Referenced by FluxDriver(), and FluxDriverPtr().

string genie::GMCJDriver::fFluxIntFileName
private

whether to save pre-generated flux tree for use in later jobs

Definition at line 146 of file GMCJDriver.h.

TFile* genie::GMCJDriver::fFluxIntProbFile
private

[input] pre-generated flux interaction probability file

Definition at line 139 of file GMCJDriver.h.

TTree* genie::GMCJDriver::fFluxIntTree
private

[computed-or-loaded] pre-computed flux interaction probabilities (expected tree name is "gFlxIntProbs")

Definition at line 140 of file GMCJDriver.h.

string genie::GMCJDriver::fFluxIntTreeName
private

name for tree holding flux probabilities

Definition at line 147 of file GMCJDriver.h.

bool genie::GMCJDriver::fForceInteraction
private

[config] force intearction?

Definition at line 137 of file GMCJDriver.h.

bool genie::GMCJDriver::fGenerateUnweighted
private

[config] force single probability scale?

Definition at line 136 of file GMCJDriver.h.

GeomAnalyzerI* genie::GMCJDriver::fGeomAnalyzer
private

[input] detector geometry analyzer

Definition at line 112 of file GMCJDriver.h.

Referenced by GeomAnalyzer(), and GeomAnalyzerPtr().

double genie::GMCJDriver::fGlobPmax
private

[computed at init] global interaction probability scale for given flux & geometry

Definition at line 128 of file GMCJDriver.h.

Referenced by GlobProbScale().

GEVGPool* genie::GMCJDriver::fGPool
private

A pool of GEVGDrivers properly configured event generation drivers / one per init state.

Definition at line 110 of file GMCJDriver.h.

bool genie::GMCJDriver::fKeepThrowingFluxNu
private

[config] keep firing flux neutrinos till one of them interacts

Definition at line 135 of file GMCJDriver.h.

PathLengthList genie::GMCJDriver::fMaxPathLengths
private

[declared by the geom driver] maximum path length list

Definition at line 116 of file GMCJDriver.h.

string genie::GMCJDriver::fMaxPlXmlFilename
private

[config] input file with max density-weighted path lengths for all materials

Definition at line 131 of file GMCJDriver.h.

double genie::GMCJDriver::fNFluxNeutrinos
private

[current] number of flux nuetrinos fired by the flux driver so far

Definition at line 122 of file GMCJDriver.h.

Referenced by NFluxNeutrinos().

PDGCodeList genie::GMCJDriver::fNuList
private

[declared by the flux driver] list of neutrino codes

Definition at line 114 of file GMCJDriver.h.

map<int,TH1D*> genie::GMCJDriver::fPmax
private

[computed at init] interaction probability scale /neutrino /energy for given geometry

Definition at line 127 of file GMCJDriver.h.

bool genie::GMCJDriver::fPmaxLogBinning
private

[config] maximum interaction probability is computed in logarithmic energy bins

Definition at line 124 of file GMCJDriver.h.

int genie::GMCJDriver::fPmaxNbins
private

[config] number of bins in energy used in the maximum interaction probability

Definition at line 125 of file GMCJDriver.h.

double genie::GMCJDriver::fPmaxSafetyFactor
private

[config] safety factor to compute the maximum interaction probability

Definition at line 126 of file GMCJDriver.h.

bool genie::GMCJDriver::fPreSelect
private

[config] set whether to pre-select events using max interaction paths

Definition at line 138 of file GMCJDriver.h.

int genie::GMCJDriver::fSelTgtPdg
private

[current] selected target material PDG code

Definition at line 120 of file GMCJDriver.h.

map<int, double> genie::GMCJDriver::fSumFluxIntProbs
private

map where the key is flux pdg code and the value is sum of fBrFluxWeight * fBrFluxIntProb for all these flux neutrinos

Definition at line 148 of file GMCJDriver.h.

Referenced by SumFluxIntProbs().

PDGCodeList genie::GMCJDriver::fTgtList
private

[declared by the geom driver] list of target codes

Definition at line 115 of file GMCJDriver.h.

TBits* genie::GMCJDriver::fUnphysEventMask
private

[config] controls whether unphysical events are returned (what used to be the $GUNPHYSMASK setting)

Definition at line 130 of file GMCJDriver.h.

bool genie::GMCJDriver::fUseExtMaxPl
private

[config] using external max path length estimate?

Definition at line 132 of file GMCJDriver.h.

bool genie::GMCJDriver::fUseLogE
private

[config] build splines = f(logE) (rather than f(E)) ?

Definition at line 134 of file GMCJDriver.h.

bool genie::GMCJDriver::fUseSplines
private

[config] compute all needed & not-loaded splines at init

Definition at line 133 of file GMCJDriver.h.

int genie::GMCJDriver::fXSecSplineNbins
private

[config] number of bins in energy used in the xsec splines

Definition at line 123 of file GMCJDriver.h.


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