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

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

#include <DMDISKinematicsGenerator.h>

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

Public Member Functions

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

Private Member Functions

void LoadConfig (void)
 
double ComputeMaxXSec (const Interaction *interaction) 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 DIS DM interaction events. Is a concrete implementation of the EventRecordVisitorI interface.

      Part of its implementation, related with the caching and retrieval of
      previously computed values, is inherited from the KineGeneratorWithCache
      abstract class.
Author
Joshua Berger <jberger physics.wisc.edu> University of Wisconsin-Madison Costas Andreopoulos <c.andreopoulos cern.ch> University of Liverpool
Created:
September 1, 2017
License:
Copyright (c) 2003-2024, The GENIE Collaboration For the full text of the license visit http://copyright.genie-mc.org

Definition at line 35 of file DMDISKinematicsGenerator.h.

Constructor & Destructor Documentation

DMDISKinematicsGenerator::DMDISKinematicsGenerator ( )

Definition at line 41 of file DMDISKinematicsGenerator.cxx.

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

Definition at line 47 of file DMDISKinematicsGenerator.cxx.

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

Definition at line 53 of file DMDISKinematicsGenerator.cxx.

54 {
55 
56 }

Member Function Documentation

double DMDISKinematicsGenerator::ComputeMaxXSec ( const Interaction interaction) const
privatevirtual

Implements genie::KineGeneratorWithCache.

Definition at line 241 of file DMDISKinematicsGenerator.cxx.

References genie::Interaction::AsString(), genie::KineGeneratorWithCache::fXSecModel, genie::Interaction::InitState(), genie::Interaction::KinePtr(), genie::kKVx, genie::kKVy, genie::kPSxyfE, genie::KPhaseSpace::Limits(), LOG, genie::Range1D_t::max, genie::Range1D_t::min, pDEBUG, genie::Interaction::PhaseSpace(), pINFO, genie::Kinematics::Setx(), genie::Kinematics::Sety(), SLOG, genie::InitialState::Tgt(), genie::utils::kinematics::UpdateWQ2FromXY(), and genie::XSecAlgorithmI::XSec().

243 {
244 // Computes the maximum differential cross section in the requested phase
245 // space. This method overloads KineGeneratorWithCache::ComputeMaxXSec
246 // method and the value is cached at a circular cache branch for retrieval
247 // during subsequent event generation.
248 // The computed max differential cross section does not need to be the exact
249 // maximum. The number used in the rejection method will be scaled up by a
250 // safety factor. But this needs to be fast - do not use a very fine grid.
251 
252 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
253  LOG("DMDISKinematics", pDEBUG)<< "Computing max xsec in allowed phase space";
254 #endif
255  double max_xsec = 0.0;
256 
257  const InitialState & init_state = interaction->InitState();
258  //double Ev = init_state.ProbeE(kRfHitNucRest);
259  //const ProcessInfo & proc = interaction->ProcInfo();
260  const Target & tgt = init_state.Tgt();
261 
262  int Ny = 20;
263  int Nx = 40;
264  int Nxb = 3;
265 
266  const KPhaseSpace & kps = interaction->PhaseSpace();
267  Range1D_t xl = kps.Limits(kKVx);
268  Range1D_t yl = kps.Limits(kKVy);
269 
270  // Look over full region or risk not finding valid phase space
271  double xmin = TMath::Max(xl.min, 5E-3);
272  double xmax = xl.max;
273  double ymin = TMath::Max(yl.min, 2E-3);
274  double ymax = yl.max;
275  double dx = (xmax-xmin)/(Nx-1);
276  double dy = (ymax-ymin)/(Ny-1);
277 
278 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
279  LOG("DMDISKinematics", pDEBUG)
280  << "Searching max. in x [" << xmin << ", " << xmax << "], y [" << ymin << ", " << ymax << "]";
281 #endif
282  double xseclast_y = -1;
283  bool increasing_y;
284 
285  for(int i=0; i<Ny; i++) {
286  double gy = ymin + i*dy;
287  //double gy = TMath::Power(10., logymin + i*dlogy);
288  interaction->KinePtr()->Sety(gy);
289 
290 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
291  LOG("DMDISKinematics", pDEBUG) << "y = " << gy;
292 #endif
293  double xseclast_x = -1;
294  bool increasing_x;
295 
296  for(int j=0; j<Nx; j++) {
297  double gx = xmin + j*dx;
298  //double gx = TMath::Power(10., logxmin + j*dlogx);
299  interaction->KinePtr()->Setx(gx);
300  kinematics::UpdateWQ2FromXY(interaction);
301 
302  double xsec = fXSecModel->XSec(interaction, kPSxyfE);
303 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
304  LOG("DMDISKinematics", pINFO)
305  << "xsec(y=" << gy << ", x=" << gx << ") = " << xsec;
306 #endif
307  // update maximum xsec
308  max_xsec = TMath::Max(xsec, max_xsec);
309 
310  increasing_x = xsec-xseclast_x>=0;
311  xseclast_x = xsec;
312 
313  // once the cross section stops increasing, I reduce the step size and
314  // step backwards a little bit to handle cases that the max cross section
315  // is grossly underestimated (very peaky distribution & large step)
316  if(!increasing_x) {
317 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
318  LOG("DMDISKinematics", pDEBUG)
319  << "d2xsec/dxdy|x stopped increasing. Stepping back & exiting x loop";
320 #endif
321  //double dlogxn = dlogx/(Nxb+1);
322  double dxn = dx/(Nxb+1);
323  for(int ik=0; ik<Nxb; ik++) {
324  //gx = TMath::Exp(TMath::Log(gx) - dlogxn);
325  gx = gx - dxn;
326  interaction->KinePtr()->Setx(gx);
327  kinematics::UpdateWQ2FromXY(interaction);
328  xsec = fXSecModel->XSec(interaction, kPSxyfE);
329 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
330  LOG("DMDISKinematics", pINFO)
331  << "xsec(y=" << gy << ", x=" << gx << ") = " << xsec;
332 #endif
333  }
334  break;
335  } // stepping back within last bin
336  } // x
337  increasing_y = max_xsec-xseclast_y>=0;
338  xseclast_y = max_xsec;
339  if(!increasing_y) {
340 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
341  LOG("DMDISKinematics", pDEBUG)
342  << "d2xsec/dxdy stopped increasing. Exiting y loop";
343 #endif
344  break;
345  }
346  }// y
347 
348  // Apply safety factor, since value retrieved from the cache might
349  // correspond to a slightly different energy
350  // max_xsec *= fSafetyFactor;
351  //max_xsec *= ( (Ev<3.0) ? 2.5 : fSafetyFactor);
352  max_xsec *= 3;
353 
354 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
355  SLOG("DMDISKinematics", pDEBUG) << interaction->AsString();
356  SLOG("DMDISKinematics", pDEBUG) << "Max xsec in phase space = " << max_xsec;
357  SLOG("DMDISKinematics", pDEBUG) << "Computed using alg = " << *fXSecModel;
358 #endif
359 
360  return max_xsec;
361 }
const KPhaseSpace & PhaseSpace(void) const
Definition: Interaction.h:73
Kinematics * KinePtr(void) const
Definition: Interaction.h:76
A simple [min,max] interval for doubles.
Definition: Range1.h:42
Range1D_t Limits(KineVar_t kvar) const
Return the kinematical variable limits.
string AsString(void) const
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
Kinematical phase space.
Definition: KPhaseSpace.h:33
#define pINFO
Definition: Messenger.h:62
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
void UpdateWQ2FromXY(const Interaction *in)
Definition: KineUtils.cxx:1290
double max
Definition: Range1.h:53
void Sety(double y, bool selected=false)
Definition: Kinematics.cxx:243
const InitialState & InitState(void) const
Definition: Interaction.h:69
double min
Definition: Range1.h:52
const Target & Tgt(void) const
Definition: InitialState.h:66
#define SLOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a short string (using the FUNCTION and...
Definition: Messenger.h:84
Initial State information.
Definition: InitialState.h:48
#define pDEBUG
Definition: Messenger.h:63
void DMDISKinematicsGenerator::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 206 of file DMDISKinematicsGenerator.cxx.

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

207 {
208  Algorithm::Configure(config);
209  this->LoadConfig();
210 }
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
void DMDISKinematicsGenerator::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 212 of file DMDISKinematicsGenerator.cxx.

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

213 {
214  Algorithm::Configure(config);
215  this->LoadConfig();
216 }
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
void DMDISKinematicsGenerator::LoadConfig ( void  )
private

Definition at line 218 of file DMDISKinematicsGenerator.cxx.

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

Referenced by Configure().

219 {
220 // Reads its configuration data from its configuration Registry and loads them
221 // in private data members to avoid looking up at the Registry all the time.
222 
223  //-- Safety factor for the maximum differential cross section
224  GetParamDef( "MaxXSec-SafetyFactor", fSafetyFactor, 1.25 ) ;
225 
226  //-- Minimum energy for which max xsec would be cached, forcing explicit
227  // calculation for lower eneries
228  GetParamDef( "Cache-MinEnergy", fEMin, 0.8 ) ;
229 
230  //-- Maximum allowed fractional cross section deviation from maxim cross
231  // section used in rejection method
232  GetParamDef( "MaxXSec-DiffTolerance", fMaxXSecDiffTolerance, 999999. ) ;
233  assert(fMaxXSecDiffTolerance>=0);
234 
235  //-- Generate kinematics uniformly over allowed phase space and compute
236  // an event weight?
237  GetParamDef( "UniformOverPhaseSpace", fGenerateUniformly, false ) ;
238 
239 }
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 DMDISKinematicsGenerator::ProcessEventRecord ( GHepRecord event_rec) const
virtual

Implements genie::EventRecordVisitorI.

Definition at line 58 of file DMDISKinematicsGenerator.cxx.

References genie::KineGeneratorWithCache::AssertXSecLimits(), genie::Kinematics::ClearRunningValues(), genie::EventGeneratorI::CrossSectionAlg(), genie::GHepRecord::EventFlags(), genie::KineGeneratorWithCache::fGenerateUniformly, genie::KineGeneratorWithCache::fXSecModel, gQ2, genie::Target::HitNucP4(), genie::Interaction::InitState(), genie::RunningThreadInfo::Instance(), genie::RandomGen::Instance(), genie::utils::mec::J(), genie::Interaction::KinePtr(), genie::kISkipKinematicChk, genie::kISkipProcessChk, genie::kKineGenErr, genie::kKVW, genie::kKVx, genie::kKVy, genie::kPSxyfE, 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(), pNOTICE, genie::InitialState::ProbeE(), pWARN, genie::Kinematics::Q2(), 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::utils::kinematics::UpdateWQ2FromXY(), genie::Kinematics::W(), genie::utils::kinematics::W(), genie::GHepRecord::Weight(), genie::XSecAlgorithmI::XSec(), genie::GHepRecord::XSec(), and genie::utils::kinematics::XYtoWQ2().

59 {
60  if(fGenerateUniformly) {
61  LOG("DMDISKinematics", pNOTICE)
62  << "Generating kinematics uniformly over the allowed phase space";
63  }
64 
65  //-- Get the random number generators
67 
68  //-- Access cross section algorithm for running thread
70  const EventGeneratorI * evg = rtinfo->RunningThread();
71  fXSecModel = evg->CrossSectionAlg();
72 
73  //-- Get the interaction
74  Interaction * interaction = evrec->Summary();
75  interaction->SetBit(kISkipProcessChk);
76 
77  //-- Get dark matter energy and hit 'nucleon mass'
78  const InitialState & init_state = interaction->InitState();
79  double Ev = init_state.ProbeE(kRfHitNucRest);
80  double M = init_state.Tgt().HitNucP4().M(); // can be off m-shell
81 
82  //-- Get the physical W range
83  const KPhaseSpace & kps = interaction->PhaseSpace();
84  Range1D_t W = kps.Limits(kKVW);
85  if(W.max <=0 || W.min>=W.max) {
86  LOG("DMDISKinematics", pWARN) << "No available phase space";
87  evrec->EventFlags()->SetBitNumber(kKineGenErr, true);
89  exception.SetReason("No available phase space");
90  exception.SwitchOnFastForward();
91  throw exception;
92  }
93 
94  Range1D_t xl = kps.Limits(kKVx);
95  Range1D_t yl = kps.Limits(kKVy);
96 
97  LOG("DMDISKinematics", pNOTICE) << "x: [" << xl.min << ", " << xl.max << "]";
98  LOG("DMDISKinematics", pNOTICE) << "y: [" << yl.min << ", " << yl.max << "]";
99 
100  assert(xl.min>0 && yl.min>0);
101 
102  //-- For the subsequent kinematic selection with the rejection method:
103  // Calculate the max differential cross section or retrieve it from the
104  // cache. Throw an exception and quit the evg thread if a non-positive
105  // value is found.
106  // If the kinematics are generated uniformly over the allowed phase
107  // space the max xsec is irrelevant
108  double xsec_max = (fGenerateUniformly) ? -1 : this->MaxXSec(evrec);
109 
110  //-- Try to select a valid (x,y) pair using the rejection method
111 
112  double dx = xl.max - xl.min;
113  double dy = yl.max - yl.min;
114  double gx=-1, gy=-1, gW=-1, gQ2=-1, xsec=-1;
115 
116  unsigned int iter = 0;
117  bool accept = false;
118  while(1) {
119  iter++;
120  if(iter > kRjMaxIterations) {
121  LOG("DMDISKinematics", pWARN)
122  << " Couldn't select kinematics after " << iter << " iterations";
123  evrec->EventFlags()->SetBitNumber(kKineGenErr, true);
125  exception.SetReason("Couldn't select kinematics");
126  exception.SwitchOnFastForward();
127  throw exception;
128  }
129 
130  //-- random x,y
131  gx = xl.min + dx * rnd->RndKine().Rndm();
132  gy = yl.min + dy * rnd->RndKine().Rndm();
133  interaction->KinePtr()->Setx(gx);
134  interaction->KinePtr()->Sety(gy);
135  kinematics::UpdateWQ2FromXY(interaction);
136 
137  LOG("DMDISKinematics", pNOTICE)
138  << "Trying: x = " << gx << ", y = " << gy
139  << " (W = " << interaction->KinePtr()->W() << ","
140  << " (Q2 = " << interaction->KinePtr()->Q2() << ")";
141 
142  //-- compute the cross section for current kinematics
143  xsec = fXSecModel->XSec(interaction, kPSxyfE);
144 
145  //-- decide whether to accept the current kinematics
146  if(!fGenerateUniformly) {
147  this->AssertXSecLimits(interaction, xsec, xsec_max);
148  double t = xsec_max * rnd->RndKine().Rndm();
149  double J = 1;
150 
151 #ifdef __GENIE_LOW_LEVEL_MESG_ENABLED__
152  LOG("DMDISKinematics", pDEBUG)
153  << "xsec= " << xsec << ", J= " << J << ", Rnd= " << t;
154 #endif
155  accept = (t < J*xsec);
156  }
157  else {
158  accept = (xsec>0);
159  }
160 
161  //-- If the generated kinematics are accepted, finish-up module's job
162  if(accept) {
163  LOG("DMDISKinematics", pNOTICE)
164  << "Selected: x = " << gx << ", y = " << gy
165  << " (W = " << interaction->KinePtr()->W() << ","
166  << " (Q2 = " << interaction->KinePtr()->Q2() << ")";
167 
168  // reset trust bits
169  interaction->ResetBit(kISkipProcessChk);
170  interaction->ResetBit(kISkipKinematicChk);
171 
172  // set the cross section for the selected kinematics
173  evrec->SetDiffXSec(xsec,kPSxyfE);
174 
175  // for uniform kinematics, compute an event weight as
176  // wght = (phase space volume)*(differential xsec)/(event total xsec)
177  if(fGenerateUniformly) {
178  double vol = kinematics::PhaseSpaceVolume(interaction,kPSxyfE);
179  double totxsec = evrec->XSec();
180  double wght = (vol/totxsec)*xsec;
181  LOG("DMDISKinematics", pNOTICE) << "Kinematics wght = "<< wght;
182 
183  // apply computed weight to the current event weight
184  wght *= evrec->Weight();
185  LOG("DMDISKinematics", pNOTICE) << "Current event wght = " << wght;
186  evrec->SetWeight(wght);
187  }
188 
189  // compute W,Q2 for selected x,y
190  kinematics::XYtoWQ2(Ev,M,gW,gQ2,gx,gy);
191 
192  LOG("DMDISKinematics", pNOTICE)
193  << "Selected x,y => W = " << gW << ", Q2 = " << gQ2;
194 
195  // lock selected kinematics & clear running values
196  interaction->KinePtr()->SetW (gW, true);
197  interaction->KinePtr()->SetQ2(gQ2, true);
198  interaction->KinePtr()->Setx (gx, true);
199  interaction->KinePtr()->Sety (gy, true);
200  interaction->KinePtr()->ClearRunningValues();
201  return;
202  }
203  } // iterations
204 }
double W(bool selected=false) const
Definition: Kinematics.cxx:157
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
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
A simple [min,max] interval for doubles.
Definition: Range1.h:42
virtual void AssertXSecLimits(const Interaction *in, double xsec, double xsec_max) const
Defines the EventGeneratorI interface.
virtual double MaxXSec(GHepRecord *evrec, const int nkey=0) const
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition: RandomGen.h:29
Range1D_t Limits(KineVar_t kvar) const
Return the kinematical variable limits.
double PhaseSpaceVolume(const Interaction *const i, KinePhaseSpace_t ps)
Definition: KineUtils.cxx:36
double W(const Interaction *const i)
Definition: KineUtils.cxx:1101
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 XYtoWQ2(double Ev, double M, double &W, double &Q2, double x, double y)
Definition: KineUtils.cxx:1155
Kinematical phase space.
Definition: KPhaseSpace.h:33
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
void UpdateWQ2FromXY(const Interaction *in)
Definition: KineUtils.cxx:1290
static RunningThreadInfo * Instance(void)
double max
Definition: Range1.h:53
void SetW(double W, bool selected=false)
Definition: Kinematics.cxx:279
void Sety(double y, bool selected=false)
Definition: Kinematics.cxx:243
virtual const XSecAlgorithmI * CrossSectionAlg(void) const =0
double gQ2
Definition: gtestDISSF.cxx:55
const InitialState & InitState(void) const
Definition: Interaction.h:69
double min
Definition: Range1.h:52
double Q2(bool selected=false) const
Definition: Kinematics.cxx:125
#define pNOTICE
Definition: Messenger.h:61
void ClearRunningValues(void)
Definition: Kinematics.cxx:347
const Target & Tgt(void) const
Definition: InitialState.h:66
static const unsigned int kRjMaxIterations
Definition: Controls.h:26
const EventGeneratorI * RunningThread(void)
double ProbeE(RefFrame_t rf) const
Keep info on the event generation thread currently on charge. This is used so that event generation m...
const UInt_t kISkipProcessChk
if set, skip process validity checks
Definition: Interaction.h:47
Initial State information.
Definition: InitialState.h:48
#define pDEBUG
Definition: Messenger.h:63

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