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

Kinematical phase space. More...

#include <KPhaseSpace.h>

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

Public Member Functions

 KPhaseSpace (void)
 
 KPhaseSpace (const Interaction *in)
 
 ~KPhaseSpace (void)
 
void UseInteraction (const Interaction *in)
 
double Threshold (void) const
 Energy threshold. More...
 
double Threshold_SPP_iso (void) const
 Energy limit for resonance single pion production on isoscalar nucleon. More...
 
bool IsAboveThreshold (void) const
 Checks whether the interaction is above the energy threshold. More...
 
bool IsAllowed (void) const
 Check whether the current kinematics is in the allowed phase space. More...
 
Range1D_t Limits (KineVar_t kvar) const
 Return the kinematical variable limits. More...
 
double Minimum (KineVar_t kvar) const
 
double Maximum (KineVar_t kvar) const
 
Range1D_t WLim (void) const
 W limits. More...
 
Range1D_t Q2Lim_W (void) const
 Q2 limits @ fixed W. More...
 
Range1D_t q2Lim_W (void) const
 q2 limits @ fixed W More...
 
Range1D_t Q2Lim (void) const
 Q2 limits. More...
 
Range1D_t q2Lim (void) const
 q2 limits More...
 
Range1D_t XLim (void) const
 x limits More...
 
Range1D_t YLim (void) const
 y limits More...
 
Range1D_t YLim_X (void) const
 y limits @ fixed x More...
 
Range1D_t YLim (double xsi) const
 y limits (COH) More...
 
Range1D_t YLim_X (double xsi) const
 y limits @ fixed x (COH) More...
 
Range1D_t TLim (void) const
 t limits More...
 
Range1D_t WLim_SPP (void) const
 W limits for single pion production models. More...
 
Range1D_t WLim_SPP_iso (void) const
 W limits for resonance single pion production on isoscalar nucleon. More...
 
Range1D_t Q2Lim_W_SPP (void) const
 Q2 limits @ fixed W for single pion production models. More...
 
Range1D_t Q2Lim_W_SPP_iso (void) const
 Q2 limits @ fixed W for resonance single pion production on isoscalar nucleon. More...
 

Static Public Member Functions

static double GetTMaxDFR ()
 

Private Member Functions

void Init (void)
 

Private Attributes

const InteractionfInteraction
 

Detailed Description

Kinematical phase space.

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

Definition at line 33 of file KPhaseSpace.h.

Constructor & Destructor Documentation

KPhaseSpace::KPhaseSpace ( void  )

Definition at line 42 of file KPhaseSpace.cxx.

42  :
43 TObject(), fInteraction(NULL)
44 {
45  this->UseInteraction(0);
46 }
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
void UseInteraction(const Interaction *in)
Definition: KPhaseSpace.cxx:77
KPhaseSpace::KPhaseSpace ( const Interaction in)

Definition at line 48 of file KPhaseSpace.cxx.

References UseInteraction().

48  :
49 TObject(), fInteraction(NULL)
50 {
51  this->UseInteraction(in);
52 }
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
void UseInteraction(const Interaction *in)
Definition: KPhaseSpace.cxx:77
KPhaseSpace::~KPhaseSpace ( void  )

Definition at line 54 of file KPhaseSpace.cxx.

55 {
56 
57 }

Member Function Documentation

double KPhaseSpace::GetTMaxDFR ( )
static

Definition at line 59 of file KPhaseSpace.cxx.

References genie::AlgConfigPool::CommonList(), genie::Registry::GetDouble(), and genie::AlgConfigPool::Instance().

Referenced by genie::DFRKinematicsGenerator::ProcessEventRecord(), and TLim().

60 {
61  static bool tMaxLoaded = false;
62  static double DFR_tMax = -1;
63 
64  if (!tMaxLoaded)
65  {
67  const Registry * r = confp->CommonList( "Param", "Diffractive" ) ;
68  double tmax = r->GetDouble("DFR-t-max");
69  DFR_tMax = tmax;
70  tMaxLoaded = true;
71  }
72 
73  return DFR_tMax;
74 
75 }
A singleton class holding all configuration registries built while parsing all loaded XML configurati...
Definition: AlgConfigPool.h:40
RgDbl GetDouble(RgKey key) const
Definition: Registry.cxx:474
Registry * CommonList(const string &file_id, const string &set_name) const
A registry. Provides the container for algorithm configuration parameters.
Definition: Registry.h:65
static AlgConfigPool * Instance()
void genie::KPhaseSpace::Init ( void  )
private
bool KPhaseSpace::IsAboveThreshold ( void  ) const

Checks whether the interaction is above the energy threshold.

Definition at line 284 of file KPhaseSpace.cxx.

References fInteraction, genie::Interaction::InitState(), genie::ProcessInfo::IsAMNuGamma(), genie::ProcessInfo::IsCoherentElastic(), genie::ProcessInfo::IsCoherentProduction(), genie::ProcessInfo::IsDarkMatterDeepInelastic(), genie::ProcessInfo::IsDarkMatterElastic(), genie::ProcessInfo::IsDarkMatterElectronElastic(), genie::ProcessInfo::IsDeepInelastic(), genie::ProcessInfo::IsDiffractive(), genie::ProcessInfo::IsGlashowResonance(), genie::ProcessInfo::IsIMDAnnihilation(), genie::ProcessInfo::IsInverseBetaDecay(), genie::ProcessInfo::IsInverseMuDecay(), genie::ProcessInfo::IsMEC(), genie::ProcessInfo::IsNuElectronElastic(), genie::ProcessInfo::IsPhotonCoherent(), genie::ProcessInfo::IsPhotonResonance(), genie::ProcessInfo::IsQuasiElastic(), genie::ProcessInfo::IsResonant(), genie::ProcessInfo::IsSingleKaon(), genie::kRfHitNucRest, genie::kRfLab, LOG, pDEBUG, genie::InitialState::ProbeE(), genie::Interaction::ProcInfo(), and Threshold().

Referenced by genie::COHXSecAR::Integrate(), genie::IMDXSec::Integrate(), genie::COHXSec::Integrate(), genie::DISXSec::Integrate(), genie::RESXSec::Integrate(), genie::CEvNSXSec::Integrate(), genie::QELXSec::Integrate(), genie::COHDNuXSec::Integrate(), genie::HEDISXSec::Integrate(), genie::DMDISXSec::Integrate(), genie::DMELXSec::Integrate(), genie::HELeptonXSec::Integrate(), genie::MECXSec::Integrate(), genie::AlamSimoAtharVacasSKXSec::Integrate(), genie::NuElectronXSec::Integrate(), genie::DFRXSec::Integrate(), genie::DMElectronXSec::Integrate(), genie::ReinSehgalRESXSec::Integrate(), genie::SmithMonizQELCCXSec::Integrate(), genie::ReinSehgalSPPXSec::Integrate(), genie::ReinSehgalRESXSecFast::Integrate(), genie::FermiMover::KickHitNucleon(), genie::XSecAlgorithmI::ValidKinematics(), genie::StrumiaVissaniIBDPXSec::ValidKinematics(), genie::KLVOxygenIBDPXSec::ValidKinematics(), and genie::BertuzzoDNuCOHPXSec::ValidKinematics().

285 {
286  double E = 0.;
287  double Ethr = this->Threshold();
288 
289  const ProcessInfo & pi = fInteraction->ProcInfo();
290  const InitialState & init_state = fInteraction->InitState();
291 
292  if (pi.IsCoherentElastic() ||
293  pi.IsCoherentProduction() ||
294  pi.IsInverseMuDecay() ||
295  pi.IsIMDAnnihilation() ||
296  pi.IsNuElectronElastic() ||
298  pi.IsMEC() ||
299  pi.IsPhotonCoherent() ||
300  pi.IsPhotonResonance() ||
301  pi.IsGlashowResonance())
302  {
303  E = init_state.ProbeE(kRfLab);
304  }
305 
306  if(pi.IsQuasiElastic() ||
307  pi.IsDarkMatterElastic() ||
308  pi.IsInverseBetaDecay() ||
309  pi.IsResonant() ||
310  pi.IsDeepInelastic() ||
312  pi.IsDiffractive() ||
313  pi.IsSingleKaon() ||
314  pi.IsAMNuGamma())
315  {
316  E = init_state.ProbeE(kRfHitNucRest);
317  }
318 
319  LOG("KPhaseSpace", pDEBUG) << "E = " << E << ", Ethr = " << Ethr;
320  return (E>Ethr);
321 }
bool IsResonant(void) const
Definition: ProcessInfo.cxx:99
bool IsPhotonResonance(void) const
bool IsDarkMatterElectronElastic(void) const
double Threshold(void) const
Energy threshold.
Definition: KPhaseSpace.cxx:82
bool IsInverseMuDecay(void) const
bool IsQuasiElastic(void) const
Definition: ProcessInfo.cxx:69
bool IsInverseBetaDecay(void) const
bool IsDiffractive(void) const
bool IsCoherentProduction(void) const
bool IsIMDAnnihilation(void) const
bool IsSingleKaon(void) const
Definition: ProcessInfo.cxx:84
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
bool IsCoherentElastic(void) const
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
bool IsNuElectronElastic(void) const
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
bool IsAMNuGamma(void) const
bool IsDarkMatterElastic(void) const
Definition: ProcessInfo.cxx:74
bool IsPhotonCoherent(void) const
bool IsMEC(void) const
bool IsDeepInelastic(void) const
Definition: ProcessInfo.cxx:89
const InitialState & InitState(void) const
Definition: Interaction.h:69
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
bool IsGlashowResonance(void) const
bool IsDarkMatterDeepInelastic(void) const
Definition: ProcessInfo.cxx:94
double ProbeE(RefFrame_t rf) const
Initial State information.
Definition: InitialState.h:48
#define pDEBUG
Definition: Messenger.h:63
bool KPhaseSpace::IsAllowed ( void  ) const

Check whether the current kinematics is in the allowed phase space.

Definition at line 323 of file KPhaseSpace.cxx.

References fInteraction, genie::ProcessInfo::IsCoherentElastic(), genie::ProcessInfo::IsCoherentProduction(), genie::ProcessInfo::IsDarkMatterDeepInelastic(), genie::ProcessInfo::IsDarkMatterElastic(), genie::ProcessInfo::IsDarkMatterElectronElastic(), genie::ProcessInfo::IsDeepInelastic(), genie::ProcessInfo::IsDiffractive(), genie::ProcessInfo::IsIMDAnnihilation(), genie::ProcessInfo::IsInverseBetaDecay(), genie::ProcessInfo::IsInverseMuDecay(), genie::ProcessInfo::IsMEC(), genie::ProcessInfo::IsNuElectronElastic(), genie::ProcessInfo::IsQuasiElastic(), genie::ProcessInfo::IsResonant(), genie::ProcessInfo::IsSingleKaon(), genie::utils::math::IsWithinLimits(), genie::Interaction::Kine(), LOG, genie::Range1D_t::max, genie::Range1D_t::min, pDEBUG, genie::Interaction::ProcInfo(), genie::Kinematics::Q2(), genie::utils::kinematics::Q2(), Q2Lim(), Q2Lim_W(), genie::Kinematics::t(), TLim(), genie::utils::kinematics::UpdateWQ2FromXY(), genie::Kinematics::W(), genie::utils::kinematics::W(), WLim(), genie::Kinematics::x(), XLim(), genie::Kinematics::y(), and YLim().

Referenced by genie::XSecAlgorithmI::ValidKinematics().

324 {
325  const ProcessInfo & pi = fInteraction->ProcInfo();
326  const Kinematics & kine = fInteraction->Kine();
327 
328  // ASK single kaon:
329  // XSec code returns zero when kinematics are not allowed
330  // Here just let kinematics always be allowed
331  if(pi.IsSingleKaon()) {
332  return true;
333  }
334 
335  // QEL:
336  // Check the running Q2 vs the Q2 limits
337  if(pi.IsQuasiElastic() || pi.IsInverseBetaDecay() || pi.IsDarkMatterElastic()) {
338  Range1D_t Q2l = this->Q2Lim();
339  double Q2 = kine.Q2();
340  bool in_phys = math::IsWithinLimits(Q2, Q2l);
341  bool allowed = in_phys;
342  return allowed;
343  }
344 
345  // RES
346  // Check the running W vs the W limits
347  // & the running Q2 vs Q2 limits for the given W
348  if(pi.IsResonant()) {
349  Range1D_t Wl = this->WLim();
350  Range1D_t Q2l = this->Q2Lim_W();
351  double W = kine.W();
352  double Q2 = kine.Q2();
353  bool in_phys = (math::IsWithinLimits(Q2, Q2l) && math::IsWithinLimits(W, Wl));
354  bool allowed = in_phys;
355  return allowed;
356  }
357 
358  // DIS
359  if(pi.IsDeepInelastic() || pi.IsDarkMatterDeepInelastic()) {
360  Range1D_t Wl = this->WLim();
361  Range1D_t Q2l = this->Q2Lim_W();
362  double W = kine.W();
363  double Q2 = kine.Q2();
364  bool in_phys = (math::IsWithinLimits(Q2, Q2l) && math::IsWithinLimits(W, Wl));
365  bool allowed = in_phys;
366  return allowed;
367  }
368 
369  //IMD
371  Range1D_t yl = this->YLim();
372  double y = kine.y();
373  bool in_phys = math::IsWithinLimits(y, yl);
374  bool allowed = in_phys;
375  return allowed;
376  }
377 
378  //COH
379  if (pi.IsCoherentProduction()) {
380  Range1D_t xl = this->XLim();
381  Range1D_t yl = this->YLim();
382  double x = kine.x();
383  double y = kine.y();
384  bool in_phys = (math::IsWithinLimits(x, xl) && math::IsWithinLimits(y, yl));
385  bool allowed = in_phys;
386  return allowed;
387  }
388 
389  // CEvNS
390  if (pi.IsCoherentElastic()) {
391  double Q2 = kine.Q2();
392  bool allowed (Q2 > 0);
393  return allowed;
394  }
395 
396  // DFR
397  if (pi.IsDiffractive()) {
398  // first two checks are the same as RES & DIS
399  Range1D_t Wl = this->WLim();
400  Range1D_t Q2l = this->Q2Lim_W();
401 
403  double W = kine.W();
404  double Q2 = kine.Q2();
405 
406  LOG("KPhaseSpace", pDEBUG) << " W = " << W << ", limits = [" << Wl.min << "," << Wl.max << "];";
407  LOG("KPhaseSpace", pDEBUG) << " Q2 = " << Q2 << ", limits = [" << Q2l.min << "," << Q2l.max << "];";
408  bool in_phys = math::IsWithinLimits(W, Wl);
409  in_phys = in_phys && math::IsWithinLimits(Q2, Q2l);
410 
411  // extra check: there's a t minimum.
412  // but only check if W, Q2 is reasonable
413  // (otherwise get NaNs in tmin)
414  if (in_phys)
415  {
416  double t = kine.t();
417  Range1D_t tl = this->TLim();
418  LOG("KPhaseSpace", pDEBUG) << " t = " << t << ", limits = [" << tl.min << "," << tl.max << "];";
419  in_phys = in_phys && math::IsWithinLimits(t, tl);
420  }
421  LOG("KPhaseSpace", pDEBUG) << " phase space point is " << ( in_phys ? "ALLOWED" : "NOT ALLOWED");
422 
423 
424  bool allowed = in_phys;
425  return allowed;
426  }
427 
428  // was MECTensor
429  if (pi.IsMEC()){
430  Range1D_t Q2l = this->Q2Lim();
431  double Q2 = kine.Q2();
432  bool in_phys = math::IsWithinLimits(Q2, Q2l);
433  bool allowed = in_phys;
434  return allowed;
435  }
436 
437  return false;
438 }
bool IsResonant(void) const
Definition: ProcessInfo.cxx:99
double W(bool selected=false) const
Definition: Kinematics.cxx:157
bool IsDarkMatterElectronElastic(void) const
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:1077
A simple [min,max] interval for doubles.
Definition: Range1.h:42
bool IsInverseMuDecay(void) const
bool IsQuasiElastic(void) const
Definition: ProcessInfo.cxx:69
bool IsWithinLimits(double x, Range1D_t range)
Definition: MathUtils.cxx:258
Generated/set kinematical variables for an event.
Definition: Kinematics.h:39
bool IsInverseBetaDecay(void) const
double x(bool selected=false) const
Definition: Kinematics.cxx:99
Range1D_t YLim(void) const
y limits
bool IsDiffractive(void) const
bool IsCoherentProduction(void) const
bool IsIMDAnnihilation(void) const
Range1D_t Q2Lim_W(void) const
Q2 limits @ fixed W.
Range1D_t Q2Lim(void) const
Q2 limits.
double y(bool selected=false) const
Definition: Kinematics.cxx:112
bool IsSingleKaon(void) const
Definition: ProcessInfo.cxx:84
double W(const Interaction *const i)
Definition: KineUtils.cxx:1101
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
bool IsCoherentElastic(void) const
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
bool IsNuElectronElastic(void) const
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
const Kinematics & Kine(void) const
Definition: Interaction.h:71
bool IsDarkMatterElastic(void) const
Definition: ProcessInfo.cxx:74
bool IsMEC(void) const
bool IsDeepInelastic(void) const
Definition: ProcessInfo.cxx:89
void UpdateWQ2FromXY(const Interaction *in)
Definition: KineUtils.cxx:1290
double max
Definition: Range1.h:53
Range1D_t XLim(void) const
x limits
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
double min
Definition: Range1.h:52
double t(bool selected=false) const
Definition: Kinematics.cxx:170
double Q2(bool selected=false) const
Definition: Kinematics.cxx:125
Range1D_t TLim(void) const
t limits
bool IsDarkMatterDeepInelastic(void) const
Definition: ProcessInfo.cxx:94
Range1D_t WLim(void) const
W limits.
#define pDEBUG
Definition: Messenger.h:63
Range1D_t KPhaseSpace::Limits ( KineVar_t  kvar) const

Return the kinematical variable limits.

Definition at line 250 of file KPhaseSpace.cxx.

References genie::KineVar::AsString(), fInteraction, genie::kKVQ2, genie::kKVq2, genie::kKVt, genie::kKVW, genie::kKVx, genie::kKVy, LOG, pERROR, Q2Lim(), q2Lim(), TLim(), WLim(), XLim(), and YLim().

Referenced by genie::ReinSehgalRESXSecWithCache::CacheResExcitationXSec(), genie::IBDKinematicsGenerator::ComputeMaxXSec(), genie::NuEKinematicsGenerator::ComputeMaxXSec(), genie::QELKinematicsGenerator::ComputeMaxXSec(), genie::DISKinematicsGenerator::ComputeMaxXSec(), genie::DMEKinematicsGenerator::ComputeMaxXSec(), genie::DMELKinematicsGenerator::ComputeMaxXSec(), genie::DMDISKinematicsGenerator::ComputeMaxXSec(), genie::IMDXSec::Integrate(), genie::COHXSecAR::Integrate(), genie::RESXSec::Integrate(), genie::COHXSec::Integrate(), genie::QELXSec::Integrate(), genie::DMELXSec::Integrate(), genie::DFRXSec::Integrate(), genie::NuElectronXSec::Integrate(), genie::DMElectronXSec::Integrate(), genie::ReinSehgalRESXSec::Integrate(), genie::SmithMonizQELCCXSec::Integrate(), Maximum(), Minimum(), genie::utils::kinematics::PhaseSpaceVolume(), PrintLimits(), genie::DFRKinematicsGenerator::ProcessEventRecord(), genie::NuEKinematicsGenerator::ProcessEventRecord(), genie::DMEKinematicsGenerator::ProcessEventRecord(), genie::QELKinematicsGenerator::ProcessEventRecord(), genie::RESKinematicsGenerator::ProcessEventRecord(), genie::IBDKinematicsGenerator::ProcessEventRecord(), genie::DMELKinematicsGenerator::ProcessEventRecord(), genie::DISKinematicsGenerator::ProcessEventRecord(), genie::DMDISKinematicsGenerator::ProcessEventRecord(), genie::QELKinematicsGenerator::SpectralFuncExperimentalCode(), and genie::DMELKinematicsGenerator::SpectralFuncExperimentalCode().

251 {
252  // Compute limits for the input kinematic variable irrespective of any other
253  // relevant kinematical variable
254  //
255  assert(fInteraction);
256 
257  switch(kvar) {
258  case(kKVW) : return this->WLim(); break;
259  case(kKVQ2) : return this->Q2Lim(); break;
260  case(kKVq2) : return this->q2Lim(); break;
261  case(kKVx) : return this->XLim(); break;
262  case(kKVy) : return this->YLim(); break;
263  case(kKVt) : return this->TLim(); break;
264  default:
265  LOG("KPhaseSpace", pERROR)
266  << "Couldn't compute limits for " << KineVar::AsString(kvar);
267  Range1D_t R(-1.,-1);
268  return R;
269  }
270 }
#define pERROR
Definition: Messenger.h:59
A simple [min,max] interval for doubles.
Definition: Range1.h:42
Range1D_t q2Lim(void) const
q2 limits
Range1D_t YLim(void) const
y limits
Range1D_t Q2Lim(void) const
Q2 limits.
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
static string AsString(KineVar_t kv)
Definition: KineVar.h:77
Range1D_t XLim(void) const
x limits
Range1D_t TLim(void) const
t limits
Range1D_t WLim(void) const
W limits.
double KPhaseSpace::Maximum ( KineVar_t  kvar) const

Definition at line 278 of file KPhaseSpace.cxx.

References Limits(), and genie::Range1D_t::max.

279 {
280  Range1D_t lim = this->Limits(kvar);
281  return lim.max;
282 }
A simple [min,max] interval for doubles.
Definition: Range1.h:42
Range1D_t Limits(KineVar_t kvar) const
Return the kinematical variable limits.
double max
Definition: Range1.h:53
double KPhaseSpace::Minimum ( KineVar_t  kvar) const

Definition at line 272 of file KPhaseSpace.cxx.

References Limits(), and genie::Range1D_t::min.

273 {
274  Range1D_t lim = this->Limits(kvar);
275  return lim.min;
276 }
A simple [min,max] interval for doubles.
Definition: Range1.h:42
Range1D_t Limits(KineVar_t kvar) const
Return the kinematical variable limits.
double min
Definition: Range1.h:52
Range1D_t KPhaseSpace::Q2Lim ( void  ) const

Q2 limits.

Definition at line 571 of file KPhaseSpace.cxx.

References genie::utils::kinematics::CEvNSQ2Lim(), genie::XclsTag::CharmHadronPdg(), genie::utils::kinematics::CohQ2Lim(), genie::utils::kinematics::DarkQ2Lim(), genie::utils::kinematics::DarkQ2Lim_W(), genie::Interaction::ExclTag(), genie::PDGLibrary::Find(), fInteraction, genie::Interaction::FSPrimLepton(), genie::Target::HitNucP4Ptr(), genie::utils::kinematics::InelQ2Lim(), genie::utils::kinematics::InelQ2Lim_W(), genie::Interaction::InitState(), genie::PDGLibrary::Instance(), genie::XclsTag::IsCharmEvent(), genie::ProcessInfo::IsCoherentElastic(), genie::ProcessInfo::IsCoherentProduction(), genie::ProcessInfo::IsDarkMatterDeepInelastic(), genie::ProcessInfo::IsDarkMatterElastic(), genie::ProcessInfo::IsDeepInelastic(), genie::ProcessInfo::IsEM(), genie::ProcessInfo::IsInverseBetaDecay(), genie::ProcessInfo::IsMEC(), genie::ProcessInfo::IsQuasiElastic(), genie::ProcessInfo::IsResonant(), genie::XclsTag::IsStrangeEvent(), genie::ProcessInfo::IsWeakCC(), genie::controls::kASmallNum, genie::controls::kMinQ2Limit_VLE, genie::constants::kPi0Mass, genie::constants::kPionMass, genie::kRfHitNucRest, genie::kRfLab, genie::Range1D_t::max, genie::Range1D_t::min, genie::XclsTag::NPions(), genie::InitialState::ProbeE(), genie::Interaction::ProcInfo(), genie::Interaction::RecoilNucleon(), genie::XclsTag::StrangeHadronPdg(), genie::InitialState::Tgt(), and genie::utils::kinematics::W().

Referenced by genie::DISXSec::CacheFreeNucleonXSec(), genie::DMDISXSec::CacheFreeNucleonXSec(), genie::utils::ComputeFullDMELPXSec(), genie::utils::ComputeFullQELPXSec(), genie::CEvNSEventGenerator::GenerateKinematics(), genie::CEvNSXSec::Integrate(), genie::DISXSec::Integrate(), genie::HEDISXSec::Integrate(), genie::DMDISXSec::Integrate(), IsAllowed(), Limits(), genie::COHKinematicsGenerator::MaxXSec_BergerSehgal(), genie::COHKinematicsGenerator::MaxXSec_BergerSehgalFM(), genie::HEDISKinematicsGenerator::ProcessEventRecord(), q2Lim(), and Q2Lim_W().

572 {
573  // Computes momentum transfer (Q2>0) limits irrespective of the invariant mass
574  // For QEL this is identical to Q2Lim_W (since W is fixed)
575  // For RES & DIS, the calculation proceeds as in kinematics::InelQ2Lim().
576  //
577  Range1D_t Q2l;
578  Q2l.min = -1;
579  Q2l.max = -1;
580 
581  const ProcessInfo & pi = fInteraction->ProcInfo();
582 
583  bool is_em = pi.IsEM();
584  bool is_qel = pi.IsQuasiElastic() || pi.IsInverseBetaDecay();
585  bool is_inel = pi.IsDeepInelastic() || pi.IsResonant();
586  bool is_coh = pi.IsCoherentProduction();
587  bool is_cevns = pi.IsCoherentElastic();
588  bool is_dme = pi.IsDarkMatterElastic();
589  bool is_dmdis = pi.IsDarkMatterDeepInelastic();
590 
591  if(!is_qel && !is_inel && !is_coh && !is_cevns && !is_dme && !is_dmdis) return Q2l;
592 
593  const InitialState & init_state = fInteraction->InitState();
594  double Ev = init_state.ProbeE(kRfHitNucRest);
595  double M = init_state.Tgt().HitNucP4Ptr()->M(); // can be off m/shell
596  double ml = fInteraction->FSPrimLepton()->Mass();
597 
598  if(is_cevns) {
599  double Ev_lab = init_state.ProbeE(kRfLab);
600  Q2l = kinematics::CEvNSQ2Lim(Ev_lab);
601  return Q2l;
602  }
603 
604  const XclsTag & xcls = fInteraction->ExclTag();
605 
606  if(is_coh) {
607 
608  double m_other = controls::kASmallNum ;
609  // as a default the mass of hadronic system is the mass of the photon.
610  // which is assumed to be a small number to avoid divergences
611 
612  if ( xcls.NPions() > 0 ) {
613  bool pionIsCharged = pi.IsWeakCC();
614  m_other = pionIsCharged ? kPionMass : kPi0Mass;
615  }
616 
617  Q2l = kinematics::CohQ2Lim(M, m_other, ml, Ev);
618  return Q2l;
619  }
620 
621  // quasi-elastic
622  if(is_qel) {
623  double W = fInteraction->RecoilNucleon()->Mass();
624  if(xcls.IsCharmEvent()) {
625  int charm_pdgc = xcls.CharmHadronPdg();
626  W = PDGLibrary::Instance()->Find(charm_pdgc)->Mass();
627  } else if(xcls.IsStrangeEvent()) {
628  int strange_pdgc = xcls.StrangeHadronPdg();
629  W = PDGLibrary::Instance()->Find(strange_pdgc)->Mass();
630  }
631  if (pi.IsInverseBetaDecay()) {
633  } else {
634  Q2l = is_em ? kinematics::electromagnetic::InelQ2Lim_W(Ev,ml,M,W) : kinematics::InelQ2Lim_W(Ev,M,ml,W);
635  }
636 
637  return Q2l;
638  }
639 
640  // dark mattter elastic
641  if(is_dme) {
642  double W = fInteraction->RecoilNucleon()->Mass();
643  if(xcls.IsCharmEvent()) {
644  int charm_pdgc = xcls.CharmHadronPdg();
645  W = PDGLibrary::Instance()->Find(charm_pdgc)->Mass();
646  } else if(xcls.IsStrangeEvent()) {
647  int strange_pdgc = xcls.StrangeHadronPdg();
648  W = PDGLibrary::Instance()->Find(strange_pdgc)->Mass();
649  }
650  if (pi.IsInverseBetaDecay()) {
652  } else {
653  Q2l = kinematics::DarkQ2Lim_W(Ev,M,ml,W);
654  }
655 
656 
657  return Q2l;
658  }
659 
660  // was MECTensor
661  // TODO: Q2maxConfig
662  if (pi.IsMEC()){
663  double W = fInteraction->RecoilNucleon()->Mass();
664  Q2l = is_em ? kinematics::electromagnetic::InelQ2Lim_W(Ev,ml,M,W) : kinematics::InelQ2Lim_W(Ev,M,ml,W);
665  double Q2maxConfig = 1.44; // need to pull from config file somehow?
666  if (Q2l.max > Q2maxConfig) Q2l.max = Q2maxConfig;
667  return Q2l;
668  }
669 
670  if (is_dmdis) {
671  Q2l = kinematics::DarkQ2Lim(Ev,M,ml);
672  return Q2l;
673  }
674 
675  // inelastic
676  Q2l = is_em ? kinematics::electromagnetic::InelQ2Lim(Ev,ml,M) : kinematics::InelQ2Lim(Ev,M,ml);
677  return Q2l;
678 }
bool IsResonant(void) const
Definition: ProcessInfo.cxx:99
bool IsWeakCC(void) const
static const double kMinQ2Limit_VLE
Definition: Controls.h:42
TParticlePDG * RecoilNucleon(void) const
recoil nucleon
A simple [min,max] interval for doubles.
Definition: Range1.h:42
int NPions(void) const
Definition: XclsTag.h:62
bool IsQuasiElastic(void) const
Definition: ProcessInfo.cxx:69
int CharmHadronPdg(void) const
Definition: XclsTag.h:52
bool IsStrangeEvent(void) const
Definition: XclsTag.h:53
Range1D_t DarkQ2Lim_W(double Ev, double M, double ml, double W, double Q2min_cut=controls::kMinQ2Limit)
Definition: KineUtils.cxx:915
bool IsInverseBetaDecay(void) const
Range1D_t CEvNSQ2Lim(double Ev)
Definition: KineUtils.cxx:886
bool IsCoherentProduction(void) const
Range1D_t InelQ2Lim(double Ev, double M, double ml, double Q2min_cut=controls::kMinQ2Limit)
Definition: KineUtils.cxx:429
Contains minimal information for tagging exclusive processes.
Definition: XclsTag.h:39
Range1D_t InelQ2Lim_W(double Ev, double M, double ml, double W, double Q2min_cut=controls::kMinQ2Limit)
Definition: KineUtils.cxx:379
bool IsCharmEvent(void) const
Definition: XclsTag.h:50
double W(const Interaction *const i)
Definition: KineUtils.cxx:1101
Range1D_t DarkQ2Lim(double Ev, double M, double ml, double Q2min_cut=controls::kMinQ2Limit)
Definition: KineUtils.cxx:973
int StrangeHadronPdg(void) const
Definition: XclsTag.h:55
bool IsCoherentElastic(void) const
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
static const double kASmallNum
Definition: Controls.h:40
bool IsDarkMatterElastic(void) const
Definition: ProcessInfo.cxx:74
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
bool IsMEC(void) const
bool IsEM(void) const
bool IsDeepInelastic(void) const
Definition: ProcessInfo.cxx:89
double max
Definition: Range1.h:53
TLorentzVector * HitNucP4Ptr(void) const
Definition: Target.cxx:247
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:68
const XclsTag & ExclTag(void) const
Definition: Interaction.h:72
Range1D_t CohQ2Lim(double Mn, double m_produced, double mlep, double Ev)
Definition: KineUtils.cxx:743
const InitialState & InitState(void) const
Definition: Interaction.h:69
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
double min
Definition: Range1.h:52
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:86
const Target & Tgt(void) const
Definition: InitialState.h:66
bool IsDarkMatterDeepInelastic(void) const
Definition: ProcessInfo.cxx:94
double ProbeE(RefFrame_t rf) const
Initial State information.
Definition: InitialState.h:48
Range1D_t KPhaseSpace::q2Lim ( void  ) const

q2 limits

Definition at line 680 of file KPhaseSpace.cxx.

References genie::Range1D_t::max, genie::Range1D_t::min, genie::utils::kinematics::Q2(), and Q2Lim().

Referenced by Limits().

681 {
682 // As Q2Lim(void) but with reversed sign (Q2 -> q2)
683 //
684  Range1D_t Q2 = this->Q2Lim();
685  Range1D_t q2;
686  q2.min = - Q2.max;
687  q2.max = - Q2.min;
688  return q2;
689 }
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:1077
A simple [min,max] interval for doubles.
Definition: Range1.h:42
Range1D_t Q2Lim(void) const
Q2 limits.
double max
Definition: Range1.h:53
double min
Definition: Range1.h:52
Range1D_t KPhaseSpace::Q2Lim_W ( void  ) const

Q2 limits @ fixed W.

Definition at line 510 of file KPhaseSpace.cxx.

References genie::utils::kinematics::DarkQ2Lim_W(), fInteraction, genie::Interaction::FSPrimLepton(), genie::Target::HitNucP4Ptr(), genie::utils::kinematics::InelQ2Lim_W(), genie::Interaction::InitState(), genie::ProcessInfo::IsCoherentProduction(), genie::ProcessInfo::IsDarkMatterDeepInelastic(), genie::ProcessInfo::IsDarkMatterElastic(), genie::ProcessInfo::IsDeepInelastic(), genie::ProcessInfo::IsDiffractive(), genie::ProcessInfo::IsEM(), genie::ProcessInfo::IsInverseBetaDecay(), genie::ProcessInfo::IsQuasiElastic(), genie::ProcessInfo::IsResonant(), genie::controls::kMinQ2Limit_VLE, genie::kRfHitNucRest, genie::Range1D_t::max, genie::Range1D_t::min, genie::InitialState::ProbeE(), genie::Interaction::ProcInfo(), Q2Lim(), genie::Interaction::RecoilNucleon(), genie::InitialState::Tgt(), and genie::utils::kinematics::W().

Referenced by genie::DFRKinematicsGenerator::ComputeMaxXSec(), genie::RESKinematicsGenerator::ComputeMaxXSec(), IsAllowed(), genie::utils::kinematics::PhaseSpaceVolume(), genie::RESKinematicsGenerator::ProcessEventRecord(), and q2Lim_W().

511 {
512  // Computes momentum transfer (Q2>0) limits @ the input invariant mass
513  // The calculation proceeds as in kinematics::InelQ2Lim_W().
514  // For QEL, W is set to the recoil nucleon mass
515  //
516  // TODO: For now, choosing to handle Q2 at fixed W for coherent in the
517  // same way as for the general Q2 limits... but shouldn't we just use
518  // W = m_pi? - which we do in Q2Lim() anyway... seems like there are
519  // cleanup opportunities here.
520 
521  Range1D_t Q2l;
522  Q2l.min = -1;
523  Q2l.max = -1;
524 
525  const ProcessInfo & pi = fInteraction->ProcInfo();
526 
527  bool is_em = pi.IsEM();
528  bool is_qel = pi.IsQuasiElastic() || pi.IsInverseBetaDecay();
529  bool is_inel = pi.IsDeepInelastic() || pi.IsResonant() || pi.IsDiffractive();
530  bool is_coh = pi.IsCoherentProduction();
531  bool is_dme = pi.IsDarkMatterElastic();
532  bool is_dmdis = pi.IsDarkMatterDeepInelastic();
533 
534  if(!is_qel && !is_inel && !is_coh && !is_dme && !is_dmdis) return Q2l;
535 
536  if(is_coh) {
537  return Q2Lim();
538  }
539 
540  const InitialState & init_state = fInteraction->InitState();
541  double Ev = init_state.ProbeE(kRfHitNucRest);
542  double M = init_state.Tgt().HitNucP4Ptr()->M(); // can be off m/shell
543  double ml = fInteraction->FSPrimLepton()->Mass();
544 
545  double W = 0;
546  if(is_qel || is_dme) W = fInteraction->RecoilNucleon()->Mass();
547  else W = kinematics::W(fInteraction);
548 
549  if (pi.IsInverseBetaDecay()) {
551  } else if (is_dme || is_dmdis) {
552  Q2l = kinematics::DarkQ2Lim_W(Ev,M,ml,W);
553  } else {
554  Q2l = is_em ? kinematics::electromagnetic::InelQ2Lim_W(Ev,ml,M,W) : kinematics::InelQ2Lim_W(Ev,M,ml,W);
555  }
556 
557  return Q2l;
558 }
bool IsResonant(void) const
Definition: ProcessInfo.cxx:99
static const double kMinQ2Limit_VLE
Definition: Controls.h:42
TParticlePDG * RecoilNucleon(void) const
recoil nucleon
A simple [min,max] interval for doubles.
Definition: Range1.h:42
bool IsQuasiElastic(void) const
Definition: ProcessInfo.cxx:69
Range1D_t DarkQ2Lim_W(double Ev, double M, double ml, double W, double Q2min_cut=controls::kMinQ2Limit)
Definition: KineUtils.cxx:915
bool IsInverseBetaDecay(void) const
bool IsDiffractive(void) const
bool IsCoherentProduction(void) const
Range1D_t Q2Lim(void) const
Q2 limits.
Range1D_t InelQ2Lim_W(double Ev, double M, double ml, double W, double Q2min_cut=controls::kMinQ2Limit)
Definition: KineUtils.cxx:379
double W(const Interaction *const i)
Definition: KineUtils.cxx:1101
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
bool IsDarkMatterElastic(void) const
Definition: ProcessInfo.cxx:74
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
bool IsEM(void) const
bool IsDeepInelastic(void) const
Definition: ProcessInfo.cxx:89
double max
Definition: Range1.h:53
TLorentzVector * HitNucP4Ptr(void) const
Definition: Target.cxx:247
const InitialState & InitState(void) const
Definition: Interaction.h:69
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
double min
Definition: Range1.h:52
const Target & Tgt(void) const
Definition: InitialState.h:66
bool IsDarkMatterDeepInelastic(void) const
Definition: ProcessInfo.cxx:94
double ProbeE(RefFrame_t rf) const
Initial State information.
Definition: InitialState.h:48
Range1D_t KPhaseSpace::q2Lim_W ( void  ) const

q2 limits @ fixed W

Definition at line 560 of file KPhaseSpace.cxx.

References genie::Range1D_t::max, genie::Range1D_t::min, genie::utils::kinematics::Q2(), and Q2Lim_W().

561 {
562 // As Q2Lim_W(void) but with reversed sign (Q2 -> q2)
563 //
564  Range1D_t Q2 = this->Q2Lim_W();
565  Range1D_t q2;
566  q2.min = - Q2.max;
567  q2.max = - Q2.min;
568  return q2;
569 }
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:1077
A simple [min,max] interval for doubles.
Definition: Range1.h:42
Range1D_t Q2Lim_W(void) const
Q2 limits @ fixed W.
double max
Definition: Range1.h:53
double min
Definition: Range1.h:52
Range1D_t KPhaseSpace::Q2Lim_W_SPP ( void  ) const

Q2 limits @ fixed W for single pion production models.

Definition at line 1072 of file KPhaseSpace.cxx.

References genie::InitialState::CMEnergy(), epsilon, genie::PDGLibrary::Find(), fInteraction, genie::SppChannel::FromInteraction(), genie::Interaction::FSPrimLepton(), genie::Interaction::InitState(), genie::SppChannel::InitStateNucleon(), genie::PDGLibrary::Instance(), genie::utils::res::Mass(), genie::Range1D_t::max, genie::Range1D_t::min, genie::InitialState::ProbePdg(), genie::units::s, and genie::utils::kinematics::W().

1073 {
1074  Range1D_t Q2l;
1075  const InitialState & init_state = fInteraction->InitState();
1077  PDGLibrary * pdglib = PDGLibrary::Instance();
1078  double Mi = pdglib->Find( SppChannel::InitStateNucleon(spp_channel) )->Mass();
1079  double mi = pdglib->Find( init_state.ProbePdg() )->Mass();
1080  double mf = fInteraction->FSPrimLepton()->Mass();
1081  double mi2 = mi*mi;
1082  double mf2 = mf*mf;
1083  double W = kinematics::W(fInteraction);
1084 
1085  double ECM = init_state.CMEnergy();
1086  double s = ECM*ECM;
1087 
1088  double Ei_CM = (s + mi2 - Mi*Mi)/2/ECM;
1089  double Ef_CM = (s + mf2 - W*W)/2/ECM;
1090  double Pi_CM = (Ei_CM - mi)<0?0:TMath::Sqrt(Ei_CM*Ei_CM - mi2);
1091  double Pf_CM = (Ef_CM - mf)<0?0:TMath::Sqrt(Ef_CM*Ef_CM - mf2);
1092  // kinematic Q2-limits
1093  Q2l.min = 2*(Ei_CM*Ef_CM - Pi_CM*Pf_CM) - mi2 - mf2;
1094  Q2l.max = 2*(Ei_CM*Ef_CM + Pi_CM*Pf_CM) - mi2 - mf2;
1095 
1096  if ( (Q2l.max - Q2l.min) < (Q2l.max + Q2l.min)*std::numeric_limits<double>::epsilon() )
1097  {
1098  Q2l.min = 2*Q2l.max*Q2l.min/(Q2l.max + Q2l.min);
1099  Q2l.max = Q2l.min;
1100  }
1101  else
1102  {
1105  }
1106 
1107  return Q2l;
1108 }
static SppChannel_t FromInteraction(const Interaction *interaction)
Definition: SppChannel.h:402
A simple [min,max] interval for doubles.
Definition: Range1.h:42
static constexpr double s
Definition: Units.h:95
double Mass(Resonance_t res)
resonance mass (GeV)
const double epsilon
enum genie::ESppChannel SppChannel_t
double W(const Interaction *const i)
Definition: KineUtils.cxx:1101
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
static int InitStateNucleon(SppChannel_t channel)
Definition: SppChannel.h:103
int ProbePdg(void) const
Definition: InitialState.h:64
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
double max
Definition: Range1.h:53
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:68
Singleton class to load &amp; serve a TDatabasePDG.
Definition: PDGLibrary.h:35
double CMEnergy() const
centre-of-mass energy (sqrt s)
const InitialState & InitState(void) const
Definition: Interaction.h:69
double min
Definition: Range1.h:52
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:86
Initial State information.
Definition: InitialState.h:48
Range1D_t KPhaseSpace::Q2Lim_W_SPP_iso ( void  ) const

Q2 limits @ fixed W for resonance single pion production on isoscalar nucleon.

Definition at line 1110 of file KPhaseSpace.cxx.

References epsilon, genie::PDGLibrary::Find(), fInteraction, genie::Interaction::FSPrimLepton(), genie::Interaction::InitState(), genie::PDGLibrary::Instance(), genie::kPdgNeutron, genie::kPdgProton, genie::kRfHitNucRest, genie::utils::res::Mass(), genie::Range1D_t::max, genie::Range1D_t::min, genie::InitialState::ProbeE(), genie::InitialState::ProbePdg(), genie::units::s, and genie::utils::kinematics::W().

Referenced by genie::SPPEventGenerator::ProcessEventRecord(), and genie::MKSPPPXSec2020::ValidKinematics().

1111 {
1112  Range1D_t Q2l;
1113  const InitialState & init_state = fInteraction->InitState();
1114  PDGLibrary * pdglib = PDGLibrary::Instance();
1115  // imply isospin symmetry
1116  double M = (pdglib->Find(kPdgProton)->Mass() + pdglib->Find(kPdgNeutron)->Mass())/2;
1117  double mi = pdglib->Find( init_state.ProbePdg() )->Mass();
1118  double mf = fInteraction->FSPrimLepton()->Mass();
1119  double mi2 = mi*mi;
1120  double mf2 = mf*mf;
1121  double W = kinematics::W(fInteraction);
1122 
1123  double Ei = init_state.ProbeE(kRfHitNucRest);
1124  double s = M*(M + 2*Ei) + mi2;
1125  double ECM = TMath::Sqrt(s);
1126 
1127  double Ei_CM = (s + mi2 - M*M)/2/ECM;
1128  double Ef_CM = (s + mf2 - W*W)/2/ECM;
1129  double Pi_CM = (Ei_CM - mi)<0?0:TMath::Sqrt(Ei_CM*Ei_CM - mi2);
1130  double Pf_CM = (Ef_CM - mf)<0?0:TMath::Sqrt(Ef_CM*Ef_CM - mf2);
1131  // kinematic Q2-limits
1132  Q2l.min = 2*(Ei_CM*Ef_CM - Pi_CM*Pf_CM) - mi2 - mf2;
1133  Q2l.max = 2*(Ei_CM*Ef_CM + Pi_CM*Pf_CM) - mi2 - mf2;
1134 
1135  if ( (Q2l.max - Q2l.min) < (Q2l.max + Q2l.min)*std::numeric_limits<double>::epsilon() )
1136  {
1137  Q2l.min = 2*Q2l.max*Q2l.min/(Q2l.max + Q2l.min);
1138  Q2l.max = Q2l.min;
1139  }
1140  else
1141  {
1144  }
1145 
1146  return Q2l;
1147 }
A simple [min,max] interval for doubles.
Definition: Range1.h:42
static constexpr double s
Definition: Units.h:95
double Mass(Resonance_t res)
resonance mass (GeV)
const double epsilon
double W(const Interaction *const i)
Definition: KineUtils.cxx:1101
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
int ProbePdg(void) const
Definition: InitialState.h:64
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
double max
Definition: Range1.h:53
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:68
Singleton class to load &amp; serve a TDatabasePDG.
Definition: PDGLibrary.h:35
const InitialState & InitState(void) const
Definition: Interaction.h:69
double min
Definition: Range1.h:52
const int kPdgProton
Definition: PDGCodes.h:81
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:86
double ProbeE(RefFrame_t rf) const
const int kPdgNeutron
Definition: PDGCodes.h:83
Initial State information.
Definition: InitialState.h:48
double KPhaseSpace::Threshold ( void  ) const

Energy threshold.

Definition at line 82 of file KPhaseSpace.cxx.

References genie::XclsTag::CharmHadronPdg(), genie::Interaction::ExclTag(), genie::PDGLibrary::Find(), fInteraction, genie::Interaction::FSPrimLepton(), genie::Target::HitNucIsSet(), genie::Target::HitNucP4Ptr(), genie::Target::HitNucPdg(), genie::Interaction::InitState(), genie::PDGLibrary::Instance(), genie::ProcessInfo::IsAMNuGamma(), genie::XclsTag::IsCharmEvent(), genie::ProcessInfo::IsCoherentElastic(), genie::ProcessInfo::IsCoherentProduction(), genie::ProcessInfo::IsDarkMatterDeepInelastic(), genie::ProcessInfo::IsDarkMatterElastic(), genie::ProcessInfo::IsDarkMatterElectronElastic(), genie::ProcessInfo::IsDeepInelastic(), genie::ProcessInfo::IsDiffractive(), genie::ProcessInfo::IsGlashowResonance(), genie::ProcessInfo::IsIMDAnnihilation(), genie::XclsTag::IsInclusiveCharm(), genie::ProcessInfo::IsInverseBetaDecay(), genie::ProcessInfo::IsInverseMuDecay(), genie::ProcessInfo::IsKnown(), genie::ProcessInfo::IsMEC(), genie::ProcessInfo::IsNorm(), genie::pdg::IsNuE(), genie::ProcessInfo::IsNuElectronElastic(), genie::pdg::IsNuMu(), genie::pdg::IsNuTau(), genie::ProcessInfo::IsPhotonCoherent(), genie::ProcessInfo::IsPhotonResonance(), genie::ProcessInfo::IsQuasiElastic(), genie::ProcessInfo::IsResonant(), genie::ProcessInfo::IsSingleKaon(), genie::ProcessInfo::IsSinglePion(), genie::ProcessInfo::IsWeakCC(), genie::controls::kASmallNum, genie::constants::kElectronMass, genie::constants::kElectronMass2, genie::constants::kLightestChmHad, genie::constants::kMuonMass, genie::constants::kMuonMass2, genie::constants::kMw, genie::constants::kNeutronMass, genie::constants::kNucleonMass, genie::kPdgPi0, genie::kPdgPiM, genie::kPdgPiP, genie::constants::kPhotontest, genie::constants::kPi0Mass, genie::constants::kPionMass, genie::constants::kProtonMass, genie::constants::kTauMass, genie::units::m, genie::units::m2, genie::utils::res::Mass(), genie::Target::Mass(), genie::Target::N(), genie::XclsTag::NPi0(), genie::XclsTag::NPiMinus(), genie::XclsTag::NPions(), genie::XclsTag::NPiPlus(), genie::XclsTag::NProtons(), genie::Target::Pdg(), pERROR, genie::InitialState::ProbePdg(), genie::Interaction::ProcInfo(), genie::Interaction::RecoilNucleon(), SLOG, genie::XclsTag::StrangeHadronPdg(), genie::pdg::SwitchProtonNeutron(), genie::InitialState::Tgt(), and genie::Target::Z().

Referenced by genie::DISXSec::CacheFreeNucleonXSec(), genie::DMDISXSec::CacheFreeNucleonXSec(), genie::ReinSehgalRESXSecWithCache::CacheResExcitationXSec(), genie::ReinSehgalRESXSecWithCacheFast::CacheResExcitationXSec(), genie::XSecSplineList::CreateSpline(), IsAboveThreshold(), and genie::FermiMover::KickHitNucleon().

83 {
84  const ProcessInfo & pi = fInteraction->ProcInfo();
85  const InitialState & init_state = fInteraction->InitState();
86  const XclsTag & xcls = fInteraction->ExclTag();
87  const Target & tgt = init_state.Tgt();
88 
89  double ml = fInteraction->FSPrimLepton()->Mass();
90 
91  if( ! pi.IsKnown() ) return 0;
92 
93  if (pi.IsSinglePion()) {
94  double Mi = tgt.HitNucP4Ptr()->M(); // initial nucleon mass
95  double Mf = (xcls.NProtons()==1) ? kProtonMass : kNeutronMass;
96  int pion_pdgc = kPdgPi0;
97  if ( xcls.NPiPlus() == 1 )
98  pion_pdgc = kPdgPiP;
99  else if ( xcls.NPiMinus() == 1 )
100  pion_pdgc = kPdgPiM;
101  else if ( xcls.NPi0() != 1 )
102  throw genie::exceptions::InteractionException("Can't compute threshold");
103  double mpi = PDGLibrary::Instance()->Find(pion_pdgc)->Mass();
104  double mi = PDGLibrary::Instance()->Find( init_state.ProbePdg() )->Mass();
105  double mf = ml;
106  double mtot = Mf + mf + mpi; // total mass of FS particles
107  double Ethresh = (mtot*mtot - Mi*Mi - mi*mi)/2/Mi;
108  return Ethresh;
109  }
110 
111  if (pi.IsNorm() ) return 0;
112 
113  if (pi.IsSingleKaon()) {
114  int kaon_pdgc = xcls.StrangeHadronPdg();
115  double Mi = tgt.HitNucP4Ptr()->M(); // initial nucleon mass
116  // Final nucleon can be different for K0 interaction
117  double Mf = (xcls.NProtons()==1) ? kProtonMass : kNeutronMass;
118  double mk = PDGLibrary::Instance()->Find(kaon_pdgc)->Mass();
119  double mtot = Mf + ml + mk; // total mass of FS particles
120  double Ethresh = (mtot*mtot - Mi*Mi)/2/Mi;
121  return Ethresh;
122  }
123 
124  if(pi.IsCoherentElastic()) {
125  return ml + 0.5*ml*ml/tgt.Mass();
126  }
127 
128  if (pi.IsCoherentProduction()) {
129 
130  int tgtpdgc = tgt.Pdg(); // nuclear target PDG code (10LZZZAAAI)
131  double MA = PDGLibrary::Instance()->Find(tgtpdgc)->Mass();
132 
133  double m_other = controls::kASmallNum ;
134  // as a default the mass of hadronic system is the mass of the photon.
135  // which is assumed to be a small number to avoid divergences
136 
137  if ( xcls.NPions() > 0 ) {
138  m_other = pi.IsWeakCC() ? kPionMass : kPi0Mass;
139  }
140 
141  double m = ml + m_other ;
142  double m2 = TMath::Power(m,2);
143  double Ethr = m + 0.5*m2/MA;
144 
145  return TMath::Max(0.,Ethr);
146  }
147 
148  if(pi.IsQuasiElastic() ||
149  pi.IsDarkMatterElastic() ||
150  pi.IsInverseBetaDecay() ||
151  pi.IsResonant() ||
152  pi.IsDeepInelastic() ||
154  pi.IsDiffractive())
155  {
156  assert(tgt.HitNucIsSet());
157  double Mn = tgt.HitNucP4Ptr()->M();
158  double Mn2 = TMath::Power(Mn,2);
159  double Wmin = kNucleonMass + kPionMass;
160  if ( pi.IsQuasiElastic() || pi.IsDarkMatterElastic() || pi.IsInverseBetaDecay() ) {
161  int finalNucPDG = tgt.HitNucPdg();
162  if ( pi.IsWeakCC() ) finalNucPDG = pdg::SwitchProtonNeutron( finalNucPDG );
163  Wmin = PDGLibrary::Instance()->Find( finalNucPDG )->Mass();
164  }
165  if (pi.IsResonant()) {
166  Wmin = kNucleonMass + kPhotontest;
167  }
168 
169  if(xcls.IsCharmEvent()) {
170  if(xcls.IsInclusiveCharm()) {
172  } else {
173  int cpdg = xcls.CharmHadronPdg();
174  double mchm = PDGLibrary::Instance()->Find(cpdg)->Mass();
175  if(pi.IsQuasiElastic() || pi.IsInverseBetaDecay()) {
176  Wmin = mchm + controls::kASmallNum;
177  }
178  else {
179  Wmin = kNeutronMass + mchm + controls::kASmallNum;
180  }
181  }//incl.?
182  }//charm?
183 
184  double smin = TMath::Power(Wmin+ml,2.);
185  double Ethr = 0.5*(smin-Mn2)/Mn;
186  // threshold is different for dark matter case
188  // Correction to minimum kinematic variables
189  Wmin = Mn;
190  smin = TMath::Power(Wmin+ml,2);
191  Ethr = TMath::Max(0.5*(smin-Mn2-ml*ml)/Mn,ml);
192  }
193 
194  return TMath::Max(0.,Ethr);
195  }
196 
197  if(pi.IsInverseMuDecay() || pi.IsIMDAnnihilation()) {
198  double Ethr = 0.5 * (kMuonMass2-kElectronMass2)/kElectronMass;
199  return TMath::Max(0.,Ethr);
200  }
201 
203  return 0;
204  }
205  if(pi.IsAMNuGamma()) {
206  return 0;
207  }
208  if (pi.IsMEC()) {
209  if (tgt.HitNucIsSet()) {
210  double Mn = tgt.HitNucP4Ptr()->M();
211  double Mn2 = TMath::Power(Mn,2);
212  double Wmin = fInteraction->RecoilNucleon()->Mass(); // mass of the recoil nucleon cluster
213  double smin = TMath::Power(Wmin+ml,2.);
214  double Ethr = 0.5*(smin-Mn2)/Mn;
215  return TMath::Max(0.,Ethr);
216  }
217  else {
218  // this was ... if (pi.IsMECTensor())
219  return ml;
220  }
221  }
222  if(pi.IsGlashowResonance()) {
223  double Ethr = 0.5 * (ml*ml-kElectronMass2)/kElectronMass;
224  return TMath::Max(0.,Ethr);
225  }
226  if(pi.IsPhotonResonance()) {
227  double Mn = tgt.HitNucP4Ptr()->M();
228  double Ethr = 0.5 * (ml*ml-TMath::Power(Mn,2))/Mn;
229  return TMath::Max(0.,Ethr);
230  }
231  if(pi.IsPhotonCoherent()) {
232  double ml = 0;
233  if (pdg::IsNuE (TMath::Abs(init_state.ProbePdg()))) ml = kElectronMass;
234  else if (pdg::IsNuMu (TMath::Abs(init_state.ProbePdg()))) ml = kMuonMass;
235  else if (pdg::IsNuTau(TMath::Abs(init_state.ProbePdg()))) ml = kTauMass;
236  double MA = init_state.Tgt().Z()*kProtonMass + init_state.Tgt().N()*kNeutronMass;
237  double Ethr = 0.5 * (TMath::Power(kMw+ml,2)-TMath::Power(MA,2))/MA;
238  return TMath::Max(0.,Ethr);
239  }
240 
241 
242  SLOG("KPhaseSpace", pERROR)
243  << "Can't compute threshold for \n" << *fInteraction;
244  throw genie::exceptions::InteractionException("Can't compute threshold");
245  //exit(1);
246 
247  return 99999999;
248 }
bool IsResonant(void) const
Definition: ProcessInfo.cxx:99
bool IsPhotonResonance(void) const
bool IsNuTau(int pdgc)
Definition: PDGUtils.cxx:168
bool IsWeakCC(void) const
#define pERROR
Definition: Messenger.h:59
bool IsDarkMatterElectronElastic(void) const
static const double kNucleonMass
int HitNucPdg(void) const
Definition: Target.cxx:304
TParticlePDG * RecoilNucleon(void) const
recoil nucleon
int NPions(void) const
Definition: XclsTag.h:62
bool IsInverseMuDecay(void) const
bool IsQuasiElastic(void) const
Definition: ProcessInfo.cxx:69
int CharmHadronPdg(void) const
Definition: XclsTag.h:52
int Pdg(void) const
Definition: Target.h:71
bool IsNuE(int pdgc)
Definition: PDGUtils.cxx:158
int SwitchProtonNeutron(int pdgc)
Definition: PDGUtils.cxx:356
int NPiMinus(void) const
Definition: XclsTag.h:60
bool IsInverseBetaDecay(void) const
static const double kLightestChmHad
double Mass(Resonance_t res)
resonance mass (GeV)
bool IsDiffractive(void) const
bool IsCoherentProduction(void) const
bool IsIMDAnnihilation(void) const
double Mass(void) const
Definition: Target.cxx:224
static const double kElectronMass
bool IsKnown(void) const
Contains minimal information for tagging exclusive processes.
Definition: XclsTag.h:39
bool IsCharmEvent(void) const
Definition: XclsTag.h:50
bool IsSingleKaon(void) const
Definition: ProcessInfo.cxx:84
bool IsNuMu(int pdgc)
Definition: PDGUtils.cxx:163
int StrangeHadronPdg(void) const
Definition: XclsTag.h:55
bool IsCoherentElastic(void) const
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
static const double kElectronMass2
Exception used inside Interaction classes.
bool IsNuElectronElastic(void) const
static constexpr double m2
Definition: Units.h:72
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
bool IsAMNuGamma(void) const
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
int ProbePdg(void) const
Definition: InitialState.h:64
static const double kNeutronMass
int NPiPlus(void) const
Definition: XclsTag.h:59
const int kPdgPiP
Definition: PDGCodes.h:158
const int kPdgPi0
Definition: PDGCodes.h:160
int Z(void) const
Definition: Target.h:68
static const double kASmallNum
Definition: Controls.h:40
bool IsDarkMatterElastic(void) const
Definition: ProcessInfo.cxx:74
int NPi0(void) const
Definition: XclsTag.h:58
bool IsPhotonCoherent(void) const
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
bool IsMEC(void) const
bool IsNorm(void) const
bool IsDeepInelastic(void) const
Definition: ProcessInfo.cxx:89
int N(void) const
Definition: Target.h:69
TLorentzVector * HitNucP4Ptr(void) const
Definition: Target.cxx:247
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:68
bool HitNucIsSet(void) const
Definition: Target.cxx:283
const XclsTag & ExclTag(void) const
Definition: Interaction.h:72
bool IsInclusiveCharm(void) const
Definition: XclsTag.cxx:54
const int kPdgPiM
Definition: PDGCodes.h:159
const InitialState & InitState(void) const
Definition: Interaction.h:69
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:86
const Target & Tgt(void) const
Definition: InitialState.h:66
bool IsSinglePion(void) const
Definition: ProcessInfo.cxx:79
bool IsGlashowResonance(void) const
#define SLOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a short string (using the FUNCTION and...
Definition: Messenger.h:84
bool IsDarkMatterDeepInelastic(void) const
Definition: ProcessInfo.cxx:94
static constexpr double m
Definition: Units.h:71
Initial State information.
Definition: InitialState.h:48
int NProtons(void) const
Definition: XclsTag.h:56
double KPhaseSpace::Threshold_SPP_iso ( void  ) const

Energy limit for resonance single pion production on isoscalar nucleon.

Definition at line 999 of file KPhaseSpace.cxx.

References genie::PDGLibrary::Find(), fInteraction, genie::Interaction::FSPrimLepton(), genie::Interaction::InitState(), genie::PDGLibrary::Instance(), genie::kPdgNeutron, genie::kPdgPi0, genie::kPdgPiM, genie::kPdgPiP, genie::kPdgProton, genie::utils::res::Mass(), and genie::InitialState::ProbePdg().

Referenced by genie::SPPXSecWithCache::CacheResExcitationXSec(), genie::SPPXSec::Integrate(), and genie::MKSPPPXSec2020::ValidKinematics().

1000 {
1001  const InitialState & init_state = fInteraction->InitState();
1002  PDGLibrary * pdglib = PDGLibrary::Instance();
1003 
1004  // imply isospin symmetry
1005  double mpi = (pdglib->Find(kPdgPiP)->Mass() + pdglib->Find(kPdgPi0)->Mass() + pdglib->Find(kPdgPiM)->Mass())/3;
1006  double M = (pdglib->Find(kPdgProton)->Mass() + pdglib->Find(kPdgNeutron)->Mass())/2;
1007  double mi = PDGLibrary::Instance()->Find( init_state.ProbePdg() )->Mass();
1008  double mf = fInteraction->FSPrimLepton()->Mass();
1009  double mtot = M + mf + mpi; // total mass of FS particles
1010  double Ethresh = (mtot*mtot - M*M - mi*mi)/2/M;
1011  return Ethresh;
1012 }
double Mass(Resonance_t res)
resonance mass (GeV)
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
int ProbePdg(void) const
Definition: InitialState.h:64
const int kPdgPiP
Definition: PDGCodes.h:158
const int kPdgPi0
Definition: PDGCodes.h:160
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:68
Singleton class to load &amp; serve a TDatabasePDG.
Definition: PDGLibrary.h:35
const int kPdgPiM
Definition: PDGCodes.h:159
const InitialState & InitState(void) const
Definition: Interaction.h:69
const int kPdgProton
Definition: PDGCodes.h:81
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:86
const int kPdgNeutron
Definition: PDGCodes.h:83
Initial State information.
Definition: InitialState.h:48
Range1D_t KPhaseSpace::TLim ( void  ) const

t limits

Definition at line 915 of file KPhaseSpace.cxx.

References fInteraction, GetTMaxDFR(), genie::Target::HitNucMass(), genie::Interaction::InitState(), genie::ProcessInfo::IsCoherentProduction(), genie::ProcessInfo::IsDiffractive(), genie::ProcessInfo::IsWeakCC(), genie::controls::kASmallNum, genie::Interaction::Kine(), genie::constants::kPi0Mass, genie::constants::kPionMass, genie::kRfHitNucRest, LOG, genie::Range1D_t::max, genie::Range1D_t::min, genie::XclsTag::NPions(), pERROR, genie::InitialState::ProbeE(), genie::Interaction::ProcInfo(), pWARN, genie::Kinematics::Q2(), genie::utils::kinematics::Q2(), genie::InitialState::Tgt(), genie::utils::kinematics::UpdateWQ2FromXY(), and genie::Kinematics::y().

Referenced by genie::DFRKinematicsGenerator::ComputeMaxXSec(), IsAllowed(), and Limits().

916 {
917  // t limits for Coherent pion production from
918  // Kartavtsev, Paschos, and Gounaris, PRD 74 054007, and
919  // Paschos and Schalla, PRD 80, 03305
920  // TODO: Attempt to assign t bounds for other reactions?
921  Range1D_t tl;
922  tl.min = -1;
923  tl.max = -1;
924 
925  const InitialState & init_state = fInteraction->InitState();
926  const ProcessInfo & pi = fInteraction->ProcInfo();
927  const Kinematics & kine = fInteraction->Kine();
929  double Ev = init_state.ProbeE(kRfHitNucRest);
930  double Q2 = kine.Q2();
931  double nu = Ev * kine.y();
932 
933  //COH
934  if(pi.IsCoherentProduction()) {
935 
936  double m_other = controls::kASmallNum ;
937  // as a default the mass of hadronic system is the mass of the photon.
938  // which is assumed to be a small number to avoid divergences
939 
940  const XclsTag & xcls = fInteraction -> ExclTag() ;
941 
942  if ( xcls.NPions() > 0 ) {
943  bool pionIsCharged = pi.IsWeakCC();
944  m_other = pionIsCharged ? kPionMass : kPi0Mass;
945  }
946 
947  double m_other2 = m_other * m_other ;
948 
949  tl.min = 1.0 * (Q2 + m_other2)/(2.0 * nu) * (Q2 + m_other2)/(2.0 * nu);
950  tl.max = 0.05;
951  return tl;
952  }
953  // DFR
954  else if (pi.IsDiffractive()) {
955 
956  // diffractive tmin from Nucl.Phys.B278,61 (1986), eq. 12
957 
958  bool pionIsCharged = pi.IsWeakCC();
959  double mpi = pionIsCharged ? kPionMass : kPi0Mass;
960  double mpi2 = mpi*mpi;
961 
962  double M = init_state.Tgt().HitNucMass();
963  double M2 = M*M;
964  double nuSqPlusQ2 = nu*nu + Q2;
965  double nuOverM = nu / M;
966  double mpiQ2term = mpi2 - Q2 - 2*nu*nu;
967  double A1 = 1 + 2*nuOverM + nuOverM*nuOverM - nuSqPlusQ2/M2;
968  double A2 = (1+nuOverM) * mpiQ2term + 2*nuOverM*nuSqPlusQ2;
969  double A3 = mpiQ2term*mpiQ2term - 4*nuSqPlusQ2*(nu*nu - mpi2);
970 
971  tl.min = std::abs( (A2 + sqrt(A2*A2 - A1*A3)) / A1 ); // GENIE's convention is that t is positive
972  bool tminIsNaN;
973  // use std::isnan when C++11 is around
974 #if __cplusplus >= 201103L
975  tminIsNaN = std::isnan(tl.min);
976 #else
977  // this the old-fashioned way to check for NaN:
978  // NaN's aren't equal to anything, including themselves
979  tminIsNaN = tl.min != tl.min;
980 #endif
981  if (tminIsNaN)
982  {
983  LOG("KPhaseSpace", pERROR)
984  << "tmin for diffractive scattering is NaN "
985  << "( Enu = " << Ev << ", Q2 = " << Q2 << ", nu = " << nu << ")";
986  throw genie::exceptions::InteractionException("NaN tmin for diffractive scattering");
987  }
988  tl.max = this->GetTMaxDFR();
989 
990  return tl;
991  }
992 
993  // RES+DIS
994  // IMD
995  LOG("KPhaseSpace", pWARN) << "It is not sensible to ask for t limits for events that are not coherent or diffractive.";
996  return tl;
997 }
bool IsWeakCC(void) const
#define pERROR
Definition: Messenger.h:59
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:1077
A simple [min,max] interval for doubles.
Definition: Range1.h:42
int NPions(void) const
Definition: XclsTag.h:62
double HitNucMass(void) const
Definition: Target.cxx:233
Generated/set kinematical variables for an event.
Definition: Kinematics.h:39
bool IsDiffractive(void) const
bool IsCoherentProduction(void) const
static double GetTMaxDFR()
Definition: KPhaseSpace.cxx:59
Contains minimal information for tagging exclusive processes.
Definition: XclsTag.h:39
double y(bool selected=false) const
Definition: Kinematics.cxx:112
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
Exception used inside Interaction classes.
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
const Kinematics & Kine(void) const
Definition: Interaction.h:71
static const double kASmallNum
Definition: Controls.h:40
#define pWARN
Definition: Messenger.h:60
void UpdateWQ2FromXY(const Interaction *in)
Definition: KineUtils.cxx:1290
double max
Definition: Range1.h:53
const InitialState & InitState(void) const
Definition: Interaction.h:69
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
double min
Definition: Range1.h:52
double Q2(bool selected=false) const
Definition: Kinematics.cxx:125
const Target & Tgt(void) const
Definition: InitialState.h:66
double ProbeE(RefFrame_t rf) const
Initial State information.
Definition: InitialState.h:48
void KPhaseSpace::UseInteraction ( const Interaction in)

Definition at line 77 of file KPhaseSpace.cxx.

References fInteraction.

Referenced by KPhaseSpace().

78 {
79  fInteraction = in;
80 }
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
Range1D_t KPhaseSpace::WLim ( void  ) const

W limits.

Definition at line 440 of file KPhaseSpace.cxx.

References genie::utils::kinematics::DarkWLim(), genie::Interaction::ExclTag(), fInteraction, genie::Interaction::FSPrimLepton(), genie::Target::HitNucP4Ptr(), genie::utils::kinematics::InelWLim(), genie::Interaction::InitState(), genie::XclsTag::IsCharmEvent(), genie::ProcessInfo::IsDarkMatterDeepInelastic(), genie::ProcessInfo::IsDarkMatterElastic(), genie::ProcessInfo::IsDeepInelastic(), genie::ProcessInfo::IsDiffractive(), genie::ProcessInfo::IsEM(), genie::ProcessInfo::IsInverseBetaDecay(), genie::ProcessInfo::IsQuasiElastic(), genie::ProcessInfo::IsResonant(), genie::constants::kLightestChmHad, genie::constants::kNeutronMass, genie::constants::kPionMass, genie::kRfHitNucRest, LOG, genie::Range1D_t::max, genie::Range1D_t::min, pDEBUG, genie::InitialState::ProbeE(), genie::Interaction::ProcInfo(), genie::Interaction::RecoilNucleon(), and genie::InitialState::Tgt().

Referenced by genie::DISXSec::CacheFreeNucleonXSec(), genie::DMDISXSec::CacheFreeNucleonXSec(), genie::DFRKinematicsGenerator::ComputeMaxXSec(), genie::RESKinematicsGenerator::ComputeMaxXSec(), genie::utils::gsl::d2XSecRESFast_dWQ2_E::d2XSecRESFast_dWQ2_E(), genie::DISXSec::Integrate(), genie::DMDISXSec::Integrate(), IsAllowed(), Limits(), and genie::HEDISPXSec::XSec().

441 {
442 // Computes hadronic invariant mass limits.
443 // For QEL the range reduces to the recoil nucleon mass.
444 // For DIS & RES the calculation proceeds as in kinematics::InelWLim().
445 // It is not computed for other interactions
446 //
447  Range1D_t Wl;
448  Wl.min = -1;
449  Wl.max = -1;
450 
451  const ProcessInfo & pi = fInteraction->ProcInfo();
452 
453  bool is_em = pi.IsEM();
454  bool is_qel = pi.IsQuasiElastic() || pi.IsInverseBetaDecay() || pi.IsDarkMatterElastic();
455  bool is_inel = pi.IsDeepInelastic() || pi.IsResonant() || pi.IsDiffractive();
456  bool is_dmdis = pi.IsDarkMatterDeepInelastic();
457 
458  if(is_qel) {
459  double MR = fInteraction->RecoilNucleon()->Mass();
460  Wl.min = MR;
461  Wl.max = MR;
462  return Wl;
463  }
464  if(is_inel) {
465  const InitialState & init_state = fInteraction->InitState();
466  double Ev = init_state.ProbeE(kRfHitNucRest);
467  double M = init_state.Tgt().HitNucP4Ptr()->M(); //can be off m/shell
468  double ml = fInteraction->FSPrimLepton()->Mass();
469 
470  Wl = is_em ? kinematics::electromagnetic::InelWLim(Ev,ml,M) : kinematics::InelWLim(Ev,M,ml);
471 
473  //Wl.min = TMath::Max(Wl.min, kNeutronMass+kPionMass+kLightestChmHad);
474  Wl.min = TMath::Max(Wl.min, kNeutronMass+kLightestChmHad);
475  }
476  else if (fInteraction->ProcInfo().IsDiffractive())
477  Wl.min = TMath::Max(Wl.min, kNeutronMass+kPionMass);
478  else if ( fInteraction->ProcInfo().IsDeepInelastic() ) {
479  Wl.min = TMath::Max(Wl.min, kNeutronMass + kPionMass );
480  }
481 
482  // sanity check
483  if(Wl.min>Wl.max) {Wl.min=-1; Wl.max=-1;}
484 
485  return Wl;
486  }
487  if(is_dmdis) {
488  const InitialState & init_state = fInteraction->InitState();
489  double Ev = init_state.ProbeE(kRfHitNucRest);
490  double M = init_state.Tgt().HitNucP4Ptr()->M(); //can be off m/shell
491  double ml = fInteraction->FSPrimLepton()->Mass();
492  Wl = kinematics::DarkWLim(Ev,M,ml);
494  //Wl.min = TMath::Max(Wl.min, kNeutronMass+kPionMass+kLightestChmHad);
495  Wl.min = TMath::Max(Wl.min, kNeutronMass+kLightestChmHad);
496  }
497  else if (fInteraction->ProcInfo().IsDiffractive())
498  Wl.min = TMath::Max(Wl.min, kNeutronMass+kPionMass);
499 
500  LOG("KPhaseSpace", pDEBUG) << "Found nominal limits: " << Wl.min << ", " << Wl.max;
501 
502  // sanity check
503  if(Wl.min>Wl.max) {Wl.min=-1; Wl.max=-1;}
504 
505  return Wl;
506  }
507  return Wl;
508 }
bool IsResonant(void) const
Definition: ProcessInfo.cxx:99
TParticlePDG * RecoilNucleon(void) const
recoil nucleon
Range1D_t InelWLim(double Ev, double M, double ml)
Definition: KineUtils.cxx:358
A simple [min,max] interval for doubles.
Definition: Range1.h:42
bool IsQuasiElastic(void) const
Definition: ProcessInfo.cxx:69
bool IsInverseBetaDecay(void) const
static const double kLightestChmHad
bool IsDiffractive(void) const
bool IsCharmEvent(void) const
Definition: XclsTag.h:50
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
static const double kNeutronMass
bool IsDarkMatterElastic(void) const
Definition: ProcessInfo.cxx:74
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
bool IsEM(void) const
bool IsDeepInelastic(void) const
Definition: ProcessInfo.cxx:89
double max
Definition: Range1.h:53
TLorentzVector * HitNucP4Ptr(void) const
Definition: Target.cxx:247
const XclsTag & ExclTag(void) const
Definition: Interaction.h:72
const InitialState & InitState(void) const
Definition: Interaction.h:69
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
double min
Definition: Range1.h:52
const Target & Tgt(void) const
Definition: InitialState.h:66
bool IsDarkMatterDeepInelastic(void) const
Definition: ProcessInfo.cxx:94
double ProbeE(RefFrame_t rf) const
Range1D_t DarkWLim(double Ev, double M, double ml)
Definition: KineUtils.cxx:892
Initial State information.
Definition: InitialState.h:48
#define pDEBUG
Definition: Messenger.h:63
Range1D_t KPhaseSpace::WLim_SPP ( void  ) const

W limits for single pion production models.

Definition at line 1014 of file KPhaseSpace.cxx.

References genie::InitialState::CMEnergy(), epsilon, genie::PDGLibrary::Find(), genie::SppChannel::FinStateNucleon(), genie::SppChannel::FinStatePion(), fInteraction, genie::SppChannel::FromInteraction(), genie::Interaction::FSPrimLepton(), genie::Interaction::InitState(), genie::PDGLibrary::Instance(), genie::Range1D_t::max, and genie::Range1D_t::min.

1015 {
1016  Range1D_t Wl;
1017  const InitialState & init_state = fInteraction->InitState();
1019  PDGLibrary * pdglib = PDGLibrary::Instance();
1020  double Mf = pdglib->Find( SppChannel::FinStateNucleon(spp_channel) )->Mass();
1021  double mpi = pdglib->Find( SppChannel::FinStatePion(spp_channel) )->Mass();
1022  double mf = fInteraction->FSPrimLepton()->Mass();
1023  double ECM = init_state.CMEnergy();
1024  // kinematic W-limits
1025  Wl.min = Mf + mpi;
1026  Wl.max = ECM - mf;
1027 
1028  if ( (Wl.max - Wl.min) < (Wl.max + Wl.min)*std::numeric_limits<double>::epsilon() )
1029  {
1030  Wl.min = 2*Wl.max*Wl.min/(Wl.max + Wl.min);
1031  Wl.max = Wl.min;
1032  }
1033  else
1034  {
1037  }
1038 
1039  return Wl;
1040 }
static SppChannel_t FromInteraction(const Interaction *interaction)
Definition: SppChannel.h:402
A simple [min,max] interval for doubles.
Definition: Range1.h:42
static int FinStateNucleon(SppChannel_t channel)
Definition: SppChannel.h:130
const double epsilon
enum genie::ESppChannel SppChannel_t
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
static int FinStatePion(SppChannel_t channel)
Definition: SppChannel.h:157
double max
Definition: Range1.h:53
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:68
Singleton class to load &amp; serve a TDatabasePDG.
Definition: PDGLibrary.h:35
double CMEnergy() const
centre-of-mass energy (sqrt s)
const InitialState & InitState(void) const
Definition: Interaction.h:69
double min
Definition: Range1.h:52
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:86
Initial State information.
Definition: InitialState.h:48
Range1D_t KPhaseSpace::WLim_SPP_iso ( void  ) const

W limits for resonance single pion production on isoscalar nucleon.

Definition at line 1042 of file KPhaseSpace.cxx.

References epsilon, genie::PDGLibrary::Find(), fInteraction, genie::Interaction::FSPrimLepton(), genie::Interaction::InitState(), genie::PDGLibrary::Instance(), genie::kPdgNeutron, genie::kPdgPi0, genie::kPdgPiM, genie::kPdgPiP, genie::kPdgProton, genie::kRfHitNucRest, genie::utils::res::Mass(), genie::Range1D_t::max, genie::Range1D_t::min, genie::InitialState::ProbeE(), and genie::InitialState::ProbePdg().

Referenced by genie::SPPEventGenerator::ComputeMaxXSec(), genie::utils::gsl::d3XSecMK_dWQ2CosTheta_E::d3XSecMK_dWQ2CosTheta_E(), genie::utils::gsl::d4XSecMK_dWQ2CosThetaPhi_E::d4XSecMK_dWQ2CosThetaPhi_E(), genie::SPPEventGenerator::ProcessEventRecord(), and genie::MKSPPPXSec2020::ValidKinematics().

1043 {
1044  Range1D_t Wl;
1045  const InitialState & init_state = fInteraction->InitState();
1046  PDGLibrary * pdglib = PDGLibrary::Instance();
1047  // imply isospin symmetry
1048  double M = (pdglib->Find(kPdgProton)->Mass() + pdglib->Find(kPdgNeutron)->Mass())/2;
1049  double mpi = (pdglib->Find(kPdgPiP)->Mass() + pdglib->Find(kPdgPi0)->Mass() + pdglib->Find(kPdgPiM)->Mass())/3;
1050  double mi = PDGLibrary::Instance()->Find( init_state.ProbePdg() )->Mass();
1051  double mf = fInteraction->FSPrimLepton()->Mass();
1052  double Ei = init_state.ProbeE(kRfHitNucRest);
1053  double ECM = TMath::Sqrt(M*(M + 2*Ei) + mi*mi);
1054  // kinematic W-limits
1055  Wl.min = M + mpi;
1056  Wl.max = ECM - mf;
1057 
1058  if ( (Wl.max - Wl.min) < (Wl.max + Wl.min)*std::numeric_limits<double>::epsilon() )
1059  {
1060  Wl.min = 2*Wl.max*Wl.min/(Wl.max + Wl.min);
1061  Wl.max = Wl.min;
1062  }
1063  else
1064  {
1067  }
1068 
1069  return Wl;
1070 }
A simple [min,max] interval for doubles.
Definition: Range1.h:42
double Mass(Resonance_t res)
resonance mass (GeV)
const double epsilon
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
int ProbePdg(void) const
Definition: InitialState.h:64
const int kPdgPiP
Definition: PDGCodes.h:158
const int kPdgPi0
Definition: PDGCodes.h:160
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
double max
Definition: Range1.h:53
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:68
Singleton class to load &amp; serve a TDatabasePDG.
Definition: PDGLibrary.h:35
const int kPdgPiM
Definition: PDGCodes.h:159
const InitialState & InitState(void) const
Definition: Interaction.h:69
double min
Definition: Range1.h:52
const int kPdgProton
Definition: PDGCodes.h:81
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:86
double ProbeE(RefFrame_t rf) const
const int kPdgNeutron
Definition: PDGCodes.h:83
Initial State information.
Definition: InitialState.h:48
Range1D_t KPhaseSpace::XLim ( void  ) const

x limits

Definition at line 691 of file KPhaseSpace.cxx.

References genie::utils::kinematics::CohXLim(), genie::utils::kinematics::DarkXLim(), fInteraction, genie::Interaction::FSPrimLepton(), genie::Target::HitNucP4Ptr(), genie::utils::kinematics::InelXLim(), genie::Interaction::InitState(), genie::ProcessInfo::IsCoherentProduction(), genie::ProcessInfo::IsDarkMatterDeepInelastic(), genie::ProcessInfo::IsDarkMatterElastic(), genie::ProcessInfo::IsDeepInelastic(), genie::ProcessInfo::IsDiffractive(), genie::ProcessInfo::IsEM(), genie::ProcessInfo::IsInverseBetaDecay(), genie::ProcessInfo::IsQuasiElastic(), genie::ProcessInfo::IsResonant(), genie::controls::kASmallNum, genie::kRfHitNucRest, genie::Range1D_t::max, genie::Range1D_t::min, genie::InitialState::ProbeE(), genie::Interaction::ProcInfo(), and genie::InitialState::Tgt().

Referenced by genie::DFRKinematicsGenerator::ComputeMaxXSec(), genie::HEDISXSec::Integrate(), IsAllowed(), Limits(), and genie::HEDISKinematicsGenerator::ProcessEventRecord().

692 {
693  // Computes x-limits;
694 
695  Range1D_t xl;
696  xl.min = -1;
697  xl.max = -1;
698 
699  const ProcessInfo & pi = fInteraction->ProcInfo();
700  bool is_em = pi.IsEM();
701 
702  //RES+DIS
703  bool is_inel = pi.IsDeepInelastic() || pi.IsResonant();
704  if(is_inel) {
705  const InitialState & init_state = fInteraction->InitState();
706  double Ev = init_state.ProbeE(kRfHitNucRest);
707  double M = init_state.Tgt().HitNucP4Ptr()->M(); // can be off m/shell
708  double ml = fInteraction->FSPrimLepton()->Mass();
709  xl = is_em ? kinematics::electromagnetic::InelXLim(Ev,ml,M) : kinematics::InelXLim(Ev,M,ml);
710  return xl;
711  }
712  //DMDIS
713  bool is_dmdis = pi.IsDarkMatterDeepInelastic();
714  if(is_dmdis) {
715  const InitialState & init_state = fInteraction->InitState();
716  double Ev = init_state.ProbeE(kRfHitNucRest);
717  double M = init_state.Tgt().HitNucP4Ptr()->M(); // can be off m/shell
718  double ml = fInteraction->FSPrimLepton()->Mass();
719  xl = kinematics::DarkXLim(Ev,M,ml);
720  return xl;
721  }
722  //COH
723  bool is_coh = pi.IsCoherentProduction();
724  if(is_coh) {
725  xl = kinematics::CohXLim();
726  return xl;
727  }
728  //QEL
729  bool is_qel = pi.IsQuasiElastic() || pi.IsInverseBetaDecay() || pi.IsDarkMatterElastic();
730  if(is_qel) {
731  xl.min = 1;
732  xl.max = 1;
733  return xl;
734  }
735  bool is_dfr = pi.IsDiffractive();
736  if(is_dfr) {
738  xl.max = 1. - controls::kASmallNum;
739  return xl;
740  }
741 
742  return xl;
743 }
bool IsResonant(void) const
Definition: ProcessInfo.cxx:99
A simple [min,max] interval for doubles.
Definition: Range1.h:42
bool IsQuasiElastic(void) const
Definition: ProcessInfo.cxx:69
bool IsInverseBetaDecay(void) const
bool IsDiffractive(void) const
bool IsCoherentProduction(void) const
Range1D_t InelXLim(double Ev, double M, double ml)
Definition: KineUtils.cxx:457
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
Range1D_t CohXLim(void)
Definition: KineUtils.cxx:735
static const double kASmallNum
Definition: Controls.h:40
bool IsDarkMatterElastic(void) const
Definition: ProcessInfo.cxx:74
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
bool IsEM(void) const
bool IsDeepInelastic(void) const
Definition: ProcessInfo.cxx:89
double max
Definition: Range1.h:53
TLorentzVector * HitNucP4Ptr(void) const
Definition: Target.cxx:247
const InitialState & InitState(void) const
Definition: Interaction.h:69
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
double min
Definition: Range1.h:52
const Target & Tgt(void) const
Definition: InitialState.h:66
bool IsDarkMatterDeepInelastic(void) const
Definition: ProcessInfo.cxx:94
double ProbeE(RefFrame_t rf) const
Range1D_t DarkXLim(double Ev, double M, double ml)
Definition: KineUtils.cxx:1001
Initial State information.
Definition: InitialState.h:48
Range1D_t KPhaseSpace::YLim ( void  ) const

y limits

Definition at line 745 of file KPhaseSpace.cxx.

References genie::utils::kinematics::CohYLim(), genie::utils::kinematics::DarkYLim(), fInteraction, genie::Interaction::FSPrimLepton(), genie::Target::HitNucP4Ptr(), genie::utils::kinematics::InelYLim(), genie::Interaction::InitState(), genie::ProcessInfo::IsCoherentProduction(), genie::ProcessInfo::IsDarkMatterDeepInelastic(), genie::ProcessInfo::IsDarkMatterElectronElastic(), genie::ProcessInfo::IsDeepInelastic(), genie::ProcessInfo::IsDiffractive(), genie::ProcessInfo::IsEM(), genie::ProcessInfo::IsIMDAnnihilation(), genie::ProcessInfo::IsInverseMuDecay(), genie::ProcessInfo::IsNuElectronElastic(), genie::ProcessInfo::IsResonant(), genie::controls::kASmallNum, genie::constants::kElectronMass, genie::constants::kPionMass, genie::kRfHitNucRest, genie::kRfLab, genie::Range1D_t::max, genie::Range1D_t::min, genie::InitialState::ProbeE(), genie::Interaction::ProcInfo(), and genie::InitialState::Tgt().

Referenced by genie::COHKinematicsGenerator::CalculateKin_BergerSehgal(), genie::COHKinematicsGenerator::CalculateKin_BergerSehgalFM(), genie::COHKinematicsGenerator::CalculateKin_ReinSehgal(), genie::DFRKinematicsGenerator::ComputeMaxXSec(), IsAllowed(), Limits(), genie::COHKinematicsGenerator::MaxXSec_AlvarezRuso(), genie::COHKinematicsGenerator::MaxXSec_BergerSehgal(), genie::COHKinematicsGenerator::MaxXSec_BergerSehgalFM(), genie::COHKinematicsGenerator::MaxXSec_ReinSehgal(), YLim(), and YLim_X().

746 {
747  Range1D_t yl;
748  yl.min = -1;
749  yl.max = -1;
750 
751  const ProcessInfo & pi = fInteraction->ProcInfo();
752  bool is_em = pi.IsEM();
753 
754  //RES+DIS
755  bool is_inel = pi.IsDeepInelastic() || pi.IsResonant();
756  if(is_inel) {
757  const InitialState & init_state = fInteraction->InitState();
758  double Ev = init_state.ProbeE(kRfHitNucRest);
759  double M = init_state.Tgt().HitNucP4Ptr()->M(); // can be off m/shell
760  double ml = fInteraction->FSPrimLepton()->Mass();
761  yl = is_em ? kinematics::electromagnetic::InelYLim(Ev,ml,M) : kinematics::InelYLim(Ev,M,ml);
762  return yl;
763  }
764  //DMDIS
765  bool is_dmdis = pi.IsDarkMatterDeepInelastic();
766  if(is_dmdis) {
767  const InitialState & init_state = fInteraction->InitState();
768  double Ev = init_state.ProbeE(kRfHitNucRest);
769  double M = init_state.Tgt().HitNucP4Ptr()->M(); // can be off m/shell
770  double ml = fInteraction->FSPrimLepton()->Mass();
771  yl = kinematics::DarkYLim(Ev,M,ml);
772  return yl;
773  }
774  //COH
775  bool is_coh = pi.IsCoherentProduction();
776  if(is_coh) {
777  const InitialState & init_state = fInteraction->InitState();
778  double EvL = init_state.ProbeE(kRfLab);
779  double ml = fInteraction->FSPrimLepton()->Mass();
780  yl = kinematics::CohYLim(EvL,ml);
781  return yl;
782  }
783  // IMD
784  if(pi.IsInverseMuDecay() || pi.IsIMDAnnihilation() || pi.IsNuElectronElastic()) {
785  const InitialState & init_state = fInteraction->InitState();
786  double Ev = init_state.ProbeE(kRfLab);
787  double ml = fInteraction->FSPrimLepton()->Mass();
788  double me = kElectronMass;
790  yl.max = 1 - (ml*ml + me*me)/(2*me*Ev) - controls::kASmallNum;
791  return yl;
792  }
793  // EDIT: y limits are different for massive probe
794  if(pi.IsDarkMatterElectronElastic()) {
795  const InitialState & init_state = fInteraction->InitState();
796  double Ev = init_state.ProbeE(kRfLab);
797  double ml = fInteraction->FSPrimLepton()->Mass();
798  double me = kElectronMass;
799  yl.min = (Ev*me*me + ml*ml*(Ev + 2.0*me)) / (Ev * (2.0*Ev*me + me*me + ml*ml)) + controls::kASmallNum;
800  yl.max = 1.0 - controls::kASmallNum;
801  return yl;
802  }
803  bool is_dfr = pi.IsDiffractive();
804  if(is_dfr) {
805  const InitialState & init_state = fInteraction -> InitState();
806  double Ev = init_state.ProbeE(kRfHitNucRest);
807  double ml = fInteraction->FSPrimLepton()->Mass();
809  yl.max = 1. -ml/Ev - controls::kASmallNum;
810  return yl;
811  }
812  return yl;
813 }
bool IsResonant(void) const
Definition: ProcessInfo.cxx:99
bool IsDarkMatterElectronElastic(void) const
A simple [min,max] interval for doubles.
Definition: Range1.h:42
bool IsInverseMuDecay(void) const
bool IsDiffractive(void) const
bool IsCoherentProduction(void) const
bool IsIMDAnnihilation(void) const
static const double kElectronMass
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
bool IsNuElectronElastic(void) const
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
static const double kASmallNum
Definition: Controls.h:40
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
bool IsEM(void) const
bool IsDeepInelastic(void) const
Definition: ProcessInfo.cxx:89
double max
Definition: Range1.h:53
TLorentzVector * HitNucP4Ptr(void) const
Definition: Target.cxx:247
Range1D_t CohYLim(double Mn, double m_produced, double mlep, double Ev, double Q2, double xsi)
Definition: KineUtils.cxx:840
Range1D_t DarkYLim(double Ev, double M, double ml)
Definition: KineUtils.cxx:1021
const InitialState & InitState(void) const
Definition: Interaction.h:69
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
double min
Definition: Range1.h:52
const Target & Tgt(void) const
Definition: InitialState.h:66
bool IsDarkMatterDeepInelastic(void) const
Definition: ProcessInfo.cxx:94
Range1D_t InelYLim(double Ev, double M, double ml)
Definition: KineUtils.cxx:477
double ProbeE(RefFrame_t rf) const
Initial State information.
Definition: InitialState.h:48
Range1D_t KPhaseSpace::YLim ( double  xsi) const

y limits (COH)

Definition at line 860 of file KPhaseSpace.cxx.

References genie::utils::kinematics::CohYLim(), fInteraction, genie::Interaction::FSPrimLepton(), genie::Interaction::InitState(), genie::ProcessInfo::IsCoherentProduction(), genie::ProcessInfo::IsWeakCC(), genie::controls::kASmallNum, genie::Interaction::Kine(), genie::constants::kPi0Mass, genie::constants::kPionMass, genie::kRfHitNucRest, genie::Target::Mass(), genie::Range1D_t::max, genie::Range1D_t::min, genie::XclsTag::NPions(), genie::InitialState::ProbeE(), genie::Interaction::ProcInfo(), genie::Kinematics::Q2(), genie::utils::kinematics::Q2(), genie::InitialState::Tgt(), and YLim().

861 {
862  // Paschos-Schalla xsi parameter for y-limits in COH
863  // From PRD 80, 033005 (2009)
864 
865  Range1D_t yl;
866  yl.min = -1;
867  yl.max = -1;
868 
869  const ProcessInfo & pi = fInteraction->ProcInfo();
870 
871  //COH
872  bool is_coh = pi.IsCoherentProduction();
873  if(is_coh) {
874  const InitialState & init_state = fInteraction->InitState();
875  const Kinematics & kine = fInteraction->Kine();
876  double Ev = init_state.ProbeE(kRfHitNucRest);
877  double Q2 = kine.Q2();
878  double Mn = init_state.Tgt().Mass();
879  double mlep = fInteraction->FSPrimLepton()->Mass();
880 
881  double m_other = controls::kASmallNum ;
882  // as a default the mass of hadronic system is the mass of the photon.
883  // which is assumed to be a small number to avoid divergences
884 
885  const XclsTag & xcls = fInteraction -> ExclTag() ;
886 
887  if ( xcls.NPions() > 0 ) {
888  bool pionIsCharged = pi.IsWeakCC();
889  m_other = pionIsCharged ? kPionMass : kPi0Mass;
890  }
891 
892  yl = kinematics::CohYLim(Mn, m_other, mlep, Ev, Q2, xsi);
893  return yl;
894  } else {
895  return this->YLim();
896  }
897 }
bool IsWeakCC(void) const
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:1077
A simple [min,max] interval for doubles.
Definition: Range1.h:42
int NPions(void) const
Definition: XclsTag.h:62
Generated/set kinematical variables for an event.
Definition: Kinematics.h:39
Range1D_t YLim(void) const
y limits
bool IsCoherentProduction(void) const
double Mass(void) const
Definition: Target.cxx:224
Contains minimal information for tagging exclusive processes.
Definition: XclsTag.h:39
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
const Kinematics & Kine(void) const
Definition: Interaction.h:71
static const double kASmallNum
Definition: Controls.h:40
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
double max
Definition: Range1.h:53
Range1D_t CohYLim(double Mn, double m_produced, double mlep, double Ev, double Q2, double xsi)
Definition: KineUtils.cxx:840
const InitialState & InitState(void) const
Definition: Interaction.h:69
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
double min
Definition: Range1.h:52
double Q2(bool selected=false) const
Definition: Kinematics.cxx:125
const Target & Tgt(void) const
Definition: InitialState.h:66
double ProbeE(RefFrame_t rf) const
Initial State information.
Definition: InitialState.h:48
Range1D_t KPhaseSpace::YLim_X ( void  ) const

y limits @ fixed x

Definition at line 815 of file KPhaseSpace.cxx.

References genie::utils::kinematics::CohYLim(), genie::utils::kinematics::DarkYLim_X(), fInteraction, genie::Interaction::FSPrimLepton(), genie::Target::HitNucP4Ptr(), genie::utils::kinematics::InelYLim_X(), genie::Interaction::InitState(), genie::ProcessInfo::IsCoherentProduction(), genie::ProcessInfo::IsDarkMatterDeepInelastic(), genie::ProcessInfo::IsDeepInelastic(), genie::ProcessInfo::IsEM(), genie::ProcessInfo::IsResonant(), genie::Interaction::Kine(), genie::kRfHitNucRest, genie::kRfLab, genie::Range1D_t::max, genie::Range1D_t::min, genie::InitialState::ProbeE(), genie::Interaction::ProcInfo(), genie::InitialState::Tgt(), and genie::Kinematics::x().

Referenced by YLim_X().

816 {
817 // Computes kinematical limits for y @ the input x
818 
819  Range1D_t yl;
820  yl.min = -1;
821  yl.max = -1;
822 
823  const ProcessInfo & pi = fInteraction->ProcInfo();
824  bool is_em = pi.IsEM();
825 
826  //RES+DIS
827  bool is_inel = pi.IsDeepInelastic() || pi.IsResonant();
828  if(is_inel) {
829  const InitialState & init_state = fInteraction->InitState();
830  double Ev = init_state.ProbeE(kRfHitNucRest);
831  double M = init_state.Tgt().HitNucP4Ptr()->M(); // can be off m/shell
832  double ml = fInteraction->FSPrimLepton()->Mass();
833  double x = fInteraction->Kine().x();
834  yl = is_em ? kinematics::electromagnetic::InelYLim_X(Ev,ml,M,x) : kinematics::InelYLim_X(Ev,M,ml,x);
835  return yl;
836  }
837  //DMDIS
838  bool is_dmdis = pi.IsDarkMatterDeepInelastic();
839  if(is_dmdis) {
840  const InitialState & init_state = fInteraction->InitState();
841  double Ev = init_state.ProbeE(kRfHitNucRest);
842  double M = init_state.Tgt().HitNucP4Ptr()->M(); // can be off m/shell
843  double ml = fInteraction->FSPrimLepton()->Mass();
844  double x = fInteraction->Kine().x();
845  yl = kinematics::DarkYLim_X(Ev,M,ml,x);
846  return yl;
847  }
848  //COH
849  bool is_coh = pi.IsCoherentProduction();
850  if(is_coh) {
851  const InitialState & init_state = fInteraction->InitState();
852  double EvL = init_state.ProbeE(kRfLab);
853  double ml = fInteraction->FSPrimLepton()->Mass();
854  yl = kinematics::CohYLim(EvL,ml);
855  return yl;
856  }
857  return yl;
858 }
bool IsResonant(void) const
Definition: ProcessInfo.cxx:99
A simple [min,max] interval for doubles.
Definition: Range1.h:42
double x(bool selected=false) const
Definition: Kinematics.cxx:99
bool IsCoherentProduction(void) const
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
const Kinematics & Kine(void) const
Definition: Interaction.h:71
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
bool IsEM(void) const
bool IsDeepInelastic(void) const
Definition: ProcessInfo.cxx:89
double max
Definition: Range1.h:53
TLorentzVector * HitNucP4Ptr(void) const
Definition: Target.cxx:247
Range1D_t InelYLim_X(double Ev, double M, double ml, double x)
Definition: KineUtils.cxx:512
Range1D_t CohYLim(double Mn, double m_produced, double mlep, double Ev, double Q2, double xsi)
Definition: KineUtils.cxx:840
const InitialState & InitState(void) const
Definition: Interaction.h:69
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
double min
Definition: Range1.h:52
Range1D_t DarkYLim_X(double Ev, double M, double ml, double x)
Definition: KineUtils.cxx:1037
const Target & Tgt(void) const
Definition: InitialState.h:66
bool IsDarkMatterDeepInelastic(void) const
Definition: ProcessInfo.cxx:94
double ProbeE(RefFrame_t rf) const
Initial State information.
Definition: InitialState.h:48
Range1D_t KPhaseSpace::YLim_X ( double  xsi) const

y limits @ fixed x (COH)

Definition at line 899 of file KPhaseSpace.cxx.

References fInteraction, genie::ProcessInfo::IsCoherentProduction(), genie::Interaction::ProcInfo(), YLim(), and YLim_X().

900 {
901  // Paschos-Schalla xsi parameter for y-limits in COH
902  // From PRD 80, 033005 (2009)
903 
904  const ProcessInfo & pi = fInteraction->ProcInfo();
905 
906  //COH
907  bool is_coh = pi.IsCoherentProduction();
908  if(is_coh) {
909  return this->YLim(xsi);
910  } else {
911  return this->YLim_X();
912  }
913 }
Range1D_t YLim_X(void) const
y limits @ fixed x
Range1D_t YLim(void) const
y limits
bool IsCoherentProduction(void) const
const Interaction * fInteraction
Definition: KPhaseSpace.h:78
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70

Member Data Documentation

const Interaction* genie::KPhaseSpace::fInteraction
private

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