Public Member Functions | Static Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes | Friends | List of all members
genie::Interaction Class Reference

Summary information for an interaction. More...

#include <Interaction.h>

Inheritance diagram for genie::Interaction:

Public Member Functions

 Interaction ()
 
 Interaction (const InitialState &init, const ProcessInfo &proc)
 
 Interaction (const Interaction &i)
 
 Interaction (TRootIOCtor *)
 
 ~Interaction ()
 
const InitialStateInitState (void) const
 
const ProcessInfoProcInfo (void) const
 
const KinematicsKine (void) const
 
const XclsTagExclTag (void) const
 
const KPhaseSpacePhaseSpace (void) const
 
InitialStateInitStatePtr (void) const
 
ProcessInfoProcInfoPtr (void) const
 
KinematicsKinePtr (void) const
 
XclsTagExclTagPtr (void) const
 
KPhaseSpacePhaseSpacePtr (void) const
 
void SetInitState (const InitialState &init)
 
void SetProcInfo (const ProcessInfo &proc)
 
void SetKine (const Kinematics &kine)
 
void SetExclTag (const XclsTag &xcls)
 
int FSPrimLeptonPdg (void) const
 final state primary lepton pdg More...
 
int RecoilNucleonPdg (void) const
 recoil nucleon pdg More...
 
TParticlePDG * FSPrimLepton (void) const
 final state primary lepton More...
 
TParticlePDG * RecoilNucleon (void) const
 recoil nucleon More...
 
void Reset (void)
 
void Copy (const Interaction &i)
 
string AsString (void) const
 
void Print (ostream &stream) const
 
Interactionoperator= (const Interaction &i)
 copy More...
 

Static Public Member Functions

static InteractionDISCC (int tgt, int nuc, int probe, double E=0)
 
static InteractionDISCC (int tgt, int nuc, int qrk, bool sea, int probe, double E=0)
 
static InteractionDISCC (int tgt, int nuc, int qrk, bool sea, int fqrk, int probe, double E=0)
 
static InteractionDISCC (int tgt, int nuc, int probe, const TLorentzVector &p4probe)
 
static InteractionDISCC (int tgt, int nuc, int qrk, bool sea, int probe, const TLorentzVector &p4probe)
 
static InteractionDISNC (int tgt, int nuc, int probe, double E=0)
 
static InteractionDISNC (int tgt, int nuc, int qrk, bool sea, int probe, double E=0)
 
static InteractionDISNC (int tgt, int nuc, int qrk, bool sea, int fqrk, int probe, double E=0)
 
static InteractionDISNC (int tgt, int nuc, int probe, const TLorentzVector &p4probe)
 
static InteractionDISNC (int tgt, int nuc, int qrk, bool sea, int probe, const TLorentzVector &p4probe)
 
static InteractionDISEM (int tgt, int nuc, int probe, double E=0)
 
static InteractionDISEM (int tgt, int nuc, int qrk, bool sea, int probe, double E=0)
 
static InteractionDISEM (int tgt, int nuc, int probe, const TLorentzVector &p4probe)
 
static InteractionDISEM (int tgt, int nuc, int qrk, bool sea, int probe, const TLorentzVector &p4probe)
 
static InteractionQELCC (int tgt, int nuc, int probe, double E=0)
 
static InteractionQELCC (int tgt, int nuc, int probe, const TLorentzVector &p4probe)
 
static InteractionQELNC (int tgt, int nuc, int probe, double E=0)
 
static InteractionQELNC (int tgt, int nuc, int probe, const TLorentzVector &p4probe)
 
static InteractionQELEM (int tgt, int nuc, int probe, double E=0)
 
static InteractionQELEM (int tgt, int nuc, int probe, const TLorentzVector &p4probe)
 
static InteractionIBD (int tgt, int nuc, int probe, double E=0)
 
static InteractionIBD (int tgt, int nuc, int probe, const TLorentzVector &p4probe)
 
static InteractionRESCC (int tgt, int nuc, int probe, double E=0)
 
static InteractionRESCC (int tgt, int nuc, int probe, const TLorentzVector &p4probe)
 
static InteractionRESNC (int tgt, int nuc, int probe, double E=0)
 
static InteractionRESNC (int tgt, int nuc, int probe, const TLorentzVector &p4probe)
 
static InteractionRESEM (int tgt, int nuc, int probe, double E=0)
 
static InteractionRESEM (int tgt, int nuc, int probe, const TLorentzVector &p4probe)
 
static InteractionDFRCC (int tgt, int nuc, int probe, double E=0)
 
static InteractionDFRCC (int tgt, int nuc, int probe, const TLorentzVector &p4probe)
 
static InteractionCOHCC (int tgt, int probe, unsigned int prod_pdg, double E=0)
 
static InteractionCOHCC (int tgt, int probe, unsigned int prod_pdg, const TLorentzVector &p4probe)
 
static InteractionCOHNC (int tgt, int probe, unsigned int prod_pdg, double E=0)
 
static InteractionCOHNC (int tgt, int probe, unsigned int prod_pdg, const TLorentzVector &p4probe)
 
static InteractionCEvNS (int tgt, int probe, double E=0)
 
static InteractionCEvNS (int tgt, int probe, const TLorentzVector &p4probe)
 
static InteractionIMD (int tgt, double E=0)
 
static InteractionIMD (int tgt, const TLorentzVector &p4probe)
 
static InteractionAMNuGamma (int tgt, int nuc, int probe, double E=0)
 
static InteractionAMNuGamma (int tgt, int nuc, int probe, const TLorentzVector &p4probe)
 
static InteractionMECCC (int tgt, int nuccluster, int probe, double E=0)
 
static InteractionMECCC (int tgt, int nuccluster, int probe, const TLorentzVector &p4probe)
 
static InteractionMECCC (int tgt, int probe, double E=0)
 
static InteractionMECCC (int tgt, int probe, const TLorentzVector &p4probe)
 
static InteractionMECNC (int tgt, int nuccluster, int probe, double E=0)
 
static InteractionMECNC (int tgt, int nuccluster, int probe, const TLorentzVector &p4probe)
 
static InteractionMECEM (int tgt, int nuccluster, int probe, double E=0)
 
static InteractionMECEM (int tgt, int probe, double E=0)
 
static InteractionMECEM (int tgt, int nuccluster, int probe, const TLorentzVector &p4probe)
 
static InteractionGLR (int tgt, double E=0)
 
static InteractionGLR (int tgt, const TLorentzVector &p4probe)
 
static InteractionNDecay (int tgt, int decay_mode=-1, int decayed_nucleon=0)
 
static InteractionNOsc (int tgt, int annihilation_mode=-1)
 
static InteractionNHL (double E=0, int decayed_mode=-1)
 
static InteractionASK (int tgt, int probe, double E=0)
 
static InteractionASK (int tgt, int probe, const TLorentzVector &p4probe)
 
static InteractionDME (int tgt, int nuc, int probe, double E=0)
 
static InteractionDME (int tgt, int nuc, int probe, const TLorentzVector &p4probe)
 
static InteractionDMDI (int tgt, int nuc, int probe, double E=0)
 
static InteractionDMDI (int tgt, int nuc, int qrk, bool sea, int probe, double E=0)
 
static InteractionDMDI (int tgt, int nuc, int probe, const TLorentzVector &p4probe)
 
static InteractionDMDI (int tgt, int nuc, int qrk, bool sea, int probe, const TLorentzVector &p4probe)
 

Private Member Functions

void Init (void)
 
void CleanUp (void)
 

Static Private Member Functions

static InteractionCreate (int tgt, int probe, ScatteringType_t st, InteractionType_t it)
 

Private Attributes

InitialStatefInitialState
 Initial State info. More...
 
ProcessInfofProcInfo
 Process info (scattering, weak current,...) More...
 
KinematicsfKinematics
 kinematical variables More...
 
XclsTagfExclusiveTag
 Additional info for exclusive channels. More...
 
KPhaseSpacefKinePhSp
 Kinematic phase space. More...
 

Friends

ostream & operator<< (ostream &stream, const Interaction &i)
 print More...
 

Detailed Description

Summary information for an interaction.

It is a container of an InitialState, a ProcessInfo, an XclsTag and a Kinematics object.

Author
Costas Andreopoulos <constantinos.andreopoulos cern.ch> University of Liverpool & STFC Rutherford Appleton Laboratory

Changes required to implement the GENIE Boosted Dark Matter module were installed by Josh Berger (Univ. of Wisconsin)

April 25, 2004

Copyright (c) 2003-2020, The GENIE Collaboration For the full text of the license visit http://copyright.genie-mc.org

Definition at line 56 of file Interaction.h.

Constructor & Destructor Documentation

Interaction::Interaction ( )

Definition at line 45 of file Interaction.cxx.

45  :
46 TObject()
47 {
48  this->Init();
49 }
Interaction::Interaction ( const InitialState init,
const ProcessInfo proc 
)

Definition at line 51 of file Interaction.cxx.

51  :
52 TObject()
53 {
54  this->Init();
55 
56  fInitialState -> Copy (ist);
57  fProcInfo -> Copy (prc);
58 }
InitialState * fInitialState
Initial State info.
Definition: Interaction.h:179
void Copy(const Interaction &i)
ProcessInfo * fProcInfo
Process info (scattering, weak current,...)
Definition: Interaction.h:180
Interaction::Interaction ( const Interaction i)

Definition at line 60 of file Interaction.cxx.

60  :
61 TObject()
62 {
63  this->Init();
64  this->Copy(interaction);
65 }
void Copy(const Interaction &i)
Interaction::Interaction ( TRootIOCtor *  )

Definition at line 67 of file Interaction.cxx.

67  :
68 TObject(),
69 fInitialState(0),
70 fProcInfo(0),
71 fKinematics(0),
72 fExclusiveTag(0),
73 fKinePhSp(0)
74 {
75 
76 }
InitialState * fInitialState
Initial State info.
Definition: Interaction.h:179
XclsTag * fExclusiveTag
Additional info for exclusive channels.
Definition: Interaction.h:182
Kinematics * fKinematics
kinematical variables
Definition: Interaction.h:181
ProcessInfo * fProcInfo
Process info (scattering, weak current,...)
Definition: Interaction.h:180
KPhaseSpace * fKinePhSp
Kinematic phase space.
Definition: Interaction.h:183
Interaction::~Interaction ( )

Definition at line 78 of file Interaction.cxx.

79 {
80  this->CleanUp();
81 }
void CleanUp(void)
Definition: Interaction.cxx:98

Member Function Documentation

Interaction * Interaction::AMNuGamma ( int  tgt,
int  nuc,
int  probe,
double  E = 0 
)
static

Definition at line 819 of file Interaction.cxx.

820 {
823 
824  InitialState * init_state = interaction->InitStatePtr();
825  init_state->SetProbeE(E);
826  init_state->TgtPtr()->SetHitNucPdg(nuc);
827 
828  return interaction;
829 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::AMNuGamma ( int  tgt,
int  nuc,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 831 of file Interaction.cxx.

833 {
836 
837  InitialState * init_state = interaction->InitStatePtr();
838  init_state->SetProbeP4(p4probe);
839  init_state->TgtPtr()->SetHitNucPdg(nuc);
840 
841  return interaction;
842 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::ASK ( int  tgt,
int  probe,
double  E = 0 
)
static

Definition at line 1011 of file Interaction.cxx.

1012 {
1015 
1016  InitialState * init_state = interaction->InitStatePtr();
1017  init_state->SetProbeE(E);
1018 
1019  return interaction;
1020 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::ASK ( int  tgt,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 1022 of file Interaction.cxx.

1024 {
1027 
1028  InitialState * init_state = interaction->InitStatePtr();
1029  init_state->SetProbeP4(p4probe);
1030 
1031  return interaction;
1032 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
string Interaction::AsString ( void  ) const

Definition at line 249 of file Interaction.cxx.

250 {
251 // Code-ify the interaction in a string to be used as (part of a) cache
252 // branch key.
253 // Template:
254 // nu:x;tgt:x;N:x;q:x(s/v);proc:x;xclv_tag
255 
256  const Target & tgt = fInitialState->Tgt();
257 
258  ostringstream interaction;
259 
260  // If the probe has non-zero mass, then it is DM
261  if (fInitialState->Probe()->PdgCode() == kPdgDarkMatter) {
262  interaction << "dm;";
263  }
264  else if (fInitialState->Probe()->PdgCode() == kPdgAntiDarkMatter) {
265  interaction << "dmb;";
266  }
267  else {
268  interaction << "nu:" << fInitialState->ProbePdg() << ";";
269  }
270  interaction << "tgt:" << tgt.Pdg() << ";";
271 
272  if(tgt.HitNucIsSet()) {
273  interaction << "N:" << tgt.HitNucPdg() << ";";
274  }
275  if(tgt.HitQrkIsSet()) {
276  interaction << "q:" << tgt.HitQrkPdg()
277  << (tgt.HitSeaQrk() ? "(s)" : "(v)") << ";";
278  }
279 
280  interaction << "proc:" << fProcInfo->InteractionTypeAsString()
281  << "," << fProcInfo->ScatteringTypeAsString() << ";";
282 
283  string xcls = fExclusiveTag->AsString();
284  interaction << xcls;
285  if(xcls.size()>0) interaction << ";";
286 
287  return interaction.str();
288 }
bool HitSeaQrk(void) const
Definition: Target.cxx:299
string ScatteringTypeAsString(void) const
int HitNucPdg(void) const
Definition: Target.cxx:304
int HitQrkPdg(void) const
Definition: Target.cxx:242
InitialState * fInitialState
Initial State info.
Definition: Interaction.h:179
const int kPdgDarkMatter
Definition: PDGCodes.h:218
int Pdg(void) const
Definition: Target.h:71
TParticlePDG * Probe(void) const
XclsTag * fExclusiveTag
Additional info for exclusive channels.
Definition: Interaction.h:182
const int kPdgAntiDarkMatter
Definition: PDGCodes.h:219
string AsString(void) const
pack into a string code
Definition: XclsTag.cxx:212
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
int ProbePdg(void) const
Definition: InitialState.h:64
string InteractionTypeAsString(void) const
bool HitNucIsSet(void) const
Definition: Target.cxx:283
bool HitQrkIsSet(void) const
Definition: Target.cxx:292
ProcessInfo * fProcInfo
Process info (scattering, weak current,...)
Definition: Interaction.h:180
const Target & Tgt(void) const
Definition: InitialState.h:66
Interaction * Interaction::CEvNS ( int  tgt,
int  probe,
double  E = 0 
)
static

Definition at line 774 of file Interaction.cxx.

775 {
778 
779  InitialState * init_state = interaction->InitStatePtr();
780  init_state->SetProbeE(E);
781 
782  return interaction;
783 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::CEvNS ( int  tgt,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 785 of file Interaction.cxx.

787 {
790 
791  InitialState * init_state = interaction->InitStatePtr();
792  init_state->SetProbeP4(p4probe);
793 
794  return interaction;
795 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
void Interaction::CleanUp ( void  )
private

Definition at line 98 of file Interaction.cxx.

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 }
InitialState * fInitialState
Initial State info.
Definition: Interaction.h:179
XclsTag * fExclusiveTag
Additional info for exclusive channels.
Definition: Interaction.h:182
Kinematics * fKinematics
kinematical variables
Definition: Interaction.h:181
ProcessInfo * fProcInfo
Process info (scattering, weak current,...)
Definition: Interaction.h:180
KPhaseSpace * fKinePhSp
Kinematic phase space.
Definition: Interaction.h:183
Interaction * Interaction::COHCC ( int  tgt,
int  probe,
unsigned int  prod_pdg,
double  E = 0 
)
static

Definition at line 708 of file Interaction.cxx.

709 {
712 
713  XclsTag * xcl = interaction -> ExclTagPtr() ;
714  if ( pdg::IsPion( prod_pdg ) ) {
715  if ( pdg::IsNeutrino( probe ) ) xcl -> SetNPions( 1,0,0 ) ;
716  else if ( pdg::IsAntiNeutrino( probe ) ) xcl -> SetNPions( 0,0,1 ) ;
717  }
718 
719  InitialState * init_state = interaction->InitStatePtr();
720  init_state->SetProbeE(E);
721 
722  return interaction;
723 }
bool IsPion(int pdgc)
Definition: PDGUtils.cxx:323
bool IsNeutrino(int pdgc)
Definition: PDGUtils.cxx:107
void SetProbeE(double E)
Contains minimal information for tagging exclusive processes.
Definition: XclsTag.h:39
Summary information for an interaction.
Definition: Interaction.h:56
bool IsAntiNeutrino(int pdgc)
Definition: PDGUtils.cxx:115
XclsTag * ExclTagPtr(void) const
Definition: Interaction.h:77
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::COHCC ( int  tgt,
int  probe,
unsigned int  prod_pdg,
const TLorentzVector &  p4probe 
)
static

Definition at line 725 of file Interaction.cxx.

727 {
730 
731  XclsTag * xcl = interaction -> ExclTagPtr() ;
732  if ( pdg::IsPion( prod_pdg ) ) {
733  if ( pdg::IsNeutrino( probe ) ) xcl -> SetNPions( 1,0,0 ) ;
734  else if ( pdg::IsAntiNeutrino( probe ) ) xcl -> SetNPions( 0,0,1 ) ;
735  }
736 
737  InitialState * init_state = interaction->InitStatePtr();
738  init_state->SetProbeP4(p4probe);
739 
740  return interaction;
741 }
bool IsPion(int pdgc)
Definition: PDGUtils.cxx:323
bool IsNeutrino(int pdgc)
Definition: PDGUtils.cxx:107
void SetProbeP4(const TLorentzVector &P4)
Contains minimal information for tagging exclusive processes.
Definition: XclsTag.h:39
Summary information for an interaction.
Definition: Interaction.h:56
bool IsAntiNeutrino(int pdgc)
Definition: PDGUtils.cxx:115
XclsTag * ExclTagPtr(void) const
Definition: Interaction.h:77
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::COHNC ( int  tgt,
int  probe,
unsigned int  prod_pdg,
double  E = 0 
)
static

Definition at line 743 of file Interaction.cxx.

744 {
747 
748  XclsTag * xcl = interaction -> ExclTagPtr() ;
749  if ( pdg::IsPion( prod_pdg ) ) xcl -> SetNPions( 0,1,0 ) ;
750  else if ( prod_pdg == kPdgGamma ) xcl -> SetNSingleGammas(1) ;
751 
752  InitialState * init_state = interaction->InitStatePtr();
753  init_state->SetProbeE(E);
754 
755  return interaction;
756 }
bool IsPion(int pdgc)
Definition: PDGUtils.cxx:323
void SetProbeE(double E)
Contains minimal information for tagging exclusive processes.
Definition: XclsTag.h:39
Summary information for an interaction.
Definition: Interaction.h:56
const int kPdgGamma
Definition: PDGCodes.h:189
XclsTag * ExclTagPtr(void) const
Definition: Interaction.h:77
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::COHNC ( int  tgt,
int  probe,
unsigned int  prod_pdg,
const TLorentzVector &  p4probe 
)
static

Definition at line 758 of file Interaction.cxx.

760 {
763 
764  XclsTag * xcl = interaction -> ExclTagPtr() ;
765  if ( pdg::IsPion( prod_pdg ) ) xcl -> SetNPions( 0,1,0 ) ;
766  else if ( prod_pdg == kPdgGamma ) xcl -> SetNSingleGammas(1) ;
767 
768  InitialState * init_state = interaction->InitStatePtr();
769  init_state->SetProbeP4(p4probe);
770 
771  return interaction;
772 }
bool IsPion(int pdgc)
Definition: PDGUtils.cxx:323
void SetProbeP4(const TLorentzVector &P4)
Contains minimal information for tagging exclusive processes.
Definition: XclsTag.h:39
Summary information for an interaction.
Definition: Interaction.h:56
const int kPdgGamma
Definition: PDGCodes.h:189
XclsTag * ExclTagPtr(void) const
Definition: Interaction.h:77
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
void Interaction::Copy ( const Interaction i)

Definition at line 113 of file Interaction.cxx.

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 }
InitialState * fInitialState
Initial State info.
Definition: Interaction.h:179
Generated/set kinematical variables for an event.
Definition: Kinematics.h:39
init
Definition: train.py:42
XclsTag * fExclusiveTag
Additional info for exclusive channels.
Definition: Interaction.h:182
Contains minimal information for tagging exclusive processes.
Definition: XclsTag.h:39
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
void Copy(const Interaction &i)
Kinematics * fKinematics
kinematical variables
Definition: Interaction.h:181
ProcessInfo * fProcInfo
Process info (scattering, weak current,...)
Definition: Interaction.h:180
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::Create ( int  tgt,
int  probe,
ScatteringType_t  st,
InteractionType_t  it 
)
staticprivate

Definition at line 318 of file Interaction.cxx.

320 {
321  InitialState init_state (target, probe);
322  ProcessInfo proc_info (st, it);
323 
324  Interaction * interaction = new Interaction(init_state, proc_info);
325  return interaction;
326 }
Summary information for an interaction.
Definition: Interaction.h:56
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::DFRCC ( int  tgt,
int  nuc,
int  probe,
double  E = 0 
)
static

Definition at line 683 of file Interaction.cxx.

684 {
687 
688  InitialState * init_state = interaction->InitStatePtr();
689  init_state->SetProbeE(E);
690  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
691 
692  return interaction;
693 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::DFRCC ( int  tgt,
int  nuc,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 695 of file Interaction.cxx.

697 {
700 
701  InitialState * init_state = interaction->InitStatePtr();
702  init_state->SetProbeP4(p4probe);
703  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
704 
705  return interaction;
706 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::DISCC ( int  tgt,
int  nuc,
int  probe,
double  E = 0 
)
static

Definition at line 328 of file Interaction.cxx.

329 {
332 
333  InitialState * init_state = interaction->InitStatePtr();
334  init_state->SetProbeE(E);
335  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
336 
337  return interaction;
338 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::DISCC ( int  tgt,
int  nuc,
int  qrk,
bool  sea,
int  probe,
double  E = 0 
)
static

Definition at line 340 of file Interaction.cxx.

342 {
344 
345  Target * tgt = interaction->InitStatePtr()->TgtPtr();
346  tgt -> SetHitQrkPdg (hitqrk);
347  tgt -> SetHitSeaQrk (fromsea);
348 
349  return interaction;
350 }
Summary information for an interaction.
Definition: Interaction.h:56
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
static Interaction * DISCC(int tgt, int nuc, int probe, double E=0)
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
Interaction * Interaction::DISCC ( int  tgt,
int  nuc,
int  qrk,
bool  sea,
int  fqrk,
int  probe,
double  E = 0 
)
static

Definition at line 352 of file Interaction.cxx.

354 {
356 
357  Target * tgt = interaction->InitStatePtr()->TgtPtr();
358  tgt -> SetHitQrkPdg (hitqrk);
359  tgt -> SetHitSeaQrk (fromsea);
360 
361  XclsTag * xclstag = interaction->ExclTagPtr();
362  xclstag->SetFinalQuark(fqrk);
363 
364  return interaction;
365 }
Contains minimal information for tagging exclusive processes.
Definition: XclsTag.h:39
Summary information for an interaction.
Definition: Interaction.h:56
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
XclsTag * ExclTagPtr(void) const
Definition: Interaction.h:77
static Interaction * DISCC(int tgt, int nuc, int probe, double E=0)
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
void SetFinalQuark(int finalquark_pdgc=0)
Definition: XclsTag.cxx:138
Interaction * Interaction::DISCC ( int  tgt,
int  nuc,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 367 of file Interaction.cxx.

369 {
372 
373  InitialState * init_state = interaction->InitStatePtr();
374  init_state->SetProbeP4(p4probe);
375  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
376 
377  return interaction;
378 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::DISCC ( int  tgt,
int  nuc,
int  qrk,
bool  sea,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 380 of file Interaction.cxx.

383 {
384  Interaction* interaction = Interaction::DISCC(target,hitnuc,probe,p4probe);
385 
386  Target * tgt = interaction->InitStatePtr()->TgtPtr();
387  tgt -> SetHitQrkPdg (hitqrk);
388  tgt -> SetHitSeaQrk (fromsea);
389 
390  return interaction;
391 }
Summary information for an interaction.
Definition: Interaction.h:56
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
static Interaction * DISCC(int tgt, int nuc, int probe, double E=0)
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
Interaction * Interaction::DISEM ( int  tgt,
int  nuc,
int  probe,
double  E = 0 
)
static

Definition at line 458 of file Interaction.cxx.

459 {
462 
463  InitialState * init_state = interaction->InitStatePtr();
464  init_state->SetProbeE(E);
465  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
466 
467  return interaction;
468 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::DISEM ( int  tgt,
int  nuc,
int  qrk,
bool  sea,
int  probe,
double  E = 0 
)
static

Definition at line 470 of file Interaction.cxx.

472 {
474 
475  Target * tgt = interaction->InitStatePtr()->TgtPtr();
476  tgt -> SetHitQrkPdg (hitqrk);
477  tgt -> SetHitSeaQrk (fromsea);
478 
479  return interaction;
480 }
Summary information for an interaction.
Definition: Interaction.h:56
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * DISEM(int tgt, int nuc, int probe, double E=0)
Interaction * Interaction::DISEM ( int  tgt,
int  nuc,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 482 of file Interaction.cxx.

484 {
487 
488  InitialState * init_state = interaction->InitStatePtr();
489  init_state->SetProbeP4(p4probe);
490  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
491 
492  return interaction;
493 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::DISEM ( int  tgt,
int  nuc,
int  qrk,
bool  sea,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 495 of file Interaction.cxx.

498 {
499  Interaction * interaction = Interaction::DISEM(target,hitnuc,probe,p4probe);
500 
501  Target * tgt = interaction->InitStatePtr()->TgtPtr();
502  tgt -> SetHitQrkPdg (hitqrk);
503  tgt -> SetHitSeaQrk (fromsea);
504 
505  return interaction;
506 }
Summary information for an interaction.
Definition: Interaction.h:56
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * DISEM(int tgt, int nuc, int probe, double E=0)
Interaction * Interaction::DISNC ( int  tgt,
int  nuc,
int  probe,
double  E = 0 
)
static

Definition at line 393 of file Interaction.cxx.

394 {
397 
398  InitialState * init_state = interaction->InitStatePtr();
399  init_state->SetProbeE(E);
400  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
401 
402  return interaction;
403 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::DISNC ( int  tgt,
int  nuc,
int  qrk,
bool  sea,
int  probe,
double  E = 0 
)
static

Definition at line 405 of file Interaction.cxx.

407 {
409 
410  Target * tgt = interaction->InitStatePtr()->TgtPtr();
411  tgt -> SetHitQrkPdg (hitqrk);
412  tgt -> SetHitSeaQrk (fromsea);
413 
414  return interaction;
415 }
Summary information for an interaction.
Definition: Interaction.h:56
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * DISNC(int tgt, int nuc, int probe, double E=0)
Interaction * Interaction::DISNC ( int  tgt,
int  nuc,
int  qrk,
bool  sea,
int  fqrk,
int  probe,
double  E = 0 
)
static

Definition at line 417 of file Interaction.cxx.

419 {
421 
422  Target * tgt = interaction->InitStatePtr()->TgtPtr();
423  tgt -> SetHitQrkPdg (hitqrk);
424  tgt -> SetHitSeaQrk (fromsea);
425 
426  XclsTag * xclstag = interaction->ExclTagPtr();
427  xclstag->SetFinalQuark(fqrk);
428 
429  return interaction;
430 }
Contains minimal information for tagging exclusive processes.
Definition: XclsTag.h:39
Summary information for an interaction.
Definition: Interaction.h:56
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
XclsTag * ExclTagPtr(void) const
Definition: Interaction.h:77
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
void SetFinalQuark(int finalquark_pdgc=0)
Definition: XclsTag.cxx:138
static Interaction * DISNC(int tgt, int nuc, int probe, double E=0)
Interaction * Interaction::DISNC ( int  tgt,
int  nuc,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 432 of file Interaction.cxx.

434 {
437 
438  InitialState * init_state = interaction->InitStatePtr();
439  init_state->SetProbeP4(p4probe);
440  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
441 
442  return interaction;
443 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::DISNC ( int  tgt,
int  nuc,
int  qrk,
bool  sea,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 445 of file Interaction.cxx.

448 {
449  Interaction * interaction = Interaction::DISNC(target,hitnuc,probe,p4probe);
450 
451  Target * tgt = interaction->InitStatePtr()->TgtPtr();
452  tgt -> SetHitQrkPdg (hitqrk);
453  tgt -> SetHitSeaQrk (fromsea);
454 
455  return interaction;
456 }
Summary information for an interaction.
Definition: Interaction.h:56
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * DISNC(int tgt, int nuc, int probe, double E=0)
Interaction * Interaction::DMDI ( int  tgt,
int  nuc,
int  probe,
double  E = 0 
)
static

Definition at line 1061 of file Interaction.cxx.

1062 {
1065 
1066  InitialState * init_state = interaction->InitStatePtr();
1067  init_state->SetProbeE(E);
1068  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
1069 
1070  return interaction;
1071 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::DMDI ( int  tgt,
int  nuc,
int  qrk,
bool  sea,
int  probe,
double  E = 0 
)
static

Definition at line 1073 of file Interaction.cxx.

1075 {
1077 
1078  Target * tgt = interaction->InitStatePtr()->TgtPtr();
1079  tgt -> SetHitQrkPdg (hitqrk);
1080  tgt -> SetHitSeaQrk (fromsea);
1081 
1082  return interaction;
1083 }
Summary information for an interaction.
Definition: Interaction.h:56
static Interaction * DMDI(int tgt, int nuc, int probe, double E=0)
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
Interaction * Interaction::DMDI ( int  tgt,
int  nuc,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 1085 of file Interaction.cxx.

1087 {
1090 
1091  InitialState * init_state = interaction->InitStatePtr();
1092  init_state->SetProbeP4(p4probe);
1093  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
1094 
1095  return interaction;
1096 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::DMDI ( int  tgt,
int  nuc,
int  qrk,
bool  sea,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 1098 of file Interaction.cxx.

1101 {
1102  Interaction * interaction = Interaction::DMDI(target,hitnuc,probe,p4probe);
1103 
1104  Target * tgt = interaction->InitStatePtr()->TgtPtr();
1105  tgt -> SetHitQrkPdg (hitqrk);
1106  tgt -> SetHitSeaQrk (fromsea);
1107 
1108  return interaction;
1109 }
Summary information for an interaction.
Definition: Interaction.h:56
static Interaction * DMDI(int tgt, int nuc, int probe, double E=0)
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
Interaction * Interaction::DME ( int  tgt,
int  nuc,
int  probe,
double  E = 0 
)
static

Definition at line 1034 of file Interaction.cxx.

1035 {
1036  // EDIT: need to be able to create dark matter elastic
1039 
1040  InitialState * init_state = interaction->InitStatePtr();
1041  init_state->SetProbeE(E);
1042  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
1043 
1044  return interaction;
1045 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::DME ( int  tgt,
int  nuc,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 1047 of file Interaction.cxx.

1049 {
1050  // EDIT: need to be able to create dark matter elastic
1053 
1054  InitialState * init_state = interaction->InitStatePtr();
1055  init_state->SetProbeP4(p4probe);
1056  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
1057 
1058  return interaction;
1059 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
const XclsTag& genie::Interaction::ExclTag ( void  ) const
inline

Definition at line 72 of file Interaction.h.

72 { return *fExclusiveTag; }
XclsTag * fExclusiveTag
Additional info for exclusive channels.
Definition: Interaction.h:182
XclsTag* genie::Interaction::ExclTagPtr ( void  ) const
inline

Definition at line 77 of file Interaction.h.

77 { return fExclusiveTag; }
XclsTag * fExclusiveTag
Additional info for exclusive channels.
Definition: Interaction.h:182
TParticlePDG * Interaction::FSPrimLepton ( void  ) const

final state primary lepton

Definition at line 126 of file Interaction.cxx.

127 {
128  int pdgc = this->FSPrimLeptonPdg();
129 
130  if(pdgc) return PDGLibrary::Instance()->Find(pdgc);
131  else return 0;
132 }
int FSPrimLeptonPdg(void) const
final state primary lepton pdg
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:57
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:75
int Interaction::FSPrimLeptonPdg ( void  ) const

final state primary lepton pdg

Definition at line 134 of file Interaction.cxx.

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  // vN (Weak-NC) or eN (EM)
148  if (proc_info.IsWeakNC() || proc_info.IsEM() || proc_info.IsWeakMix() || proc_info.IsDarkMatter()) return pdgc; // EDIT: DM does not change in FS
149 
150  // vN (Weak-CC)
151  else if (proc_info.IsWeakCC()) {
152  int clpdgc;
153  if (proc_info.IsIMDAnnihilation())
154  clpdgc = kPdgMuon;
155  else if (proc_info.IsGlashowResonance()) {
156  if ( pdg::IsMuon(xclstag.FinalLeptonPdg()) ) clpdgc = kPdgMuon;
157  else if ( pdg::IsTau(xclstag.FinalLeptonPdg()) ) clpdgc = kPdgTau;
158  else if ( pdg::IsElectron(xclstag.FinalLeptonPdg()) ) clpdgc = kPdgElectron;
159  else if ( pdg::IsPion(xclstag.FinalLeptonPdg()) ) clpdgc = kPdgPiP;
160  }
161  else
162  clpdgc = pdg::Neutrino2ChargedLepton(pdgc);
163  return clpdgc;
164  }
165  else if (proc_info.IsDarkNeutralCurrent()){
166  return kPdgDarkNeutrino;
167  }
168  LOG("Interaction", pWARN)
169  << "Could not figure out the final state primary lepton pdg code!!";
170 
171  return 0;
172 }
bool IsPion(int pdgc)
Definition: PDGUtils.cxx:323
bool IsWeakMix(void) const
bool IsWeakCC(void) const
int FinalLeptonPdg(void) const
Definition: XclsTag.h:74
bool IsDarkNeutralCurrent(void) const
bool IsIMDAnnihilation(void) const
const int kPdgElectron
Definition: PDGCodes.h:35
Contains minimal information for tagging exclusive processes.
Definition: XclsTag.h:39
bool IsWeakNC(void) const
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
const int kPdgTau
Definition: PDGCodes.h:39
bool IsNuElectronElastic(void) const
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
int ProbePdg(void) const
Definition: InitialState.h:64
bool IsTau(int pdgc)
Definition: PDGUtils.cxx:205
const int kPdgPiP
Definition: PDGCodes.h:158
bool IsMuon(int pdgc)
Definition: PDGUtils.cxx:195
#define pWARN
Definition: Messenger.h:60
bool IsEM(void) const
bool IsDarkMatter(void) const
const XclsTag & ExclTag(void) const
Definition: Interaction.h:72
int Neutrino2ChargedLepton(int pdgc)
Definition: PDGUtils.cxx:215
const InitialState & InitState(void) const
Definition: Interaction.h:69
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
const int kPdgMuon
Definition: PDGCodes.h:37
bool IsGlashowResonance(void) const
const int kPdgDarkNeutrino
Definition: PDGCodes.h:222
bool IsElectron(int pdgc)
Definition: PDGUtils.cxx:185
Initial State information.
Definition: InitialState.h:48
#define pDEBUG
Definition: Messenger.h:63
Interaction * Interaction::GLR ( int  tgt,
double  E = 0 
)
static

Definition at line 955 of file Interaction.cxx.

956 {
959 
960  InitialState * init_state = interaction->InitStatePtr();
961  init_state->SetProbeE(E);
962  init_state->TgtPtr()->SetHitNucPdg(0);
963 
964  return interaction;
965 }
const int kPdgAntiNuE
Definition: PDGCodes.h:29
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::GLR ( int  tgt,
const TLorentzVector &  p4probe 
)
static

Definition at line 967 of file Interaction.cxx.

968 {
971 
972  InitialState * init_state = interaction->InitStatePtr();
973  init_state->SetProbeP4(p4probe);
974  init_state->TgtPtr()->SetHitNucPdg(0);
975 
976  return interaction;
977 }
void SetProbeP4(const TLorentzVector &P4)
const int kPdgAntiNuE
Definition: PDGCodes.h:29
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::IBD ( int  tgt,
int  nuc,
int  probe,
double  E = 0 
)
static

Definition at line 583 of file Interaction.cxx.

584 {
587 
588  InitialState * init_state = interaction->InitStatePtr();
589  init_state->SetProbeE(E);
590  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
591 
592  return interaction;
593 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::IBD ( int  tgt,
int  nuc,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 595 of file Interaction.cxx.

597 {
600 
601  InitialState * init_state = interaction->InitStatePtr();
602  init_state->SetProbeP4(p4probe);
603  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
604 
605  return interaction;
606 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::IMD ( int  tgt,
double  E = 0 
)
static

Definition at line 797 of file Interaction.cxx.

798 {
801 
802  InitialState * init_state = interaction->InitStatePtr();
803  init_state->SetProbeE(E);
804 
805  return interaction;
806 }
const int kPdgNuMu
Definition: PDGCodes.h:30
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::IMD ( int  tgt,
const TLorentzVector &  p4probe 
)
static

Definition at line 808 of file Interaction.cxx.

809 {
812 
813  InitialState * init_state = interaction->InitStatePtr();
814  init_state->SetProbeP4(p4probe);
815 
816  return interaction;
817 }
void SetProbeP4(const TLorentzVector &P4)
const int kPdgNuMu
Definition: PDGCodes.h:30
Summary information for an interaction.
Definition: Interaction.h:56
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
void Interaction::Init ( void  )
private

Definition at line 89 of file Interaction.cxx.

90 {
91  fInitialState = new InitialState ();
92  fProcInfo = new ProcessInfo ();
93  fKinematics = new Kinematics ();
94  fExclusiveTag = new XclsTag ();
95  fKinePhSp = new KPhaseSpace (this);
96 }
InitialState * fInitialState
Initial State info.
Definition: Interaction.h:179
Generated/set kinematical variables for an event.
Definition: Kinematics.h:39
XclsTag * fExclusiveTag
Additional info for exclusive channels.
Definition: Interaction.h:182
Contains minimal information for tagging exclusive processes.
Definition: XclsTag.h:39
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
Kinematical phase space.
Definition: KPhaseSpace.h:33
Kinematics * fKinematics
kinematical variables
Definition: Interaction.h:181
ProcessInfo * fProcInfo
Process info (scattering, weak current,...)
Definition: Interaction.h:180
KPhaseSpace * fKinePhSp
Kinematic phase space.
Definition: Interaction.h:183
Initial State information.
Definition: InitialState.h:48
const InitialState& genie::Interaction::InitState ( void  ) const
inline

Definition at line 69 of file Interaction.h.

69 { return *fInitialState; }
InitialState * fInitialState
Initial State info.
Definition: Interaction.h:179
InitialState* genie::Interaction::InitStatePtr ( void  ) const
inline

Definition at line 74 of file Interaction.h.

74 { return fInitialState; }
InitialState * fInitialState
Initial State info.
Definition: Interaction.h:179
const Kinematics& genie::Interaction::Kine ( void  ) const
inline

Definition at line 71 of file Interaction.h.

71 { return *fKinematics; }
Kinematics * fKinematics
kinematical variables
Definition: Interaction.h:181
Kinematics* genie::Interaction::KinePtr ( void  ) const
inline

Definition at line 76 of file Interaction.h.

76 { return fKinematics; }
Kinematics * fKinematics
kinematical variables
Definition: Interaction.h:181
Interaction * Interaction::MECCC ( int  tgt,
int  nuccluster,
int  probe,
double  E = 0 
)
static

Definition at line 844 of file Interaction.cxx.

845 {
847  Interaction::Create(tgt, probe, kScMEC, kIntWeakCC);
848 
849  InitialState * init_state = interaction->InitStatePtr();
850  init_state->SetProbeE(E);
851  init_state->TgtPtr()->SetHitNucPdg(ncluster);
852 
853  return interaction;
854 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::MECCC ( int  tgt,
int  nuccluster,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 856 of file Interaction.cxx.

858 {
860  Interaction::Create(tgt, probe, kScMEC, kIntWeakCC);
861 
862  InitialState * init_state = interaction->InitStatePtr();
863  init_state->SetProbeP4(p4probe);
864  init_state->TgtPtr()->SetHitNucPdg(ncluster);
865 
866  return interaction;
867 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::MECCC ( int  tgt,
int  probe,
double  E = 0 
)
static

Definition at line 869 of file Interaction.cxx.

870 {
872  Interaction::Create(tgt, probe, kScMEC, kIntWeakCC);
873 
874  InitialState * init_state = interaction->InitStatePtr();
875  init_state->SetProbeE(E);
876 
877  return interaction;
878 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::MECCC ( int  tgt,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 880 of file Interaction.cxx.

882 {
884  Interaction::Create(tgt, probe, kScMEC, kIntWeakCC);
885 
886  InitialState * init_state = interaction->InitStatePtr();
887  init_state->SetProbeP4(p4probe);
888 
889  return interaction;
890 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::MECEM ( int  tgt,
int  nuccluster,
int  probe,
double  E = 0 
)
static

Definition at line 930 of file Interaction.cxx.

931 {
933  Interaction::Create(tgt, probe, kScMEC, kIntEM);
934 
935  InitialState * init_state = interaction->InitStatePtr();
936  init_state->SetProbeE(E);
937  init_state->TgtPtr()->SetHitNucPdg(ncluster);
938 
939  return interaction;
940 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::MECEM ( int  tgt,
int  probe,
double  E = 0 
)
static

Definition at line 918 of file Interaction.cxx.

919 {
920 
922  Interaction::Create(tgt, probe, kScMEC, kIntEM);
923 
924  InitialState * init_state = interaction->InitStatePtr();
925  init_state->SetProbeE(E);
926 
927  return interaction;
928 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::MECEM ( int  tgt,
int  nuccluster,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 942 of file Interaction.cxx.

944 {
946  Interaction::Create(tgt, probe, kScMEC, kIntEM);
947 
948  InitialState * init_state = interaction->InitStatePtr();
949  init_state->SetProbeP4(p4probe);
950  init_state->TgtPtr()->SetHitNucPdg(ncluster);
951 
952  return interaction;
953 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::MECNC ( int  tgt,
int  nuccluster,
int  probe,
double  E = 0 
)
static

Definition at line 893 of file Interaction.cxx.

894 {
896  Interaction::Create(tgt, probe, kScMEC, kIntWeakNC);
897 
898  InitialState * init_state = interaction->InitStatePtr();
899  init_state->SetProbeE(E);
900  init_state->TgtPtr()->SetHitNucPdg(ncluster);
901 
902  return interaction;
903 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::MECNC ( int  tgt,
int  nuccluster,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 905 of file Interaction.cxx.

907 {
909  Interaction::Create(tgt, probe, kScMEC, kIntWeakNC);
910 
911  InitialState * init_state = interaction->InitStatePtr();
912  init_state->SetProbeP4(p4probe);
913  init_state->TgtPtr()->SetHitNucPdg(ncluster);
914 
915  return interaction;
916 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::NDecay ( int  tgt,
int  decay_mode = -1,
int  decayed_nucleon = 0 
)
static

Definition at line 979 of file Interaction.cxx.

980 {
983  interaction->ExclTagPtr()->SetDecayMode(decay_mode);
984 
985  InitialState * init_state = interaction->InitStatePtr();
986  init_state->TgtPtr()->SetHitNucPdg(decayed_nucleon);
987 
988  return interaction;
989 }
Summary information for an interaction.
Definition: Interaction.h:56
void SetDecayMode(int decay_mode)
Definition: XclsTag.cxx:133
XclsTag * ExclTagPtr(void) const
Definition: Interaction.h:77
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::NHL ( double  E = 0,
int  decayed_mode = -1 
)
static

Definition at line 999 of file Interaction.cxx.

1000 {
1003  interaction->ExclTagPtr()->SetDecayMode(decayed_mode);
1004 
1005  InitialState * init_state = interaction->InitStatePtr();
1006  init_state->SetProbeE(E);
1007 
1008  return interaction;
1009 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
void SetDecayMode(int decay_mode)
Definition: XclsTag.cxx:133
XclsTag * ExclTagPtr(void) const
Definition: Interaction.h:77
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::NOsc ( int  tgt,
int  annihilation_mode = -1 
)
static

Definition at line 991 of file Interaction.cxx.

992 {
995  interaction->ExclTagPtr()->SetDecayMode(annihilation_mode);
996  return interaction;
997 }
Summary information for an interaction.
Definition: Interaction.h:56
void SetDecayMode(int decay_mode)
Definition: XclsTag.cxx:133
XclsTag * ExclTagPtr(void) const
Definition: Interaction.h:77
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Interaction & Interaction::operator= ( const Interaction i)

copy

Definition at line 308 of file Interaction.cxx.

309 {
310  this->Copy(interaction);
311  return (*this);
312 }
void Copy(const Interaction &i)
const KPhaseSpace& genie::Interaction::PhaseSpace ( void  ) const
inline

Definition at line 73 of file Interaction.h.

73 { return *fKinePhSp; }
KPhaseSpace * fKinePhSp
Kinematic phase space.
Definition: Interaction.h:183
KPhaseSpace* genie::Interaction::PhaseSpacePtr ( void  ) const
inline

Definition at line 78 of file Interaction.h.

78 { return fKinePhSp; }
KPhaseSpace * fKinePhSp
Kinematic phase space.
Definition: Interaction.h:183
void Interaction::Print ( ostream &  stream) const

Definition at line 290 of file Interaction.cxx.

291 {
292  const string line(110, '-');
293 
294  stream << endl;
295  stream << line << endl;
296 
297  stream << "GENIE Interaction Summary" << endl;
298  stream << line << endl;
299 
300  stream << *fInitialState << endl; // print initial state
301  stream << *fProcInfo; // print process info
302  stream << *fKinematics; // print scattering parameters
303  stream << *fExclusiveTag; // print exclusive process tag
304 
305  stream << line << endl;
306 }
InitialState * fInitialState
Initial State info.
Definition: Interaction.h:179
XclsTag * fExclusiveTag
Additional info for exclusive channels.
Definition: Interaction.h:182
Kinematics * fKinematics
kinematical variables
Definition: Interaction.h:181
void line(double t, double *p, double &x, double &y, double &z)
ProcessInfo * fProcInfo
Process info (scattering, weak current,...)
Definition: Interaction.h:180
QTextStream & endl(QTextStream &s)
const ProcessInfo& genie::Interaction::ProcInfo ( void  ) const
inline

Definition at line 70 of file Interaction.h.

70 { return *fProcInfo; }
ProcessInfo * fProcInfo
Process info (scattering, weak current,...)
Definition: Interaction.h:180
ProcessInfo* genie::Interaction::ProcInfoPtr ( void  ) const
inline

Definition at line 75 of file Interaction.h.

75 { return fProcInfo; }
ProcessInfo * fProcInfo
Process info (scattering, weak current,...)
Definition: Interaction.h:180
Interaction * Interaction::QELCC ( int  tgt,
int  nuc,
int  probe,
double  E = 0 
)
static

Definition at line 508 of file Interaction.cxx.

509 {
512 
513  InitialState * init_state = interaction->InitStatePtr();
514  init_state->SetProbeE(E);
515  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
516 
517  return interaction;
518 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::QELCC ( int  tgt,
int  nuc,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 520 of file Interaction.cxx.

522 {
525 
526  InitialState * init_state = interaction->InitStatePtr();
527  init_state->SetProbeP4(p4probe);
528  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
529 
530  return interaction;
531 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::QELEM ( int  tgt,
int  nuc,
int  probe,
double  E = 0 
)
static

Definition at line 558 of file Interaction.cxx.

559 {
562 
563  InitialState * init_state = interaction->InitStatePtr();
564  init_state->SetProbeE(E);
565  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
566 
567  return interaction;
568 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::QELEM ( int  tgt,
int  nuc,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 570 of file Interaction.cxx.

572 {
575 
576  InitialState * init_state = interaction->InitStatePtr();
577  init_state->SetProbeP4(p4probe);
578  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
579 
580  return interaction;
581 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::QELNC ( int  tgt,
int  nuc,
int  probe,
double  E = 0 
)
static

Definition at line 533 of file Interaction.cxx.

534 {
537 
538  InitialState * init_state = interaction->InitStatePtr();
539  init_state->SetProbeE(E);
540  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
541 
542  return interaction;
543 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::QELNC ( int  tgt,
int  nuc,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 545 of file Interaction.cxx.

547 {
550 
551  InitialState * init_state = interaction->InitStatePtr();
552  init_state->SetProbeP4(p4probe);
553  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
554 
555  return interaction;
556 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
TParticlePDG * Interaction::RecoilNucleon ( void  ) const

recoil nucleon

Definition at line 174 of file Interaction.cxx.

175 {
176  int rnuc = this->RecoilNucleonPdg();
177 
178  if(rnuc) return PDGLibrary::Instance()->Find(rnuc);
179  else return 0;
180 }
int RecoilNucleonPdg(void) const
recoil nucleon pdg
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:57
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:75
int Interaction::RecoilNucleonPdg ( void  ) const

recoil nucleon pdg

Definition at line 182 of file Interaction.cxx.

183 {
184 // Determine the recoil nucleon PDG code
185 
186  const Target & target = fInitialState->Tgt();
187 
188  int recoil_nuc = 0;
189  int struck_nuc = target.HitNucPdg();
190 
192  bool struck_is_nuc = pdg::IsNucleon(struck_nuc);
193  bool is_weak = fProcInfo->IsWeak();
194  bool is_em = fProcInfo->IsEM();
195  bool is_dm = fProcInfo->IsDarkMatter();
196  assert(struck_is_nuc && (is_weak || is_em || is_dm));
197  if(fProcInfo->IsWeakCC()) {
198  recoil_nuc = pdg::SwitchProtonNeutron(struck_nuc); // CC
199  } else {
200  recoil_nuc = struck_nuc; // NC, EM
201  }
202  }
203 
204  if (fProcInfo->IsMEC()) {
205  bool struck_is_2nuc_cluster = pdg::Is2NucleonCluster(struck_nuc);
206  bool is_weak = fProcInfo->IsWeak();
207  bool is_em = fProcInfo->IsEM();
208  assert(struck_is_2nuc_cluster && (is_weak || is_em));
209  if(fProcInfo->IsWeakCC()) {
210  bool isnu = pdg::IsNeutrino(fInitialState->ProbePdg());
211  // nucleon cluster charge should be incremented by +1 for
212  // neutrino CC and by -1 for antineutrino CC
213  int dQ = (isnu) ? +1 : -1;
214  recoil_nuc = pdg::ModifyNucleonCluster(struck_nuc,dQ); // CC
215  }
216  else {
217  recoil_nuc = struck_nuc; // NC, EM
218  }
219  }
220 
221  LOG("Interaction", pDEBUG) << "Recoil nucleon PDG = " << recoil_nuc;
222  return recoil_nuc;
223 }
bool IsWeak(void) const
bool IsWeakCC(void) const
bool IsNeutrino(int pdgc)
Definition: PDGUtils.cxx:107
int HitNucPdg(void) const
Definition: Target.cxx:304
bool IsNucleon(int pdgc)
Definition: PDGUtils.cxx:343
InitialState * fInitialState
Initial State info.
Definition: Interaction.h:179
bool IsQuasiElastic(void) const
Definition: ProcessInfo.cxx:69
int SwitchProtonNeutron(int pdgc)
Definition: PDGUtils.cxx:353
bool IsInverseBetaDecay(void) const
int ModifyNucleonCluster(int pdgc, int dQ)
Definition: PDGUtils.cxx:361
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
int ProbePdg(void) const
Definition: InitialState.h:64
bool IsDarkMatterElastic(void) const
Definition: ProcessInfo.cxx:74
bool IsMEC(void) const
bool IsEM(void) const
bool Is2NucleonCluster(int pdgc)
Definition: PDGUtils.cxx:399
bool IsDarkMatter(void) const
ProcessInfo * fProcInfo
Process info (scattering, weak current,...)
Definition: Interaction.h:180
const Target & Tgt(void) const
Definition: InitialState.h:66
#define pDEBUG
Definition: Messenger.h:63
Interaction * Interaction::RESCC ( int  tgt,
int  nuc,
int  probe,
double  E = 0 
)
static

Definition at line 608 of file Interaction.cxx.

609 {
612 
613  InitialState * init_state = interaction->InitStatePtr();
614  init_state->SetProbeE(E);
615  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
616 
617  return interaction;
618 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::RESCC ( int  tgt,
int  nuc,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 620 of file Interaction.cxx.

622 {
625 
626  InitialState * init_state = interaction->InitStatePtr();
627  init_state->SetProbeP4(p4probe);
628  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
629 
630  return interaction;
631 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::RESEM ( int  tgt,
int  nuc,
int  probe,
double  E = 0 
)
static

Definition at line 658 of file Interaction.cxx.

659 {
662 
663  InitialState * init_state = interaction->InitStatePtr();
664  init_state->SetProbeE(E);
665  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
666 
667  return interaction;
668 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::RESEM ( int  tgt,
int  nuc,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 670 of file Interaction.cxx.

672 {
675 
676  InitialState * init_state = interaction->InitStatePtr();
677  init_state->SetProbeP4(p4probe);
678  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
679 
680  return interaction;
681 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
void Interaction::Reset ( void  )

Definition at line 83 of file Interaction.cxx.

84 {
85  this->CleanUp();
86  this->Init();
87 }
void CleanUp(void)
Definition: Interaction.cxx:98
Interaction * Interaction::RESNC ( int  tgt,
int  nuc,
int  probe,
double  E = 0 
)
static

Definition at line 633 of file Interaction.cxx.

634 {
637 
638  InitialState * init_state = interaction->InitStatePtr();
639  init_state->SetProbeE(E);
640  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
641 
642  return interaction;
643 }
void SetProbeE(double E)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
Interaction * Interaction::RESNC ( int  tgt,
int  nuc,
int  probe,
const TLorentzVector &  p4probe 
)
static

Definition at line 645 of file Interaction.cxx.

647 {
650 
651  InitialState * init_state = interaction->InitStatePtr();
652  init_state->SetProbeP4(p4probe);
653  init_state->TgtPtr()->SetHitNucPdg(hitnuc);
654 
655  return interaction;
656 }
void SetProbeP4(const TLorentzVector &P4)
Summary information for an interaction.
Definition: Interaction.h:56
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
Target * TgtPtr(void) const
Definition: InitialState.h:67
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
static Interaction * Create(int tgt, int probe, ScatteringType_t st, InteractionType_t it)
Initial State information.
Definition: InitialState.h:48
void Interaction::SetExclTag ( const XclsTag xcls)

Definition at line 243 of file Interaction.cxx.

244 {
245  if (!fExclusiveTag) fExclusiveTag = new XclsTag();
246  fExclusiveTag->Copy(xcls_tag);
247 }
void Copy(const XclsTag &xcls)
copy input XclsTag object
Definition: XclsTag.cxx:173
XclsTag * fExclusiveTag
Additional info for exclusive channels.
Definition: Interaction.h:182
Contains minimal information for tagging exclusive processes.
Definition: XclsTag.h:39
void Interaction::SetInitState ( const InitialState init)

Definition at line 225 of file Interaction.cxx.

226 {
228  fInitialState->Copy(init_state);
229 }
InitialState * fInitialState
Initial State info.
Definition: Interaction.h:179
void Copy(const InitialState &init_state)
Initial State information.
Definition: InitialState.h:48
void Interaction::SetKine ( const Kinematics kine)

Definition at line 237 of file Interaction.cxx.

238 {
239  if (!fKinematics) fKinematics = new Kinematics();
240  fKinematics->Copy(kinematics);
241 }
Generated/set kinematical variables for an event.
Definition: Kinematics.h:39
void Copy(const Kinematics &kine)
Definition: Kinematics.cxx:83
Kinematics * fKinematics
kinematical variables
Definition: Interaction.h:181
void Interaction::SetProcInfo ( const ProcessInfo proc)

Definition at line 231 of file Interaction.cxx.

232 {
233  if (!fProcInfo) fProcInfo = new ProcessInfo();
234  fProcInfo->Copy(proc_info);
235 }
A class encapsulating an enumeration of interaction types (EM, Weak-CC, Weak-NC) and scattering types...
Definition: ProcessInfo.h:46
ProcessInfo * fProcInfo
Process info (scattering, weak current,...)
Definition: Interaction.h:180
void Copy(const ProcessInfo &proc)

Friends And Related Function Documentation

ostream& operator<< ( ostream &  stream,
const Interaction i 
)
friend

print

Member Data Documentation

XclsTag* genie::Interaction::fExclusiveTag
private

Additional info for exclusive channels.

Definition at line 182 of file Interaction.h.

InitialState* genie::Interaction::fInitialState
private

Initial State info.

Definition at line 179 of file Interaction.h.

Kinematics* genie::Interaction::fKinematics
private

kinematical variables

Definition at line 181 of file Interaction.h.

KPhaseSpace* genie::Interaction::fKinePhSp
private

Kinematic phase space.

Definition at line 183 of file Interaction.h.

ProcessInfo* genie::Interaction::fProcInfo
private

Process info (scattering, weak current,...)

Definition at line 180 of file Interaction.h.


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