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

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

#include <QELKinematicsGenerator.h>

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

Public Member Functions

 QELKinematicsGenerator ()
 
 QELKinematicsGenerator (string config)
 
 ~QELKinematicsGenerator ()
 
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 SpectralFuncExperimentalCode (GHepRecord *event_rec) const
 
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 QEL neutrino interaction events. Is a concrete implementation of the EventRecordVisitorI interface.

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

Definition at line 28 of file QELKinematicsGenerator.h.

Constructor & Destructor Documentation

QELKinematicsGenerator::QELKinematicsGenerator ( )

Definition at line 37 of file QELKinematicsGenerator.cxx.

37  :
38 KineGeneratorWithCache("genie::QELKinematicsGenerator")
39 {
40 
41 }
QELKinematicsGenerator::QELKinematicsGenerator ( string  config)

Definition at line 43 of file QELKinematicsGenerator.cxx.

43  :
44 KineGeneratorWithCache("genie::QELKinematicsGenerator", config)
45 {
46 
47 }
QELKinematicsGenerator::~QELKinematicsGenerator ( )

Definition at line 49 of file QELKinematicsGenerator.cxx.

50 {
51 
52 }

Member Function Documentation

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

Implements genie::KineGeneratorWithCache.

Definition at line 465 of file QELKinematicsGenerator.cxx.

References genie::Interaction::AsString(), genie::KineGeneratorWithCache::fSafetyFactor, genie::KineGeneratorWithCache::fXSecModel, genie::controls::kASmallNum, 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().

467 {
468 // Computes the maximum differential cross section in the requested phase
469 // space. This method overloads KineGeneratorWithCache::ComputeMaxXSec
470 // method and the value is cached at a circular cache branch for retrieval
471 // during subsequent event generation.
472 // The computed max differential cross section does not need to be the exact
473 // maximum. The number used in the rejection method will be scaled up by a
474 // safety factor. But it needs to be fast - do not use a very small dQ2 step.
475 
476  double max_xsec = 0.0;
477 
478  const KPhaseSpace & kps = interaction->PhaseSpace();
479  Range1D_t rQ2 = kps.Limits(kKVQ2);
480  if(rQ2.min <=0 || rQ2.max <= rQ2.min) return 0.;
481 
482  const double logQ2min = TMath::Log(rQ2.min + kASmallNum);
483  const double logQ2max = TMath::Log(rQ2.max - kASmallNum);
484 
485  const int N = 15;
486  const int Nb = 10;
487 
488  double dlogQ2 = (logQ2max - logQ2min) /(N-1);
489  double xseclast = -1;
490  bool increasing;
491 
492  for(int i=0; i<N; i++) {
493  double Q2 = TMath::Exp(logQ2min + i * dlogQ2);
494  interaction->KinePtr()->SetQ2(Q2);
495  double xsec = fXSecModel->XSec(interaction, kPSQ2fE);
496 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
497  LOG("QELKinematics", pDEBUG) << "xsec(Q2= " << Q2 << ") = " << xsec;
498 #endif
499  max_xsec = TMath::Max(xsec, max_xsec);
500  increasing = xsec-xseclast>=0;
501  xseclast = xsec;
502 
503  // once the cross section stops increasing, I reduce the step size and
504  // step backwards a little bit to handle cases that the max cross section
505  // is grossly underestimated (very peaky distribution & large step)
506  if(!increasing) {
507  dlogQ2/=(Nb+1);
508  for(int ib=0; ib<Nb; ib++) {
509  Q2 = TMath::Exp(TMath::Log(Q2) - dlogQ2);
510  if(Q2 < rQ2.min) continue;
511  interaction->KinePtr()->SetQ2(Q2);
512  xsec = fXSecModel->XSec(interaction, kPSQ2fE);
513 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
514  LOG("QELKinematics", pDEBUG) << "xsec(Q2= " << Q2 << ") = " << xsec;
515 #endif
516  max_xsec = TMath::Max(xsec, max_xsec);
517  }
518  break;
519  }
520  }//Q^2
521 
522  // Apply safety factor, since value retrieved from the cache might
523  // correspond to a slightly different energy
524  max_xsec *= fSafetyFactor;
525 
526 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
527  SLOG("QELKinematics", pDEBUG) << interaction->AsString();
528  SLOG("QELKinematics", pDEBUG) << "Max xsec in phase space = " << max_xsec;
529  SLOG("QELKinematics", pDEBUG) << "Computed using alg = " << *fXSecModel;
530 #endif
531 
532  return max_xsec;
533 }
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
static const double kASmallNum
Definition: Controls.h:40
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 QELKinematicsGenerator::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 430 of file QELKinematicsGenerator.cxx.

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

431 {
432  Algorithm::Configure(config);
433  this->LoadConfig();
434 }
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
void QELKinematicsGenerator::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 436 of file QELKinematicsGenerator.cxx.

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

437 {
438  Algorithm::Configure(config);
439  this->LoadConfig();
440 }
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
void QELKinematicsGenerator::LoadConfig ( void  )
private

Definition at line 442 of file QELKinematicsGenerator.cxx.

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

Referenced by Configure().

443 {
444 // Load sub-algorithms and config data to reduce the number of registry
445 // lookups
446 
447  //-- Safety factor for the maximum differential cross section
448  GetParamDef( "MaxXSec-SafetyFactor", fSafetyFactor , 1.25 ) ;
449 
450  //-- Minimum energy for which max xsec would be cached, forcing explicit
451  // calculation for lower eneries
452  GetParamDef( "Cache-MinEnergy", fEMin, 1.00 ) ;
453 
454  //-- Maximum allowed fractional cross section deviation from maxim cross
455  // section used in rejection method
456  GetParamDef( "MaxXSec-DiffTolerance", fMaxXSecDiffTolerance, 999999. ) ;
457  assert(fMaxXSecDiffTolerance>=0);
458 
459  //-- Generate kinematics uniformly over allowed phase space and compute
460  // an event weight?
461  GetParamDef( "UniformOverPhaseSpace", fGenerateUniformly, false ) ;
462 
463 }
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 QELKinematicsGenerator::ProcessEventRecord ( GHepRecord event_rec) const
virtual

Implements genie::EventRecordVisitorI.

Definition at line 54 of file QELKinematicsGenerator.cxx.

References genie::KineGeneratorWithCache::AssertXSecLimits(), genie::XclsTag::CharmHadronPdg(), genie::Kinematics::ClearRunningValues(), genie::EventGeneratorI::CrossSectionAlg(), genie::GHepRecord::EventFlags(), genie::Interaction::ExclTag(), genie::KineGeneratorWithCache::fGenerateUniformly, genie::PDGLibrary::Find(), genie::KineGeneratorWithCache::fXSecModel, gQ2, genie::GHepRecord::HitNucleon(), genie::Target::HitNucP4(), genie::Interaction::InitState(), genie::Interaction::InitStatePtr(), genie::RunningThreadInfo::Instance(), genie::RandomGen::Instance(), genie::PDGLibrary::Instance(), genie::XclsTag::IsCharmEvent(), genie::XclsTag::IsStrangeEvent(), genie::utils::mec::J(), genie::controls::kASmallNum, 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::Target::SetHitNucPosition(), genie::Kinematics::SetQ2(), genie::exceptions::EVGThreadException::SetReason(), genie::Kinematics::SetW(), genie::GHepRecord::SetWeight(), genie::Kinematics::Setx(), genie::Kinematics::Sety(), genie::XclsTag::StrangeHadronPdg(), genie::GHepRecord::Summary(), genie::exceptions::EVGThreadException::SwitchOnFastForward(), genie::InitialState::Tgt(), genie::InitialState::TgtPtr(), genie::GHepRecord::Weight(), genie::utils::kinematics::WQ2toXY(), genie::GHepParticle::X4(), genie::XSecAlgorithmI::XSec(), and genie::GHepRecord::XSec().

55 {
56  if(fGenerateUniformly) {
57  LOG("QELKinematics", pNOTICE)
58  << "Generating kinematics uniformly over the allowed phase space";
59  }
60 
61  //-- Get the random number generators
63 
64  //-- Access cross section algorithm for running thread
66  const EventGeneratorI * evg = rtinfo->RunningThread();
67  fXSecModel = evg->CrossSectionAlg();
68 
69  //-- Get the interaction and set the 'trust' bits
70  Interaction * interaction = evrec->Summary();
71  interaction->SetBit(kISkipProcessChk);
72  interaction->SetBit(kISkipKinematicChk);
73 
74  // store the struck nucleon position for use by the xsec method
75  double hitNucPos = evrec->HitNucleon()->X4()->Vect().Mag();
76  interaction->InitStatePtr()->TgtPtr()->SetHitNucPosition(hitNucPos);
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  Range1D_t Q2 = kps.Limits(kKVQ2);
90 
91  if(Q2.max <=0 || Q2.min>=Q2.max) {
92  LOG("QELKinematics", 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  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  double Q2min = Q2.min+kASmallNum;
112  double Q2max = Q2.max-kASmallNum;
113 //double QD2min = utils::kinematics::Q2toQD2(Q2min);
114 //double QD2max = utils::kinematics::Q2toQD2(Q2max);
115  double xsec = -1.;
116  double gQ2 = 0.;
117 
118  unsigned int iter = 0;
119  bool accept = false;
120  while(1) {
121  iter++;
122  if(iter > kRjMaxIterations) {
123  LOG("QELKinematics", pWARN)
124  << "Couldn't select a valid Q^2 after " << iter << " iterations";
125  evrec->EventFlags()->SetBitNumber(kKineGenErr, true);
127  exception.SetReason("Couldn't select kinematics");
128  exception.SwitchOnFastForward();
129  throw exception;
130  }
131 
132  //-- Generate a Q2 value within the allowed phase space
133 /*
134  if(fGenerateUniformly) {
135  gQ2 = Q2min + (Q2max-Q2min) * rnd->RndKine().Rndm();
136  } else {
137  // In unweighted mode - use transform that takes out the dipole form
138  double gQD2 = QD2min + (QD2max-QD2min) * rnd->RndKine().Rndm();
139  gQ2 = utils::kinematics::QD2toQ2(gQD2);
140  }
141 */
142  gQ2 = Q2min + (Q2max-Q2min) * rnd->RndKine().Rndm();
143  interaction->KinePtr()->SetQ2(gQ2);
144  LOG("QELKinematics", pINFO) << "Trying: Q^2 = " << gQ2;
145 
146  //-- Computing cross section for the current kinematics
147  xsec = fXSecModel->XSec(interaction, kPSQ2fE);
148 
149  //-- Decide whether to accept the current kinematics
150  if(!fGenerateUniformly) {
151  this->AssertXSecLimits(interaction, xsec, xsec_max);
152 
153  double t = xsec_max * rnd->RndKine().Rndm();
154  //double J = kinematics::Jacobian(interaction,kPSQ2fE,kPSQD2fE);
155  double J = 1.;
156 
157 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
158  LOG("QELKinematics", pDEBUG)
159  << "xsec= " << xsec << ", J= " << J << ", Rnd= " << t;
160 #endif
161  accept = (t < J*xsec);
162  } else {
163  accept = (xsec>0);
164  }
165 
166  //-- If the generated kinematics are accepted, finish-up module's job
167  if(accept) {
168  LOG("QELKinematics", pINFO) << "Selected: Q^2 = " << gQ2;
169 
170  // reset bits
171  interaction->ResetBit(kISkipProcessChk);
172  interaction->ResetBit(kISkipKinematicChk);
173  interaction->ResetBit(kIAssumeFreeNucleon);
174 
175  // compute the rest of the kinematical variables
176 
177  // get neutrino energy at struck nucleon rest frame and the
178  // struck nucleon mass (can be off the mass shell)
179  const InitialState & init_state = interaction->InitState();
180  double E = init_state.ProbeE(kRfHitNucRest);
181  double M = init_state.Tgt().HitNucP4().M();
182 
183  LOG("QELKinematics", pNOTICE) << "E = " << E << ", M = "<< M;
184 
185  // The hadronic inv. mass is equal to the recoil nucleon on-shell mass.
186  // For QEL/Charm events it is set to be equal to the on-shell mass of
187  // the generated charm baryon (Lamda_c+, Sigma_c+ or Sigma_c++)
188  // Similarly for strange baryons
189  //
190  const XclsTag & xcls = interaction->ExclTag();
191  int rpdgc = 0;
192  if(xcls.IsCharmEvent()) { rpdgc = xcls.CharmHadronPdg(); }
193  else if(xcls.IsStrangeEvent()) { rpdgc = xcls.StrangeHadronPdg(); }
194  else { rpdgc = interaction->RecoilNucleonPdg(); }
195  assert(rpdgc);
196  double gW = PDGLibrary::Instance()->Find(rpdgc)->Mass();
197 
198  LOG("QELKinematics", pNOTICE) << "Selected: W = "<< gW;
199 
200  // (W,Q2) -> (x,y)
201  double gx=0, gy=0;
202  kinematics::WQ2toXY(E,M,gW,gQ2,gx,gy);
203 
204  // set the cross section for the selected kinematics
205  evrec->SetDiffXSec(xsec,kPSQ2fE);
206 
207  // for uniform kinematics, compute an event weight as
208  // wght = (phase space volume)*(differential xsec)/(event total xsec)
209  if(fGenerateUniformly) {
210  double vol = kinematics::PhaseSpaceVolume(interaction,kPSQ2fE);
211  double totxsec = evrec->XSec();
212  double wght = (vol/totxsec)*xsec;
213  LOG("QELKinematics", pNOTICE) << "Kinematics wght = "<< wght;
214 
215  // apply computed weight to the current event weight
216  wght *= evrec->Weight();
217  LOG("QELKinematics", pNOTICE) << "Current event wght = " << wght;
218  evrec->SetWeight(wght);
219  }
220 
221  // lock selected kinematics & clear running values
222  interaction->KinePtr()->SetQ2(gQ2, true);
223  interaction->KinePtr()->SetW (gW, true);
224  interaction->KinePtr()->Setx (gx, true);
225  interaction->KinePtr()->Sety (gy, true);
226  interaction->KinePtr()->ClearRunningValues();
227 
228  return;
229  }
230  }// iterations
231 }
const KPhaseSpace & PhaseSpace(void) const
Definition: Interaction.h:73
bool fGenerateUniformly
uniform over allowed phase space + event weight?
double J(double q0, double q3, double Enu, double ml)
Definition: MECUtils.cxx:147
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
int CharmHadronPdg(void) const
Definition: XclsTag.h:52
bool IsStrangeEvent(void) const
Definition: XclsTag.h:53
Defines the EventGeneratorI interface.
virtual double MaxXSec(GHepRecord *evrec, const int nkey=0) const
void SetHitNucPosition(double r)
Definition: Target.cxx:210
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.
Contains minimal information for tagging exclusive processes.
Definition: XclsTag.h:39
double PhaseSpaceVolume(const Interaction *const i, KinePhaseSpace_t ps)
Definition: KineUtils.cxx:36
bool IsCharmEvent(void) const
Definition: XclsTag.h:50
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...
int StrangeHadronPdg(void) const
Definition: XclsTag.h:55
#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
static const double kASmallNum
Definition: Controls.h:40
#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
const XclsTag & ExclTag(void) const
Definition: Interaction.h:72
Target * TgtPtr(void) const
Definition: InitialState.h:67
virtual const XSecAlgorithmI * CrossSectionAlg(void) const =0
double gQ2
Definition: gtestDISSF.cxx:55
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
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
void QELKinematicsGenerator::SpectralFuncExperimentalCode ( GHepRecord event_rec) const
private

Definition at line 233 of file QELKinematicsGenerator.cxx.

References genie::KineGeneratorWithCache::AssertXSecLimits(), genie::XclsTag::CharmHadronPdg(), genie::Kinematics::ClearRunningValues(), genie::EventGeneratorI::CrossSectionAlg(), genie::GHepRecord::EventFlags(), genie::Interaction::ExclTag(), genie::PDGLibrary::Find(), genie::KineGeneratorWithCache::fXSecModel, gQ2, genie::GHepRecord::HitNucleon(), genie::Target::HitNucMass(), genie::Target::HitNucP4(), genie::Target::HitNucP4Ptr(), genie::Interaction::InitState(), genie::Interaction::InitStatePtr(), genie::RunningThreadInfo::Instance(), genie::RandomGen::Instance(), genie::PDGLibrary::Instance(), genie::XclsTag::IsCharmEvent(), genie::controls::kASmallNum, 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(), pNOTICE, genie::InitialState::ProbeE(), pWARN, genie::utils::kinematics::Q2(), genie::Interaction::RecoilNucleonPdg(), genie::RandomGen::RndKine(), genie::RunningThreadInfo::RunningThread(), genie::GHepRecord::SetDiffXSec(), genie::Target::SetHitNucPosition(), genie::Kinematics::SetQ2(), genie::exceptions::EVGThreadException::SetReason(), genie::Kinematics::SetW(), genie::Kinematics::Setx(), genie::Kinematics::Sety(), genie::GHepRecord::Summary(), genie::exceptions::EVGThreadException::SwitchOnFastForward(), genie::InitialState::Tgt(), genie::InitialState::TgtPtr(), genie::utils::kinematics::WQ2toXY(), genie::GHepParticle::X4(), and genie::XSecAlgorithmI::XSec().

235 {
236  //-- Get the random number generators
237  RandomGen * rnd = RandomGen::Instance();
238 
239  //-- Access cross section algorithm for running thread
241  const EventGeneratorI * evg = rtinfo->RunningThread();
242  fXSecModel = evg->CrossSectionAlg();
243 
244  //-- Get the interaction and set the 'trust' bits
245  Interaction * interaction = new Interaction(*evrec->Summary());
246  interaction->SetBit(kISkipProcessChk);
247  interaction->SetBit(kISkipKinematicChk);
248 
249  // store the struck nucleon position for use by the xsec method
250  double hitNucPos = evrec->HitNucleon()->X4()->Vect().Mag();
251  interaction->InitStatePtr()->TgtPtr()->SetHitNucPosition(hitNucPos);
252 
253  //-- Note: The kinematic generator would be using the free nucleon cross
254  // section (even for nuclear targets) so as not to double-count nuclear
255  // suppression. This assumes that a) the nuclear suppression was turned
256  // on when computing the cross sections for selecting the current event
257  // and that b) if the event turns out to be unphysical (Pauli-blocked)
258  // the next attempted event will be forced to QEL again.
259  // (discussion with Hugh - GENIE/NeuGEN integration workshop - 07APR2006
260  interaction->SetBit(kIAssumeFreeNucleon);
261 
262  //-- Assume scattering off a nucleon on the mass shell (PWIA prescription)
263  double Mn = interaction->InitState().Tgt().HitNucMass(); // PDG mass, take it to be on-shell
264  double pxn = interaction->InitState().Tgt().HitNucP4().Px();
265  double pyn = interaction->InitState().Tgt().HitNucP4().Py();
266  double pzn = interaction->InitState().Tgt().HitNucP4().Pz();
267  double En = interaction->InitState().Tgt().HitNucP4().Energy();
268  double En0 = TMath::Sqrt(pxn*pxn + pyn*pyn + pzn*pzn + Mn*Mn);
269  double Eb = En0 - En;
270  interaction->InitStatePtr()->TgtPtr()->HitNucP4Ptr()->SetE(En0);
271 
272  //-- Get the limits for the generated Q2
273  const KPhaseSpace & kps = interaction->PhaseSpace();
274  Range1D_t Q2 = kps.Limits(kKVQ2);
275 
276  if(Q2.max <=0 || Q2.min>=Q2.max) {
277  LOG("QELKinematics", pWARN) << "No available phase space";
278  evrec->EventFlags()->SetBitNumber(kKineGenErr, true);
280  exception.SetReason("No available phase space");
281  exception.SwitchOnFastForward();
282  throw exception;
283  }
284 
285  //-- For the subsequent kinematic selection with the rejection method:
286  // Calculate the max differential cross section or retrieve it from the
287  // cache. Throw an exception and quit the evg thread if a non-positive
288  // value is found.
289  // If the kinematics are generated uniformly over the allowed phase
290  // space the max xsec is irrelevant
291 // double xsec_max = (fGenerateUniformly) ? -1 : this->MaxXSec(evrec);
292  double xsec_max = this->MaxXSec(evrec);
293 
294  // get neutrino energy at struck nucleon rest frame and the
295  // struck nucleon mass (can be off the mass shell)
296  const InitialState & init_state = interaction->InitState();
297  double E = init_state.ProbeE(kRfHitNucRest);
298 
299  LOG("QELKinematics", pNOTICE) << "E = " << E << ", M = "<< Mn;
300 
301  //-- Try to select a valid Q2 using the rejection method
302 
303  // kinematical limits
304  double Q2min = Q2.min+kASmallNum;
305  double Q2max = Q2.max-kASmallNum;
306  double xsec = -1.;
307  double gQ2 = 0.;
308  double gW = 0.;
309  double gx = 0.;
310  double gy = 0.;
311 
312  unsigned int iter = 0;
313  bool accept = false;
314  while(1) {
315  iter++;
316  if(iter > kRjMaxIterations) {
317  LOG("QELKinematics", pWARN)
318  << "Couldn't select a valid Q^2 after " << iter << " iterations";
319  evrec->EventFlags()->SetBitNumber(kKineGenErr, true);
321  exception.SetReason("Couldn't select kinematics");
322  exception.SwitchOnFastForward();
323  throw exception;
324  }
325 
326  //-- Generate a Q2 value within the allowed phase space
327  gQ2 = Q2min + (Q2max-Q2min) * rnd->RndKine().Rndm();
328  LOG("QELKinematics", pNOTICE) << "Trying: Q^2 = " << gQ2;
329 
330  // The hadronic inv. mass is equal to the recoil nucleon on-shell mass.
331  // For QEL/Charm events it is set to be equal to the on-shell mass of
332  // the generated charm baryon (Lamda_c+, Sigma_c+ or Sigma_c++)
333  //
334  const XclsTag & xcls = interaction->ExclTag();
335  int rpdgc = 0;
336  if(xcls.IsCharmEvent()) { rpdgc = xcls.CharmHadronPdg(); }
337  else { rpdgc = interaction->RecoilNucleonPdg(); }
338  assert(rpdgc);
339  gW = PDGLibrary::Instance()->Find(rpdgc)->Mass();
340 
341  // (W,Q2) -> (x,y)
342  kinematics::WQ2toXY(E,Mn,gW,gQ2,gx,gy);
343 
344  LOG("QELKinematics", pNOTICE) << "W = "<< gW;
345  LOG("QELKinematics", pNOTICE) << "x = "<< gx;
346  LOG("QELKinematics", pNOTICE) << "y = "<< gy;
347 
348  // v
349  double gv = gy * E;
350  double gv2 = gv*gv;
351 
352  LOG("QELKinematics", pNOTICE) << "v = "<< gv;
353 
354  // v -> v~
355  double gvtilde = gv + Mn - Eb - TMath::Sqrt(Mn*Mn+pxn*pxn+pyn*pyn+pzn*pzn);
356  double gvtilde2 = gvtilde*gvtilde;
357 
358  LOG("QELKinematics", pNOTICE) << "v~ = "<< gvtilde;
359 
360  // Q~^2
361  double gQ2tilde = gQ2 - gv2 + gvtilde2;
362 
363  LOG("QELKinematics", pNOTICE) << "Q~^2 = "<< gQ2tilde;
364 
365  // Set updated Q2
366  interaction->KinePtr()->SetQ2(gQ2tilde);
367 
368  //-- Computing cross section for the current kinematics
369  xsec = fXSecModel->XSec(interaction, kPSQ2fE);
370 
371  //-- Decide whether to accept the current kinematics
372 // if(!fGenerateUniformly) {
373  this->AssertXSecLimits(interaction, xsec, xsec_max);
374 
375  double t = xsec_max * rnd->RndKine().Rndm();
376 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
377  LOG("QELKinematics", pDEBUG)
378  << "xsec= " << xsec << ", Rnd= " << t;
379 #endif
380  accept = (t < xsec);
381 // } else {
382 // accept = (xsec>0);
383 // }
384 
385  //-- If the generated kinematics are accepted, finish-up module's job
386  if(accept) {
387  LOG("QELKinematics", pNOTICE) << "Selected: Q^2 = " << gQ2;
388 
389  // reset bits
390 // interaction->ResetBit(kISkipProcessChk);
391 // interaction->ResetBit(kISkipKinematicChk);
392 // interaction->ResetBit(kIAssumeFreeNucleon);
393 
394  // set the cross section for the selected kinematics
395  evrec->SetDiffXSec(xsec,kPSQ2fE);
396 
397  // for uniform kinematics, compute an event weight as
398  // wght = (phase space volume)*(differential xsec)/(event total xsec)
399 // if(fGenerateUniformly) {
400 // double vol = kinematics::PhaseSpaceVolume(interaction,kPSQ2fE);
401 // double totxsec = evrec->XSec();
402 // double wght = (vol/totxsec)*xsec;
403 // LOG("QELKinematics", pNOTICE) << "Kinematics wght = "<< wght;
404 
405  // apply computed weight to the current event weight
406 // wght *= evrec->Weight();
407 // LOG("QELKinematics", pNOTICE) << "Current event wght = " << wght;
408 // evrec->SetWeight(wght);
409 // }
410 
411  // lock selected kinematics & clear running values
412 // interaction->KinePtr()->SetQ2(gQ2, true);
413 // interaction->KinePtr()->SetW (gW, true);
414 // interaction->KinePtr()->Setx (gx, true);
415 // interaction->KinePtr()->Sety (gy, true);
416 // interaction->KinePtr()->ClearRunningValues();
417 
418  evrec->Summary()->KinePtr()->SetQ2(gQ2, true);
419  evrec->Summary()->KinePtr()->SetW (gW, true);
420  evrec->Summary()->KinePtr()->Setx (gx, true);
421  evrec->Summary()->KinePtr()->Sety (gy, true);
422  evrec->Summary()->KinePtr()->ClearRunningValues();
423  delete interaction;
424 
425  return;
426  }
427  }// iterations
428 }
const KPhaseSpace & PhaseSpace(void) const
Definition: Interaction.h:73
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
double HitNucMass(void) const
Definition: Target.cxx:233
int CharmHadronPdg(void) const
Definition: XclsTag.h:52
Defines the EventGeneratorI interface.
virtual double MaxXSec(GHepRecord *evrec, const int nkey=0) const
void SetHitNucPosition(double r)
Definition: Target.cxx:210
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.
Contains minimal information for tagging exclusive processes.
Definition: XclsTag.h:39
bool IsCharmEvent(void) const
Definition: XclsTag.h:50
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
static const double kASmallNum
Definition: Controls.h:40
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.
static RunningThreadInfo * Instance(void)
double max
Definition: Range1.h:53
TLorentzVector * HitNucP4Ptr(void) const
Definition: Target.cxx:247
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:68
const UInt_t kIAssumeFreeNucleon
Definition: Interaction.h:49
const XclsTag & ExclTag(void) const
Definition: Interaction.h:72
Target * TgtPtr(void) const
Definition: InitialState.h:67
virtual const XSecAlgorithmI * CrossSectionAlg(void) const =0
double gQ2
Definition: gtestDISSF.cxx:55
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
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
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: