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

Generates values for the kinematic variables describing IBD neutrino interaction events. Is a concrete implementation of the EventRecordVisitorI interface. More...

#include <IBDKinematicsGenerator.h>

Inheritance diagram for genie::IBDKinematicsGenerator:
Inheritance graph
[legend]
Collaboration diagram for genie::IBDKinematicsGenerator:
Collaboration graph
[legend]

Public Member Functions

 IBDKinematicsGenerator ()
 
 IBDKinematicsGenerator (string config)
 
virtual ~IBDKinematicsGenerator ()
 
void ProcessEventRecord (GHepRecord *event_rec) const
 
void Configure (const Registry &config)
 
void Configure (string config)
 
- Public Member Functions inherited from genie::EventRecordVisitorI
virtual ~EventRecordVisitorI ()
 
- Public Member Functions inherited from genie::Algorithm
virtual ~Algorithm ()
 
virtual void FindConfig (void)
 
virtual const RegistryGetConfig (void) const
 
RegistryGetOwnedConfig (void)
 
virtual const AlgIdId (void) const
 Get algorithm ID. More...
 
virtual AlgStatus_t GetStatus (void) const
 Get algorithm status. More...
 
virtual bool AllowReconfig (void) const
 
virtual AlgCmp_t Compare (const Algorithm *alg) const
 Compare with input algorithm. More...
 
virtual void SetId (const AlgId &id)
 Set algorithm ID. More...
 
virtual void SetId (string name, string config)
 
const AlgorithmSubAlg (const RgKey &registry_key) const
 
void AdoptConfig (void)
 
void AdoptSubstructure (void)
 
virtual void Print (ostream &stream) const
 Print algorithm info. More...
 

Private Member Functions

void LoadConfig (void)
 
double ComputeMaxXSec (const Interaction *in) const
 

Additional Inherited Members

- Static Public Member Functions inherited from genie::Algorithm
static string BuildParamVectKey (const std::string &comm_name, unsigned int i)
 
static string BuildParamVectSizeKey (const std::string &comm_name)
 
static string BuildParamMatKey (const std::string &comm_name, unsigned int i, unsigned int j)
 
static string BuildParamMatRowSizeKey (const std::string &comm_name)
 
static string BuildParamMatColSizeKey (const std::string &comm_name)
 
- Protected Member Functions inherited from genie::KineGeneratorWithCache
 KineGeneratorWithCache ()
 
 KineGeneratorWithCache (string name)
 
 KineGeneratorWithCache (string name, string config)
 
 ~KineGeneratorWithCache ()
 
virtual double ComputeMaxXSec (const Interaction *in, const int nkey) const
 
virtual double MaxXSec (GHepRecord *evrec, const int nkey=0) const
 
virtual double FindMaxXSec (const Interaction *in, const int nkey=0) const
 
virtual void CacheMaxXSec (const Interaction *in, double xsec, const int nkey=0) const
 
virtual double Energy (const Interaction *in) const
 
virtual CacheBranchFxAccessCacheBranch (const Interaction *in, const int nkey=0) const
 
virtual void AssertXSecLimits (const Interaction *in, double xsec, double xsec_max) const
 
- Protected Member Functions inherited from genie::EventRecordVisitorI
 EventRecordVisitorI ()
 
 EventRecordVisitorI (string name)
 
 EventRecordVisitorI (string name, string config)
 
- Protected Member Functions inherited from genie::Algorithm
 Algorithm ()
 
 Algorithm (string name)
 
 Algorithm (string name, string config)
 
void Initialize (void)
 
void DeleteConfig (void)
 
void DeleteSubstructure (void)
 
RegistryExtractLocalConfig (const Registry &in) const
 
RegistryExtractLowerConfig (const Registry &in, const string &alg_key) const
 Split an incoming configuration Registry into a block valid for the sub-algo identified by alg_key. More...
 
template<class T >
bool GetParam (const RgKey &name, T &p, bool is_top_call=true) const
 
template<class T >
bool GetParamDef (const RgKey &name, T &p, const T &def) const
 
template<class T >
int GetParamVect (const std::string &comm_name, std::vector< T > &v, bool is_top_call=true) const
 Handle to load vectors of parameters. More...
 
int GetParamVectKeys (const std::string &comm_name, std::vector< RgKey > &k, bool is_top_call=true) const
 
template<class T >
int GetParamMat (const std::string &comm_name, TMatrixT< T > &mat, bool is_top_call=true) const
 Handle to load matrix of parameters. More...
 
template<class T >
int GetParamMatSym (const std::string &comm_name, TMatrixTSym< T > &mat, bool is_top_call=true) const
 
int GetParamMatKeys (const std::string &comm_name, std::vector< RgKey > &k, bool is_top_call=true) const
 
int AddTopRegistry (Registry *rp, bool owns=true)
 add registry with top priority, also update ownership More...
 
int AddLowRegistry (Registry *rp, bool owns=true)
 add registry with lowest priority, also update ownership More...
 
int MergeTopRegistry (const Registry &r)
 
int AddTopRegisties (const vector< Registry * > &rs, bool owns=false)
 Add registries with top priority, also udated Ownerships. More...
 
- Protected Attributes inherited from genie::KineGeneratorWithCache
const XSecAlgorithmIfXSecModel
 
double fSafetyFactor
 ComputeMaxXSec -> ComputeMaxXSec * fSafetyFactor. More...
 
std::vector< double > vSafetyFactors
 MaxXSec -> MaxXSec * fSafetyFactors[nkey]. More...
 
int fNumOfSafetyFactors
 Number of given safety factors. More...
 
std::vector< string > vInterpolatorTypes
 Type of interpolator for each key in a branch. More...
 
int fNumOfInterpolatorTypes
 Number of given interpolators types. More...
 
double fMaxXSecDiffTolerance
 max{100*(xsec-maxxsec)/.5*(xsec+maxxsec)} if xsec>maxxsec More...
 
double fEMin
 min E for which maxxsec is cached - forcing explicit calc. More...
 
bool fGenerateUniformly
 uniform over allowed phase space + event weight? More...
 
- Protected Attributes inherited from genie::Algorithm
bool fAllowReconfig
 
bool fOwnsSubstruc
 true if it owns its substructure (sub-algs,...) More...
 
AlgId fID
 algorithm name and configuration set More...
 
vector< Registry * > fConfVect
 
vector< bool > fOwnerships
 ownership for every registry in fConfVect More...
 
AlgStatus_t fStatus
 algorithm execution status More...
 
AlgMapfOwnedSubAlgMp
 local pool for owned sub-algs (taken out of the factory pool) More...
 

Detailed Description

Generates values for the kinematic variables describing IBD neutrino interaction events. Is a concrete implementation of the EventRecordVisitorI interface.

Author
Corey Reed <cjreed nikhef.nl> - October 29, 2009 using code from the QELKinematicGenerator written by Costas Andreopoulos <c.andreopoulos cern.ch> University of Liverpool
Created:
October 29, 2009
License:
Copyright (c) 2003-2024, The GENIE Collaboration For the full text of the license visit http://copyright.genie-mc.org

Definition at line 31 of file IBDKinematicsGenerator.h.

Constructor & Destructor Documentation

IBDKinematicsGenerator::IBDKinematicsGenerator ( )

Definition at line 41 of file IBDKinematicsGenerator.cxx.

41  :
42 KineGeneratorWithCache("genie::IBDKinematicsGenerator")
43 {
44 
45 }
IBDKinematicsGenerator::IBDKinematicsGenerator ( string  config)

Definition at line 47 of file IBDKinematicsGenerator.cxx.

47  :
48 KineGeneratorWithCache("genie::IBDKinematicsGenerator", config)
49 {
50 
51 }
IBDKinematicsGenerator::~IBDKinematicsGenerator ( )
virtual

Definition at line 53 of file IBDKinematicsGenerator.cxx.

54 {
55 
56 }

Member Function Documentation

double IBDKinematicsGenerator::ComputeMaxXSec ( const Interaction in) const
privatevirtual

Implements genie::KineGeneratorWithCache.

Definition at line 242 of file IBDKinematicsGenerator.cxx.

References genie::Interaction::AsString(), genie::KineGeneratorWithCache::fSafetyFactor, genie::KineGeneratorWithCache::fXSecModel, genie::Interaction::KinePtr(), genie::kKVQ2, genie::kPSQ2fE, genie::KPhaseSpace::Limits(), LOG, genie::Range1D_t::max, genie::Range1D_t::min, pDEBUG, genie::Interaction::PhaseSpace(), genie::utils::kinematics::Q2(), genie::Kinematics::SetQ2(), SLOG, and genie::XSecAlgorithmI::XSec().

244 {
245 // Computes the maximum differential cross section in the requested phase
246 // space. This method overloads KineGeneratorWithCache::ComputeMaxXSec
247 // method and the value is cached at a circular cache branch for retrieval
248 // during subsequent event generation.
249 // The computed max differential cross section does not need to be the exact
250 // maximum. The number used in the rejection method will be scaled up by a
251 // safety factor. But it needs to be fast - do not use a very small dQ2 step.
252 
253  double max_xsec = 0.0;
254 
255  const KPhaseSpace & kps = interaction->PhaseSpace();
256  Range1D_t rQ2 = kps.Limits(kKVQ2);
257  if(rQ2.min <=0 || rQ2.max <= rQ2.min) return 0.;
258 
259  //const double logQ2min = TMath::Log(rQ2.min + kASmallNum);
260  //const double logQ2max = TMath::Log(rQ2.max - kASmallNum);
261  const double logQ2min = TMath::Log(rQ2.min);
262  const double logQ2max = TMath::Log(rQ2.max);
263 
264  const int N = 15;
265  const int Nb = 10;
266 
267  double dlogQ2 = (logQ2max - logQ2min) /(N-1);
268  double xseclast = -1;
269  bool increasing;
270 
271  for(int i=0; i<N; i++) {
272  double Q2 = TMath::Exp(logQ2min + i * dlogQ2);
273  interaction->KinePtr()->SetQ2(Q2);
274  double xsec = fXSecModel->XSec(interaction, kPSQ2fE);
275 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
276  LOG("IBD", pDEBUG) << "xsec(Q2= " << Q2 << ") = " << xsec;
277 #endif
278  max_xsec = TMath::Max(xsec, max_xsec);
279  increasing = xsec-xseclast>=0;
280  xseclast = xsec;
281 
282  // once the cross section stops increasing, I reduce the step size and
283  // step backwards a little bit to handle cases that the max cross section
284  // is grossly underestimated (very peaky distribution & large step)
285  if(!increasing) {
286  dlogQ2/=(Nb+1);
287  for(int ib=0; ib<Nb; ib++) {
288  Q2 = TMath::Exp(TMath::Log(Q2) - dlogQ2);
289  if(Q2 < rQ2.min) continue;
290  interaction->KinePtr()->SetQ2(Q2);
291  xsec = fXSecModel->XSec(interaction, kPSQ2fE);
292 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
293  LOG("IBD", pDEBUG) << "xsec(Q2= " << Q2 << ") = " << xsec;
294 #endif
295  max_xsec = TMath::Max(xsec, max_xsec);
296  }
297  break;
298  }
299  }//Q^2
300 
301  // Apply safety factor, since value retrieved from the cache might
302  // correspond to a slightly different energy
303  max_xsec *= fSafetyFactor;
304 
305 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
306  SLOG("IBD", pDEBUG) << interaction->AsString();
307  SLOG("IBD", pDEBUG) << "Max xsec in phase space = " << max_xsec;
308  SLOG("IBD", pDEBUG) << "Computed using alg = " << *fXSecModel;
309 #endif
310 
311  return max_xsec;
312 }
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:1077
double fSafetyFactor
ComputeMaxXSec -&gt; ComputeMaxXSec * fSafetyFactor.
A simple [min,max] interval for doubles.
Definition: Range1.h:42
Range1D_t Limits(KineVar_t kvar) const
Return the kinematical variable limits.
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
Kinematical phase space.
Definition: KPhaseSpace.h:33
virtual double XSec(const Interaction *i, KinePhaseSpace_t k=kPSfE) const =0
Compute the cross section for the input interaction.
double max
Definition: Range1.h:53
double min
Definition: Range1.h:52
#define SLOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a short string (using the FUNCTION and...
Definition: Messenger.h:84
#define pDEBUG
Definition: Messenger.h:63
void IBDKinematicsGenerator::Configure ( const Registry config)
virtual

Configure the algorithm with an external registry The registry is merged with the top level registry if it is owned, Otherwise a copy of it is added with the highest priority

Reimplemented from genie::Algorithm.

Definition at line 210 of file IBDKinematicsGenerator.cxx.

References genie::Algorithm::Configure(), and LoadConfig().

211 {
212  Algorithm::Configure(config);
213  this->LoadConfig();
214 }
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
void IBDKinematicsGenerator::Configure ( string  config)
virtual

Configure the algorithm from the AlgoConfigPool based on param_set string given in input An algorithm contains a vector of registries coming from different xml configuration files, which are loaded according a very precise prioriy This methods will load a number registries in order of priority: 1) "Tunable" parameter set from CommonParametes. This is loaded with the highest prioriry and it is designed to be used for tuning procedure Usage not expected from the user. 2) For every string defined in "CommonParame" the corresponding parameter set will be loaded from CommonParameter.xml 3) parameter set specified by the config string and defined in the xml file of the algorithm 4) if config is not "Default" also the Default parameter set from the same xml file will be loaded Effectively this avoids the repetion of a parameter when it is not changed in the requested configuration

Reimplemented from genie::Algorithm.

Definition at line 216 of file IBDKinematicsGenerator.cxx.

References genie::Algorithm::Configure(), and LoadConfig().

217 {
218  Algorithm::Configure(config);
219  this->LoadConfig();
220 }
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
void IBDKinematicsGenerator::LoadConfig ( void  )
private

Definition at line 222 of file IBDKinematicsGenerator.cxx.

References genie::KineGeneratorWithCache::fEMin, genie::KineGeneratorWithCache::fGenerateUniformly, genie::KineGeneratorWithCache::fMaxXSecDiffTolerance, genie::KineGeneratorWithCache::fSafetyFactor, and genie::Algorithm::GetParamDef().

Referenced by Configure().

223 {
224  //-- Safety factor for the maximum differential cross section
225  GetParamDef( "MaxXSec-SafetyFactor", fSafetyFactor, 1.25 ) ;
226 
227  //-- Minimum energy for which max xsec would be cached, forcing explicit
228  // calculation for lower eneries
229  GetParamDef( "Cache-MinEnergy", fEMin, 1.00 ) ;
230 
231  //-- Maximum allowed fractional cross section deviation from maxim cross
232  // section used in rejection method
233  GetParamDef( "MaxXSec-DiffTolerance", fMaxXSecDiffTolerance, 999999. ) ;
234  assert(fMaxXSecDiffTolerance>=0);
235 
236  //-- Generate kinematics uniformly over allowed phase space and compute
237  // an event weight?
238  GetParamDef( "UniformOverPhaseSpace", fGenerateUniformly, false ) ;
239 
240 }
bool fGenerateUniformly
uniform over allowed phase space + event weight?
double fSafetyFactor
ComputeMaxXSec -&gt; ComputeMaxXSec * fSafetyFactor.
double fMaxXSecDiffTolerance
max{100*(xsec-maxxsec)/.5*(xsec+maxxsec)} if xsec&gt;maxxsec
double fEMin
min E for which maxxsec is cached - forcing explicit calc.
bool GetParamDef(const RgKey &name, T &p, const T &def) const
void IBDKinematicsGenerator::ProcessEventRecord ( GHepRecord event_rec) const
virtual

Implements genie::EventRecordVisitorI.

Definition at line 58 of file IBDKinematicsGenerator.cxx.

References genie::KineGeneratorWithCache::AssertXSecLimits(), genie::Kinematics::ClearRunningValues(), genie::EventGeneratorI::CrossSectionAlg(), genie::GHepRecord::EventFlags(), genie::KineGeneratorWithCache::fGenerateUniformly, genie::PDGLibrary::Find(), genie::KineGeneratorWithCache::fXSecModel, gQ2, genie::Target::HitNucP4(), genie::Interaction::InitState(), genie::RunningThreadInfo::Instance(), genie::RandomGen::Instance(), genie::PDGLibrary::Instance(), genie::kIAssumeFreeNucleon, genie::Interaction::KinePtr(), genie::kISkipKinematicChk, genie::kISkipProcessChk, genie::kKineGenErr, genie::kKVQ2, genie::kPSQ2fE, genie::kRfHitNucRest, genie::controls::kRjMaxIterations, genie::KPhaseSpace::Limits(), LOG, genie::Range1D_t::max, genie::KineGeneratorWithCache::MaxXSec(), genie::Range1D_t::min, pDEBUG, genie::Interaction::PhaseSpace(), genie::utils::kinematics::PhaseSpaceVolume(), pINFO, pNOTICE, genie::InitialState::ProbeE(), pWARN, genie::utils::kinematics::Q2(), genie::Interaction::RecoilNucleonPdg(), genie::RandomGen::RndKine(), genie::RunningThreadInfo::RunningThread(), genie::GHepRecord::SetDiffXSec(), genie::Kinematics::SetQ2(), genie::exceptions::EVGThreadException::SetReason(), genie::Kinematics::SetW(), genie::GHepRecord::SetWeight(), genie::Kinematics::Setx(), genie::Kinematics::Sety(), genie::GHepRecord::Summary(), genie::exceptions::EVGThreadException::SwitchOnFastForward(), genie::InitialState::Tgt(), genie::GHepRecord::Weight(), genie::utils::kinematics::WQ2toXY(), genie::XSecAlgorithmI::XSec(), and genie::GHepRecord::XSec().

59 {
60  if(fGenerateUniformly) {
61  LOG("IBD", pNOTICE)
62  << "Generating kinematics uniformly over the allowed phase space";
63  }
64 
65  //-- Get the random number generators
67 
68  //-- Access cross section algorithm for running thread
70  const EventGeneratorI * evg = rtinfo->RunningThread();
71  fXSecModel = evg->CrossSectionAlg();
72 
73  //-- Get the interaction and set the 'trust' bits
74  Interaction * interaction = evrec->Summary();
75  interaction->SetBit(kISkipProcessChk);
76  interaction->SetBit(kISkipKinematicChk);
77 
78  //-- Note: The kinematic generator would be using the free nucleon cross
79  // section (even for nuclear targets) so as not to double-count nuclear
80  // suppression. This assumes that a) the nuclear suppression was turned
81  // on when computing the cross sections for selecting the current event
82  // and that b) if the event turns out to be unphysical (Pauli-blocked)
83  // the next attempted event will be forced to QEL again.
84  // (discussion with Hugh - GENIE/NeuGEN integration workshop - 07APR2006
85  interaction->SetBit(kIAssumeFreeNucleon);
86 
87  //-- Get the limits for the generated Q2
88  const KPhaseSpace & kps = interaction->PhaseSpace();
89  const Range1D_t Q2 = kps.Limits(kKVQ2);
90 
91  if(Q2.max <=0 || Q2.min>=Q2.max) {
92  LOG("IBD", pWARN) << "No available phase space";
93  evrec->EventFlags()->SetBitNumber(kKineGenErr, true);
95  exception.SetReason("No available phase space");
96  exception.SwitchOnFastForward();
97  throw exception;
98  }
99 
100  //-- For the subsequent kinematic selection with the rejection method:
101  // Calculate the max differential cross section or retrieve it from the
102  // cache. Throw an exception and quit the evg thread if a non-positive
103  // value is found.
104  // If the kinematics are generated uniformly over the allowed phase
105  // space the max xsec is irrelevant
106  const double xsec_max = (fGenerateUniformly) ? -1 : this->MaxXSec(evrec);
107 
108  //-- Try to select a valid Q2 using the rejection method
109 
110  // kinematical limits
111  const double Q2min = Q2.min;
112  const double Q2max = Q2.max;
113  double xsec = -1.;
114  double gQ2 = 0.;
115 
116  unsigned int iter = 0;
117  bool accept = false;
118  while(1) {
119  iter++;
120  if(iter > kRjMaxIterations) {
121  LOG("IBD", pWARN)
122  << "Couldn't select a valid Q^2 after " << iter << " iterations";
123  evrec->EventFlags()->SetBitNumber(kKineGenErr, true);
125  exception.SetReason("Couldn't select kinematics");
126  exception.SwitchOnFastForward();
127  throw exception;
128  }
129 
130  //-- Generate a Q2 value within the allowed phase space
131  gQ2 = Q2min + (Q2max-Q2min) * rnd->RndKine().Rndm();
132  interaction->KinePtr()->SetQ2(gQ2);
133  LOG("IBD", pINFO) << "Trying: Q^2 = " << gQ2;
134 
135  //-- Computing cross section for the current kinematics
136  xsec = fXSecModel->XSec(interaction, kPSQ2fE);
137 
138  //-- Decide whether to accept the current kinematics
139  if(!fGenerateUniformly) {
140  this->AssertXSecLimits(interaction, xsec, xsec_max);
141  const double t = xsec_max * rnd->RndKine().Rndm();
142 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
143  LOG("IBD", pDEBUG)
144  << "dxsec/dQ2 = " << xsec << ", rnd = " << t;
145 #endif
146  accept = (t < xsec);
147  } else {
148  accept = (xsec>0);
149  }
150 
151  //-- If the generated kinematics are accepted, finish-up module's job
152  if(accept) {
153  LOG("IBD", pINFO) << "Selected: Q^2 = " << gQ2;
154 
155  // reset bits
156  interaction->ResetBit(kISkipProcessChk);
157  interaction->ResetBit(kISkipKinematicChk);
158  interaction->ResetBit(kIAssumeFreeNucleon);
159 
160  // compute the rest of the kinematical variables
161 
162  // get neutrino energy at struck nucleon rest frame and the
163  // struck nucleon mass (can be off the mass shell)
164  const InitialState & init_state = interaction->InitState();
165  const double E = init_state.ProbeE(kRfHitNucRest);
166  const double M = init_state.Tgt().HitNucP4().M();
167 
168  LOG("IBD", pNOTICE) << "E = " << E << ", M = "<< M;
169 
170  // The hadronic inv. mass is equal to the recoil nucleon on-shell mass.
171  const int rpdgc = interaction->RecoilNucleonPdg();
172  assert(rpdgc);
173  const double gW = PDGLibrary::Instance()->Find(rpdgc)->Mass();
174 
175  LOG("IBD", pNOTICE) << "Selected: W = "<< gW;
176 
177  // (W,Q2) -> (x,y)
178  double gx=0, gy=0;
179  kinematics::WQ2toXY(E,M,gW,gQ2,gx,gy);
180 
181  // set the cross section for the selected kinematics
182  evrec->SetDiffXSec(xsec,kPSQ2fE);
183 
184  // for uniform kinematics, compute an event weight as
185  // wght = (phase space volume)*(differential xsec)/(event total xsec)
186  if(fGenerateUniformly) {
187  const double vol = kinematics::PhaseSpaceVolume(interaction,kPSQ2fE);
188  const double totxsec = evrec->XSec();
189  double wght = (vol/totxsec)*xsec;
190  LOG("IBD", pNOTICE) << "Kinematics wght = "<< wght;
191 
192  // apply computed weight to the current event weight
193  wght *= evrec->Weight();
194  LOG("IBD", pNOTICE) << "Current event wght = " << wght;
195  evrec->SetWeight(wght);
196  }
197 
198  // lock selected kinematics & clear running values
199  interaction->KinePtr()->SetQ2(gQ2, true);
200  interaction->KinePtr()->SetW (gW, true);
201  interaction->KinePtr()->Setx (gx, true);
202  interaction->KinePtr()->Sety (gy, true);
203  interaction->KinePtr()->ClearRunningValues();
204 
205  return;
206  }
207  }// iterations
208 }
const KPhaseSpace & PhaseSpace(void) const
Definition: Interaction.h:73
bool fGenerateUniformly
uniform over allowed phase space + event weight?
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:1077
static RandomGen * Instance()
Access instance.
Definition: RandomGen.cxx:71
void SetQ2(double Q2, bool selected=false)
Definition: Kinematics.cxx:255
Kinematics * KinePtr(void) const
Definition: Interaction.h:76
int RecoilNucleonPdg(void) const
recoil nucleon pdg
A simple [min,max] interval for doubles.
Definition: Range1.h:42
virtual void AssertXSecLimits(const Interaction *in, double xsec, double xsec_max) const
Defines the EventGeneratorI interface.
virtual double MaxXSec(GHepRecord *evrec, const int nkey=0) const
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition: RandomGen.h:29
Range1D_t Limits(KineVar_t kvar) const
Return the kinematical variable limits.
double PhaseSpaceVolume(const Interaction *const i, KinePhaseSpace_t ps)
Definition: KineUtils.cxx:36
Summary information for an interaction.
Definition: Interaction.h:56
const TLorentzVector & HitNucP4(void) const
Definition: Target.h:91
An exception thrown by EventRecordVisitorI when the normal processing sequence has to be disrupted (f...
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
void WQ2toXY(double Ev, double M, double W, double Q2, double &x, double &y)
Definition: KineUtils.cxx:1132
Kinematical phase space.
Definition: KPhaseSpace.h:33
#define pINFO
Definition: Messenger.h:62
TRandom3 & RndKine(void) const
rnd number generator used by kinematics generators
Definition: RandomGen.h:50
#define pWARN
Definition: Messenger.h:60
const UInt_t kISkipKinematicChk
if set, skip kinematic validity checks
Definition: Interaction.h:48
virtual double XSec(const Interaction *i, KinePhaseSpace_t k=kPSfE) const =0
Compute the cross section for the input interaction.
void Setx(double x, bool selected=false)
Definition: Kinematics.cxx:231
static RunningThreadInfo * Instance(void)
double max
Definition: Range1.h:53
void SetW(double W, bool selected=false)
Definition: Kinematics.cxx:279
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:68
void Sety(double y, bool selected=false)
Definition: Kinematics.cxx:243
const UInt_t kIAssumeFreeNucleon
Definition: Interaction.h:49
virtual const XSecAlgorithmI * CrossSectionAlg(void) const =0
double gQ2
Definition: gtestDISSF.cxx:55
const InitialState & InitState(void) const
Definition: Interaction.h:69
double min
Definition: Range1.h:52
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:86
#define pNOTICE
Definition: Messenger.h:61
void ClearRunningValues(void)
Definition: Kinematics.cxx:347
const Target & Tgt(void) const
Definition: InitialState.h:66
static const unsigned int kRjMaxIterations
Definition: Controls.h:26
const EventGeneratorI * RunningThread(void)
double ProbeE(RefFrame_t rf) const
Keep info on the event generation thread currently on charge. This is used so that event generation m...
const UInt_t kISkipProcessChk
if set, skip process validity checks
Definition: Interaction.h:47
Initial State information.
Definition: InitialState.h:48
#define pDEBUG
Definition: Messenger.h:63

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