GENIEGenerator
Loading...
Searching...
No Matches
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:
[legend]
Collaboration diagram for genie::QELKinematicsGenerator:
[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.
virtual AlgStatus_t GetStatus (void) const
 Get algorithm status.
virtual bool AllowReconfig (void) const
virtual AlgCmp_t Compare (const Algorithm *alg) const
 Compare with input algorithm.
virtual void SetId (const AlgId &id)
 Set algorithm ID.
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.

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.
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.
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.
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
int AddLowRegistry (Registry *rp, bool owns=true)
 add registry with lowest priority, also update ownership
int MergeTopRegistry (const Registry &r)
int AddTopRegisties (const vector< Registry * > &rs, bool owns=false)
 Add registries with top priority, also udated Ownerships.
Protected Attributes inherited from genie::KineGeneratorWithCache
const XSecAlgorithmIfXSecModel
double fSafetyFactor
 ComputeMaxXSec -> ComputeMaxXSec * fSafetyFactor.
std::vector< double > vSafetyFactors
 MaxXSec -> MaxXSec * fSafetyFactors[nkey].
int fNumOfSafetyFactors
 Number of given safety factors.
std::vector< string > vInterpolatorTypes
 Type of interpolator for each key in a branch.
int fNumOfInterpolatorTypes
 Number of given interpolators types.
double fMaxXSecDiffTolerance
 max{100*(xsec-maxxsec)/.5*(xsec+maxxsec)} if xsec>maxxsec
double fEMin
 min E for which maxxsec is cached - forcing explicit calc.
bool fGenerateUniformly
 uniform over allowed phase space + event weight?
Protected Attributes inherited from genie::Algorithm
bool fAllowReconfig
bool fOwnsSubstruc
 true if it owns its substructure (sub-algs,...)
AlgId fID
 algorithm name and configuration set
vector< Registry * > fConfVect
vector< bool > fOwnerships
 ownership for every registry in fConfVect
AlgStatus_t fStatus
 algorithm execution status
AlgMapfOwnedSubAlgMp
 local pool for owned sub-algs (taken out of the factory pool)

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 \at cern.ch> University of Liverpool
Created:\n October 03, 2004
License:\n Copyright (c) 2003-2025, 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() [1/2]

QELKinematicsGenerator::QELKinematicsGenerator ( )

Definition at line 37 of file QELKinematicsGenerator.cxx.

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

References genie::KineGeneratorWithCache::KineGeneratorWithCache().

◆ QELKinematicsGenerator() [2/2]

QELKinematicsGenerator::QELKinematicsGenerator ( string config)

Definition at line 43 of file QELKinematicsGenerator.cxx.

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

References genie::KineGeneratorWithCache::KineGeneratorWithCache().

◆ ~QELKinematicsGenerator()

QELKinematicsGenerator::~QELKinematicsGenerator ( )

Definition at line 49 of file QELKinematicsGenerator.cxx.

50{
51
52}

Member Function Documentation

◆ ComputeMaxXSec()

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

Implements genie::KineGeneratorWithCache.

Definition at line 465 of file QELKinematicsGenerator.cxx.

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}
#define pDEBUG
Definition Messenger.h:63
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE,...
Definition Messenger.h:96
#define SLOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a short string (using the FUNCTION and...
Definition Messenger.h:84
Range1D_t Limits(KineVar_t kvar) const
Return the kinematical variable limits.
double fSafetyFactor
ComputeMaxXSec -> ComputeMaxXSec * fSafetyFactor.
static const double kASmallNum
Definition Controls.h:40
double Q2(const Interaction *const i)
@ kKVQ2
Definition KineVar.h:33

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::Kinematics::SetQ2(), and SLOG.

◆ Configure() [1/2]

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.

431{
432 Algorithm::Configure(config);
433 this->LoadConfig();
434}
virtual void Configure(const Registry &config)
Definition Algorithm.cxx:62

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

◆ Configure() [2/2]

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.

437{
438 Algorithm::Configure(config);
439 this->LoadConfig();
440}

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

◆ LoadConfig()

void QELKinematicsGenerator::LoadConfig ( void )
private

Definition at line 442 of file QELKinematicsGenerator.cxx.

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 GetParamDef(const RgKey &name, T &p, const T &def) const
bool fGenerateUniformly
uniform over allowed phase space + event weight?
double fMaxXSecDiffTolerance
max{100*(xsec-maxxsec)/.5*(xsec+maxxsec)} if xsec>maxxsec
double fEMin
min E for which maxxsec is cached - forcing explicit calc.

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

Referenced by Configure(), and Configure().

◆ ProcessEventRecord()

void QELKinematicsGenerator::ProcessEventRecord ( GHepRecord * event_rec) const
virtual

Implements genie::EventRecordVisitorI.

Definition at line 54 of file QELKinematicsGenerator.cxx.

55{
57 LOG("QELKinematics", pNOTICE)
58 << "Generating kinematics uniformly over the allowed phase space";
59 }
60
61 //-- Get the random number generators
62 RandomGen * rnd = RandomGen::Instance();
63
64 //-- Access cross section algorithm for running thread
65 RunningThreadInfo * rtinfo = RunningThreadInfo::Instance();
66 const EventGeneratorI * evg = rtinfo->RunningThread();
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);
94 genie::exceptions::EVGThreadException exception;
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);
126 genie::exceptions::EVGThreadException exception;
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)
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}
#define pNOTICE
Definition Messenger.h:61
#define pINFO
Definition Messenger.h:62
#define pWARN
Definition Messenger.h:60
virtual const XSecAlgorithmI * CrossSectionAlg(void) const =0
const Target & Tgt(void) const
double ProbeE(RefFrame_t rf) const
Target * TgtPtr(void) const
const XclsTag & ExclTag(void) const
Definition Interaction.h:72
InitialState * InitStatePtr(void) const
Definition Interaction.h:74
int RecoilNucleonPdg(void) const
recoil nucleon pdg
const KPhaseSpace & PhaseSpace(void) const
Definition Interaction.h:73
const InitialState & InitState(void) const
Definition Interaction.h:69
Kinematics * KinePtr(void) const
Definition Interaction.h:76
virtual void AssertXSecLimits(const Interaction *in, double xsec, double xsec_max) const
virtual double MaxXSec(GHepRecord *evrec, const int nkey=0) const
void Setx(double x, bool selected=false)
void SetQ2(double Q2, bool selected=false)
void ClearRunningValues(void)
void Sety(double y, bool selected=false)
void SetW(double W, bool selected=false)
static PDGLibrary * Instance(void)
TParticlePDG * Find(int pdgc, bool must_exist=true)
static RandomGen * Instance()
Access instance.
Definition RandomGen.cxx:74
TRandom3 & RndKine(void) const
rnd number generator used by kinematics generators
Definition RandomGen.h:50
static RunningThreadInfo * Instance(void)
const EventGeneratorI * RunningThread(void)
void SetHitNucPosition(double r)
Definition Target.cxx:210
const TLorentzVector & HitNucP4(void) const
Definition Target.h:91
bool IsStrangeEvent(void) const
Definition XclsTag.h:53
bool IsCharmEvent(void) const
Definition XclsTag.h:50
int StrangeHadronPdg(void) const
Definition XclsTag.h:55
int CharmHadronPdg(void) const
Definition XclsTag.h:52
double gQ2
static const unsigned int kRjMaxIterations
Definition Controls.h:26
void WQ2toXY(double Ev, double M, double W, double Q2, double &x, double &y)
double PhaseSpaceVolume(const Interaction *const i, KinePhaseSpace_t ps)
Definition KineUtils.cxx:36
double J(double q0, double q3, double Enu, double ml)
Definition MECUtils.cxx:147
const UInt_t kISkipKinematicChk
if set, skip kinematic validity checks
Definition Interaction.h:48
@ kRfHitNucRest
Definition RefFrame.h:30
const UInt_t kISkipProcessChk
if set, skip process validity checks
Definition Interaction.h:47
@ kKineGenErr
Definition GHepFlags.h:31
const UInt_t kIAssumeFreeNucleon
Definition Interaction.h:49

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::PDGLibrary::Instance(), genie::RandomGen::Instance(), genie::RunningThreadInfo::Instance(), genie::XclsTag::IsCharmEvent(), genie::XclsTag::IsStrangeEvent(), 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::KineGeneratorWithCache::MaxXSec(), pDEBUG, genie::Interaction::PhaseSpace(), genie::utils::kinematics::PhaseSpaceVolume(), pINFO, pNOTICE, genie::InitialState::ProbeE(), pWARN, 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(), and genie::GHepRecord::XSec().

◆ SpectralFuncExperimentalCode()

void QELKinematicsGenerator::SpectralFuncExperimentalCode ( GHepRecord * event_rec) const
private

Definition at line 233 of file QELKinematicsGenerator.cxx.

235{
236 //-- Get the random number generators
237 RandomGen * rnd = RandomGen::Instance();
238
239 //-- Access cross section algorithm for running thread
240 RunningThreadInfo * rtinfo = RunningThreadInfo::Instance();
241 const EventGeneratorI * evg = rtinfo->RunningThread();
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);
279 genie::exceptions::EVGThreadException exception;
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);
320 genie::exceptions::EVGThreadException exception;
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}
TLorentzVector * HitNucP4Ptr(void) const
Definition Target.cxx:247
double HitNucMass(void) const
Definition Target.cxx:233

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::PDGLibrary::Instance(), genie::RandomGen::Instance(), genie::RunningThreadInfo::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::KineGeneratorWithCache::MaxXSec(), pDEBUG, genie::Interaction::PhaseSpace(), pNOTICE, genie::InitialState::ProbeE(), pWARN, 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(), and genie::GHepParticle::X4().


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