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

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

#include <QELEventGeneratorSM.h>

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

Public Member Functions

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

Private Member Functions

void LoadConfig (void)
 
double ComputeMaxXSec (const Interaction *in) const
 
double ComputeMaxXSec (const Interaction *in, const int nkey) const
 
void AddTargetNucleusRemnant (GHepRecord *evrec) const
 add a recoiled nucleus remnant More...
 

Private Attributes

SmithMonizUtilssm_utils
 
KinePhaseSpace_t fkps
 
bool fGenerateNucleonInNucleus
 generate struck nucleon in nucleus More...
 
double fQ2Min
 Q2-threshold for seeking the second maximum. More...
 

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 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 for Smith-Moniz model. Is a concrete implementation of the EventRecordVisitorI interface.

References:
[1] R.A.Smith and E.J.Moniz, Nuclear Physics B43, (1972) 605-622
[2] K.S. Kuzmin, V.V. Lyubushkin, V.A.Naumov Eur. Phys. J. C54, (2008) 517-538
Author
Igor Kakorin kakor.nosp@m.in@j.nosp@m.inr.r.nosp@m.u
Joint Institute for Nuclear Research
adapted from fortran code provided by:

Konstantin Kuzmin kkuzm.nosp@m.in@t.nosp@m.heor..nosp@m.jinr.nosp@m..ru,
Joint Institute for Nuclear Research, Institute for Theoretical and Experimental Physics
Vadim Naumov vnaum.nosp@m.ov@t.nosp@m.heor..nosp@m.jinr.nosp@m..ru,
Joint Institute for Nuclear Research
based on code of: Costas Andreopoulos <c.andreopoulos cern.ch> University of Liverpool

Created:
May 05, 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 47 of file QELEventGeneratorSM.h.

Constructor & Destructor Documentation

QELEventGeneratorSM::QELEventGeneratorSM ( )

Definition at line 65 of file QELEventGeneratorSM.cxx.

65  :
66 KineGeneratorWithCache("genie::QELEventGeneratorSM")
67 {
68 
69 }
QELEventGeneratorSM::QELEventGeneratorSM ( string  config)

Definition at line 71 of file QELEventGeneratorSM.cxx.

71  :
72 KineGeneratorWithCache("genie::QELEventGeneratorSM", config)
73 {
74 
75 }
QELEventGeneratorSM::~QELEventGeneratorSM ( )

Definition at line 77 of file QELEventGeneratorSM.cxx.

78 {
79 
80 }

Member Function Documentation

void QELEventGeneratorSM::AddTargetNucleusRemnant ( GHepRecord evrec) const
private

add a recoiled nucleus remnant

Definition at line 339 of file QELEventGeneratorSM.cxx.

References genie::units::A, genie::GHepParticle::A(), genie::GHepRecord::AddParticle(), genie::GHepParticle::E(), genie::PDGLibrary::Find(), genie::GHepParticle::FirstDaughter(), genie::PDGLibrary::Instance(), genie::pdg::IonPdgCode(), genie::pdg::IsNeutron(), genie::pdg::IsProton(), genie::kIStStableFinalState, genie::GHepParticle::LastDaughter(), LOG, genie::GHepParticle::Mass(), genie::GHepRecord::Particle(), genie::GHepParticle::Pdg(), pFATAL, pINFO, genie::GHepParticle::Px(), genie::GHepParticle::Py(), genie::GHepParticle::Pz(), genie::GHepRecord::TargetNucleus(), genie::GHepRecord::TargetNucleusPosition(), and genie::GHepParticle::Z().

Referenced by ProcessEventRecord().

340 {
341 // add the remnant nuclear target at the GHEP record
342 
343  LOG("QELEvent", pINFO) << "Adding final state nucleus";
344 
345  double Px = 0;
346  double Py = 0;
347  double Pz = 0;
348  double E = 0;
349 
350  GHepParticle * nucleus = evrec->TargetNucleus();
351  int A = nucleus->A();
352  int Z = nucleus->Z();
353 
354  int fd = nucleus->FirstDaughter();
355  int ld = nucleus->LastDaughter();
356 
357  for(int id = fd; id <= ld; id++)
358  {
359 
360  // compute A,Z for final state nucleus & get its PDG code and its mass
361  GHepParticle * particle = evrec->Particle(id);
362  assert(particle);
363  int pdgc = particle->Pdg();
364  bool is_p = pdg::IsProton (pdgc);
365  bool is_n = pdg::IsNeutron(pdgc);
366 
367  if (is_p) Z--;
368  if (is_p || is_n) A--;
369 
370  Px += particle->Px();
371  Py += particle->Py();
372  Pz += particle->Pz();
373  E += particle->E();
374 
375  }//daughters
376 
377  TParticlePDG * remn = 0;
378  int ipdgc = pdg::IonPdgCode(A, Z);
379  remn = PDGLibrary::Instance()->Find(ipdgc);
380  if(!remn)
381  {
382  LOG("HadronicVtx", pFATAL)
383  << "No particle with [A = " << A << ", Z = " << Z
384  << ", pdgc = " << ipdgc << "] in PDGLibrary!";
385  assert(remn);
386  }
387 
388  double Mi = nucleus->Mass();
389  Px *= -1;
390  Py *= -1;
391  Pz *= -1;
392  E = Mi-E;
393 
394  // Add the nucleus to the event record
395  LOG("QELEvent", pINFO)
396  << "Adding nucleus [A = " << A << ", Z = " << Z
397  << ", pdgc = " << ipdgc << "]";
398 
399  int imom = evrec->TargetNucleusPosition();
400  evrec->AddParticle(
401  ipdgc,kIStStableFinalState, imom,-1,-1,-1, Px,Py,Pz,E, 0,0,0,0);
402 
403  LOG("QELEvent", pINFO) << "Done";
404  LOG("QELEvent", pINFO) << *evrec;
405 }
int Z(void) const
virtual GHepParticle * Particle(int position) const
Definition: GHepRecord.cxx:104
double E(void) const
Get energy.
Definition: GHepParticle.h:91
int FirstDaughter(void) const
Definition: GHepParticle.h:68
#define pFATAL
Definition: Messenger.h:56
double Mass(void) const
Mass that corresponds to the PDG code.
double Pz(void) const
Get Pz.
Definition: GHepParticle.h:90
double Px(void) const
Get Px.
Definition: GHepParticle.h:88
int Pdg(void) const
Definition: GHepParticle.h:63
bool IsNeutron(int pdgc)
Definition: PDGUtils.cxx:341
int LastDaughter(void) const
Definition: GHepParticle.h:69
bool IsProton(int pdgc)
Definition: PDGUtils.cxx:336
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
static constexpr double A
Definition: Units.h:74
virtual GHepParticle * TargetNucleus(void) const
Definition: GHepRecord.cxx:293
#define pINFO
Definition: Messenger.h:62
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:68
int IonPdgCode(int A, int Z)
Definition: PDGUtils.cxx:71
virtual void AddParticle(const GHepParticle &p)
Definition: GHepRecord.cxx:499
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:86
int A(void) const
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
virtual int TargetNucleusPosition(void) const
Definition: GHepRecord.cxx:370
double Py(void) const
Get Py.
Definition: GHepParticle.h:89
double QELEventGeneratorSM::ComputeMaxXSec ( const Interaction in) const
privatevirtual

Implements genie::KineGeneratorWithCache.

Definition at line 456 of file QELEventGeneratorSM.cxx.

References genie::Target::A(), e, fkps, fQ2Min, genie::KineGeneratorWithCache::fXSecModel, genie::SmithMonizUtils::GetFermiMomentum(), genie::Interaction::KinePtr(), genie::kKVPn, genie::kKVQ2, genie::kKVv, genie::Range1D_t::max, genie::Range1D_t::min, genie::utils::kinematics::Q2(), genie::SmithMonizUtils::Q2QES_SM_lim(), genie::Kinematics::SetKV(), sm_utils, genie::InitialState::TgtPtr(), genie::SmithMonizUtils::vQES_SM_lim(), genie::SmithMonizUtils::vQES_SM_max(), genie::SmithMonizUtils::vQES_SM_min(), and genie::XSecAlgorithmI::XSec().

Referenced by ComputeMaxXSec().

457 {
458  double xsec_max = -1;
459  double tmp_xsec_max = -1;
460  double Q20, v0;
461  const int N_Q2 = 32;
462  const InitialState & init_state = interaction -> InitState();
463  Target * tgt = init_state.TgtPtr();
464  bool isHeavyNucleus = tgt->A()>3;
465  int N_v = isHeavyNucleus?32:0;
466 
468  const double logQ2min = TMath::Log(TMath::Max(rQ2.min, eps));
469  const double logQ2max = TMath::Log(TMath::Min(rQ2.max, fQ2Min));
470  Kinematics * kinematics = interaction->KinePtr();
471  const double pFmax = sm_utils->GetFermiMomentum();
472  // Now scan through kinematical variables Q2,v,kF
473  for (int Q2_n=0; Q2_n <= N_Q2; Q2_n++)
474  {
475  // Scan around Q2
476  double Q2 = TMath::Exp(Q2_n*(logQ2max-logQ2min)/N_Q2 + logQ2min);
477  kinematics->SetKV(kKVQ2, Q2);
478  Range1D_t rv = sm_utils->vQES_SM_lim(Q2);
479  const double logvmin = TMath::Log(TMath::Max(rv.min, eps));
480  const double logvmax = TMath::Log(TMath::Max(rv.max, TMath::Max(rv.min, eps)));
481  for (int v_n=0; v_n <= N_v; v_n++)
482  {
483  // Scan around v
484  double v = TMath::Exp(v_n*(logvmax-logvmin)/N_v + logvmin);
485  kinematics->SetKV(kKVv, v);
486  kinematics->SetKV(kKVPn, pFmax);
487  // Compute the QE cross section for the current kinematics
488  double xs = fXSecModel->XSec(interaction, fkps);
489  if (xs > tmp_xsec_max)
490  {
491  tmp_xsec_max = xs;
492  Q20 = Q2;
493  v0 = v;
494  }
495  } // Done with v scan
496  }// Done with Q2 scan
497 
498  const double Q2min = rQ2.min;
499  const double Q2max = TMath::Min(rQ2.max, fQ2Min);
500  const double vmin = sm_utils->vQES_SM_min(Q2min, Q2max);
501  const double vmax = sm_utils->vQES_SM_max(Q2min, Q2max);
502 
503  ROOT::Math::Minimizer * min = ROOT::Math::Factory::CreateMinimizer("Minuit", "Minimize");
504  ROOT::Math::IBaseFunctionMultiDim * f = isHeavyNucleus?static_cast<ROOT::Math::IBaseFunctionMultiDim*>(new d3XSecSM_dQ2dvdkF_E(fXSecModel, interaction, pFmax)):
505  static_cast<ROOT::Math::IBaseFunctionMultiDim*>(new d1XSecSM_dQ2_E(fXSecModel, interaction));
506  min->SetFunction( *f );
507  min->SetMaxFunctionCalls(10000); // for Minuit/Minuit2
508  min->SetMaxIterations(10000); // for GSL
509  min->SetTolerance(0.001);
510  min->SetPrintLevel(0);
511  double step = 1e-7;
512  min->SetVariable(0, "Q2", Q20, step);
513  min->SetVariableLimits(0, Q2min, Q2max);
514  if (isHeavyNucleus)
515  {
516  min->SetVariable(1, "v", v0, step);
517  min->SetVariableLimits(1, vmin, vmax);
518  }
519  min->Minimize();
520  xsec_max = -min->MinValue();
521  if (tmp_xsec_max > xsec_max)
522  {
523  xsec_max = tmp_xsec_max;
524  }
525 
526  return xsec_max;
527 
528 }
Range1D_t Q2QES_SM_lim(void) const
double vQES_SM_max(double Q2min, double Q2max) const
double fQ2Min
Q2-threshold for seeking the second maximum.
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:1077
int A(void) const
Definition: Target.h:70
A simple [min,max] interval for doubles.
Definition: Range1.h:42
Generated/set kinematical variables for an event.
Definition: Kinematics.h:39
Range1D_t vQES_SM_lim(double Q2) const
const double e
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
double vQES_SM_min(double Q2min, double Q2max) const
virtual double XSec(const Interaction *i, KinePhaseSpace_t k=kPSfE) const =0
Compute the cross section for the input interaction.
void SetKV(KineVar_t kv, double value)
Definition: Kinematics.cxx:335
double max
Definition: Range1.h:53
Target * TgtPtr(void) const
Definition: InitialState.h:67
double min
Definition: Range1.h:52
double GetFermiMomentum(void) const
Initial State information.
Definition: InitialState.h:48
double QELEventGeneratorSM::ComputeMaxXSec ( const Interaction in,
const int  nkey 
) const
privatevirtual

Reimplemented from genie::KineGeneratorWithCache.

Definition at line 530 of file QELEventGeneratorSM.cxx.

References genie::Target::A(), ComputeMaxXSec(), e, fkps, fQ2Min, genie::KineGeneratorWithCache::fXSecModel, genie::SmithMonizUtils::GetFermiMomentum(), genie::Interaction::KinePtr(), genie::kKVPn, genie::kKVQ2, genie::kKVv, genie::Range1D_t::max, genie::Range1D_t::min, genie::utils::kinematics::Q2(), genie::SmithMonizUtils::Q2QES_SM_lim(), genie::Kinematics::SetKV(), sm_utils, genie::InitialState::TgtPtr(), genie::SmithMonizUtils::vQES_SM_lim(), genie::SmithMonizUtils::vQES_SM_max(), genie::SmithMonizUtils::vQES_SM_min(), and genie::XSecAlgorithmI::XSec().

531 {
532  switch (nkey){
533  case 0:
534  return this->ComputeMaxXSec(interaction);
535 
536  case 1:
537  {
539  if (rQ2.max<fQ2Min)
540  {
541  return -1.;
542  }
543  double xsec_max = -1;
544  double tmp_xsec_max = -1;
545  double Q20, v0;
546  const int N_Q2 = 32;
547  const InitialState & init_state = interaction -> InitState();
548  Target * tgt = init_state.TgtPtr();
549  bool isHeavyNucleus = tgt->A()>3;
550  int N_v = isHeavyNucleus?32:0;
551 
552  const double logQ2min = TMath::Log(fQ2Min);
553  const double logQ2max = TMath::Log(rQ2.max);
554  Kinematics * kinematics = interaction->KinePtr();
555  const double pFmax = sm_utils->GetFermiMomentum();
556  // Now scan through kinematical variables Q2,v,kF
557  for (int Q2_n=0; Q2_n <= N_Q2; Q2_n++)
558  {
559  // Scan around Q2
560  double Q2 = TMath::Exp(Q2_n*(logQ2max-logQ2min)/N_Q2 + logQ2min);
561  kinematics->SetKV(kKVQ2, Q2);
562  Range1D_t rv = sm_utils->vQES_SM_lim(Q2);
563  const double logvmin = TMath::Log(TMath::Max(rv.min, eps));
564  const double logvmax = TMath::Log(TMath::Max(rv.max, TMath::Max(rv.min, eps)));
565  for (int v_n=0; v_n <= N_v; v_n++)
566  {
567  // Scan around v
568  double v = TMath::Exp(v_n*(logvmax-logvmin)/N_v + logvmin);
569  kinematics->SetKV(kKVv, v);
570  kinematics->SetKV(kKVPn, pFmax);
571  // Compute the QE cross section for the current kinematics
572  double xs = fXSecModel->XSec(interaction, fkps);
573  if (xs > tmp_xsec_max)
574  {
575  tmp_xsec_max = xs;
576  Q20 = Q2;
577  v0 = v;
578  }
579  } // Done with v scan
580  }// Done with Q2 scan
581 
582  const double Q2min = fQ2Min;
583  const double Q2max = rQ2.max;
584  const double vmin = sm_utils->vQES_SM_min(Q2min, Q2max);
585  const double vmax = sm_utils->vQES_SM_max(Q2min, Q2max);
586  ROOT::Math::Minimizer * min = ROOT::Math::Factory::CreateMinimizer("Minuit", "Minimize");
587  ROOT::Math::IBaseFunctionMultiDim * f = isHeavyNucleus?static_cast<ROOT::Math::IBaseFunctionMultiDim*>(new d3XSecSM_dQ2dvdkF_E(fXSecModel, interaction, pFmax)):
588  static_cast<ROOT::Math::IBaseFunctionMultiDim*>(new d1XSecSM_dQ2_E(fXSecModel, interaction));
589  min->SetFunction( *f );
590  min->SetMaxFunctionCalls(10000); // for Minuit/Minuit2
591  min->SetMaxIterations(10000); // for GSL
592  min->SetTolerance(0.001);
593  min->SetPrintLevel(0);
594  double step = 1e-7;
595  min->SetVariable(0, "Q2", Q20, step);
596  min->SetVariableLimits(0, Q2min, Q2max);
597  if (isHeavyNucleus)
598  {
599  min->SetVariable(1, "v", v0, step);
600  min->SetVariableLimits(1, vmin, vmax);
601  }
602  min->Minimize();
603  xsec_max = -min->MinValue();
604  if (tmp_xsec_max > xsec_max)
605  {
606  xsec_max = tmp_xsec_max;
607  }
608  return xsec_max;
609  }
610  case 2:
611  {
612  double diffv_max = -1;
613  double tmp_diffv_max = -1;
614  const int N_Q2 = 100;
615  double Q20;
617  for (int Q2_n = 0; Q2_n<=N_Q2; Q2_n++) // Scan around Q2
618  {
619  double Q2 = rQ2.min + 1.*Q2_n * (rQ2.max-rQ2.min)/N_Q2;
620  Range1D_t rv = sm_utils->vQES_SM_lim(Q2);
621  if (rv.max-rv.min > tmp_diffv_max)
622  {
623  tmp_diffv_max = rv.max-rv.min;
624  Q20 = Q2;
625  }
626  } // Done with Q2 scan
627 
628  ROOT::Math::Minimizer * min = ROOT::Math::Factory::CreateMinimizer("Minuit", "Minimize");
629  ROOT::Math::IBaseFunctionMultiDim * f = new dv_dQ2_E(interaction);
630  min->SetFunction( *f );
631  min->SetMaxFunctionCalls(10000); // for Minuit/Minuit2
632  min->SetMaxIterations(10000); // for GSL
633  min->SetTolerance(0.001);
634  min->SetPrintLevel(0);
635  double step = 1e-7;
636  min->SetVariable(0, "Q2", Q20, step);
637  min->SetVariableLimits(0, rQ2.min, rQ2.max);
638  min->Minimize();
639  diffv_max = -min->MinValue();
640 
641  if (tmp_diffv_max > diffv_max)
642  {
643  diffv_max = tmp_diffv_max;
644  }
645  return diffv_max;
646  }
647  default:
648  return -1.;
649  }
650 }
Range1D_t Q2QES_SM_lim(void) const
double vQES_SM_max(double Q2min, double Q2max) const
double fQ2Min
Q2-threshold for seeking the second maximum.
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:1077
int A(void) const
Definition: Target.h:70
A simple [min,max] interval for doubles.
Definition: Range1.h:42
Generated/set kinematical variables for an event.
Definition: Kinematics.h:39
double ComputeMaxXSec(const Interaction *in) const
Range1D_t vQES_SM_lim(double Q2) const
const double e
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
double vQES_SM_min(double Q2min, double Q2max) const
virtual double XSec(const Interaction *i, KinePhaseSpace_t k=kPSfE) const =0
Compute the cross section for the input interaction.
void SetKV(KineVar_t kv, double value)
Definition: Kinematics.cxx:335
double max
Definition: Range1.h:53
Target * TgtPtr(void) const
Definition: InitialState.h:67
double min
Definition: Range1.h:52
double GetFermiMomentum(void) const
Initial State information.
Definition: InitialState.h:48
void QELEventGeneratorSM::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 407 of file QELEventGeneratorSM.cxx.

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

408 {
409  Algorithm::Configure(config);
410  this->LoadConfig();
411 }
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
void QELEventGeneratorSM::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 413 of file QELEventGeneratorSM.cxx.

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

414 {
415  Algorithm::Configure(config);
416 
417  Registry r( "QELEventGeneratorSM_specific", false ) ;
418  r.Set( "sm_utils_algo", RgAlg("genie::SmithMonizUtils","Default") ) ;
419 
421 
422  this->LoadConfig();
423 }
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
A registry. Provides the container for algorithm configuration parameters.
Definition: Registry.h:65
void QELEventGeneratorSM::LoadConfig ( void  )
private

Definition at line 425 of file QELEventGeneratorSM.cxx.

References genie::KineGeneratorWithCache::fEMin, fGenerateNucleonInNucleus, genie::KineGeneratorWithCache::fGenerateUniformly, genie::KineGeneratorWithCache::fMaxXSecDiffTolerance, genie::KineGeneratorWithCache::fNumOfInterpolatorTypes, genie::KineGeneratorWithCache::fNumOfSafetyFactors, fQ2Min, genie::Algorithm::GetParamDef(), genie::Algorithm::GetParamVect(), sm_utils, genie::Algorithm::SubAlg(), genie::KineGeneratorWithCache::vInterpolatorTypes, and genie::KineGeneratorWithCache::vSafetyFactors.

Referenced by Configure().

426 {
427 
428  // Safety factors for the maximum differential cross section
429  fNumOfSafetyFactors = GetParamVect("Safety-Factors", vSafetyFactors, false);
430 
431  // Interpolator types for the maximum differential cross section
432  fNumOfInterpolatorTypes = GetParamVect("Interpolator-Types", vInterpolatorTypes, false);
433 
434 
435  // Minimum energy for which max xsec would be cached, forcing explicit
436  // calculation for lower eneries
437  GetParamDef( "Cache-MinEnergy", fEMin, 1.00) ;
438  GetParamDef( "Threshold-Q2", fQ2Min, 2.00);
439 
440  // Maximum allowed fractional cross section deviation from maxim cross
441  // section used in rejection method
442  GetParamDef( "MaxXSec-DiffTolerance", fMaxXSecDiffTolerance, 999999.);
443  assert(fMaxXSecDiffTolerance>=0);
444 
445  //-- Generate kinematics uniformly over allowed phase space and compute
446  // an event weight?
447  GetParamDef( "UniformOverPhaseSpace", fGenerateUniformly, false);
448 
449  // Generate nucleon in nucleus?
450  GetParamDef( "IsNucleonInNucleus", fGenerateNucleonInNucleus, true);
451 
452 
453  sm_utils = const_cast<genie::SmithMonizUtils *>(dynamic_cast<const genie::SmithMonizUtils *>( this -> SubAlg("sm_utils_algo") ) ) ;
454 }
double fQ2Min
Q2-threshold for seeking the second maximum.
bool fGenerateUniformly
uniform over allowed phase space + event weight?
int GetParamVect(const std::string &comm_name, std::vector< T > &v, bool is_top_call=true) const
Handle to load vectors of parameters.
double fMaxXSecDiffTolerance
max{100*(xsec-maxxsec)/.5*(xsec+maxxsec)} if xsec&gt;maxxsec
int fNumOfSafetyFactors
Number of given safety factors.
int fNumOfInterpolatorTypes
Number of given interpolators types.
std::vector< double > vSafetyFactors
MaxXSec -&gt; MaxXSec * fSafetyFactors[nkey].
bool fGenerateNucleonInNucleus
generate struck nucleon in nucleus
Contains auxiliary functions for Smith-Moniz model. .
double fEMin
min E for which maxxsec is cached - forcing explicit calc.
bool GetParamDef(const RgKey &name, T &p, const T &def) const
std::vector< string > vInterpolatorTypes
Type of interpolator for each key in a branch.
const Algorithm * SubAlg(const RgKey &registry_key) const
Definition: Algorithm.cxx:345
void QELEventGeneratorSM::ProcessEventRecord ( GHepRecord event_rec) const
virtual

Implements genie::EventRecordVisitorI.

Definition at line 82 of file QELEventGeneratorSM.cxx.

References genie::Target::A(), genie::GHepRecord::AddParticle(), AddTargetNucleusRemnant(), genie::KineGeneratorWithCache::AssertXSecLimits(), genie::Kinematics::ClearRunningValues(), genie::EventGeneratorI::CrossSectionAlg(), genie::GHepRecord::EventFlags(), fGenerateNucleonInNucleus, genie::KineGeneratorWithCache::fGenerateUniformly, genie::PDGLibrary::Find(), fkps, fQ2Min, genie::Interaction::FSPrimLeptonPdg(), genie::KineGeneratorWithCache::fXSecModel, genie::gAbortingInErr, genie::SmithMonizUtils::GetBindingEnergy(), genie::SmithMonizUtils::GetFermiMomentum(), genie::GHepParticle::GetP4(), genie::GHepRecord::HitNucleon(), genie::GHepRecord::HitNucleonPosition(), genie::Target::HitNucP4(), genie::Target::HitNucP4Ptr(), genie::Interaction::InitState(), genie::RunningThreadInfo::Instance(), genie::RandomGen::Instance(), genie::PDGLibrary::Instance(), genie::Target::IsNucleus(), genie::SmithMonizUtils::kFQES_SM_lim(), genie::Interaction::KinePtr(), genie::kISkipKinematicChk, genie::kISkipProcessChk, genie::kIStHadronInTheNucleus, genie::kIStStableFinalState, genie::kKineGenErr, genie::kKVPn, genie::kKVQ2, genie::kKVv, genie::kPSQ2fE, genie::kPSQ2vpfE, genie::kRfHitNucRest, genie::controls::kRjMaxIterations, LOG, genie::Range1D_t::max, genie::KineGeneratorWithCache::MaxXSec(), genie::Range1D_t::min, pFATAL, genie::SmithMonizUtils::PhaseSpaceVolume(), pINFO, pNOTICE, genie::GHepRecord::Probe(), genie::InitialState::ProbeE(), genie::GHepRecord::ProbePosition(), pWARN, genie::utils::kinematics::Q2(), genie::SmithMonizUtils::Q2QES_SM_lim(), genie::Interaction::RecoilNucleonPdg(), genie::RandomGen::RndGen(), genie::RandomGen::RndKine(), genie::RunningThreadInfo::RunningThread(), genie::GHepRecord::SetDiffXSec(), genie::Target::SetHitNucPosition(), genie::SmithMonizUtils::SetInteraction(), genie::Kinematics::SetKV(), genie::GHepParticle::SetMomentum(), genie::utils::SetPrimaryLeptonPolarization(), genie::Kinematics::SetQ2(), genie::exceptions::EVGThreadException::SetReason(), genie::GHepParticle::SetRemovalEnergy(), genie::Kinematics::SetW(), genie::GHepRecord::SetWeight(), genie::Kinematics::Setx(), genie::Kinematics::Sety(), sm_utils, genie::GHepRecord::Summary(), genie::exceptions::EVGThreadException::SwitchOnFastForward(), genie::InitialState::Tgt(), genie::InitialState::TgtPtr(), genie::SmithMonizUtils::vQES_SM_lim(), genie::utils::kinematics::W(), genie::GHepRecord::Weight(), genie::utils::kinematics::WQ2toXY(), genie::GHepParticle::X4(), genie::XSecAlgorithmI::XSec(), and genie::GHepRecord::XSec().

83 {
84  LOG("QELEvent", pINFO) << "Generating QE event kinematics...";
85 
86  if(fGenerateUniformly) {
87  LOG("QELEvent", pNOTICE)
88  << "Generating kinematics uniformly over the allowed phase space";
89  }
90  // Get the interaction and set the 'trust' bits
91  Interaction * interaction = evrec->Summary();
92  const InitialState & init_state = interaction -> InitState();
93  interaction->SetBit(kISkipProcessChk);
94  interaction->SetBit(kISkipKinematicChk);
95 
96  // Skip if no hit nucleon is set
97  if(! evrec->HitNucleon())
98  {
99  LOG("QELEvent", pFATAL) << "No hit nucleon was set";
100  gAbortingInErr = true;
101  exit(1);
102  }
103 
104  // Access the target from the interaction summary
105  Target * tgt = init_state.TgtPtr();
106  GHepParticle * nucleon = evrec->HitNucleon();
107  // Store position of nucleon
108  double hitNucPos = nucleon->X4()->Vect().Mag();
109  tgt->SetHitNucPosition( hitNucPos );
110 
111  // Get the random number generators
112  RandomGen * rnd = RandomGen::Instance();
113 
114  // Access cross section algorithm for running thread
116  const EventGeneratorI * evg = rtinfo->RunningThread();
117  fXSecModel = evg->CrossSectionAlg();
118 
119  // heavy nucleus is nucleus that heavier than tritium or 3He.
120  bool isHeavyNucleus = tgt->A()>3;
121 
122  sm_utils->SetInteraction(interaction);
123  // phase space for heavy nucleus is different from light one
124  fkps = isHeavyNucleus?kPSQ2vpfE:kPSQ2fE;
126  // Try to calculate the maximum cross-section in kinematical limits
127  // if not pre-computed already
128  double xsec_max1 = fGenerateUniformly ? -1 : this->MaxXSec(evrec);
129  // this make correct calculation of probability
130  double xsec_max2 = fGenerateUniformly ? -1 : (rQ2.max<fQ2Min)? 0:this->MaxXSec(evrec, 1);
131  double dvmax= isHeavyNucleus ? this->MaxXSec(evrec, 2) : 0.;
132 
133 
134  // generate Q2, v, pF
135  double Q2, v, kF, xsec;
136  unsigned int iter = 0;
137  bool accept = false;
138  TLorentzVector q;
139  while(1)
140  {
141  LOG("QELEvent", pINFO) << "Attempt #: " << iter;
142  if(iter > 100*kRjMaxIterations)
143  {
144  LOG("QELEvent", pWARN)
145  << "Couldn't select a valid kinematics after " << iter << " iterations";
146  evrec->EventFlags()->SetBitNumber(kKineGenErr, true);
148  exception.SetReason("Couldn't select kinematics");
149  exception.SwitchOnFastForward();
150  throw exception;
151  }
152 
153  // Pick Q2, v and pF
154  double xsec_max = 0.;
155  double pth = rnd->RndKine().Rndm();
156  //pth < prob1/(prob1+prob2), where prob1,prob2 - probabilities to generate event in area1 (Q2<fQ2Min) and area2 (Q2>fQ2Min) which are not normalized
157  if (pth <= xsec_max1*(TMath::Min(rQ2.max, fQ2Min)-rQ2.min)/(xsec_max1*(TMath::Min(rQ2.max, fQ2Min)-rQ2.min)+xsec_max2*(rQ2.max-fQ2Min)))
158  {
159  xsec_max = xsec_max1;
160  Q2 = (rnd->RndKine().Rndm() * (TMath::Min(rQ2.max, fQ2Min)-rQ2.min)) + rQ2.min;
161  }
162  else
163  {
164  xsec_max = xsec_max2;
165  Q2 = (rnd->RndKine().Rndm() * (rQ2.max-fQ2Min)) + fQ2Min;
166  }
167  Range1D_t rv = sm_utils->vQES_SM_lim(Q2);
168  // for nuclei heavier than deuterium generate energy transfer in defined energy interval
169  double dv = 0.;
170  if (isHeavyNucleus)
171  {
172  dv = dvmax * rnd->RndKine().Rndm();
173  if (dv > (rv.max-rv.min))
174  {
175  continue;
176  }
177  }
178  v = rv.min + dv;
179 
180  Range1D_t rkF = sm_utils->kFQES_SM_lim(Q2, v);
181  // rkF.max = Fermi momentum
182  kF = rnd->RndKine().Rndm()*sm_utils->GetFermiMomentum();
183  if (kF < rkF.min)
184  {
185  continue;
186  }
187 
188  Kinematics * kinematics = interaction->KinePtr();
189  kinematics->SetKV(kKVQ2, Q2);
190  kinematics->SetKV(kKVv, v);
191  kinematics->SetKV(kKVPn, kF);
192  xsec = fXSecModel->XSec(interaction, fkps);
193  //-- Decide whether to accept the current kinematics
194  if(!fGenerateUniformly)
195  {
196  this->AssertXSecLimits(interaction, xsec, xsec_max);
197  double t = xsec_max * rnd->RndKine().Rndm();
198 
199  accept = (t < xsec);
200  }
201  else
202  {
203  accept = (xsec>0);
204  }
205 
206  // If the generated kinematics are accepted, finish-up module's job
207  if(accept)
208  {
209  interaction->ResetBit(kISkipProcessChk);
210  interaction->ResetBit(kISkipKinematicChk);
211  break;
212  }
213  iter++;
214  }
215 
216  // Z-frame is frame where momentum transfer is (v,0,0,qv)
217  double qv = TMath::Sqrt(v*v+Q2);
218  TLorentzVector transferMom(0, 0, qv, v);
219 
220  // Momentum of initial neutrino in LAB frame
221  TLorentzVector * tempTLV = evrec->Probe()->GetP4();
222  TLorentzVector neutrinoMom = *tempTLV;
223  delete tempTLV;
224 
225  // define all angles in Z frame
226  double theta = neutrinoMom.Vect().Theta();
227  double phi = neutrinoMom.Vect().Phi();
228  double theta_k = sm_utils-> GetTheta_k(v, qv);
229  double costheta_k = TMath::Cos(theta_k);
230  double sintheta_k = TMath::Sin(theta_k);
231  double E_p; //energy of initial nucleon
232  double theta_p = sm_utils-> GetTheta_p(kF, v, qv, E_p);
233  double costheta_p = TMath::Cos(theta_p);
234  double sintheta_p = TMath::Sin(theta_p);
235  double fi_p = 2 * TMath::Pi() * rnd->RndGen().Rndm();
236  double cosfi_p = TMath::Cos(fi_p);
237  double sinfi_p = TMath::Sin(fi_p);
238  double psi = 2 * TMath::Pi() * rnd->RndGen().Rndm();
239 
240  // 4-momentum of initial neutrino in Z-frame
241  TLorentzVector neutrinoMomZ(neutrinoMom.P()*sintheta_k, 0, neutrinoMom.P()*costheta_k, neutrinoMom.E());
242 
243  // 4-momentum of final lepton in Z-frame
244  TLorentzVector outLeptonMom = neutrinoMomZ-transferMom;
245 
246  // 4-momentum of initial nucleon in Z-frame
247  TLorentzVector inNucleonMom(kF*sintheta_p*cosfi_p, kF*sintheta_p*sinfi_p, kF*costheta_p, E_p);
248 
249  // 4-momentum of final nucleon in Z-frame
250  TLorentzVector outNucleonMom = inNucleonMom+transferMom;
251 
252  // Rotate all vectors from Z frame to LAB frame
253  TVector3 yvec (0.0, 1.0, 0.0);
254  TVector3 zvec (0.0, 0.0, 1.0);
255  TVector3 unit_nudir = neutrinoMom.Vect().Unit();
256 
257  outLeptonMom.Rotate(theta-theta_k, yvec);
258  outLeptonMom.Rotate(phi, zvec);
259 
260  inNucleonMom.Rotate(theta-theta_k, yvec);
261  inNucleonMom.Rotate(phi, zvec);
262 
263  outNucleonMom.Rotate(theta-theta_k, yvec);
264  outNucleonMom.Rotate(phi, zvec);
265 
266  outLeptonMom.Rotate(psi, unit_nudir);
267  inNucleonMom.Rotate(psi, unit_nudir);
268  outNucleonMom.Rotate(psi, unit_nudir);
269 
270  // set 4-momentum of struck nucleon
271  TLorentzVector * p4 = tgt->HitNucP4Ptr();
272  p4->SetPx( inNucleonMom.Px() );
273  p4->SetPy( inNucleonMom.Py() );
274  p4->SetPz( inNucleonMom.Pz() );
275  p4->SetE ( inNucleonMom.E() );
276 
277  int rpdgc = interaction->RecoilNucleonPdg();
278  assert(rpdgc);
279  double W = PDGLibrary::Instance()->Find(rpdgc)->Mass();
280  LOG("QELEvent", pNOTICE) << "Selected: W = "<< W;
281  double M = init_state.Tgt().HitNucP4().M();
282  double E = init_state.ProbeE(kRfHitNucRest);
283 
284  // (W,Q2) -> (x,y)
285  double x=0, y=0;
286  kinematics::WQ2toXY(E,M,W,Q2,x,y);
287 
288  // lock selected kinematics & clear running values
289  interaction->KinePtr()->SetQ2(Q2, true);
290  interaction->KinePtr()->SetW (W, true);
291  interaction->KinePtr()->Setx (x, true);
292  interaction->KinePtr()->Sety (y, true);
293  interaction->KinePtr()->ClearRunningValues();
294 
295  // set the cross section for the selected kinematics
296  evrec->SetDiffXSec(xsec,fkps);
298  {
299  double vol = sm_utils->PhaseSpaceVolume(fkps);
300  double totxsec = evrec->XSec();
301  double wght = (vol/totxsec)*xsec;
302  LOG("QELEvent", pNOTICE) << "Kinematics wght = "<< wght;
303 
304  // apply computed weight to the current event weight
305  wght *= evrec->Weight();
306  LOG("QELEvent", pNOTICE) << "Current event wght = " << wght;
307  evrec->SetWeight(wght);
308  }
309  TLorentzVector x4l(*(evrec->Probe())->X4());
310 
311  // Add the final-state lepton to the event record
312  evrec->AddParticle(interaction->FSPrimLeptonPdg(), kIStStableFinalState, evrec->ProbePosition(),-1,-1,-1, outLeptonMom, x4l);
313 
314  // Set its polarization
316 
317  GHepStatus_t ist;
319  ist = kIStStableFinalState;
320  else
321  ist = (tgt->IsNucleus() && isHeavyNucleus) ? kIStHadronInTheNucleus : kIStStableFinalState;
322 
323  GHepParticle outNucleon(interaction->RecoilNucleonPdg(), ist, evrec->HitNucleonPosition(),-1,-1,-1, outNucleonMom , x4l);
324  evrec->AddParticle(outNucleon);
325 
326  // Store struck nucleon momentum
327  LOG("QELEvent",pNOTICE) << "pn: " << inNucleonMom.X() << ", " <<inNucleonMom.Y() << ", " <<inNucleonMom.Z() << ", " <<inNucleonMom.E();
328  nucleon->SetMomentum(inNucleonMom);
330 
331  // skip if not a nuclear target
332  if(evrec->Summary()->InitState().Tgt().IsNucleus())
333  // add a recoiled nucleus remnant
334  this->AddTargetNucleusRemnant(evrec);
335 
336  LOG("QELEvent", pINFO) << "Done generating QE event kinematics!";
337 }
void SetInteraction(const Interaction *i)
Range1D_t Q2QES_SM_lim(void) const
double fQ2Min
Q2-threshold for seeking the second maximum.
bool fGenerateUniformly
uniform over allowed phase space + event weight?
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:1077
static RandomGen * Instance()
Access instance.
Definition: RandomGen.cxx:71
void SetQ2(double Q2, bool selected=false)
Definition: Kinematics.cxx:255
Kinematics * KinePtr(void) const
Definition: Interaction.h:76
int RecoilNucleonPdg(void) const
recoil nucleon pdg
int A(void) const
Definition: Target.h:70
A simple [min,max] interval for doubles.
Definition: Range1.h:42
virtual void AssertXSecLimits(const Interaction *in, double xsec, double xsec_max) const
#define pFATAL
Definition: Messenger.h:56
bool IsNucleus(void) const
Definition: Target.cxx:272
Defines the EventGeneratorI interface.
virtual double MaxXSec(GHepRecord *evrec, const int nkey=0) const
Generated/set kinematical variables for an event.
Definition: Kinematics.h:39
void SetHitNucPosition(double r)
Definition: Target.cxx:210
void SetMomentum(const TLorentzVector &p4)
Range1D_t vQES_SM_lim(double Q2) const
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition: RandomGen.h:29
int FSPrimLeptonPdg(void) const
final state primary lepton pdg
double GetBindingEnergy(void) const
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
Range1D_t kFQES_SM_lim(double nu, double Q2) const
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
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
bool fGenerateNucleonInNucleus
generate struck nucleon in nucleus
#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
void SetKV(KineVar_t kv, double value)
Definition: Kinematics.cxx:335
void SetRemovalEnergy(double Erm)
static RunningThreadInfo * Instance(void)
double max
Definition: Range1.h:53
void SetW(double W, bool selected=false)
Definition: Kinematics.cxx:279
TLorentzVector * HitNucP4Ptr(void) const
Definition: Target.cxx:247
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:68
void AddTargetNucleusRemnant(GHepRecord *evrec) const
add a recoiled nucleus remnant
double PhaseSpaceVolume(KinePhaseSpace_t ps) const
const TLorentzVector * X4(void) const
Definition: GHepParticle.h:79
void Sety(double y, bool selected=false)
Definition: Kinematics.cxx:243
Target * TgtPtr(void) const
Definition: InitialState.h:67
TRandom3 & RndGen(void) const
rnd number generator for generic usage
Definition: RandomGen.h:80
virtual const XSecAlgorithmI * CrossSectionAlg(void) const =0
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)
bool gAbortingInErr
Definition: Messenger.cxx:34
void SetPrimaryLeptonPolarization(GHepRecord *ev)
double ProbeE(RefFrame_t rf) const
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
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
enum genie::EGHepStatus GHepStatus_t
double GetFermiMomentum(void) const
Initial State information.
Definition: InitialState.h:48

Member Data Documentation

bool genie::QELEventGeneratorSM::fGenerateNucleonInNucleus
private

generate struck nucleon in nucleus

Definition at line 76 of file QELEventGeneratorSM.h.

Referenced by LoadConfig(), and ProcessEventRecord().

KinePhaseSpace_t genie::QELEventGeneratorSM::fkps
mutableprivate

Definition at line 73 of file QELEventGeneratorSM.h.

Referenced by ComputeMaxXSec(), and ProcessEventRecord().

double genie::QELEventGeneratorSM::fQ2Min
private

Q2-threshold for seeking the second maximum.

Definition at line 77 of file QELEventGeneratorSM.h.

Referenced by ComputeMaxXSec(), LoadConfig(), and ProcessEventRecord().

SmithMonizUtils* genie::QELEventGeneratorSM::sm_utils
mutableprivate

Definition at line 64 of file QELEventGeneratorSM.h.

Referenced by ComputeMaxXSec(), LoadConfig(), and ProcessEventRecord().


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