GENIEGenerator
Loading...
Searching...
No Matches
KPhaseSpace.cxx
Go to the documentation of this file.
1//____________________________________________________________________________
2/*
3 Copyright (c) 2003-2025, The GENIE Collaboration
4 For the full text of the license visit http://copyright.genie-mc.org
5
6 Costas Andreopoulos <c.andreopoulos \at cern.ch>
7 University of Liverpool
8
9 Changes required to implement the GENIE Boosted Dark Matter module
10 were installed by Josh Berger (Univ. of Wisconsin)
11*/
12//____________________________________________________________________________
13
14#include <cmath>
15#include <cstdlib>
16
17#include <TMath.h>
18
20
34
35using namespace genie;
36using namespace genie::utils;
37using namespace genie::constants;
38
40
41//____________________________________________________________________________
43TObject(), fInteraction(NULL)
44{
45 this->UseInteraction(0);
46}
47//___________________________________________________________________________
49TObject(), fInteraction(NULL)
50{
51 this->UseInteraction(in);
52}
53//___________________________________________________________________________
58//___________________________________________________________________________
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}
76//___________________________________________________________________________
78{
79 fInteraction = in;
80}
81//___________________________________________________________________________
82double KPhaseSpace::Threshold(void) const
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}
249//___________________________________________________________________________
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}
271//____________________________________________________________________________
273{
274 Range1D_t lim = this->Limits(kvar);
275 return lim.min;
276}
277//___________________________________________________________________________
279{
280 Range1D_t lim = this->Limits(kvar);
281 return lim.max;
282}
283//___________________________________________________________________________
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() ||
294 pi.IsInverseMuDecay() ||
295 pi.IsIMDAnnihilation() ||
296 pi.IsNuElectronElastic() ||
298 pi.IsMEC() ||
299 pi.IsPhotonCoherent() ||
300 pi.IsPhotonResonance() ||
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.IsSinglePion() ||
315 pi.IsAMNuGamma())
316 {
317 E = init_state.ProbeE(kRfHitNucRest);
318 }
319
320 LOG("KPhaseSpace", pDEBUG) << "E = " << E << ", Ethr = " << Ethr;
321 return (E>Ethr);
322}
323//___________________________________________________________________________
325{
326 const ProcessInfo & pi = fInteraction->ProcInfo();
327 const Kinematics & kine = fInteraction->Kine();
328
329 // ASK single kaon:
330 // XSec code returns zero when kinematics are not allowed
331 // Here just let kinematics always be allowed
332 if(pi.IsSingleKaon()) {
333 return true;
334 }
335
336 // QEL:
337 // Check the running Q2 vs the Q2 limits
338 if(pi.IsQuasiElastic() || pi.IsInverseBetaDecay() || pi.IsDarkMatterElastic()) {
339 Range1D_t Q2l = this->Q2Lim();
340 double Q2 = kine.Q2();
341 bool in_phys = math::IsWithinLimits(Q2, Q2l);
342 bool allowed = in_phys;
343 return allowed;
344 }
345
346 // RES
347 // Check the running W vs the W limits
348 // & the running Q2 vs Q2 limits for the given W
349 if(pi.IsResonant()) {
350 Range1D_t Wl = this->WLim();
351 Range1D_t Q2l = this->Q2Lim_W();
352 double W = kine.W();
353 double Q2 = kine.Q2();
354 bool in_phys = (math::IsWithinLimits(Q2, Q2l) && math::IsWithinLimits(W, Wl));
355 bool allowed = in_phys;
356 return allowed;
357 }
358
359 // DIS
361 Range1D_t Wl = this->WLim();
362 Range1D_t Q2l = this->Q2Lim_W();
363 double W = kine.W();
364 double Q2 = kine.Q2();
365 bool in_phys = (math::IsWithinLimits(Q2, Q2l) && math::IsWithinLimits(W, Wl));
366 bool allowed = in_phys;
367 return allowed;
368 }
369
370 //IMD
372 Range1D_t yl = this->YLim();
373 double y = kine.y();
374 bool in_phys = math::IsWithinLimits(y, yl);
375 bool allowed = in_phys;
376 return allowed;
377 }
378
379 //COH
380 if (pi.IsCoherentProduction()) {
381 Range1D_t xl = this->XLim();
382 Range1D_t yl = this->YLim();
383 double x = kine.x();
384 double y = kine.y();
385 bool in_phys = (math::IsWithinLimits(x, xl) && math::IsWithinLimits(y, yl));
386 bool allowed = in_phys;
387 return allowed;
388 }
389
390 // CEvNS
391 if (pi.IsCoherentElastic()) {
392 double Q2 = kine.Q2();
393 bool allowed (Q2 > 0);
394 return allowed;
395 }
396
397 // DFR
398 if (pi.IsDiffractive()) {
399 // first two checks are the same as RES & DIS
400 Range1D_t Wl = this->WLim();
401 Range1D_t Q2l = this->Q2Lim_W();
402
404 double W = kine.W();
405 double Q2 = kine.Q2();
406
407 LOG("KPhaseSpace", pDEBUG) << " W = " << W << ", limits = [" << Wl.min << "," << Wl.max << "];";
408 LOG("KPhaseSpace", pDEBUG) << " Q2 = " << Q2 << ", limits = [" << Q2l.min << "," << Q2l.max << "];";
409 bool in_phys = math::IsWithinLimits(W, Wl);
410 in_phys = in_phys && math::IsWithinLimits(Q2, Q2l);
411
412 // extra check: there's a t minimum.
413 // but only check if W, Q2 is reasonable
414 // (otherwise get NaNs in tmin)
415 if (in_phys)
416 {
417 double t = kine.t();
418 Range1D_t tl = this->TLim();
419 LOG("KPhaseSpace", pDEBUG) << " t = " << t << ", limits = [" << tl.min << "," << tl.max << "];";
420 in_phys = in_phys && math::IsWithinLimits(t, tl);
421 }
422 LOG("KPhaseSpace", pDEBUG) << " phase space point is " << ( in_phys ? "ALLOWED" : "NOT ALLOWED");
423
424
425 bool allowed = in_phys;
426 return allowed;
427 }
428
429 // was MECTensor
430 if (pi.IsMEC()){
431 Range1D_t Q2l = this->Q2Lim();
432 double Q2 = kine.Q2();
433 bool in_phys = math::IsWithinLimits(Q2, Q2l);
434 bool allowed = in_phys;
435 return allowed;
436 }
437
438 return false;
439}
440//___________________________________________________________________________
442{
443// Computes hadronic invariant mass limits.
444// For QEL the range reduces to the recoil nucleon mass.
445// For DIS & RES the calculation proceeds as in kinematics::InelWLim().
446// It is not computed for other interactions
447//
448 Range1D_t Wl;
449 Wl.min = -1;
450 Wl.max = -1;
451
452 const ProcessInfo & pi = fInteraction->ProcInfo();
453
454 bool is_em = pi.IsEM();
455 bool is_qel = pi.IsQuasiElastic() || pi.IsInverseBetaDecay() || pi.IsDarkMatterElastic();
456 bool is_inel = pi.IsDeepInelastic() || pi.IsResonant() || pi.IsDiffractive();
457 bool is_dmdis = pi.IsDarkMatterDeepInelastic();
458
459 if(is_qel) {
460 double MR = fInteraction->RecoilNucleon()->Mass();
461 Wl.min = MR;
462 Wl.max = MR;
463 return Wl;
464 }
465 if(is_inel) {
466 const InitialState & init_state = fInteraction->InitState();
467 double Ev = init_state.ProbeE(kRfHitNucRest);
468 double M = init_state.Tgt().HitNucP4Ptr()->M(); //can be off m/shell
469 double ml = fInteraction->FSPrimLepton()->Mass();
470
471 Wl = is_em ? kinematics::electromagnetic::InelWLim(Ev,ml,M) : kinematics::InelWLim(Ev,M,ml);
472
473 if(fInteraction->ExclTag().IsCharmEvent()) {
474 //Wl.min = TMath::Max(Wl.min, kNeutronMass+kPionMass+kLightestChmHad);
475 Wl.min = TMath::Max(Wl.min, kNeutronMass+kLightestChmHad);
476 }
477 else if (fInteraction->ProcInfo().IsDiffractive())
478 Wl.min = TMath::Max(Wl.min, kNeutronMass+kPionMass);
479 else if ( fInteraction->ProcInfo().IsDeepInelastic() ) {
480 Wl.min = TMath::Max(Wl.min, kNeutronMass + kPionMass );
481 }
482
483 // sanity check
484 if(Wl.min>Wl.max) {Wl.min=-1; Wl.max=-1;}
485
486 return Wl;
487 }
488 if(is_dmdis) {
489 const InitialState & init_state = fInteraction->InitState();
490 double Ev = init_state.ProbeE(kRfHitNucRest);
491 double M = init_state.Tgt().HitNucP4Ptr()->M(); //can be off m/shell
492 double ml = fInteraction->FSPrimLepton()->Mass();
493 Wl = kinematics::DarkWLim(Ev,M,ml);
494 if(fInteraction->ExclTag().IsCharmEvent()) {
495 //Wl.min = TMath::Max(Wl.min, kNeutronMass+kPionMass+kLightestChmHad);
496 Wl.min = TMath::Max(Wl.min, kNeutronMass+kLightestChmHad);
497 }
498 else if (fInteraction->ProcInfo().IsDiffractive())
499 Wl.min = TMath::Max(Wl.min, kNeutronMass+kPionMass);
500
501 LOG("KPhaseSpace", pDEBUG) << "Found nominal limits: " << Wl.min << ", " << Wl.max;
502
503 // sanity check
504 if(Wl.min>Wl.max) {Wl.min=-1; Wl.max=-1;}
505
506 return Wl;
507 }
508 return Wl;
509}
510//____________________________________________________________________________
512{
513 // Computes momentum transfer (Q2>0) limits @ the input invariant mass
514 // The calculation proceeds as in kinematics::InelQ2Lim_W().
515 // For QEL, W is set to the recoil nucleon mass
516 //
517 // TODO: For now, choosing to handle Q2 at fixed W for coherent in the
518 // same way as for the general Q2 limits... but shouldn't we just use
519 // W = m_pi? - which we do in Q2Lim() anyway... seems like there are
520 // cleanup opportunities here.
521
522 Range1D_t Q2l;
523 Q2l.min = -1;
524 Q2l.max = -1;
525
526 const ProcessInfo & pi = fInteraction->ProcInfo();
527
528 bool is_em = pi.IsEM();
529 bool is_qel = pi.IsQuasiElastic() || pi.IsInverseBetaDecay();
530 bool is_inel = pi.IsDeepInelastic() || pi.IsResonant() || pi.IsDiffractive();
531 bool is_coh = pi.IsCoherentProduction();
532 bool is_dme = pi.IsDarkMatterElastic();
533 bool is_dmdis = pi.IsDarkMatterDeepInelastic();
534
535 if(!is_qel && !is_inel && !is_coh && !is_dme && !is_dmdis) return Q2l;
536
537 if(is_coh) {
538 return Q2Lim();
539 }
540
541 const InitialState & init_state = fInteraction->InitState();
542 double Ev = init_state.ProbeE(kRfHitNucRest);
543 double M = init_state.Tgt().HitNucP4Ptr()->M(); // can be off m/shell
544 double ml = fInteraction->FSPrimLepton()->Mass();
545
546 double W = 0;
547 if(is_qel || is_dme) W = fInteraction->RecoilNucleon()->Mass();
548 else W = kinematics::W(fInteraction);
549
550 if (pi.IsInverseBetaDecay()) {
552 } else if (is_dme || is_dmdis) {
553 Q2l = kinematics::DarkQ2Lim_W(Ev,M,ml,W);
554 } else {
555 Q2l = is_em ? kinematics::electromagnetic::InelQ2Lim_W(Ev,ml,M,W) : kinematics::InelQ2Lim_W(Ev,M,ml,W);
556 }
557
558 return Q2l;
559}
560//____________________________________________________________________________
562{
563// As Q2Lim_W(void) but with reversed sign (Q2 -> q2)
564//
565 Range1D_t Q2 = this->Q2Lim_W();
566 Range1D_t q2;
567 q2.min = - Q2.max;
568 q2.max = - Q2.min;
569 return q2;
570}
571//____________________________________________________________________________
573{
574 // Computes momentum transfer (Q2>0) limits irrespective of the invariant mass
575 // For QEL this is identical to Q2Lim_W (since W is fixed)
576 // For RES & DIS, the calculation proceeds as in kinematics::InelQ2Lim().
577 //
578 Range1D_t Q2l;
579 Q2l.min = -1;
580 Q2l.max = -1;
581
582 const ProcessInfo & pi = fInteraction->ProcInfo();
583
584 bool is_em = pi.IsEM();
585 bool is_qel = pi.IsQuasiElastic() || pi.IsInverseBetaDecay();
586 bool is_inel = pi.IsDeepInelastic() || pi.IsResonant();
587 bool is_coh = pi.IsCoherentProduction();
588 bool is_cevns = pi.IsCoherentElastic();
589 bool is_dme = pi.IsDarkMatterElastic();
590 bool is_dmdis = pi.IsDarkMatterDeepInelastic();
591
592 if(!is_qel && !is_inel && !is_coh && !is_cevns && !is_dme && !is_dmdis) return Q2l;
593
594 const InitialState & init_state = fInteraction->InitState();
595 double Ev = init_state.ProbeE(kRfHitNucRest);
596 double M = init_state.Tgt().HitNucP4Ptr()->M(); // can be off m/shell
597 double ml = fInteraction->FSPrimLepton()->Mass();
598
599 if(is_cevns) {
600 double Ev_lab = init_state.ProbeE(kRfLab);
601 Q2l = kinematics::CEvNSQ2Lim(Ev_lab);
602 return Q2l;
603 }
604
605 const XclsTag & xcls = fInteraction->ExclTag();
606
607 if(is_coh) {
608
609 double m_other = controls::kASmallNum ;
610 // as a default the mass of hadronic system is the mass of the photon.
611 // which is assumed to be a small number to avoid divergences
612
613 if ( xcls.NPions() > 0 ) {
614 bool pionIsCharged = pi.IsWeakCC();
615 m_other = pionIsCharged ? kPionMass : kPi0Mass;
616 }
617
618 Q2l = kinematics::CohQ2Lim(M, m_other, ml, Ev);
619 return Q2l;
620 }
621
622 // quasi-elastic
623 if(is_qel) {
624 double W = fInteraction->RecoilNucleon()->Mass();
625 if(xcls.IsCharmEvent()) {
626 int charm_pdgc = xcls.CharmHadronPdg();
627 W = PDGLibrary::Instance()->Find(charm_pdgc)->Mass();
628 } else if(xcls.IsStrangeEvent()) {
629 int strange_pdgc = xcls.StrangeHadronPdg();
630 W = PDGLibrary::Instance()->Find(strange_pdgc)->Mass();
631 }
632 if (pi.IsInverseBetaDecay()) {
634 } else {
635 Q2l = is_em ? kinematics::electromagnetic::InelQ2Lim_W(Ev,ml,M,W) : kinematics::InelQ2Lim_W(Ev,M,ml,W);
636 }
637
638 return Q2l;
639 }
640
641 // dark mattter elastic
642 if(is_dme) {
643 double W = fInteraction->RecoilNucleon()->Mass();
644 if(xcls.IsCharmEvent()) {
645 int charm_pdgc = xcls.CharmHadronPdg();
646 W = PDGLibrary::Instance()->Find(charm_pdgc)->Mass();
647 } else if(xcls.IsStrangeEvent()) {
648 int strange_pdgc = xcls.StrangeHadronPdg();
649 W = PDGLibrary::Instance()->Find(strange_pdgc)->Mass();
650 }
651 if (pi.IsInverseBetaDecay()) {
653 } else {
654 Q2l = kinematics::DarkQ2Lim_W(Ev,M,ml,W);
655 }
656
657
658 return Q2l;
659 }
660
661 // was MECTensor
662 // TODO: Q2maxConfig
663 if (pi.IsMEC()){
664 double W = fInteraction->RecoilNucleon()->Mass();
665 Q2l = is_em ? kinematics::electromagnetic::InelQ2Lim_W(Ev,ml,M,W) : kinematics::InelQ2Lim_W(Ev,M,ml,W);
666 double Q2maxConfig = 1.44; // need to pull from config file somehow?
667 if (Q2l.max > Q2maxConfig) Q2l.max = Q2maxConfig;
668 return Q2l;
669 }
670
671 if (is_dmdis) {
672 Q2l = kinematics::DarkQ2Lim(Ev,M,ml);
673 return Q2l;
674 }
675
676 // inelastic
677 Q2l = is_em ? kinematics::electromagnetic::InelQ2Lim(Ev,ml,M) : kinematics::InelQ2Lim(Ev,M,ml);
678 return Q2l;
679}
680//____________________________________________________________________________
682{
683// As Q2Lim(void) but with reversed sign (Q2 -> q2)
684//
685 Range1D_t Q2 = this->Q2Lim();
686 Range1D_t q2;
687 q2.min = - Q2.max;
688 q2.max = - Q2.min;
689 return q2;
690}
691//____________________________________________________________________________
693{
694 // Computes x-limits;
695
696 Range1D_t xl;
697 xl.min = -1;
698 xl.max = -1;
699
700 const ProcessInfo & pi = fInteraction->ProcInfo();
701 bool is_em = pi.IsEM();
702
703 //RES+DIS
704 bool is_inel = pi.IsDeepInelastic() || pi.IsResonant();
705 if(is_inel) {
706 const InitialState & init_state = fInteraction->InitState();
707 double Ev = init_state.ProbeE(kRfHitNucRest);
708 double M = init_state.Tgt().HitNucP4Ptr()->M(); // can be off m/shell
709 double ml = fInteraction->FSPrimLepton()->Mass();
710 xl = is_em ? kinematics::electromagnetic::InelXLim(Ev,ml,M) : kinematics::InelXLim(Ev,M,ml);
711 return xl;
712 }
713 //DMDIS
714 bool is_dmdis = pi.IsDarkMatterDeepInelastic();
715 if(is_dmdis) {
716 const InitialState & init_state = fInteraction->InitState();
717 double Ev = init_state.ProbeE(kRfHitNucRest);
718 double M = init_state.Tgt().HitNucP4Ptr()->M(); // can be off m/shell
719 double ml = fInteraction->FSPrimLepton()->Mass();
720 xl = kinematics::DarkXLim(Ev,M,ml);
721 return xl;
722 }
723 //COH
724 bool is_coh = pi.IsCoherentProduction();
725 if(is_coh) {
726 xl = kinematics::CohXLim();
727 return xl;
728 }
729 //QEL
730 bool is_qel = pi.IsQuasiElastic() || pi.IsInverseBetaDecay() || pi.IsDarkMatterElastic();
731 if(is_qel) {
732 xl.min = 1;
733 xl.max = 1;
734 return xl;
735 }
736 bool is_dfr = pi.IsDiffractive();
737 if(is_dfr) {
739 xl.max = 1. - controls::kASmallNum;
740 return xl;
741 }
742
743 return xl;
744}
745//____________________________________________________________________________
747{
748 Range1D_t yl;
749 yl.min = -1;
750 yl.max = -1;
751
752 const ProcessInfo & pi = fInteraction->ProcInfo();
753 bool is_em = pi.IsEM();
754
755 //RES+DIS
756 bool is_inel = pi.IsDeepInelastic() || pi.IsResonant();
757 if(is_inel) {
758 const InitialState & init_state = fInteraction->InitState();
759 double Ev = init_state.ProbeE(kRfHitNucRest);
760 double M = init_state.Tgt().HitNucP4Ptr()->M(); // can be off m/shell
761 double ml = fInteraction->FSPrimLepton()->Mass();
762 yl = is_em ? kinematics::electromagnetic::InelYLim(Ev,ml,M) : kinematics::InelYLim(Ev,M,ml);
763 return yl;
764 }
765 //DMDIS
766 bool is_dmdis = pi.IsDarkMatterDeepInelastic();
767 if(is_dmdis) {
768 const InitialState & init_state = fInteraction->InitState();
769 double Ev = init_state.ProbeE(kRfHitNucRest);
770 double M = init_state.Tgt().HitNucP4Ptr()->M(); // can be off m/shell
771 double ml = fInteraction->FSPrimLepton()->Mass();
772 yl = kinematics::DarkYLim(Ev,M,ml);
773 return yl;
774 }
775 //COH
776 bool is_coh = pi.IsCoherentProduction();
777 if(is_coh) {
778 const InitialState & init_state = fInteraction->InitState();
779 double EvL = init_state.ProbeE(kRfLab);
780 double ml = fInteraction->FSPrimLepton()->Mass();
781 yl = kinematics::CohYLim(EvL,ml);
782 return yl;
783 }
784 // IMD
786 const InitialState & init_state = fInteraction->InitState();
787 double Ev = init_state.ProbeE(kRfLab);
788 double ml = fInteraction->FSPrimLepton()->Mass();
789 double me = kElectronMass;
791 yl.max = 1 - (ml*ml + me*me)/(2*me*Ev) - controls::kASmallNum;
792 return yl;
793 }
794 // EDIT: y limits are different for massive probe
796 const InitialState & init_state = fInteraction->InitState();
797 double Ev = init_state.ProbeE(kRfLab);
798 double ml = fInteraction->FSPrimLepton()->Mass();
799 double me = kElectronMass;
800 yl.min = (Ev*me*me + ml*ml*(Ev + 2.0*me)) / (Ev * (2.0*Ev*me + me*me + ml*ml)) + controls::kASmallNum;
801 yl.max = 1.0 - controls::kASmallNum;
802 return yl;
803 }
804 bool is_dfr = pi.IsDiffractive();
805 if(is_dfr) {
806 const InitialState & init_state = fInteraction -> InitState();
807 double Ev = init_state.ProbeE(kRfHitNucRest);
808 double ml = fInteraction->FSPrimLepton()->Mass();
810 yl.max = 1. -ml/Ev - controls::kASmallNum;
811 return yl;
812 }
813 return yl;
814}
815//____________________________________________________________________________
817{
818// Computes kinematical limits for y @ the input x
819
820 Range1D_t yl;
821 yl.min = -1;
822 yl.max = -1;
823
824 const ProcessInfo & pi = fInteraction->ProcInfo();
825 bool is_em = pi.IsEM();
826
827 //RES+DIS
828 bool is_inel = pi.IsDeepInelastic() || pi.IsResonant();
829 if(is_inel) {
830 const InitialState & init_state = fInteraction->InitState();
831 double Ev = init_state.ProbeE(kRfHitNucRest);
832 double M = init_state.Tgt().HitNucP4Ptr()->M(); // can be off m/shell
833 double ml = fInteraction->FSPrimLepton()->Mass();
834 double x = fInteraction->Kine().x();
835 yl = is_em ? kinematics::electromagnetic::InelYLim_X(Ev,ml,M,x) : kinematics::InelYLim_X(Ev,M,ml,x);
836 return yl;
837 }
838 //DMDIS
839 bool is_dmdis = pi.IsDarkMatterDeepInelastic();
840 if(is_dmdis) {
841 const InitialState & init_state = fInteraction->InitState();
842 double Ev = init_state.ProbeE(kRfHitNucRest);
843 double M = init_state.Tgt().HitNucP4Ptr()->M(); // can be off m/shell
844 double ml = fInteraction->FSPrimLepton()->Mass();
845 double x = fInteraction->Kine().x();
846 yl = kinematics::DarkYLim_X(Ev,M,ml,x);
847 return yl;
848 }
849 //COH
850 bool is_coh = pi.IsCoherentProduction();
851 if(is_coh) {
852 const InitialState & init_state = fInteraction->InitState();
853 double EvL = init_state.ProbeE(kRfLab);
854 double ml = fInteraction->FSPrimLepton()->Mass();
855 yl = kinematics::CohYLim(EvL,ml);
856 return yl;
857 }
858 return yl;
859}
860//____________________________________________________________________________
862{
863 // Paschos-Schalla xsi parameter for y-limits in COH
864 // From PRD 80, 033005 (2009)
865
866 Range1D_t yl;
867 yl.min = -1;
868 yl.max = -1;
869
870 const ProcessInfo & pi = fInteraction->ProcInfo();
871
872 //COH
873 bool is_coh = pi.IsCoherentProduction();
874 if(is_coh) {
875 const InitialState & init_state = fInteraction->InitState();
876 const Kinematics & kine = fInteraction->Kine();
877 double Ev = init_state.ProbeE(kRfHitNucRest);
878 double Q2 = kine.Q2();
879 double Mn = init_state.Tgt().Mass();
880 double mlep = fInteraction->FSPrimLepton()->Mass();
881
882 double m_other = controls::kASmallNum ;
883 // as a default the mass of hadronic system is the mass of the photon.
884 // which is assumed to be a small number to avoid divergences
885
886 const XclsTag & xcls = fInteraction -> ExclTag() ;
887
888 if ( xcls.NPions() > 0 ) {
889 bool pionIsCharged = pi.IsWeakCC();
890 m_other = pionIsCharged ? kPionMass : kPi0Mass;
891 }
892
893 yl = kinematics::CohYLim(Mn, m_other, mlep, Ev, Q2, xsi);
894 return yl;
895 } else {
896 return this->YLim();
897 }
898}
899//____________________________________________________________________________
901{
902 // Paschos-Schalla xsi parameter for y-limits in COH
903 // From PRD 80, 033005 (2009)
904
905 const ProcessInfo & pi = fInteraction->ProcInfo();
906
907 //COH
908 bool is_coh = pi.IsCoherentProduction();
909 if(is_coh) {
910 return this->YLim(xsi);
911 } else {
912 return this->YLim_X();
913 }
914}
915//____________________________________________________________________________
917{
918 // t limits for Coherent pion production from
919 // Kartavtsev, Paschos, and Gounaris, PRD 74 054007, and
920 // Paschos and Schalla, PRD 80, 03305
921 // TODO: Attempt to assign t bounds for other reactions?
922 Range1D_t tl;
923 tl.min = -1;
924 tl.max = -1;
925
926 const InitialState & init_state = fInteraction->InitState();
927 const ProcessInfo & pi = fInteraction->ProcInfo();
928 const Kinematics & kine = fInteraction->Kine();
930 double Ev = init_state.ProbeE(kRfHitNucRest);
931 double Q2 = kine.Q2();
932 double nu = Ev * kine.y();
933
934 //COH
935 if(pi.IsCoherentProduction()) {
936
937 double m_other = controls::kASmallNum ;
938 // as a default the mass of hadronic system is the mass of the photon.
939 // which is assumed to be a small number to avoid divergences
940
941 const XclsTag & xcls = fInteraction -> ExclTag() ;
942
943 if ( xcls.NPions() > 0 ) {
944 bool pionIsCharged = pi.IsWeakCC();
945 m_other = pionIsCharged ? kPionMass : kPi0Mass;
946 }
947
948 double m_other2 = m_other * m_other ;
949
950 tl.min = 1.0 * (Q2 + m_other2)/(2.0 * nu) * (Q2 + m_other2)/(2.0 * nu);
951 tl.max = 0.05;
952 return tl;
953 }
954 // DFR
955 else if (pi.IsDiffractive()) {
956
957 // diffractive tmin from Nucl.Phys.B278,61 (1986), eq. 12
958
959 bool pionIsCharged = pi.IsWeakCC();
960 double mpi = pionIsCharged ? kPionMass : kPi0Mass;
961 double mpi2 = mpi*mpi;
962
963 double M = init_state.Tgt().HitNucMass();
964 double M2 = M*M;
965 double nuSqPlusQ2 = nu*nu + Q2;
966 double nuOverM = nu / M;
967 double mpiQ2term = mpi2 - Q2 - 2*nu*nu;
968 double A1 = 1 + 2*nuOverM + nuOverM*nuOverM - nuSqPlusQ2/M2;
969 double A2 = (1+nuOverM) * mpiQ2term + 2*nuOverM*nuSqPlusQ2;
970 double A3 = mpiQ2term*mpiQ2term - 4*nuSqPlusQ2*(nu*nu - mpi2);
971
972 tl.min = std::abs( (A2 + sqrt(A2*A2 - A1*A3)) / A1 ); // GENIE's convention is that t is positive
973 bool tminIsNaN;
974 // use std::isnan when C++11 is around
975#if __cplusplus >= 201103L
976 tminIsNaN = std::isnan(tl.min);
977#else
978 // this the old-fashioned way to check for NaN:
979 // NaN's aren't equal to anything, including themselves
980 tminIsNaN = tl.min != tl.min;
981#endif
982 if (tminIsNaN)
983 {
984 LOG("KPhaseSpace", pERROR)
985 << "tmin for diffractive scattering is NaN "
986 << "( Enu = " << Ev << ", Q2 = " << Q2 << ", nu = " << nu << ")";
987 throw genie::exceptions::InteractionException("NaN tmin for diffractive scattering");
988 }
989 tl.max = this->GetTMaxDFR();
990
991 return tl;
992 }
993
994 // RES+DIS
995 // IMD
996 LOG("KPhaseSpace", pWARN) << "It is not sensible to ask for t limits for events that are not coherent or diffractive.";
997 return tl;
998}
999//____________________________________________________________________________
1001{
1002 const InitialState & init_state = fInteraction->InitState();
1003 PDGLibrary * pdglib = PDGLibrary::Instance();
1004
1005 // imply isospin symmetry
1006 double mpi = (pdglib->Find(kPdgPiP)->Mass() + pdglib->Find(kPdgPi0)->Mass() + pdglib->Find(kPdgPiM)->Mass())/3;
1007 double M = (pdglib->Find(kPdgProton)->Mass() + pdglib->Find(kPdgNeutron)->Mass())/2;
1008 double mi = PDGLibrary::Instance()->Find( init_state.ProbePdg() )->Mass();
1009 double mf = fInteraction->FSPrimLepton()->Mass();
1010 double mtot = M + mf + mpi; // total mass of FS particles
1011 double Ethresh = (mtot*mtot - M*M - mi*mi)/2/M;
1012 return Ethresh;
1013}
1014//____________________________________________________________________________
1016{
1017 Range1D_t Wl;
1018 const InitialState & init_state = fInteraction->InitState();
1020 PDGLibrary * pdglib = PDGLibrary::Instance();
1021 double Mf = pdglib->Find( SppChannel::FinStateNucleon(spp_channel) )->Mass();
1022 double mpi = pdglib->Find( SppChannel::FinStatePion(spp_channel) )->Mass();
1023 double mf = fInteraction->FSPrimLepton()->Mass();
1024 double ECM = init_state.CMEnergy();
1025 // kinematic W-limits
1026 Wl.min = Mf + mpi;
1027 Wl.max = ECM - mf;
1028
1029 if ( (Wl.max - Wl.min) < (Wl.max + Wl.min)*std::numeric_limits<double>::epsilon() )
1030 {
1031 Wl.min = 2*Wl.max*Wl.min/(Wl.max + Wl.min);
1032 Wl.max = Wl.min;
1033 }
1034 else
1035 {
1036 Wl.min *= 1. + std::numeric_limits<double>::epsilon();
1037 Wl.max *= 1. - std::numeric_limits<double>::epsilon();
1038 }
1039
1040 return Wl;
1041}
1042//____________________________________________________________________________
1044{
1045 Range1D_t Wl;
1046 const InitialState & init_state = fInteraction->InitState();
1047 PDGLibrary * pdglib = PDGLibrary::Instance();
1048 // imply isospin symmetry
1049 double M = (pdglib->Find(kPdgProton)->Mass() + pdglib->Find(kPdgNeutron)->Mass())/2;
1050 double mpi = (pdglib->Find(kPdgPiP)->Mass() + pdglib->Find(kPdgPi0)->Mass() + pdglib->Find(kPdgPiM)->Mass())/3;
1051 double mi = PDGLibrary::Instance()->Find( init_state.ProbePdg() )->Mass();
1052 double mf = fInteraction->FSPrimLepton()->Mass();
1053 double Ei = init_state.ProbeE(kRfHitNucRest);
1054 double ECM = TMath::Sqrt(M*(M + 2*Ei) + mi*mi);
1055 // kinematic W-limits
1056 Wl.min = M + mpi;
1057 Wl.max = ECM - mf;
1058
1059 if ( (Wl.max - Wl.min) < (Wl.max + Wl.min)*std::numeric_limits<double>::epsilon() )
1060 {
1061 Wl.min = 2*Wl.max*Wl.min/(Wl.max + Wl.min);
1062 Wl.max = Wl.min;
1063 }
1064 else
1065 {
1066 Wl.min *= 1. + std::numeric_limits<double>::epsilon();
1067 Wl.max *= 1. - std::numeric_limits<double>::epsilon();
1068 }
1069
1070 return Wl;
1071}
1072//____________________________________________________________________________
1074{
1075 Range1D_t Q2l;
1076 const InitialState & init_state = fInteraction->InitState();
1078 PDGLibrary * pdglib = PDGLibrary::Instance();
1079 double Mi = pdglib->Find( SppChannel::InitStateNucleon(spp_channel) )->Mass();
1080 double mi = pdglib->Find( init_state.ProbePdg() )->Mass();
1081 double mf = fInteraction->FSPrimLepton()->Mass();
1082 double mi2 = mi*mi;
1083 double mf2 = mf*mf;
1084 double W = kinematics::W(fInteraction);
1085
1086 double ECM = init_state.CMEnergy();
1087 double s = ECM*ECM;
1088
1089 double Ei_CM = (s + mi2 - Mi*Mi)/2/ECM;
1090 double Ef_CM = (s + mf2 - W*W)/2/ECM;
1091 double Pi_CM = (Ei_CM - mi)<0?0:TMath::Sqrt(Ei_CM*Ei_CM - mi2);
1092 double Pf_CM = (Ef_CM - mf)<0?0:TMath::Sqrt(Ef_CM*Ef_CM - mf2);
1093 // kinematic Q2-limits
1094 Q2l.min = 2*(Ei_CM*Ef_CM - Pi_CM*Pf_CM) - mi2 - mf2;
1095 Q2l.max = 2*(Ei_CM*Ef_CM + Pi_CM*Pf_CM) - mi2 - mf2;
1096
1097 if ( (Q2l.max - Q2l.min) < (Q2l.max + Q2l.min)*std::numeric_limits<double>::epsilon() )
1098 {
1099 Q2l.min = 2*Q2l.max*Q2l.min/(Q2l.max + Q2l.min);
1100 Q2l.max = Q2l.min;
1101 }
1102 else
1103 {
1104 Q2l.min *= 1. + std::numeric_limits<double>::epsilon();
1105 Q2l.max *= 1. - std::numeric_limits<double>::epsilon();
1106 }
1107
1108 return Q2l;
1109}
1110//____________________________________________________________________________
1112{
1113 Range1D_t Q2l;
1114 const InitialState & init_state = fInteraction->InitState();
1115 PDGLibrary * pdglib = PDGLibrary::Instance();
1116 // imply isospin symmetry
1117 double M = (pdglib->Find(kPdgProton)->Mass() + pdglib->Find(kPdgNeutron)->Mass())/2;
1118 double mi = pdglib->Find( init_state.ProbePdg() )->Mass();
1119 double mf = fInteraction->FSPrimLepton()->Mass();
1120 double mi2 = mi*mi;
1121 double mf2 = mf*mf;
1122 double W = kinematics::W(fInteraction);
1123
1124 double Ei = init_state.ProbeE(kRfHitNucRest);
1125 double s = M*(M + 2*Ei) + mi2;
1126 double ECM = TMath::Sqrt(s);
1127
1128 double Ei_CM = (s + mi2 - M*M)/2/ECM;
1129 double Ef_CM = (s + mf2 - W*W)/2/ECM;
1130 double Pi_CM = (Ei_CM - mi)<0?0:TMath::Sqrt(Ei_CM*Ei_CM - mi2);
1131 double Pf_CM = (Ef_CM - mf)<0?0:TMath::Sqrt(Ef_CM*Ef_CM - mf2);
1132 // kinematic Q2-limits
1133 Q2l.min = 2*(Ei_CM*Ef_CM - Pi_CM*Pf_CM) - mi2 - mf2;
1134 Q2l.max = 2*(Ei_CM*Ef_CM + Pi_CM*Pf_CM) - mi2 - mf2;
1135
1136 if ( (Q2l.max - Q2l.min) < (Q2l.max + Q2l.min)*std::numeric_limits<double>::epsilon() )
1137 {
1138 Q2l.min = 2*Q2l.max*Q2l.min/(Q2l.max + Q2l.min);
1139 Q2l.max = Q2l.min;
1140 }
1141 else
1142 {
1143 Q2l.min *= 1. + std::numeric_limits<double>::epsilon();
1144 Q2l.max *= 1. - std::numeric_limits<double>::epsilon();
1145 }
1146
1147 return Q2l;
1148}
1149//____________________________________________________________________________
1150
ClassImp(KPhaseSpace) KPhaseSpace
#define pERROR
Definition Messenger.h:59
#define pDEBUG
Definition Messenger.h:63
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE,...
Definition Messenger.h:96
#define pWARN
Definition Messenger.h:60
#define SLOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a short string (using the FUNCTION and...
Definition Messenger.h:84
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils.
A singleton class holding all configuration registries built while parsing all loaded XML configurati...
static AlgConfigPool * Instance()
Registry * CommonList(const string &file_id, const string &set_name) const
Initial State information.
const Target & Tgt(void) const
int ProbePdg(void) const
double CMEnergy() const
centre-of-mass energy (sqrt s)
double ProbeE(RefFrame_t rf) const
Summary information for an interaction.
Definition Interaction.h:56
Kinematical phase space.
Definition KPhaseSpace.h:33
double Threshold_SPP_iso(void) const
Energy limit for resonance single pion production on isoscalar nucleon.
Range1D_t XLim(void) const
x limits
Range1D_t WLim_SPP_iso(void) const
W limits for resonance single pion production on isoscalar nucleon.
double Minimum(KineVar_t kvar) const
Range1D_t q2Lim_W(void) const
q2 limits @ fixed W
Range1D_t YLim_X(void) const
y limits @ fixed x
bool IsAllowed(void) const
Check whether the current kinematics is in the allowed phase space.
Range1D_t Q2Lim_W(void) const
Q2 limits @ fixed W.
Range1D_t Q2Lim_W_SPP(void) const
Q2 limits @ fixed W for single pion production models.
void UseInteraction(const Interaction *in)
double Threshold(void) const
Energy threshold.
Range1D_t TLim(void) const
t limits
bool IsAboveThreshold(void) const
Checks whether the interaction is above the energy threshold.
Range1D_t YLim(void) const
y limits
Range1D_t q2Lim(void) const
q2 limits
Range1D_t WLim(void) const
W limits.
const Interaction * fInteraction
Definition KPhaseSpace.h:78
Range1D_t Q2Lim_W_SPP_iso(void) const
Q2 limits @ fixed W for resonance single pion production on isoscalar nucleon.
Range1D_t Limits(KineVar_t kvar) const
Return the kinematical variable limits.
static double GetTMaxDFR()
Range1D_t Q2Lim(void) const
Q2 limits.
double Maximum(KineVar_t kvar) const
Range1D_t WLim_SPP(void) const
W limits for single pion production models.
static string AsString(KineVar_t kv)
Definition KineVar.h:77
Generated/set kinematical variables for an event.
Definition Kinematics.h:39
double Q2(bool selected=false) const
double t(bool selected=false) const
double y(bool selected=false) const
double W(bool selected=false) const
double x(bool selected=false) const
Singleton class to load & serve a TDatabasePDG.
Definition PDGLibrary.h:36
static PDGLibrary * Instance(void)
TParticlePDG * Find(int pdgc, bool must_exist=true)
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition ProcessInfo.h:46
bool IsPhotonResonance(void) const
bool IsAMNuGamma(void) const
bool IsNuElectronElastic(void) const
bool IsDiffractive(void) const
bool IsNorm(void) const
bool IsDeepInelastic(void) const
bool IsInverseMuDecay(void) const
bool IsCoherentElastic(void) const
bool IsDarkMatterElastic(void) const
bool IsPhotonCoherent(void) const
bool IsWeakCC(void) const
bool IsCoherentProduction(void) const
bool IsQuasiElastic(void) const
bool IsDarkMatterElectronElastic(void) const
bool IsIMDAnnihilation(void) const
bool IsEM(void) const
bool IsGlashowResonance(void) const
bool IsKnown(void) const
bool IsMEC(void) const
bool IsSinglePion(void) const
bool IsResonant(void) const
bool IsInverseBetaDecay(void) const
bool IsDarkMatterDeepInelastic(void) const
bool IsSingleKaon(void) const
A simple [min,max] interval for doubles.
Definition Range1.h:43
A registry. Provides the container for algorithm configuration parameters.
Definition Registry.h:65
RgDbl GetDouble(RgKey key) const
Definition Registry.cxx:474
static int FinStatePion(SppChannel_t channel)
Definition SppChannel.h:157
static SppChannel_t FromInteraction(const Interaction *interaction)
Definition SppChannel.h:402
static int FinStateNucleon(SppChannel_t channel)
Definition SppChannel.h:130
static int InitStateNucleon(SppChannel_t channel)
Definition SppChannel.h:103
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition Target.h:40
int HitNucPdg(void) const
Definition Target.cxx:304
int N(void) const
Definition Target.h:69
int Z(void) const
Definition Target.h:68
TLorentzVector * HitNucP4Ptr(void) const
Definition Target.cxx:247
double Mass(void) const
Definition Target.cxx:224
int Pdg(void) const
Definition Target.h:71
double HitNucMass(void) const
Definition Target.cxx:233
bool HitNucIsSet(void) const
Definition Target.cxx:283
Contains minimal information for tagging exclusive processes.
Definition XclsTag.h:39
int NPi0(void) const
Definition XclsTag.h:58
bool IsInclusiveCharm(void) const
Definition XclsTag.cxx:54
int NPiMinus(void) const
Definition XclsTag.h:60
bool IsStrangeEvent(void) const
Definition XclsTag.h:53
int NProtons(void) const
Definition XclsTag.h:56
int NPions(void) const
Definition XclsTag.h:62
int NPiPlus(void) const
Definition XclsTag.h:59
bool IsCharmEvent(void) const
Definition XclsTag.h:50
int StrangeHadronPdg(void) const
Definition XclsTag.h:55
int CharmHadronPdg(void) const
Definition XclsTag.h:52
Exception used inside Interaction classes.
Basic constants.
static const double kLightestChmHad
static const double kMinQ2Limit_VLE
Definition Controls.h:42
static const double kASmallNum
Definition Controls.h:40
bool IsNuE(int pdgc)
Definition PDGUtils.cxx:158
int SwitchProtonNeutron(int pdgc)
Definition PDGUtils.cxx:356
bool IsNuMu(int pdgc)
Definition PDGUtils.cxx:163
bool IsNuTau(int pdgc)
Definition PDGUtils.cxx:168
Range1D_t InelXLim(double El, double ml, double M)
Range1D_t InelWLim(double El, double ml, double M)
Range1D_t InelQ2Lim_W(double El, double ml, double M, double W)
Range1D_t InelYLim(double El, double ml, double M)
Range1D_t InelQ2Lim(double El, double ml, double M)
Range1D_t InelYLim_X(double El, double ml, double M, double x)
Range1D_t InelQ2Lim(double Ev, double M, double ml, double Q2min_cut=controls::kMinQ2Limit)
Range1D_t CohXLim(void)
Range1D_t DarkWLim(double Ev, double M, double ml)
void UpdateWQ2FromXY(const Interaction *in)
double W(const Interaction *const i)
Range1D_t DarkQ2Lim_W(double Ev, double M, double ml, double W, double Q2min_cut=controls::kMinQ2Limit)
Range1D_t CEvNSQ2Lim(double Ev)
Range1D_t DarkYLim_X(double Ev, double M, double ml, double x)
Range1D_t CohYLim(double Mn, double m_produced, double mlep, double Ev, double Q2, double xsi)
Range1D_t DarkYLim(double Ev, double M, double ml)
Range1D_t InelWLim(double Ev, double M, double ml)
Range1D_t DarkQ2Lim(double Ev, double M, double ml, double Q2min_cut=controls::kMinQ2Limit)
Range1D_t CohQ2Lim(double Mn, double m_produced, double mlep, double Ev)
Range1D_t DarkXLim(double Ev, double M, double ml)
Range1D_t InelYLim_X(double Ev, double M, double ml, double x)
Range1D_t InelXLim(double Ev, double M, double ml)
Range1D_t InelQ2Lim_W(double Ev, double M, double ml, double W, double Q2min_cut=controls::kMinQ2Limit)
Range1D_t InelYLim(double Ev, double M, double ml)
bool IsWithinLimits(double x, Range1D_t range)
Root of GENIE utility namespaces.
THE MAIN GENIE PROJECT NAMESPACE
Definition AlgCmp.h:25
const int kPdgPiM
Definition PDGCodes.h:159
const int kPdgProton
Definition PDGCodes.h:81
const int kPdgPi0
Definition PDGCodes.h:160
const int kPdgNeutron
Definition PDGCodes.h:83
@ 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
enum genie::EKineVar KineVar_t
enum genie::ESppChannel SppChannel_t
@ kRfHitNucRest
Definition RefFrame.h:30
@ kRfLab
Definition RefFrame.h:26
const int kPdgPiP
Definition PDGCodes.h:158