GENIEGenerator
Loading...
Searching...
No Matches
genie::RESKinematicsGenerator Class Reference

Generates resonance event (v+N->l+Resonance) kinematics. Is a concrete implementation of the EventRecordVisitorI interface. More...

#include <RESKinematicsGenerator.h>

Inheritance diagram for genie::RESKinematicsGenerator:
[legend]
Collaboration diagram for genie::RESKinematicsGenerator:
[legend]

Public Member Functions

 RESKinematicsGenerator ()
 RESKinematicsGenerator (string config)
 ~RESKinematicsGenerator ()
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 LoadConfig (void)
double ComputeMaxXSec (const Interaction *interaction) const

Private Attributes

TF2 * fEnvelope
 2-D envelope used for importance sampling
double fWcut
 Wcut parameter in DIS/RES join scheme.

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 resonance event (v+N->l+Resonance) kinematics. Is a concrete implementation of the EventRecordVisitorI interface.

Author
Costas Andreopoulos <c.andreopoulos \at cern.ch> University of Liverpool
Created:\n November 18, 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 29 of file RESKinematicsGenerator.h.

Constructor & Destructor Documentation

◆ RESKinematicsGenerator() [1/2]

RESKinematicsGenerator::RESKinematicsGenerator ( )

Definition at line 38 of file RESKinematicsGenerator.cxx.

38 :
39KineGeneratorWithCache("genie::RESKinematicsGenerator")
40{
41 fEnvelope = 0;
42}
TF2 * fEnvelope
2-D envelope used for importance sampling

References fEnvelope, and genie::KineGeneratorWithCache::KineGeneratorWithCache().

◆ RESKinematicsGenerator() [2/2]

RESKinematicsGenerator::RESKinematicsGenerator ( string config)

Definition at line 44 of file RESKinematicsGenerator.cxx.

44 :
45KineGeneratorWithCache("genie::RESKinematicsGenerator", config)
46{
47 fEnvelope = 0;
48}

References fEnvelope, and genie::KineGeneratorWithCache::KineGeneratorWithCache().

◆ ~RESKinematicsGenerator()

RESKinematicsGenerator::~RESKinematicsGenerator ( )

Definition at line 50 of file RESKinematicsGenerator.cxx.

51{
52 if(fEnvelope) delete fEnvelope;
53}

References fEnvelope.

Member Function Documentation

◆ ComputeMaxXSec()

double RESKinematicsGenerator::ComputeMaxXSec ( const Interaction * interaction) const
privatevirtual

Implements genie::KineGeneratorWithCache.

Definition at line 273 of file RESKinematicsGenerator.cxx.

275{
276// Computes the maximum differential cross section in the requested phase
277// space. This method overloads KineGeneratorWithCache::ComputeMaxXSec
278// method and the value is cached at a circular cache branch for retrieval
279// during subsequent event generation.
280// The computed max differential cross section does not need to be the exact
281// maximum. The number used in the rejection method will be scaled up by a
282// safety factor. But this needs to be fast - do not use a very fine grid.
283
284 double max_xsec = 0.;
285
286 const InitialState & init_state = interaction -> InitState();
287 double E = init_state.ProbeE(kRfHitNucRest);
288 bool is_em = interaction->ProcInfo().IsEM();
290
291 double md;
292 if(!interaction->ExclTag().KnownResonance()) md=1.23;
293 else {
294 Resonance_t res = interaction->ExclTag().Resonance();
295 md=res::Mass(res);
296 }
297
298 // ** 2-D Scan
299 const KPhaseSpace & kps = interaction->PhaseSpace();
300 Range1D_t rW = kps.WLim();
301
302 int NW = 20;
303 double Wmin = rW.min + kASmallNum;
304 double Wmax = rW.max - kASmallNum;
305
306 Wmax = TMath::Min(Wmax,fWcut);
307
308 Wmin = TMath::Max(Wmin, md-.3);
309 Wmax = TMath::Min(Wmax, md+.3);
310
311 if(Wmax-Wmin<0.05) { NW=1; Wmin=Wmax; }
312
313 double dW = (NW>1) ? (Wmax-Wmin)/(NW-1) : 0.;
314
315 for(int iw=0; iw<NW; iw++)
316 {
317 double W = Wmin + iw*dW;
318 interaction->KinePtr()->SetW(W);
319
320 int NQ2 = 25;
321 int NQ2b = 4;
322
323 Range1D_t rQ2 = kps.Q2Lim_W();
324 if( rQ2.max < Q2Thres || rQ2.min <=0 ) continue;
325 if( rQ2.max-rQ2.min<0.02 ) {NQ2=5; NQ2b=3;}
326
327 double logQ2min = TMath::Log(rQ2.min+kASmallNum);
328 double logQ2max = TMath::Log(rQ2.max-kASmallNum);
329 double dlogQ2 = (logQ2max - logQ2min) /(NQ2-1);
330 double xseclast = -1;
331 bool increasing = true;
332
333 for(int iq2=0; iq2<NQ2; iq2++)
334 {
335 double Q2 = TMath::Exp(logQ2min + iq2 * dlogQ2);
336 interaction->KinePtr()->SetQ2(Q2);
337 double xsec = fXSecModel->XSec(interaction, kPSWQD2fE);
338 LOG("RESKinematics", pDEBUG)
339 << "xsec(W= " << W << ", Q2= " << Q2 << ") = " << xsec;
340 max_xsec = TMath::Max(xsec, max_xsec);
341 increasing = xsec-xseclast>=0;
342 xseclast=xsec;
343
344 // once the cross section stops increasing, I reduce the step size and
345 // step backwards a little bit to handle cases that the max cross section
346 // is grossly underestimated (very peaky distribution & large step)
347 if(!increasing)
348 {
349 dlogQ2/=NQ2b;
350 for(int iq2b=0; iq2b<NQ2b; iq2b++)
351 {
352 Q2 = TMath::Exp(TMath::Log(Q2) - dlogQ2);
353 if(Q2 < rQ2.min) continue;
354 interaction->KinePtr()->SetQ2(Q2);
355 xsec = fXSecModel->XSec(interaction, kPSWQD2fE);
356 LOG("RESKinematics", pDEBUG)
357 << "xsec(W= " << W << ", Q2= " << Q2 << ") = " << xsec;
358 max_xsec = TMath::Max(xsec, max_xsec);
359 }
360 break;
361 }
362 } // Q2
363 }//W
364
365 // Apply safety factor, since value retrieved from the cache might
366 // correspond to a slightly different energy
367 // Apply larger safety factor for smaller energies.
368 max_xsec *= ( (E<md) ? 2. : fSafetyFactor);
369
370 return max_xsec;
371}
#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
double ProbeE(RefFrame_t rf) const
const XclsTag & ExclTag(void) const
Definition Interaction.h:72
const ProcessInfo & ProcInfo(void) const
Definition Interaction.h:70
const KPhaseSpace & PhaseSpace(void) const
Definition Interaction.h:73
Kinematics * KinePtr(void) const
Definition Interaction.h:76
Range1D_t Q2Lim_W(void) const
Q2 limits @ fixed W.
Range1D_t WLim(void) const
W limits.
double fSafetyFactor
ComputeMaxXSec -> ComputeMaxXSec * fSafetyFactor.
void SetQ2(double Q2, bool selected=false)
void SetW(double W, bool selected=false)
bool IsEM(void) const
double fWcut
Wcut parameter in DIS/RES join scheme.
Resonance_t Resonance(void) const
Definition XclsTag.h:69
bool KnownResonance(void) const
Definition XclsTag.h:68
static const double kMinQ2Limit
Definition Controls.h:41
static const double kASmallNum
Definition Controls.h:40
double W(const Interaction *const i)
double Q2(const Interaction *const i)
double Mass(Resonance_t res)
resonance mass (GeV)
enum genie::EResonance Resonance_t
@ kRfHitNucRest
Definition RefFrame.h:30

References genie::Interaction::ExclTag(), genie::KineGeneratorWithCache::fSafetyFactor, fWcut, genie::KineGeneratorWithCache::fXSecModel, genie::ProcessInfo::IsEM(), genie::controls::kASmallNum, genie::Interaction::KinePtr(), genie::controls::kMinQ2Limit, genie::utils::kinematics::electromagnetic::kMinQ2Limit, genie::XclsTag::KnownResonance(), genie::kPSWQD2fE, genie::kRfHitNucRest, LOG, genie::utils::res::Mass(), genie::Range1D_t::max, genie::Range1D_t::min, pDEBUG, genie::Interaction::PhaseSpace(), genie::InitialState::ProbeE(), genie::Interaction::ProcInfo(), genie::KPhaseSpace::Q2Lim_W(), genie::XclsTag::Resonance(), genie::Kinematics::SetQ2(), genie::Kinematics::SetW(), and genie::KPhaseSpace::WLim().

◆ Configure() [1/2]

void RESKinematicsGenerator::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 231 of file RESKinematicsGenerator.cxx.

232{
233 Algorithm::Configure(config);
234 this->LoadConfig();
235}
virtual void Configure(const Registry &config)
Definition Algorithm.cxx:62

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

◆ Configure() [2/2]

void RESKinematicsGenerator::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 237 of file RESKinematicsGenerator.cxx.

238{
239 Algorithm::Configure(config);
240 this->LoadConfig();
241}

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

◆ LoadConfig()

void RESKinematicsGenerator::LoadConfig ( void )
private

Definition at line 243 of file RESKinematicsGenerator.cxx.

244{
245 // Safety factor for the maximum differential cross section
246 this->GetParamDef("MaxXSec-SafetyFactor", fSafetyFactor, 1.25);
247
248 // Minimum energy for which max xsec would be cached, forcing explicit
249 // calculation for lower eneries
250 this->GetParamDef("Cache-MinEnergy", fEMin, 0.5);
251
252 // Load Wcut used in DIS/RES join scheme
253 this->GetParam("Wcut", fWcut);
254
255 // Maximum allowed fractional cross section deviation from maxim cross
256 // section used in rejection method
257 this->GetParamDef("MaxXSec-DiffTolerance", fMaxXSecDiffTolerance, 999999.);
258 assert(fMaxXSecDiffTolerance>=0);
259
260 // Generate kinematics uniformly over allowed phase space and compute
261 // an event weight?
262 this->GetParamDef("UniformOverPhaseSpace", fGenerateUniformly, false);
263
264 // Envelope employed when importance sampling is used
265 // (initialize with dummy range)
266 if(fEnvelope) delete fEnvelope;
267 fEnvelope = new TF2("res-envelope",
269 // stop ROOT from deleting this object of its own volition
270 gROOT->GetListOfFunctions()->Remove(fEnvelope);
271}
bool GetParam(const RgKey &name, T &p, bool is_top_call=true) const
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.
double RESImportanceSamplingEnvelope(double *x, double *par)

References genie::KineGeneratorWithCache::fEMin, fEnvelope, genie::KineGeneratorWithCache::fGenerateUniformly, genie::KineGeneratorWithCache::fMaxXSecDiffTolerance, genie::KineGeneratorWithCache::fSafetyFactor, fWcut, genie::Algorithm::GetParam(), genie::Algorithm::GetParamDef(), and genie::utils::kinematics::RESImportanceSamplingEnvelope().

Referenced by Configure(), and Configure().

◆ ProcessEventRecord()

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

Implements genie::EventRecordVisitorI.

Definition at line 55 of file RESKinematicsGenerator.cxx.

56{
58 LOG("RESKinematics", pNOTICE)
59 << "Generating kinematics uniformly over the allowed phase space";
60 }
61
62 //-- Get the random number generators
63 RandomGen * rnd = RandomGen::Instance();
64
65 //-- Access cross section algorithm for running thread
66 RunningThreadInfo * rtinfo = RunningThreadInfo::Instance();
67 const EventGeneratorI * evg = rtinfo->RunningThread();
69
70 //-- Get the interaction from the GHEP record
71 Interaction * interaction = evrec->Summary();
72 interaction->SetBit(kISkipProcessChk);
73
74 //-- EM or CC/NC process (different importance sampling methods)
75 bool is_em = interaction->ProcInfo().IsEM();
76
77 //-- Compute the W limits
78 // (the physically allowed W's, unless an external cut is imposed)
79 const KPhaseSpace & kps = interaction->PhaseSpace();
80 Range1D_t W = kps.Limits(kKVW);
81
82 if(W.max <=0 || W.min>=W.max) {
83 LOG("RESKinematics", pWARN) << "No available phase space";
84 evrec->EventFlags()->SetBitNumber(kKineGenErr, true);
85 genie::exceptions::EVGThreadException exception;
86 exception.SetReason("No available phase space");
87 exception.SwitchOnFastForward();
88 throw exception;
89 }
90
91 const InitialState & init_state = interaction -> InitState();
92 double E = init_state.ProbeE(kRfHitNucRest);
93
94 //-- For the subsequent kinematic selection with the rejection method:
95 // Calculate the max differential cross section or retrieve it from the
96 // cache. Throw an exception and quit the evg thread if a non-positive
97 // value is found.
98 // If the kinematics are generated uniformly over the allowed phase
99 // space the max xsec is irrelevant
100 double xsec_max = (fGenerateUniformly) ? -1 : this->MaxXSec(evrec);
101
102 //-- Try to select a valid W, Q2 pair using the rejection method
103 double dW = W.max - W.min;
104 double xsec = -1;
105
106 unsigned int iter = 0;
107 bool accept = false;
108 while(1) {
109 iter++;
110 if(iter > kRjMaxIterations) {
111 LOG("RESKinematics", pWARN)
112 << "*** Could not select a valid (W,Q^2) pair after "
113 << iter << " iterations";
114 evrec->EventFlags()->SetBitNumber(kKineGenErr, true);
115 genie::exceptions::EVGThreadException exception;
116 exception.SetReason("Couldn't select kinematics");
117 exception.SwitchOnFastForward();
118 throw exception;
119 }
120
121 double gW = 0; // current hadronic invariant mass
122 double gQ2 = 0; // current momentum transfer
123 double gQD2 = 0; // tranformed Q2 to take out dipole form
124
126 {
127 //-- Generate a W uniformly in the kinematically allowed range.
128 // For the generated W, compute the Q2 range and generate a value
129 // uniformly over that range
130 gW = W.min + dW * rnd->RndKine().Rndm();
131 interaction->KinePtr()->SetW(gW);
132 Range1D_t Q2 = kps.Q2Lim_W();
133 if(Q2.max<=0. || Q2.min>=Q2.max) continue;
134 gQ2 = Q2.min + (Q2.max-Q2.min) * rnd->RndKine().Rndm();
135 interaction->SetBit(kISkipKinematicChk);
136 }
137 else
138 {
139 // neutrino scattering
140 // Selecting unweighted event kinematics using an importance sampling
141 // method. Q2 with be transformed to QD2 to take out the dipole form.
142 interaction->KinePtr()->SetW(W.min);
143 Range1D_t Q2 = kps.Q2Lim_W();
144 double Q2min = -99.;
145 if (is_em)
146 Q2min = Q2.min + kASmallNum;
147 else
148 Q2min = 0 + kASmallNum;
149 double Q2max = Q2.max - kASmallNum;
150
151 // In unweighted mode - use transform that takes out the dipole form
152 double QD2min = utils::kinematics::Q2toQD2(Q2max);
153 double QD2max = utils::kinematics::Q2toQD2(Q2min);
154
155 gW = W.min + dW * rnd->RndKine().Rndm();
156 gQD2 = QD2min + (QD2max - QD2min) * rnd->RndKine().Rndm();
157
158 // QD2 -> Q2
160 } // uniformly over phase space?
161
162 LOG("RESKinematics", pINFO) << "Trying: W = " << gW << ", Q2 = " << gQ2;
163
164 //-- Set kinematics for current trial
165 interaction->KinePtr()->SetW(gW);
166 interaction->KinePtr()->SetQ2(gQ2);
167
168 //-- Computing cross section for the current kinematics
169 xsec = fXSecModel->XSec(interaction, kPSWQD2fE);
170 //-- Decide whether to accept the current kinematics
172 {
173 // unified neutrino / electron scattering
174 double t = xsec_max * rnd->RndKine().Rndm();
175 this->AssertXSecLimits(interaction, xsec, xsec_max);
176 accept = (t < xsec);
177 } // charged lepton or neutrino scattering?
178 else
179 {
180 accept = (xsec>0);
181 } // uniformly over phase space
182
183 //-- If the generated kinematics are accepted, finish-up module's job
184 if(accept) {
185 LOG("RESKinematics", pINFO)
186 << "Selected: W = " << gW << ", Q2 = " << gQ2;
187 // reset 'trust' bits
188 interaction->ResetBit(kISkipProcessChk);
189 interaction->ResetBit(kISkipKinematicChk);
190
191 // compute x,y for selected W,Q2
192 // note: hit nucleon can be off the mass-shell
193 double gx=-1, gy=-1;
194 double M = init_state.Tgt().HitNucP4().M();
195 kinematics::WQ2toXY(E,M,gW,gQ2,gx,gy);
196
197 // set the cross section for the selected kinematics.
198 // we're converting here to the more familiar "W*Q2" space
199 // rather than the "W*Q2D" (precomputed dipole) space that's used above
200 // for generation efficiency in the accept-reject loop
201 double J = kinematics::Jacobian(interaction, kPSWQD2fE, kPSWQ2fE);
202 xsec *= J;
203 evrec->SetDiffXSec(xsec, kPSWQ2fE);
204
205 // for uniform kinematics, compute an event weight as
206 // wght = (phase space volume)*(differential xsec)/(event total xsec)
208 double vol = kinematics::PhaseSpaceVolume(interaction,kPSWQ2fE);
209 double totxsec = evrec->XSec();
210 double wght = (vol/totxsec)*xsec;
211 LOG("RESKinematics", pNOTICE) << "Kinematics wght = "<< wght;
212
213 // apply computed weight to the current event weight
214 wght *= evrec->Weight();
215 LOG("RESKinematics", pNOTICE) << "Current event wght = " << wght;
216 evrec->SetWeight(wght);
217 }
218
219 // lock selected kinematics & clear running values
220 interaction->KinePtr()->SetQ2(gQ2, true);
221 interaction->KinePtr()->SetW (gW, true);
222 interaction->KinePtr()->Setx (gx, true);
223 interaction->KinePtr()->Sety (gy, true);
224 interaction->KinePtr()->ClearRunningValues();
225
226 return;
227 } // accept
228 } // iterations
229}
#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
Range1D_t Limits(KineVar_t kvar) const
Return the kinematical variable limits.
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 ClearRunningValues(void)
void Sety(double y, bool selected=false)
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)
const TLorentzVector & HitNucP4(void) const
Definition Target.h:91
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 QD2toQ2(double QD2)
double PhaseSpaceVolume(const Interaction *const i, KinePhaseSpace_t ps)
Definition KineUtils.cxx:36
double Jacobian(const Interaction *const i, KinePhaseSpace_t f, KinePhaseSpace_t t)
double Q2toQD2(double Q2)
double J(double q0, double q3, double Enu, double ml)
Definition MECUtils.cxx:147
@ kKVW
Definition KineVar.h:35
const UInt_t kISkipKinematicChk
if set, skip kinematic validity checks
Definition Interaction.h:48
const UInt_t kISkipProcessChk
if set, skip process validity checks
Definition Interaction.h:47
@ kKineGenErr
Definition GHepFlags.h:31

References genie::KineGeneratorWithCache::AssertXSecLimits(), genie::Kinematics::ClearRunningValues(), genie::EventGeneratorI::CrossSectionAlg(), genie::GHepRecord::EventFlags(), genie::KineGeneratorWithCache::fGenerateUniformly, genie::KineGeneratorWithCache::fXSecModel, gQ2, genie::Target::HitNucP4(), genie::RandomGen::Instance(), genie::RunningThreadInfo::Instance(), genie::ProcessInfo::IsEM(), genie::utils::kinematics::Jacobian(), genie::controls::kASmallNum, genie::Interaction::KinePtr(), genie::kISkipKinematicChk, genie::kISkipProcessChk, genie::kKineGenErr, genie::kKVW, genie::kPSWQ2fE, genie::kPSWQD2fE, genie::kRfHitNucRest, genie::controls::kRjMaxIterations, genie::KPhaseSpace::Limits(), LOG, genie::KineGeneratorWithCache::MaxXSec(), genie::Interaction::PhaseSpace(), genie::utils::kinematics::PhaseSpaceVolume(), pINFO, pNOTICE, genie::InitialState::ProbeE(), genie::Interaction::ProcInfo(), pWARN, genie::KPhaseSpace::Q2Lim_W(), genie::utils::kinematics::Q2toQD2(), genie::utils::kinematics::QD2toQ2(), 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(), and genie::GHepRecord::XSec().

Member Data Documentation

◆ fEnvelope

TF2* genie::RESKinematicsGenerator::fEnvelope
mutableprivate

2-D envelope used for importance sampling

Definition at line 48 of file RESKinematicsGenerator.h.

Referenced by LoadConfig(), RESKinematicsGenerator(), RESKinematicsGenerator(), and ~RESKinematicsGenerator().

◆ fWcut

double genie::RESKinematicsGenerator::fWcut
private

Wcut parameter in DIS/RES join scheme.

Definition at line 49 of file RESKinematicsGenerator.h.

Referenced by ComputeMaxXSec(), and LoadConfig().


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