GENIEGenerator
Loading...
Searching...
No Matches
Interaction.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 Changes required to implement the GENIE Dark Neutrino module
13 were installed by Iker de Icaza (Univ. of Sussex)
14*/
15//____________________________________________________________________________
16
17#include <sstream>
18
19#include <TRootIOCtor.h>
20
27
28using namespace genie;
29using namespace genie::constants;
30
31using std::endl;
32using std::ostringstream;
33
35
36//____________________________________________________________________________
37namespace genie {
38 ostream & operator<< (ostream& stream, const Interaction & interaction)
39 {
40 interaction.Print(stream);
41 return stream;
42 }
43}
44//___________________________________________________________________________
46TObject()
47{
48 this->Init();
49}
50//___________________________________________________________________________
52TObject()
53{
54 this->Init();
55
56 fInitialState -> Copy (ist);
57 fProcInfo -> Copy (prc);
58}
59//___________________________________________________________________________
61TObject()
62{
63 this->Init();
64 this->Copy(interaction);
65}
66//___________________________________________________________________________
68TObject(),
70fProcInfo(0),
73fKinePhSp(0)
74{
75
76}
77//___________________________________________________________________________
79{
80 this->CleanUp();
81}
82//___________________________________________________________________________
84{
85 this->CleanUp();
86 this->Init();
87}
88//___________________________________________________________________________
90{
92 fProcInfo = new ProcessInfo ();
93 fKinematics = new Kinematics ();
94 fExclusiveTag = new XclsTag ();
95 fKinePhSp = new KPhaseSpace (this);
96}
97//___________________________________________________________________________
99{
100 if ( fInitialState ) delete fInitialState;
101 if ( fProcInfo ) delete fProcInfo;
102 if ( fKinematics ) delete fKinematics;
103 if ( fExclusiveTag ) delete fExclusiveTag;
104 if ( fKinePhSp ) delete fKinePhSp;
105
106 fInitialState = 0;
107 fProcInfo = 0;
108 fKinematics = 0;
109 fExclusiveTag = 0;
110 fKinePhSp = 0;
111}
112//___________________________________________________________________________
113void Interaction::Copy(const Interaction & interaction)
114{
115 const InitialState & init = *interaction.fInitialState;
116 const ProcessInfo & proc = *interaction.fProcInfo;
117 const Kinematics & kine = *interaction.fKinematics;
118 const XclsTag & xcls = *interaction.fExclusiveTag;
119
120 fInitialState -> Copy (init);
121 fProcInfo -> Copy (proc);
122 fKinematics -> Copy (kine);
123 fExclusiveTag -> Copy (xcls);
124}
125//___________________________________________________________________________
126TParticlePDG * Interaction::FSPrimLepton(void) const
127{
128 int pdgc = this->FSPrimLeptonPdg();
129
130 if(pdgc) return PDGLibrary::Instance()->Find(pdgc);
131 else return 0;
132}
133//___________________________________________________________________________
135{
136 const ProcessInfo & proc_info = this -> ProcInfo();
137 const InitialState & init_state = this -> InitState();
138 const XclsTag & xclstag = this -> ExclTag();
139
140 int pdgc = init_state.ProbePdg();
141
142 LOG("Interaction", pDEBUG) << "Probe PDG code: " << pdgc;
143
144 if (proc_info.IsNuElectronElastic())
145 return kPdgElectron;
146
147 if (proc_info.IsGlashowResonance() || proc_info.IsPhotonResonance())
148 return xclstag.FinalLeptonPdg();
149
150 // vN (Weak-NC) or eN (EM)
151 if (proc_info.IsWeakNC() || proc_info.IsEM() || proc_info.IsWeakMix() || proc_info.IsDarkMatter()) return pdgc; // EDIT: DM does not change in FS
152
153 // vN (Weak-CC)
154 else if (proc_info.IsWeakCC()) {
155 int clpdgc;
156 if (proc_info.IsIMDAnnihilation())
157 clpdgc = kPdgMuon;
158 else
159 clpdgc = pdg::Neutrino2ChargedLepton(pdgc);
160 return clpdgc;
161 }
162 else if (proc_info.IsDarkNeutralCurrent()){
163 return kPdgDarkNeutrino;
164 }
165 LOG("Interaction", pWARN)
166 << "Could not figure out the final state primary lepton pdg code!!";
167
168 return 0;
169}
170//___________________________________________________________________________
171TParticlePDG * Interaction::RecoilNucleon(void) const
172{
173 int rnuc = this->RecoilNucleonPdg();
174
175 if(rnuc) return PDGLibrary::Instance()->Find(rnuc);
176 else return 0;
177}
178//___________________________________________________________________________
180{
181// Determine the recoil nucleon PDG code
182
183 const Target & target = fInitialState->Tgt();
184
185 int recoil_nuc = 0;
186 int struck_nuc = target.HitNucPdg();
187
188 if (fProcInfo->IsQuasiElastic() || fProcInfo->IsInverseBetaDecay() || fProcInfo->IsDarkMatterElastic()) {
189 bool struck_is_nuc = pdg::IsNucleon(struck_nuc);
190 bool is_weak = fProcInfo->IsWeak();
191 bool is_em = fProcInfo->IsEM();
192 bool is_dm = fProcInfo->IsDarkMatter();
193 assert(struck_is_nuc && (is_weak || is_em || is_dm));
194 if(fProcInfo->IsWeakCC()) {
195 recoil_nuc = pdg::SwitchProtonNeutron(struck_nuc); // CC
196 } else {
197 recoil_nuc = struck_nuc; // NC, EM
198 }
199 }
200
201 if (fProcInfo->IsMEC()) {
202 bool struck_is_2nuc_cluster = pdg::Is2NucleonCluster(struck_nuc);
203 bool is_weak = fProcInfo->IsWeak();
204 bool is_em = fProcInfo->IsEM();
205 assert(struck_is_2nuc_cluster && (is_weak || is_em));
206 if(fProcInfo->IsWeakCC()) {
207 bool isnu = pdg::IsNeutrino(fInitialState->ProbePdg());
208 // nucleon cluster charge should be incremented by +1 for
209 // neutrino CC and by -1 for antineutrino CC
210 int dQ = (isnu) ? +1 : -1;
211 recoil_nuc = pdg::ModifyNucleonCluster(struck_nuc,dQ); // CC
212 }
213 else {
214 recoil_nuc = struck_nuc; // NC, EM
215 }
216 }
217
218 LOG("Interaction", pDEBUG) << "Recoil nucleon PDG = " << recoil_nuc;
219 return recoil_nuc;
220}
221//___________________________________________________________________________
223{
225 fInitialState->Copy(init_state);
226}
227//___________________________________________________________________________
229{
230 if (!fProcInfo) fProcInfo = new ProcessInfo();
231 fProcInfo->Copy(proc_info);
232}
233//___________________________________________________________________________
234void Interaction::SetKine(const Kinematics & kinematics)
235{
236 if (!fKinematics) fKinematics = new Kinematics();
237 fKinematics->Copy(kinematics);
238}
239//___________________________________________________________________________
240void Interaction::SetExclTag(const XclsTag & xcls_tag)
241{
242 if (!fExclusiveTag) fExclusiveTag = new XclsTag();
243 fExclusiveTag->Copy(xcls_tag);
244}
245//___________________________________________________________________________
246string Interaction::AsString(void) const
247{
248// Code-ify the interaction in a string to be used as (part of a) cache
249// branch key.
250// Template:
251// nu:x;tgt:x;N:x;q:x(s/v);proc:x;xclv_tag
252
253 const Target & tgt = fInitialState->Tgt();
254
255 ostringstream interaction;
256
257 // If the probe has non-zero mass, then it is DM
258 if (fInitialState->Probe()->PdgCode() == kPdgDarkMatter) {
259 interaction << "dm;";
260 }
261 else if (fInitialState->Probe()->PdgCode() == kPdgAntiDarkMatter) {
262 interaction << "dmb;";
263 }
264 else {
265 interaction << "nu:" << fInitialState->ProbePdg() << ";";
266 }
267 interaction << "tgt:" << tgt.Pdg() << ";";
268
269 if(tgt.HitNucIsSet()) {
270 interaction << "N:" << tgt.HitNucPdg() << ";";
271 }
272 if(tgt.HitQrkIsSet()) {
273 interaction << "q:" << tgt.HitQrkPdg()
274 << (tgt.HitSeaQrk() ? "(s)" : "(v)") << ";";
275 }
276
277 interaction << "proc:" << fProcInfo->InteractionTypeAsString()
278 << "," << fProcInfo->ScatteringTypeAsString() << ";";
279
280 string xcls = fExclusiveTag->AsString();
281 interaction << xcls;
282 if(xcls.size()>0) interaction << ";";
283
284 return interaction.str();
285}
286//___________________________________________________________________________
287void Interaction::Print(ostream & stream) const
288{
289 const string line(110, '-');
290
291 stream << endl;
292 stream << line << endl;
293
294 stream << "GENIE Interaction Summary" << endl;
295 stream << line << endl;
296
297 stream << *fInitialState << endl; // print initial state
298 stream << *fProcInfo; // print process info
299 stream << *fKinematics; // print scattering parameters
300 stream << *fExclusiveTag; // print exclusive process tag
301
302 stream << line << endl;
303}
304//___________________________________________________________________________
306{
307 this->Copy(interaction);
308 return (*this);
309}
310//___________________________________________________________________________
311//
312// **** Methods using the "named constructor" C++ idiom ****
313//
314//___________________________________________________________________________
316 int target, int probe, ScatteringType_t st, InteractionType_t it)
317{
318 InitialState init_state (target, probe);
319 ProcessInfo proc_info (st, it);
320
321 Interaction * interaction = new Interaction(init_state, proc_info);
322 return interaction;
323}
324//___________________________________________________________________________
325Interaction * Interaction::DISCC(int target, int hitnuc, int probe, double E)
326{
327 Interaction * interaction =
329
330 InitialState * init_state = interaction->InitStatePtr();
331 init_state->SetProbeE(E);
332 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
333
334 return interaction;
335}
336//___________________________________________________________________________
338 int target, int hitnuc, int hitqrk, bool fromsea, int probe, double E)
339{
340 Interaction* interaction = Interaction::DISCC(target,hitnuc,probe,E);
341
342 Target * tgt = interaction->InitStatePtr()->TgtPtr();
343 tgt -> SetHitQrkPdg (hitqrk);
344 tgt -> SetHitSeaQrk (fromsea);
345
346 return interaction;
347}
348//___________________________________________________________________________
350 int target, int hitnuc, int hitqrk, bool fromsea, int fqrk, int probe, double E)
351{
352 Interaction* interaction = Interaction::DISCC(target,hitnuc,probe,E);
353
354 Target * tgt = interaction->InitStatePtr()->TgtPtr();
355 tgt -> SetHitQrkPdg (hitqrk);
356 tgt -> SetHitSeaQrk (fromsea);
357
358 XclsTag * xclstag = interaction->ExclTagPtr();
359 xclstag->SetFinalQuark(fqrk);
360
361 return interaction;
362}
363//___________________________________________________________________________
365 int target, int hitnuc, int probe, const TLorentzVector & p4probe)
366{
367 Interaction * interaction =
369
370 InitialState * init_state = interaction->InitStatePtr();
371 init_state->SetProbeP4(p4probe);
372 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
373
374 return interaction;
375}
376//___________________________________________________________________________
378 int target, int hitnuc, int hitqrk, bool fromsea, int probe,
379 const TLorentzVector & p4probe)
380{
381 Interaction* interaction = Interaction::DISCC(target,hitnuc,probe,p4probe);
382
383 Target * tgt = interaction->InitStatePtr()->TgtPtr();
384 tgt -> SetHitQrkPdg (hitqrk);
385 tgt -> SetHitSeaQrk (fromsea);
386
387 return interaction;
388}
389//___________________________________________________________________________
390Interaction * Interaction::DISNC(int target, int hitnuc, int probe, double E)
391{
392 Interaction * interaction =
394
395 InitialState * init_state = interaction->InitStatePtr();
396 init_state->SetProbeE(E);
397 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
398
399 return interaction;
400}
401//___________________________________________________________________________
403 int target, int hitnuc, int hitqrk, bool fromsea, int probe, double E)
404{
405 Interaction* interaction = Interaction::DISNC(target,hitnuc,probe,E);
406
407 Target * tgt = interaction->InitStatePtr()->TgtPtr();
408 tgt -> SetHitQrkPdg (hitqrk);
409 tgt -> SetHitSeaQrk (fromsea);
410
411 return interaction;
412}
413//___________________________________________________________________________
415 int target, int hitnuc, int hitqrk, bool fromsea, int fqrk, int probe, double E)
416{
417 Interaction* interaction = Interaction::DISNC(target,hitnuc,probe,E);
418
419 Target * tgt = interaction->InitStatePtr()->TgtPtr();
420 tgt -> SetHitQrkPdg (hitqrk);
421 tgt -> SetHitSeaQrk (fromsea);
422
423 XclsTag * xclstag = interaction->ExclTagPtr();
424 xclstag->SetFinalQuark(fqrk);
425
426 return interaction;
427}
428//___________________________________________________________________________
430 int target, int hitnuc, int probe, const TLorentzVector & p4probe)
431{
432 Interaction * interaction =
434
435 InitialState * init_state = interaction->InitStatePtr();
436 init_state->SetProbeP4(p4probe);
437 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
438
439 return interaction;
440}
441//___________________________________________________________________________
443 int target, int hitnuc, int hitqrk, bool fromsea, int probe,
444 const TLorentzVector & p4probe)
445{
446 Interaction * interaction = Interaction::DISNC(target,hitnuc,probe,p4probe);
447
448 Target * tgt = interaction->InitStatePtr()->TgtPtr();
449 tgt -> SetHitQrkPdg (hitqrk);
450 tgt -> SetHitSeaQrk (fromsea);
451
452 return interaction;
453}
454//___________________________________________________________________________
455Interaction * Interaction::DISEM(int target, int hitnuc, int probe, double E)
456{
457 Interaction * interaction =
459
460 InitialState * init_state = interaction->InitStatePtr();
461 init_state->SetProbeE(E);
462 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
463
464 return interaction;
465}
466//___________________________________________________________________________
468 int target, int hitnuc, int hitqrk, bool fromsea, int probe, double E)
469{
470 Interaction* interaction = Interaction::DISEM(target,hitnuc,probe,E);
471
472 Target * tgt = interaction->InitStatePtr()->TgtPtr();
473 tgt -> SetHitQrkPdg (hitqrk);
474 tgt -> SetHitSeaQrk (fromsea);
475
476 return interaction;
477}
478//___________________________________________________________________________
480 int target, int hitnuc, int probe, const TLorentzVector & p4probe)
481{
482 Interaction * interaction =
484
485 InitialState * init_state = interaction->InitStatePtr();
486 init_state->SetProbeP4(p4probe);
487 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
488
489 return interaction;
490}
491//___________________________________________________________________________
493 int target, int hitnuc, int hitqrk, bool fromsea, int probe,
494 const TLorentzVector & p4probe)
495{
496 Interaction * interaction = Interaction::DISEM(target,hitnuc,probe,p4probe);
497
498 Target * tgt = interaction->InitStatePtr()->TgtPtr();
499 tgt -> SetHitQrkPdg (hitqrk);
500 tgt -> SetHitSeaQrk (fromsea);
501
502 return interaction;
503}
504//___________________________________________________________________________
505Interaction * Interaction::QELCC(int target, int hitnuc, int probe, double E)
506{
507 Interaction * interaction =
509
510 InitialState * init_state = interaction->InitStatePtr();
511 init_state->SetProbeE(E);
512 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
513
514 return interaction;
515}
516//___________________________________________________________________________
518 int target, int hitnuc, int probe, const TLorentzVector & p4probe)
519{
520 Interaction * interaction =
522
523 InitialState * init_state = interaction->InitStatePtr();
524 init_state->SetProbeP4(p4probe);
525 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
526
527 return interaction;
528}
529//___________________________________________________________________________
530Interaction * Interaction::QELNC(int target, int hitnuc, int probe, double E)
531{
532 Interaction * interaction =
534
535 InitialState * init_state = interaction->InitStatePtr();
536 init_state->SetProbeE(E);
537 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
538
539 return interaction;
540}
541//___________________________________________________________________________
543 int target, int hitnuc, int probe, const TLorentzVector & p4probe)
544{
545 Interaction * interaction =
547
548 InitialState * init_state = interaction->InitStatePtr();
549 init_state->SetProbeP4(p4probe);
550 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
551
552 return interaction;
553}
554//___________________________________________________________________________
555Interaction * Interaction::QELEM(int target, int hitnuc, int probe, double E)
556{
557 Interaction * interaction =
559
560 InitialState * init_state = interaction->InitStatePtr();
561 init_state->SetProbeE(E);
562 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
563
564 return interaction;
565}
566//___________________________________________________________________________
568 int target, int hitnuc, int probe, const TLorentzVector & p4probe)
569{
570 Interaction * interaction =
572
573 InitialState * init_state = interaction->InitStatePtr();
574 init_state->SetProbeP4(p4probe);
575 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
576
577 return interaction;
578}
579//___________________________________________________________________________
580Interaction * Interaction::IBD(int target, int hitnuc, int probe, double E)
581{
582 Interaction * interaction =
584
585 InitialState * init_state = interaction->InitStatePtr();
586 init_state->SetProbeE(E);
587 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
588
589 return interaction;
590}
591//___________________________________________________________________________
593 int target, int hitnuc, int probe, const TLorentzVector & p4probe)
594{
595 Interaction * interaction =
597
598 InitialState * init_state = interaction->InitStatePtr();
599 init_state->SetProbeP4(p4probe);
600 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
601
602 return interaction;
603}
604//___________________________________________________________________________
605Interaction * Interaction::RESCC(int target, int hitnuc, int probe, double E)
606{
607 Interaction * interaction =
609
610 InitialState * init_state = interaction->InitStatePtr();
611 init_state->SetProbeE(E);
612 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
613
614 return interaction;
615}
616//___________________________________________________________________________
618 int target, int hitnuc, int probe, const TLorentzVector & p4probe)
619{
620 Interaction * interaction =
622
623 InitialState * init_state = interaction->InitStatePtr();
624 init_state->SetProbeP4(p4probe);
625 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
626
627 return interaction;
628}
629//___________________________________________________________________________
630Interaction * Interaction::RESNC(int target, int hitnuc, int probe, double E)
631{
632 Interaction * interaction =
634
635 InitialState * init_state = interaction->InitStatePtr();
636 init_state->SetProbeE(E);
637 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
638
639 return interaction;
640}
641//___________________________________________________________________________
643 int target, int hitnuc, int probe, const TLorentzVector & p4probe)
644{
645 Interaction * interaction =
647
648 InitialState * init_state = interaction->InitStatePtr();
649 init_state->SetProbeP4(p4probe);
650 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
651
652 return interaction;
653}
654//___________________________________________________________________________
655Interaction * Interaction::RESEM(int target, int hitnuc, int probe, double E)
656{
657 Interaction * interaction =
659
660 InitialState * init_state = interaction->InitStatePtr();
661 init_state->SetProbeE(E);
662 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
663
664 return interaction;
665}
666//___________________________________________________________________________
668 int target, int hitnuc, int probe, const TLorentzVector & p4probe)
669{
670 Interaction * interaction =
672
673 InitialState * init_state = interaction->InitStatePtr();
674 init_state->SetProbeP4(p4probe);
675 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
676
677 return interaction;
678}
679//___________________________________________________________________________
680Interaction * Interaction::DFRCC(int tgt,int hitnuc, int probe, double E)
681{
682 Interaction * interaction =
684
685 InitialState * init_state = interaction->InitStatePtr();
686 init_state->SetProbeE(E);
687 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
688
689 return interaction;
690}
691//___________________________________________________________________________
693 int tgt, int hitnuc, int probe, const TLorentzVector & p4probe)
694{
695 Interaction * interaction =
697
698 InitialState * init_state = interaction->InitStatePtr();
699 init_state->SetProbeP4(p4probe);
700 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
701
702 return interaction;
703}
704//___________________________________________________________________________
705Interaction * Interaction::COHCC(int tgt, int probe, unsigned int prod_pdg, double E)
706{
707 Interaction * interaction =
709
710 XclsTag * xcl = interaction -> ExclTagPtr() ;
711 if ( pdg::IsPion( prod_pdg ) ) {
712 if ( pdg::IsNeutrino( probe ) ) xcl -> SetNPions( 1,0,0 ) ;
713 else if ( pdg::IsAntiNeutrino( probe ) ) xcl -> SetNPions( 0,0,1 ) ;
714 }
715
716 InitialState * init_state = interaction->InitStatePtr();
717 init_state->SetProbeE(E);
718
719 return interaction;
720}
721//___________________________________________________________________________
723 int tgt, int probe, unsigned int prod_pdg, const TLorentzVector & p4probe)
724{
725 Interaction * interaction =
727
728 XclsTag * xcl = interaction -> ExclTagPtr() ;
729 if ( pdg::IsPion( prod_pdg ) ) {
730 if ( pdg::IsNeutrino( probe ) ) xcl -> SetNPions( 1,0,0 ) ;
731 else if ( pdg::IsAntiNeutrino( probe ) ) xcl -> SetNPions( 0,0,1 ) ;
732 }
733
734 InitialState * init_state = interaction->InitStatePtr();
735 init_state->SetProbeP4(p4probe);
736
737 return interaction;
738}
739//___________________________________________________________________________
740Interaction * Interaction::COHNC(int tgt, int probe, unsigned int prod_pdg, double E)
741{
742 Interaction * interaction =
744
745 XclsTag * xcl = interaction -> ExclTagPtr() ;
746 if ( pdg::IsPion( prod_pdg ) ) xcl -> SetNPions( 0,1,0 ) ;
747 else if ( prod_pdg == kPdgGamma ) xcl -> SetNSingleGammas(1) ;
748
749 InitialState * init_state = interaction->InitStatePtr();
750 init_state->SetProbeE(E);
751
752 return interaction;
753}
754//___________________________________________________________________________
756 int tgt, int probe, unsigned int prod_pdg, const TLorentzVector & p4probe)
757{
758 Interaction * interaction =
760
761 XclsTag * xcl = interaction -> ExclTagPtr() ;
762 if ( pdg::IsPion( prod_pdg ) ) xcl -> SetNPions( 0,1,0 ) ;
763 else if ( prod_pdg == kPdgGamma ) xcl -> SetNSingleGammas(1) ;
764
765 InitialState * init_state = interaction->InitStatePtr();
766 init_state->SetProbeP4(p4probe);
767
768 return interaction;
769}
770//___________________________________________________________________________
771Interaction * Interaction::CEvNS(int tgt, int probe, double E)
772{
773 Interaction * interaction =
775
776 InitialState * init_state = interaction->InitStatePtr();
777 init_state->SetProbeE(E);
778
779 return interaction;
780}
781//___________________________________________________________________________
783 int tgt, int probe, const TLorentzVector & p4probe)
784{
785 Interaction * interaction =
787
788 InitialState * init_state = interaction->InitStatePtr();
789 init_state->SetProbeP4(p4probe);
790
791 return interaction;
792}
793//___________________________________________________________________________
794Interaction * Interaction::IMD(int target, double E)
795{
796 Interaction * interaction =
798
799 InitialState * init_state = interaction->InitStatePtr();
800 init_state->SetProbeE(E);
801
802 return interaction;
803}
804//___________________________________________________________________________
805Interaction * Interaction::IMD(int target, const TLorentzVector & p4probe)
806{
807 Interaction * interaction =
809
810 InitialState * init_state = interaction->InitStatePtr();
811 init_state->SetProbeP4(p4probe);
812
813 return interaction;
814}
815//___________________________________________________________________________
816Interaction * Interaction::AMNuGamma(int tgt, int nuc, int probe, double E)
817{
818 Interaction * interaction =
820
821 InitialState * init_state = interaction->InitStatePtr();
822 init_state->SetProbeE(E);
823 init_state->TgtPtr()->SetHitNucPdg(nuc);
824
825 return interaction;
826}
827//___________________________________________________________________________
829 int tgt, int nuc, int probe, const TLorentzVector & p4probe)
830{
831 Interaction * interaction =
833
834 InitialState * init_state = interaction->InitStatePtr();
835 init_state->SetProbeP4(p4probe);
836 init_state->TgtPtr()->SetHitNucPdg(nuc);
837
838 return interaction;
839}
840//___________________________________________________________________________
841Interaction * Interaction::MECCC(int tgt, int ncluster, int probe, double E)
842{
843 Interaction * interaction =
845
846 InitialState * init_state = interaction->InitStatePtr();
847 init_state->SetProbeE(E);
848 init_state->TgtPtr()->SetHitNucPdg(ncluster);
849
850 return interaction;
851}
852//___________________________________________________________________________
854 int tgt, int ncluster, int probe, const TLorentzVector & p4probe)
855{
856 Interaction * interaction =
858
859 InitialState * init_state = interaction->InitStatePtr();
860 init_state->SetProbeP4(p4probe);
861 init_state->TgtPtr()->SetHitNucPdg(ncluster);
862
863 return interaction;
864}
865//___________________________________________________________________________
866Interaction * Interaction::MECCC(int tgt, int probe, double E)
867{
868 Interaction * interaction =
870
871 InitialState * init_state = interaction->InitStatePtr();
872 init_state->SetProbeE(E);
873
874 return interaction;
875}
876//___________________________________________________________________________
878 int tgt, int probe, const TLorentzVector & p4probe)
879{
880 Interaction * interaction =
882
883 InitialState * init_state = interaction->InitStatePtr();
884 init_state->SetProbeP4(p4probe);
885
886 return interaction;
887}
888
889//___________________________________________________________________________
890Interaction * Interaction::MECNC(int tgt, int ncluster, int probe, double E)
891{
892 Interaction * interaction =
894
895 InitialState * init_state = interaction->InitStatePtr();
896 init_state->SetProbeE(E);
897 init_state->TgtPtr()->SetHitNucPdg(ncluster);
898
899 return interaction;
900}
901//___________________________________________________________________________
903 int tgt, int ncluster, int probe, const TLorentzVector & p4probe)
904{
905 Interaction * interaction =
907
908 InitialState * init_state = interaction->InitStatePtr();
909 init_state->SetProbeP4(p4probe);
910 init_state->TgtPtr()->SetHitNucPdg(ncluster);
911
912 return interaction;
913}
914//___________________________________________________________________________
915Interaction * Interaction::MECEM(int tgt, int probe, double E)
916{
917
918 Interaction * interaction =
919 Interaction::Create(tgt, probe, kScMEC, kIntEM);
920
921 InitialState * init_state = interaction->InitStatePtr();
922 init_state->SetProbeE(E);
923
924 return interaction;
925}
926//___________________________________________________________________________
927Interaction * Interaction::MECEM(int tgt, int ncluster, int probe, double E)
928{
929 Interaction * interaction =
930 Interaction::Create(tgt, probe, kScMEC, kIntEM);
931
932 InitialState * init_state = interaction->InitStatePtr();
933 init_state->SetProbeE(E);
934 init_state->TgtPtr()->SetHitNucPdg(ncluster);
935
936 return interaction;
937}
938//___________________________________________________________________________
940 int tgt, int ncluster, int probe, const TLorentzVector & p4probe)
941{
942 Interaction * interaction =
943 Interaction::Create(tgt, probe, kScMEC, kIntEM);
944
945 InitialState * init_state = interaction->InitStatePtr();
946 init_state->SetProbeP4(p4probe);
947 init_state->TgtPtr()->SetHitNucPdg(ncluster);
948
949 return interaction;
950}
951//___________________________________________________________________________
952Interaction * Interaction::GLR(int tgt, double E)
953{
954 Interaction * interaction =
956
957 InitialState * init_state = interaction->InitStatePtr();
958 init_state->SetProbeE(E);
959 init_state->TgtPtr()->SetHitNucPdg(0);
960
961 return interaction;
962}
963//___________________________________________________________________________
964Interaction * Interaction::GLR(int tgt, const TLorentzVector & p4probe)
965{
966 Interaction * interaction =
968
969 InitialState * init_state = interaction->InitStatePtr();
970 init_state->SetProbeP4(p4probe);
971 init_state->TgtPtr()->SetHitNucPdg(0);
972
973 return interaction;
974}
975//___________________________________________________________________________
976Interaction * Interaction::NDecay(int tgt, int decay_mode, int decayed_nucleon)
977{
978 Interaction * interaction =
980 interaction->ExclTagPtr()->SetDecayMode(decay_mode);
981
982 InitialState * init_state = interaction->InitStatePtr();
983 init_state->TgtPtr()->SetHitNucPdg(decayed_nucleon);
984
985 return interaction;
986}
987//___________________________________________________________________________
988Interaction * Interaction::NOsc(int tgt, int annihilation_mode)
989{
990 Interaction * interaction =
992 interaction->ExclTagPtr()->SetDecayMode(annihilation_mode);
993 return interaction;
994}
995//___________________________________________________________________________
996Interaction * Interaction::ASK(int tgt, int probe, double E)
997{
998 Interaction * interaction =
1000
1001 InitialState * init_state = interaction->InitStatePtr();
1002 init_state->SetProbeE(E);
1003
1004 return interaction;
1005}
1006//___________________________________________________________________________
1008 int tgt, int probe, const TLorentzVector & p4probe)
1009{
1010 Interaction * interaction =
1012
1013 InitialState * init_state = interaction->InitStatePtr();
1014 init_state->SetProbeP4(p4probe);
1015
1016 return interaction;
1017}
1018//___________________________________________________________________________
1019Interaction * Interaction::DME(int target, int hitnuc, int probe, double E)
1020{
1021 // EDIT: need to be able to create dark matter elastic
1022 Interaction * interaction =
1024
1025 InitialState * init_state = interaction->InitStatePtr();
1026 init_state->SetProbeE(E);
1027 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
1028
1029 return interaction;
1030}
1031//___________________________________________________________________________
1033 int target, int hitnuc, int probe, const TLorentzVector & p4probe)
1034{
1035 // EDIT: need to be able to create dark matter elastic
1036 Interaction * interaction =
1038
1039 InitialState * init_state = interaction->InitStatePtr();
1040 init_state->SetProbeP4(p4probe);
1041 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
1042
1043 return interaction;
1044}
1045//___________________________________________________________________________
1046Interaction * Interaction::DMDI(int target, int hitnuc, int probe, double E)
1047{
1048 Interaction * interaction =
1050
1051 InitialState * init_state = interaction->InitStatePtr();
1052 init_state->SetProbeE(E);
1053 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
1054
1055 return interaction;
1056}
1057//___________________________________________________________________________
1059 int target, int hitnuc, int hitqrk, bool fromsea, int probe, double E)
1060{
1061 Interaction* interaction = Interaction::DMDI(target,hitnuc,probe,E);
1062
1063 Target * tgt = interaction->InitStatePtr()->TgtPtr();
1064 tgt -> SetHitQrkPdg (hitqrk);
1065 tgt -> SetHitSeaQrk (fromsea);
1066
1067 return interaction;
1068}
1069//___________________________________________________________________________
1071 int target, int hitnuc, int probe, const TLorentzVector & p4probe)
1072{
1073 Interaction * interaction =
1075
1076 InitialState * init_state = interaction->InitStatePtr();
1077 init_state->SetProbeP4(p4probe);
1078 init_state->TgtPtr()->SetHitNucPdg(hitnuc);
1079
1080 return interaction;
1081}
1082//___________________________________________________________________________
1084 int target, int hitnuc, int hitqrk, bool fromsea, int probe,
1085 const TLorentzVector & p4probe)
1086{
1087 Interaction * interaction = Interaction::DMDI(target,hitnuc,probe,p4probe);
1088
1089 Target * tgt = interaction->InitStatePtr()->TgtPtr();
1090 tgt -> SetHitQrkPdg (hitqrk);
1091 tgt -> SetHitSeaQrk (fromsea);
1092
1093 return interaction;
1094}
1095//___________________________________________________________________________
1096Interaction * Interaction::HNL(int probe, double E, int decayed_mode)
1097{
1098 Interaction * interaction =
1100 interaction->ExclTagPtr()->SetDecayMode(decayed_mode);
1101
1102 InitialState * init_state = interaction->InitStatePtr();
1103 init_state->SetProbeE(E);
1104
1105 return interaction;
1106}
ClassImp(Interaction) namespace genie
#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
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils.
Initial State information.
void SetProbeE(double E)
void SetProbeP4(const TLorentzVector &P4)
int ProbePdg(void) const
Target * TgtPtr(void) const
Summary information for an interaction.
Definition Interaction.h:56
string AsString(void) const
void Copy(const Interaction &i)
TParticlePDG * RecoilNucleon(void) const
recoil nucleon
const XclsTag & ExclTag(void) const
Definition Interaction.h:72
static Interaction * MECCC(int tgt, int nuccluster, int probe, double E=0)
static Interaction * MECNC(int tgt, int nuccluster, int probe, double E=0)
static Interaction * DISEM(int tgt, int nuc, int probe, double E=0)
static Interaction * DMDI(int tgt, int nuc, int probe, double E=0)
void SetInitState(const InitialState &init)
ProcessInfo * fProcInfo
Process info (scattering, weak current,...)
InitialState * InitStatePtr(void) const
Definition Interaction.h:74
static Interaction * DISNC(int tgt, int nuc, int probe, double E=0)
static Interaction * ASK(int tgt, int probe, double E=0)
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
static Interaction * DME(int tgt, int nuc, int probe, double E=0)
void SetProcInfo(const ProcessInfo &proc)
static Interaction * NOsc(int tgt, int annihilation_mode=-1)
int RecoilNucleonPdg(void) const
recoil nucleon pdg
int FSPrimLeptonPdg(void) const
final state primary lepton pdg
static Interaction * RESCC(int tgt, int nuc, int probe, double E=0)
XclsTag * fExclusiveTag
Additional info for exclusive channels.
static Interaction * MECEM(int tgt, int nuccluster, int probe, double E=0)
static Interaction * CEvNS(int tgt, int probe, double E=0)
void SetKine(const Kinematics &kine)
static Interaction * QELCC(int tgt, int nuc, int probe, double E=0)
Interaction & operator=(const Interaction &i)
copy
XclsTag * ExclTagPtr(void) const
Definition Interaction.h:77
Kinematics * fKinematics
kinematical variables
KPhaseSpace * fKinePhSp
Kinematic phase space.
static Interaction * QELEM(int tgt, int nuc, int probe, double E=0)
static Interaction * RESNC(int tgt, int nuc, int probe, double E=0)
static Interaction * DISCC(int tgt, int nuc, int probe, double E=0)
static Interaction * COHNC(int tgt, int probe, unsigned int prod_pdg, double E=0)
static Interaction * IMD(int tgt, double E=0)
static Interaction * NDecay(int tgt, int decay_mode=-1, int decayed_nucleon=0)
const ProcessInfo & ProcInfo(void) const
Definition Interaction.h:70
static Interaction * DFRCC(int tgt, int nuc, int probe, double E=0)
static Interaction * QELNC(int tgt, int nuc, int probe, double E=0)
static Interaction * IBD(int tgt, int nuc, int probe, double E=0)
static Interaction * AMNuGamma(int tgt, int nuc, int probe, double E=0)
static Interaction * RESEM(int tgt, int nuc, int probe, double E=0)
static Interaction * HNL(int probe, double E=0, int decayed_mode=-1)
InitialState * fInitialState
Initial State info.
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
void SetExclTag(const XclsTag &xcls)
void Print(ostream &stream) const
const InitialState & InitState(void) const
Definition Interaction.h:69
static Interaction * GLR(int tgt, double E=0)
static Interaction * COHCC(int tgt, int probe, unsigned int prod_pdg, double E=0)
Kinematical phase space.
Definition KPhaseSpace.h:33
Generated/set kinematical variables for an event.
Definition Kinematics.h:39
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 IsWeakNC(void) const
bool IsNuElectronElastic(void) const
bool IsWeakMix(void) const
bool IsWeakCC(void) const
bool IsDarkMatter(void) const
bool IsIMDAnnihilation(void) const
bool IsEM(void) const
bool IsGlashowResonance(void) const
bool IsDarkNeutralCurrent(void) const
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
void SetHitNucPdg(int pdgc)
Definition Target.cxx:171
int HitQrkPdg(void) const
Definition Target.cxx:242
bool HitSeaQrk(void) const
Definition Target.cxx:299
int Pdg(void) const
Definition Target.h:71
bool HitQrkIsSet(void) const
Definition Target.cxx:292
bool HitNucIsSet(void) const
Definition Target.cxx:283
Contains minimal information for tagging exclusive processes.
Definition XclsTag.h:39
void SetFinalQuark(int finalquark_pdgc=0)
Definition XclsTag.cxx:138
void SetDecayMode(int decay_mode)
Definition XclsTag.cxx:133
int FinalLeptonPdg(void) const
Definition XclsTag.h:74
Basic constants.
bool Is2NucleonCluster(int pdgc)
Definition PDGUtils.cxx:402
bool IsNeutrino(int pdgc)
Definition PDGUtils.cxx:110
int SwitchProtonNeutron(int pdgc)
Definition PDGUtils.cxx:356
bool IsNucleon(int pdgc)
Definition PDGUtils.cxx:346
int Neutrino2ChargedLepton(int pdgc)
Definition PDGUtils.cxx:218
bool IsPion(int pdgc)
Definition PDGUtils.cxx:326
bool IsAntiNeutrino(int pdgc)
Definition PDGUtils.cxx:118
int ModifyNucleonCluster(int pdgc, int dQ)
Definition PDGUtils.cxx:364
THE MAIN GENIE PROJECT NAMESPACE
Definition AlgCmp.h:25
const int kPdgDarkMatter
Definition PDGCodes.h:218
const int kPdgDarkNeutrino
Definition PDGCodes.h:221
enum genie::EInteractionType InteractionType_t
const int kPdgAntiNuE
Definition PDGCodes.h:29
enum genie::EScatteringType ScatteringType_t
@ kScCoherentProduction
@ kScInverseMuDecay
@ kScDarkMatterDeepInelastic
@ kScQuasiElastic
@ kScInverseBetaDecay
@ kScDarkMatterElastic
@ kScDeepInelastic
@ kScCoherentElastic
@ kScGlashowResonance
@ kScDiffractive
@ kScSingleKaon
ostream & operator<<(ostream &stream, const AlgConfigPool &config_pool)
const int kPdgMuon
Definition PDGCodes.h:37
const int kPdgNuMu
Definition PDGCodes.h:30
const int kPdgGamma
Definition PDGCodes.h:189
const int kPdgAntiDarkMatter
Definition PDGCodes.h:219
const int kPdgElectron
Definition PDGCodes.h:35