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

Generated/set kinematical variables for an event. More...

#include <Kinematics.h>

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

Public Member Functions

 Kinematics ()
 Kinematics (const Kinematics &kv)
 Kinematics (TRootIOCtor *)
 ~Kinematics ()
double x (bool selected=false) const
double y (bool selected=false) const
double Q2 (bool selected=false) const
double q2 (bool selected=false) const
double W (bool selected=false) const
double t (bool selected=false) const
double Logx (bool selected=false) const
double Logy (bool selected=false) const
double LogQ2 (bool selected=false) const
double LogW (bool selected=false) const
double Log10x (bool selected=false) const
double Log10y (bool selected=false) const
double Log10Q2 (bool selected=false) const
double Log10W (bool selected=false) const
const TLorentzVector & FSLeptonP4 (void) const
const TLorentzVector & HadSystP4 (void) const
void Setx (double x, bool selected=false)
void Sety (double y, bool selected=false)
void SetQ2 (double Q2, bool selected=false)
void Setq2 (double q2, bool selected=false)
void SetW (double W, bool selected=false)
void Sett (double t, bool selected=false)
void SetFSLeptonP4 (const TLorentzVector &p4)
void SetFSLeptonP4 (double px, double py, double pz, double E)
void SetHadSystP4 (const TLorentzVector &p4)
void SetHadSystP4 (double px, double py, double pz, double E)
bool KVSet (KineVar_t kv) const
double GetKV (KineVar_t kv) const
void SetKV (KineVar_t kv, double value)
void ClearRunningValues (void)
void UseSelectedKinematics (void)
void Reset (void)
void Copy (const Kinematics &kine)
void Print (ostream &stream) const
Kinematicsoperator= (const Kinematics &kine)

Private Member Functions

void Init (void)
 initialize
void CleanUp (void)
 clean-up

Private Attributes

map< KineVar_t, double > fKV
 selected kinematics
TLorentzVector * fP4Fsl
 generated final state primary lepton 4-p (LAB)
TLorentzVector * fP4HadSyst
 generated final state hadronic system 4-p (LAB)

Friends

ostream & operator<< (ostream &stream, const Kinematics &kine)

Detailed Description

Generated/set kinematical variables for an event.

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

Definition at line 39 of file Kinematics.h.

Constructor & Destructor Documentation

◆ Kinematics() [1/3]

Kinematics::Kinematics ( )

Definition at line 33 of file Kinematics.cxx.

33 :
34TObject()
35{
36 this->Init();
37}
void Init(void)
initialize

References Init().

Referenced by Copy(), Kinematics(), operator<<, and operator=().

◆ Kinematics() [2/3]

Kinematics::Kinematics ( const Kinematics & kv)

Definition at line 39 of file Kinematics.cxx.

39 :
40TObject()
41{
42 this->Init();
43 this->Copy(kinematics);
44}
void Copy(const Kinematics &kine)

References Copy(), Init(), and Kinematics().

◆ Kinematics() [3/3]

Kinematics::Kinematics ( TRootIOCtor * )

Definition at line 46 of file Kinematics.cxx.

46 :
47TObject(),
48fP4Fsl(0),
50{
51
52}
TLorentzVector * fP4Fsl
generated final state primary lepton 4-p (LAB)
Definition Kinematics.h:103
TLorentzVector * fP4HadSyst
generated final state hadronic system 4-p (LAB)
Definition Kinematics.h:104

References fP4Fsl, and fP4HadSyst.

◆ ~Kinematics()

Kinematics::~Kinematics ( )

Definition at line 54 of file Kinematics.cxx.

55{
56 this->CleanUp();
57}
void CleanUp(void)
clean-up

References CleanUp().

Member Function Documentation

◆ CleanUp()

void Kinematics::CleanUp ( void )
private

clean-up

Definition at line 67 of file Kinematics.cxx.

68{
69 fKV.clear();
70
71 delete fP4Fsl;
72 delete fP4HadSyst;
73}
map< KineVar_t, double > fKV
selected kinematics
Definition Kinematics.h:102

References fKV, fP4Fsl, and fP4HadSyst.

Referenced by ~Kinematics().

◆ ClearRunningValues()

void Kinematics::ClearRunningValues ( void )

Definition at line 347 of file Kinematics.cxx.

348{
349// clear the running values (leave the selected ones)
350//
351 fKV.erase( kKVx );
352 fKV.erase( kKVy );
353 fKV.erase( kKVQ2 );
354 fKV.erase( kKVq2 );
355 fKV.erase( kKVW );
356 fKV.erase( kKVt );
357}
@ kKVQ2
Definition KineVar.h:33
@ kKVx
Definition KineVar.h:31
@ kKVy
Definition KineVar.h:32
@ kKVW
Definition KineVar.h:35
@ kKVq2
Definition KineVar.h:34
@ kKVt
Definition KineVar.h:36

References fKV, genie::kKVQ2, genie::kKVq2, genie::kKVt, genie::kKVW, genie::kKVx, and genie::kKVy.

Referenced by genie::COHKinematicsGenerator::CalculateKin_AlvarezRuso(), genie::SKKinematicsGenerator::CalculateKin_AtharSingleKaon(), genie::COHKinematicsGenerator::CalculateKin_BergerSehgal(), genie::COHKinematicsGenerator::CalculateKin_BergerSehgalFM(), genie::COHKinematicsGenerator::CalculateKin_ReinSehgal(), genie::CEvNSEventGenerator::GenerateKinematics(), genie::COHDNuEventGenerator::GenerateKinematics(), main(), genie::DFRKinematicsGenerator::ProcessEventRecord(), genie::DISKinematicsGenerator::ProcessEventRecord(), genie::DMDISKinematicsGenerator::ProcessEventRecord(), genie::DMEKinematicsGenerator::ProcessEventRecord(), genie::DMELEventGenerator::ProcessEventRecord(), genie::DMELKinematicsGenerator::ProcessEventRecord(), genie::HEDISKinematicsGenerator::ProcessEventRecord(), genie::HELeptonKinematicsGenerator::ProcessEventRecord(), genie::IBDKinematicsGenerator::ProcessEventRecord(), genie::NuEKinematicsGenerator::ProcessEventRecord(), genie::QELEventGenerator::ProcessEventRecord(), genie::QELEventGeneratorSM::ProcessEventRecord(), genie::QELKinematicsGenerator::ProcessEventRecord(), genie::RESKinematicsGenerator::ProcessEventRecord(), genie::SPPEventGenerator::ProcessEventRecord(), genie::MECGenerator::SelectEmpiricalKinematics(), genie::RSPPResonanceSelector::SelectResonance(), genie::DMELKinematicsGenerator::SpectralFuncExperimentalCode(), and genie::QELKinematicsGenerator::SpectralFuncExperimentalCode().

◆ Copy()

void Kinematics::Copy ( const Kinematics & kine)

Definition at line 83 of file Kinematics.cxx.

84{
85 this->Reset();
86
87 map<KineVar_t, double>::const_iterator iter;
88
89 for(iter = kinematics.fKV.begin(); iter != kinematics.fKV.end(); ++iter) {
90 KineVar_t kv = iter->first;
91 double val = iter->second;
92 this->SetKV(kv,val);
93 }
94
95 this->SetFSLeptonP4 (*kinematics.fP4Fsl);
96 this->SetHadSystP4 (*kinematics.fP4HadSyst);
97}
void SetHadSystP4(const TLorentzVector &p4)
void SetKV(KineVar_t kv, double value)
void SetFSLeptonP4(const TLorentzVector &p4)
enum genie::EKineVar KineVar_t

References fKV, fP4Fsl, fP4HadSyst, Kinematics(), Reset(), SetFSLeptonP4(), SetHadSystP4(), and SetKV().

Referenced by Kinematics(), and operator=().

◆ FSLeptonP4()

◆ GetKV()

double Kinematics::GetKV ( KineVar_t kv) const

Definition at line 323 of file Kinematics.cxx.

324{
325 if(this->KVSet(kv)) {
326 map<KineVar_t, double>::const_iterator iter = fKV.find(kv);
327 return iter->second;
328 } else {
329 LOG("Interaction", pWARN)
330 << "Kinematic variable: " << KineVar::AsString(kv) << " was not set";
331 }
332 return -99999;
333}
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE,...
Definition Messenger.h:96
#define pWARN
Definition Messenger.h:60
static string AsString(KineVar_t kv)
Definition KineVar.h:77
bool KVSet(KineVar_t kv) const

References genie::KineVar::AsString(), fKV, KVSet(), LOG, and pWARN.

Referenced by genie::SKHadronicSystemGenerator::CalculateHadronicSystem_AtharSingleKaon(), genie::LabFrameHadronTensorI::contraction(), genie::TabulatedLabFrameHadronTensor::dSigma_dT_dCosTheta(), genie::TabulatedLabFrameHadronTensor::dSigma_dT_dCosTheta_rosenbluth(), genie::MECScaleVsW::GetScaling(), genie::utils::kinematics::Jacobian(), genie::SuSAv2MECPXSec::PairRatio(), genie::HENuElGenerator::ProcessEventRecord(), Q2(), q2(), t(), W(), genie::GLRESWdecPythia6::Wdecay(), genie::GLRESWdecPythia8::Wdecay(), genie::PhotonCOHWdecPythia6::Wdecay(), genie::PhotonCOHWdecPythia8::Wdecay(), genie::PhotonRESWdecPythia6::Wdecay(), genie::PhotonRESWdecPythia8::Wdecay(), x(), genie::EmpiricalMECPXSec2015::XSec(), genie::GLRESPXSec::XSec(), genie::HENuElPXSec::XSec(), genie::MKSPPPXSec2020::XSec(), genie::NievesSimoVacasMECPXSec2016::XSec(), genie::PhotonCOHPXSec::XSec(), genie::PhotonRESPXSec::XSec(), genie::SuSAv2MECPXSec::XSec(), genie::SuSAv2QELPXSec::XSec(), and y().

◆ HadSystP4()

◆ Init()

void Kinematics::Init ( void )
private

initialize

Definition at line 59 of file Kinematics.cxx.

60{
61 fKV.clear();
62
63 fP4Fsl = new TLorentzVector;
64 fP4HadSyst = new TLorentzVector;
65}

References fKV, fP4Fsl, and fP4HadSyst.

Referenced by Kinematics(), and Kinematics().

◆ KVSet()

◆ Log10Q2()

double Kinematics::Log10Q2 ( bool selected = false) const

Definition at line 219 of file Kinematics.cxx.

220{
221 double Q2s = this->Q2(selected);
222 return (Q2s>0) ? TMath::Log10(Q2s) : -99999;
223}
double Q2(bool selected=false) const

References Q2().

◆ Log10W()

double Kinematics::Log10W ( bool selected = false) const

Definition at line 225 of file Kinematics.cxx.

226{
227 double Ws = this->W(selected);
228 return (Ws>0) ? TMath::Log10(Ws) : -99999;
229}
double W(bool selected=false) const

References W().

◆ Log10x()

double Kinematics::Log10x ( bool selected = false) const

Definition at line 207 of file Kinematics.cxx.

208{
209 double xs = this->x(selected);
210 return (xs>0) ? TMath::Log10(xs) : -99999;
211}
double x(bool selected=false) const

References x().

◆ Log10y()

double Kinematics::Log10y ( bool selected = false) const

Definition at line 213 of file Kinematics.cxx.

214{
215 double ys = this->y(selected);
216 return (ys>0) ? TMath::Log10(ys) : -99999;
217}
double y(bool selected=false) const

References y().

◆ LogQ2()

double Kinematics::LogQ2 ( bool selected = false) const

Definition at line 195 of file Kinematics.cxx.

196{
197 double Q2s = this->Q2(selected);
198 return (Q2s>0) ? TMath::Log(Q2s) : -99999;
199}

References Q2().

◆ LogW()

double Kinematics::LogW ( bool selected = false) const

Definition at line 201 of file Kinematics.cxx.

202{
203 double Ws = this->W(selected);
204 return (Ws>0) ? TMath::Log(Ws) : -99999;
205}

References W().

◆ Logx()

double Kinematics::Logx ( bool selected = false) const

Definition at line 183 of file Kinematics.cxx.

184{
185 double xs = this->x(selected);
186 return (xs>0) ? TMath::Log(xs) : -99999;
187}

References x().

◆ Logy()

double Kinematics::Logy ( bool selected = false) const

Definition at line 189 of file Kinematics.cxx.

190{
191 double ys = this->y(selected);
192 return (ys>0) ? TMath::Log(ys) : -99999;
193}

References y().

◆ operator=()

Kinematics & Kinematics::operator= ( const Kinematics & kine)

Definition at line 391 of file Kinematics.cxx.

392{
393 this->Copy(kinematics);
394 return (*this);
395}

References Copy(), and Kinematics().

◆ Print()

void Kinematics::Print ( ostream & stream) const

Definition at line 378 of file Kinematics.cxx.

379{
380 stream << "[-] [Kinematics]" << endl;
381
382 map<KineVar_t, double>::const_iterator iter;
383
384 for(iter = fKV.begin(); iter != fKV.end(); ++iter) {
385 KineVar_t kv = iter->first;
386 double val = iter->second;
387 stream << " |--> " << KineVar::AsString(kv) << " = " << val << endl;
388 }
389}

References genie::KineVar::AsString(), and fKV.

Referenced by ClassImp().

◆ Q2()

double Kinematics::Q2 ( bool selected = false) const

Definition at line 125 of file Kinematics.cxx.

126{
127// returns the running or selected value of momentum transfer Q2 (>0)
128
129 if(selected) {
130 if (this->KVSet(kKVSelQ2) ) { return this->GetKV(kKVSelQ2); }
131 else if (this->KVSet(kKVSelq2) ) { return -1* this->GetKV(kKVSelq2); }
132 } else {
133 if (this->KVSet(kKVQ2) ) { return this->GetKV(kKVQ2); }
134 else if (this->KVSet(kKVq2) ) { return -1* this->GetKV(kKVq2); }
135 }
136
137 LOG("Interaction", pWARN) << "Kinematic variable Q2 was not set";
138 return -99999;
139}
double GetKV(KineVar_t kv) const
@ kKVSelq2
Definition KineVar.h:44
@ kKVSelQ2
Definition KineVar.h:43

References GetKV(), genie::kKVQ2, genie::kKVq2, genie::kKVSelQ2, genie::kKVSelq2, KVSet(), LOG, and pWARN.

Referenced by genie::MECGenerator::AddFinalStateLepton(), genie::HEDISGenerator::AddPrimaryLepton(), genie::COHHadronicSystemGenerator::CalculateHadronicSystem_BergerSehgalFM(), genie::DFRKinematicsGenerator::ComputeMaxXSec(), ConvertToGST(), genie::KovalenkoQELCharmPXSec::DR(), GenerateEvent(), genie::TransverseEnhancementFFModel::GetTransEnhMagFF(), genie::KPhaseSpace::IsAllowed(), genie::utils::kinematics::Jacobian(), Log10Q2(), LogQ2(), genie::utils::nuclear::NuclQELXSecSuppression(), genie::DISKinematicsGenerator::ProcessEventRecord(), genie::DMDISKinematicsGenerator::ProcessEventRecord(), genie::DMELEventGenerator::ProcessEventRecord(), genie::OutgoingDarkGenerator::ProcessEventRecord(), genie::PrimaryLeptonGenerator::ProcessEventRecord(), genie::QELEventGenerator::ProcessEventRecord(), genie::QPMDISStrucFuncBase::Q2(), genie::QPMDMDISStrucFuncBase::Q2(), genie::KPhaseSpace::TLim(), genie::utils::kinematics::UpdateWYFromXQ2(), genie::utils::kinematics::UpdateXFromQ2Y(), genie::utils::kinematics::UpdateXYFromWQ2(), genie::MKSPPPXSec2020::ValidKinematics(), genie::EmpiricalMECPXSec2015::XSec(), genie::HEDISPXSec::XSec(), genie::KovalenkoQELCharmPXSec::XSec(), genie::MKSPPPXSec2020::XSec(), genie::P33PaschosLalakulichPXSec::XSec(), and genie::KPhaseSpace::YLim().

◆ q2()

double Kinematics::q2 ( bool selected = false) const

Definition at line 141 of file Kinematics.cxx.

142{
143// returns the running or selected value of momentum transfer q2 (<0)
144
145 if(selected) {
146 if (this->KVSet(kKVSelQ2) ) { return -1* this->GetKV(kKVSelQ2); }
147 else if (this->KVSet(kKVSelq2) ) { return this->GetKV(kKVSelq2); }
148 } else {
149 if (this->KVSet(kKVQ2) ) { return -1* this->GetKV(kKVQ2); }
150 else if (this->KVSet(kKVq2) ) { return this->GetKV(kKVq2); }
151 }
152
153 LOG("Interaction", pWARN) << "Kinematic variable q2 was not set";
154 return -99999;
155}

References GetKV(), genie::kKVQ2, genie::kKVq2, genie::kKVSelQ2, genie::kKVSelq2, KVSet(), LOG, and pWARN.

Referenced by CalculateFormFactor(), genie::DipoleAxialFormFactorModel::FA(), genie::MArunAxialFormFactorModel::FA(), genie::ZExpAxialFormFactorModel::FA(), genie::LwlynSmithFF::Fp(), genie::LwlynSmithFFNC::Fp(), genie::BBA03ELFormFactorsModel::Gen(), genie::BBA07ELFormFactorsModel::Gen(), genie::GalsterELFormFactorsModel::Gen(), genie::ZExpELFormFactorModel::Gen(), genie::BBA03ELFormFactorsModel::Gep(), genie::BBA07ELFormFactorsModel::Gep(), genie::DipoleELFormFactorsModel::Gep(), genie::GalsterELFormFactorsModel::Gep(), genie::ZExpELFormFactorModel::Gep(), genie::BBA03ELFormFactorsModel::Gmn(), genie::BBA07ELFormFactorsModel::Gmn(), genie::DipoleELFormFactorsModel::Gmn(), genie::ZExpELFormFactorModel::Gmn(), genie::BBA03ELFormFactorsModel::Gmp(), genie::BBA07ELFormFactorsModel::Gmp(), genie::DipoleELFormFactorsModel::Gmp(), genie::ZExpELFormFactorModel::Gmp(), genie::utils::nuclear::NuclQELXSecSuppression(), genie::BBA05ELFormFactorsModel::tau(), genie::LwlynSmithFF::tau(), genie::MKFFCC::tau(), genie::MKFFEM::tau(), genie::BSKLNBaseRESPXSec2014::XSec(), genie::ReinSehgalRESPXSec::XSec(), and genie::StrumiaVissaniIBDPXSec::XSec().

◆ Reset()

void Kinematics::Reset ( void )

Definition at line 75 of file Kinematics.cxx.

76{
77 fKV.clear();
78
79 this->SetFSLeptonP4 (0,0,0,0);
80 this->SetHadSystP4 (0,0,0,0);
81}

References fKV, SetFSLeptonP4(), and SetHadSystP4().

Referenced by Copy().

◆ SetFSLeptonP4() [1/2]

◆ SetFSLeptonP4() [2/2]

void Kinematics::SetFSLeptonP4 ( double px,
double py,
double pz,
double E )

Definition at line 302 of file Kinematics.cxx.

303{
304 fP4Fsl->SetPxPyPzE(px,py,pz,E);
305}

References fP4Fsl.

◆ SetHadSystP4() [1/2]

◆ SetHadSystP4() [2/2]

void Kinematics::SetHadSystP4 ( double px,
double py,
double pz,
double E )

Definition at line 312 of file Kinematics.cxx.

313{
314 fP4HadSyst->SetPxPyPzE(px,py,pz,E);
315}

References fP4HadSyst.

◆ SetKV()

void Kinematics::SetKV ( KineVar_t kv,
double value )

Definition at line 335 of file Kinematics.cxx.

336{
337 LOG("Interaction", pDEBUG)
338 << "Setting " << KineVar::AsString(kv) << " to " << value;
339
340 if(this->KVSet(kv)) {
341 fKV[kv] = value;
342 } else {
343 fKV.insert( map<KineVar_t, double>::value_type(kv,value) );
344 }
345}
#define pDEBUG
Definition Messenger.h:63

References genie::KineVar::AsString(), fKV, KVSet(), LOG, and pDEBUG.

Referenced by genie::SKKinematicsGenerator::CalculateKin_AtharSingleKaon(), genie::HELeptonKinematicsGenerator::ComputeMaxXSec(), genie::SKKinematicsGenerator::ComputeMaxXSec(), Copy(), genie::utils::mec::GetMaxXSecTlctl(), genie::HELeptonKinematicsGenerator::ProcessEventRecord(), SetQ2(), Setq2(), Sett(), SetW(), Setx(), Sety(), and genie::NievesSimoVacasMECPXSec2016::XSec().

◆ SetQ2()

void Kinematics::SetQ2 ( double Q2,
bool selected = false )

Definition at line 255 of file Kinematics.cxx.

256{
257// sets the running or selected value of momentum transfer Q2 (>0)
258
259 if(Qsqrd<0) {
260 LOG("Interaction", pWARN)
261 << "Setting unphysical value for Q2 (Q2 = " << Qsqrd << ")";
262 }
263 KineVar_t kvar = (selected) ? kKVSelQ2 : kKVQ2;
264 this->SetKV(kvar, Qsqrd);
265}

References genie::kKVQ2, genie::kKVSelQ2, LOG, pWARN, and SetKV().

Referenced by BuildStdNtuple(), CalculateFormFactor(), genie::COHKinematicsGenerator::CalculateKin_AlvarezRuso(), genie::SKKinematicsGenerator::CalculateKin_AtharSingleKaon(), genie::COHKinematicsGenerator::CalculateKin_BergerSehgal(), genie::COHKinematicsGenerator::CalculateKin_BergerSehgalFM(), genie::COHKinematicsGenerator::CalculateKin_ReinSehgal(), genie::NievesQELCCPXSec::CompareNievesTensors(), genie::utils::ComputeFullDMELPXSec(), genie::utils::ComputeFullQELPXSec(), genie::DMELKinematicsGenerator::ComputeMaxXSec(), genie::IBDKinematicsGenerator::ComputeMaxXSec(), genie::QELKinematicsGenerator::ComputeMaxXSec(), genie::RESKinematicsGenerator::ComputeMaxXSec(), genie::evtlib::EventLibraryInterface::FillKinematics(), genie::LwlynSmithQELCCPXSec::FullDifferentialXSec(), genie::CEvNSEventGenerator::GenerateKinematics(), genie::COHDNuEventGenerator::GenerateKinematics(), GetCrossSection(), genie::ZExpELFormFactorModel::LoadConfig(), main(), main(), MakePlots(), genie::COHKinematicsGenerator::MaxXSec_BergerSehgal(), genie::COHKinematicsGenerator::MaxXSec_BergerSehgalFM(), genie::DFRKinematicsGenerator::ProcessEventRecord(), genie::DISKinematicsGenerator::ProcessEventRecord(), genie::DMDISKinematicsGenerator::ProcessEventRecord(), genie::DMELEventGenerator::ProcessEventRecord(), genie::DMELKinematicsGenerator::ProcessEventRecord(), genie::HEDISKinematicsGenerator::ProcessEventRecord(), genie::IBDKinematicsGenerator::ProcessEventRecord(), genie::QELEventGenerator::ProcessEventRecord(), genie::QELEventGeneratorSM::ProcessEventRecord(), genie::QELKinematicsGenerator::ProcessEventRecord(), genie::RESKinematicsGenerator::ProcessEventRecord(), genie::SPPEventGenerator::ProcessEventRecord(), genie::HEDISKinematicsGenerator::Scan(), genie::MECGenerator::SelectEmpiricalKinematics(), genie::QELEventGeneratorSuSA::SelectLeptonKinematics(), genie::MECGenerator::SelectNSVLeptonKinematics(), genie::MECGenerator::SelectSuSALeptonKinematics(), genie::DMELKinematicsGenerator::SpectralFuncExperimentalCode(), genie::QELKinematicsGenerator::SpectralFuncExperimentalCode(), genie::hnl::Decayer::UpdateEventRecord(), genie::utils::kinematics::UpdateWQ2FromXY(), UseSelectedKinematics(), VerticalSlice(), genie::EmpiricalMECPXSec2015::XSec(), and genie::NievesQELCCPXSec::XSec().

◆ Setq2()

void Kinematics::Setq2 ( double q2,
bool selected = false )

Definition at line 267 of file Kinematics.cxx.

268{
269// sets the running or selected value of momentum transfer q2 (<0)
270
271 if(qsqrd>0) {
272 LOG("Interaction", pWARN)
273 << "Setting unphysical value for q2 (q2 = " << qsqrd << ")";
274 }
275 KineVar_t kvar = (selected) ? kKVSelq2 : kKVq2;
276 this->SetKV(kvar, qsqrd);
277}

References genie::kKVq2, genie::kKVSelq2, LOG, pWARN, and SetKV().

Referenced by UseSelectedKinematics().

◆ Sett()

◆ SetW()

void Kinematics::SetW ( double W,
bool selected = false )

Definition at line 279 of file Kinematics.cxx.

280{
281// sets the running or selected value of invariant hadronic mass W
282
283 if(hadr_mass_W<0) {
284 LOG("Interaction", pWARN)
285 << "Setting unphysical value for W (W = " << hadr_mass_W << ")";
286 }
287 KineVar_t kvar = (selected) ? kKVSelW : kKVW;
288 this->SetKV(kvar, hadr_mass_W);
289}
@ kKVSelW
Definition KineVar.h:45

References genie::kKVSelW, genie::kKVW, LOG, pWARN, and SetKV().

Referenced by genie::COHKinematicsGenerator::CalculateKin_AlvarezRuso(), genie::SKKinematicsGenerator::CalculateKin_AtharSingleKaon(), genie::COHKinematicsGenerator::CalculateKin_BergerSehgal(), genie::COHKinematicsGenerator::CalculateKin_BergerSehgalFM(), genie::COHKinematicsGenerator::CalculateKin_ReinSehgal(), genie::RESKinematicsGenerator::ComputeMaxXSec(), genie::KNOTunedQPMDISPXSec::DISRESJoinSuppressionFactor(), genie::evtlib::EventLibraryInterface::FillKinematics(), GetCrossSection(), genie::LeptoHadPythia6::Hadronize(), genie::LeptoHadPythia8::Hadronize(), main(), main(), genie::utils::kinematics::PhaseSpaceVolume(), genie::DFRKinematicsGenerator::ProcessEventRecord(), genie::DISHadronicSystemGenerator::ProcessEventRecord(), genie::DISKinematicsGenerator::ProcessEventRecord(), genie::DMDISKinematicsGenerator::ProcessEventRecord(), genie::DMELEventGenerator::ProcessEventRecord(), genie::DMELKinematicsGenerator::ProcessEventRecord(), genie::HEDISKinematicsGenerator::ProcessEventRecord(), genie::IBDKinematicsGenerator::ProcessEventRecord(), genie::QELEventGenerator::ProcessEventRecord(), genie::QELEventGeneratorSM::ProcessEventRecord(), genie::QELKinematicsGenerator::ProcessEventRecord(), genie::RESKinematicsGenerator::ProcessEventRecord(), genie::SPPEventGenerator::ProcessEventRecord(), genie::MECGenerator::SelectEmpiricalKinematics(), genie::QELEventGeneratorSuSA::SelectLeptonKinematics(), genie::MECGenerator::SelectNSVLeptonKinematics(), genie::MECGenerator::SelectSuSALeptonKinematics(), genie::DMELKinematicsGenerator::SpectralFuncExperimentalCode(), genie::QELKinematicsGenerator::SpectralFuncExperimentalCode(), genie::hnl::Decayer::UpdateEventRecord(), genie::utils::kinematics::UpdateWQ2FromXY(), genie::utils::kinematics::UpdateWYFromXQ2(), UseSelectedKinematics(), and genie::EmpiricalMECPXSec2015::XSec().

◆ Setx()

void Kinematics::Setx ( double x,
bool selected = false )

Definition at line 231 of file Kinematics.cxx.

232{
233// sets the running or selected value of Bjorken scaling variable x
234
235 if(xbj<0 || xbj>1) {
236 LOG("Interaction", pWARN)
237 << "Setting unphysical value for x (x = " << xbj << ")";
238 }
239 KineVar_t kvar = (selected) ? kKVSelx : kKVx;
240 this->SetKV(kvar, xbj);
241}
@ kKVSelx
Definition KineVar.h:41

References genie::kKVSelx, genie::kKVx, LOG, pWARN, and SetKV().

Referenced by BuildStdNtuple(), genie::COHKinematicsGenerator::CalculateKin_AlvarezRuso(), genie::SKKinematicsGenerator::CalculateKin_AtharSingleKaon(), genie::COHKinematicsGenerator::CalculateKin_BergerSehgal(), genie::COHKinematicsGenerator::CalculateKin_ReinSehgal(), genie::DFRKinematicsGenerator::ComputeMaxXSec(), genie::DISKinematicsGenerator::ComputeMaxXSec(), genie::DMDISKinematicsGenerator::ComputeMaxXSec(), genie::evtlib::EventLibraryInterface::FillKinematics(), GetCrossSection(), main(), main(), MakePlots(), genie::COHKinematicsGenerator::MaxXSec_ReinSehgal(), genie::DFRKinematicsGenerator::ProcessEventRecord(), genie::DISKinematicsGenerator::ProcessEventRecord(), genie::DMDISKinematicsGenerator::ProcessEventRecord(), genie::DMELEventGenerator::ProcessEventRecord(), genie::DMELKinematicsGenerator::ProcessEventRecord(), genie::HEDISKinematicsGenerator::ProcessEventRecord(), genie::IBDKinematicsGenerator::ProcessEventRecord(), genie::QELEventGenerator::ProcessEventRecord(), genie::QELEventGeneratorSM::ProcessEventRecord(), genie::QELKinematicsGenerator::ProcessEventRecord(), genie::RESKinematicsGenerator::ProcessEventRecord(), genie::SPPEventGenerator::ProcessEventRecord(), genie::HEDISKinematicsGenerator::Scan(), genie::MECGenerator::SelectEmpiricalKinematics(), genie::QELEventGeneratorSuSA::SelectLeptonKinematics(), genie::MECGenerator::SelectNSVLeptonKinematics(), genie::MECGenerator::SelectSuSALeptonKinematics(), genie::DMELKinematicsGenerator::SpectralFuncExperimentalCode(), genie::QELKinematicsGenerator::SpectralFuncExperimentalCode(), genie::utils::kinematics::UpdateXFromQ2Y(), genie::utils::kinematics::UpdateXYFromWQ2(), UseSelectedKinematics(), and VerticalSlice().

◆ Sety()

void Kinematics::Sety ( double y,
bool selected = false )

Definition at line 243 of file Kinematics.cxx.

244{
245// sets the running or selected value of inelasticity y
246
247 if(inel_y<0 || inel_y>1) {
248 LOG("Interaction", pWARN)
249 << "Setting unphysical value for y (y = " << inel_y << ")";
250 }
251 KineVar_t kvar = (selected) ? kKVSely : kKVy;
252 this->SetKV(kvar, inel_y);
253}
@ kKVSely
Definition KineVar.h:42

References genie::kKVSely, genie::kKVy, LOG, pWARN, and SetKV().

Referenced by genie::COHKinematicsGenerator::CalculateKin_AlvarezRuso(), genie::SKKinematicsGenerator::CalculateKin_AtharSingleKaon(), genie::COHKinematicsGenerator::CalculateKin_BergerSehgal(), genie::COHKinematicsGenerator::CalculateKin_BergerSehgalFM(), genie::COHKinematicsGenerator::CalculateKin_ReinSehgal(), genie::DFRKinematicsGenerator::ComputeMaxXSec(), genie::DISKinematicsGenerator::ComputeMaxXSec(), genie::DMDISKinematicsGenerator::ComputeMaxXSec(), genie::DMEKinematicsGenerator::ComputeMaxXSec(), genie::NuEKinematicsGenerator::ComputeMaxXSec(), genie::evtlib::EventLibraryInterface::FillKinematics(), GetCrossSection(), main(), main(), genie::COHKinematicsGenerator::MaxXSec_BergerSehgal(), genie::COHKinematicsGenerator::MaxXSec_BergerSehgalFM(), genie::COHKinematicsGenerator::MaxXSec_ReinSehgal(), genie::DFRKinematicsGenerator::ProcessEventRecord(), genie::DISKinematicsGenerator::ProcessEventRecord(), genie::DMDISKinematicsGenerator::ProcessEventRecord(), genie::DMEKinematicsGenerator::ProcessEventRecord(), genie::DMELEventGenerator::ProcessEventRecord(), genie::DMELKinematicsGenerator::ProcessEventRecord(), genie::HEDISKinematicsGenerator::ProcessEventRecord(), genie::IBDKinematicsGenerator::ProcessEventRecord(), genie::NuEKinematicsGenerator::ProcessEventRecord(), genie::QELEventGenerator::ProcessEventRecord(), genie::QELEventGeneratorSM::ProcessEventRecord(), genie::QELKinematicsGenerator::ProcessEventRecord(), genie::RESKinematicsGenerator::ProcessEventRecord(), genie::SPPEventGenerator::ProcessEventRecord(), genie::MECGenerator::SelectEmpiricalKinematics(), genie::QELEventGeneratorSuSA::SelectLeptonKinematics(), genie::MECGenerator::SelectNSVLeptonKinematics(), genie::MECGenerator::SelectSuSALeptonKinematics(), genie::DMELKinematicsGenerator::SpectralFuncExperimentalCode(), genie::QELKinematicsGenerator::SpectralFuncExperimentalCode(), genie::utils::kinematics::UpdateWYFromXQ2(), genie::utils::kinematics::UpdateXYFromWQ2(), and UseSelectedKinematics().

◆ t()

double Kinematics::t ( bool selected = false) const

Definition at line 170 of file Kinematics.cxx.

171{
172// returns the running or selected value of invariant hadronic mass W
173
174 KineVar_t kvar = (selected) ? kKVSelt : kKVt;
175
176 if(this->KVSet(kvar)) { return this->GetKV(kvar); }
177 else {
178 LOG("Interaction", pWARN) << "Kinematic variable t was not set";
179 }
180 return -99999;
181}

References GetKV(), genie::kKVSelt, genie::kKVt, KVSet(), LOG, and pWARN.

Referenced by genie::COHHadronicSystemGenerator::CalculateHadronicSystem_BergerSehgalFM(), genie::COHHadronicSystemGenerator::CalculateHadronicSystem_ReinSehgal(), ConvertToGST(), genie::KPhaseSpace::IsAllowed(), and genie::DFRHadronicSystemGenerator::ProcessEventRecord().

◆ UseSelectedKinematics()

void Kinematics::UseSelectedKinematics ( void )

Definition at line 359 of file Kinematics.cxx.

360{
361// copy the selected kinematics into the running ones
362//
363 map<KineVar_t, double>::const_iterator iter;
364 iter = fKV.find(kKVSelx);
365 if(iter != fKV.end()) this->Setx(iter->second);
366 iter = fKV.find(kKVSely);
367 if(iter != fKV.end()) this->Sety(iter->second);
368 iter = fKV.find(kKVSelQ2);
369 if(iter != fKV.end()) this->SetQ2(iter->second);
370 iter = fKV.find(kKVSelq2);
371 if(iter != fKV.end()) this->Setq2(iter->second);
372 iter = fKV.find(kKVSelW);
373 if(iter != fKV.end()) this->SetW(iter->second);
374 iter = fKV.find(kKVSelt);
375 if(iter != fKV.end()) this->Sett(iter->second);
376}
void Setx(double x, bool selected=false)
void SetQ2(double Q2, bool selected=false)
void Sett(double t, bool selected=false)
void Setq2(double q2, bool selected=false)
void Sety(double y, bool selected=false)
void SetW(double W, bool selected=false)

References fKV, genie::kKVSelQ2, genie::kKVSelq2, genie::kKVSelt, genie::kKVSelW, genie::kKVSelx, genie::kKVSely, SetQ2(), Setq2(), Sett(), SetW(), Setx(), and Sety().

Referenced by main(), and genie::RSPPResonanceSelector::SelectResonance().

◆ W()

double Kinematics::W ( bool selected = false) const

◆ x()

double Kinematics::x ( bool selected = false) const

Definition at line 99 of file Kinematics.cxx.

100{
101// returns the running or selected value of Bjorken scaling variable x
102
103 KineVar_t kvar = (selected) ? kKVSelx : kKVx;
104
105 if(this->KVSet(kvar)) { return this->GetKV(kvar); }
106 else {
107 LOG("Interaction", pWARN) << "Kinematic variable x was not set";
108 }
109 return -99999;
110}

References GetKV(), genie::kKVSelx, genie::kKVx, KVSet(), LOG, and pWARN.

Referenced by genie::QPMDISStrucFuncBase::Calculate(), genie::QPMDMDISStrucFuncBase::Calculate(), genie::COHHadronicSystemGenerator::CalculateHadronicSystem_ReinSehgal(), genie::COHKinematicsGenerator::CalculateKin_BergerSehgal(), ConvertToGHad(), ConvertToGST(), genie::KNOTunedQPMDISPXSec::DISRESJoinSuppressionFactor(), GenerateEvent(), genie::KPhaseSpace::IsAllowed(), genie::utils::kinematics::Jacobian(), Log10x(), Logx(), genie::QPMDISStrucFuncBase::NuclMod(), genie::QPMDMDISStrucFuncBase::NuclMod(), genie::DFRHadronicSystemGenerator::ProcessEventRecord(), genie::QPMDISStrucFuncBase::Q2(), genie::QPMDMDISStrucFuncBase::Q2(), genie::QPMDISStrucFuncBase::R(), genie::QPMDMDISStrucFuncBase::R(), genie::BYStrucFunc::ScalingVar(), genie::DMBYStrucFunc::ScalingVar(), genie::QPMDISStrucFuncBase::ScalingVar(), genie::QPMDMDISStrucFuncBase::ScalingVar(), genie::utils::kinematics::UpdateWQ2FromXY(), genie::utils::kinematics::UpdateWYFromXQ2(), genie::AivazisCharmPXSecLO::XSec(), genie::HEDISPXSec::XSec(), genie::QPMDISPXSec::XSec(), genie::QPMDMDISPXSec::XSec(), and genie::SlowRsclCharmDISPXSecLO::XSec().

◆ y()

double Kinematics::y ( bool selected = false) const

Definition at line 112 of file Kinematics.cxx.

113{
114// returns the running or selected value of inelasticity y
115
116 KineVar_t kvar = (selected) ? kKVSely : kKVy;
117
118 if(this->KVSet(kvar)) { return this->GetKV(kvar); }
119 else {
120 LOG("Interaction", pWARN) << "Kinematic variable y was not set";
121 }
122 return -99999;
123}

References GetKV(), genie::kKVSely, genie::kKVy, KVSet(), LOG, and pWARN.

Referenced by genie::MECGenerator::AddFinalStateLepton(), genie::HEDISGenerator::AddPrimaryLepton(), genie::COHHadronicSystemGenerator::CalculateHadronicSystem_BergerSehgalFM(), genie::COHHadronicSystemGenerator::CalculateHadronicSystem_ReinSehgal(), ConvertToGHad(), ConvertToGST(), genie::KNOTunedQPMDISPXSec::DISRESJoinSuppressionFactor(), GenerateEvent(), genie::KPhaseSpace::IsAllowed(), genie::utils::kinematics::Jacobian(), Log10y(), Logy(), genie::DFRHadronicSystemGenerator::ProcessEventRecord(), genie::DMEOutgoingDarkGenerator::ProcessEventRecord(), genie::HEDISKinematicsGenerator::ProcessEventRecord(), genie::NuEPrimaryLeptonGenerator::ProcessEventRecord(), genie::OutgoingDarkGenerator::ProcessEventRecord(), genie::PrimaryLeptonGenerator::ProcessEventRecord(), genie::QPMDISStrucFuncBase::Q2(), genie::QPMDMDISStrucFuncBase::Q2(), genie::KPhaseSpace::TLim(), genie::utils::kinematics::UpdateWQ2FromXY(), genie::utils::kinematics::UpdateXFromQ2Y(), genie::AivazisCharmPXSecLO::XSec(), genie::BardinIMDRadCorPXSec::XSec(), genie::DMElectronPXSec::XSec(), genie::HEDISPXSec::XSec(), genie::IMDAnnihilationPXSec::XSec(), genie::NuElectronPXSec::XSec(), genie::QPMDISPXSec::XSec(), genie::QPMDMDISPXSec::XSec(), and genie::SlowRsclCharmDISPXSecLO::XSec().

◆ operator<<

ostream & operator<< ( ostream & stream,
const Kinematics & kine )
friend

References Kinematics().

Member Data Documentation

◆ fKV

map<KineVar_t, double> genie::Kinematics::fKV
private

selected kinematics

Definition at line 102 of file Kinematics.h.

Referenced by CleanUp(), ClearRunningValues(), Copy(), GetKV(), Init(), KVSet(), Print(), Reset(), SetKV(), and UseSelectedKinematics().

◆ fP4Fsl

TLorentzVector* genie::Kinematics::fP4Fsl
private

generated final state primary lepton 4-p (LAB)

Definition at line 103 of file Kinematics.h.

Referenced by CleanUp(), Copy(), FSLeptonP4(), Init(), Kinematics(), SetFSLeptonP4(), and SetFSLeptonP4().

◆ fP4HadSyst

TLorentzVector* genie::Kinematics::fP4HadSyst
private

generated final state hadronic system 4-p (LAB)

Definition at line 104 of file Kinematics.h.

Referenced by CleanUp(), Copy(), HadSystP4(), Init(), Kinematics(), SetHadSystP4(), and SetHadSystP4().


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