LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
MCNeutrino.h
Go to the documentation of this file.
1 #ifndef SIMB_MCNEUTRINO_H
9 #define SIMB_MCNEUTRINO_H
10 
12 
13 namespace simb {
14 
15  //......................................................................
16 
18  class MCNeutrino {
19  public:
20 
21  MCNeutrino();
22 
23  private:
24 
27  int fMode;
29  int fCCNC;
30  int fTarget;
31  int fHitNuc;
32  int fHitQuark;
33  double fW;
34  double fX;
35  double fY;
36  double fQSqr;
37 
38 #ifndef __GCCXML__
39  public:
40 
42  simb::MCParticle &lep,
43  int CCNC,
44  int mode,
45  int interactionType,
46  int target,
47  int nucleon,
48  int quark,
49  double w,
50  double x,
51  double y,
52  double qsqr);
53 
54  const simb::MCParticle& Nu() const;
55  const simb::MCParticle& Lepton() const;
56  int CCNC() const;
57  int Mode() const;
58  int InteractionType() const;
59  int Target() const;
60  int HitNuc() const;
61  int HitQuark() const;
62  double W() const;
63  double X() const;
64  double Y() const;
65  double QSqr() const;
66  double Pt() const;
67  double Theta() const;
68  friend std::ostream& operator<< (std::ostream& output, const simb::MCNeutrino &mcnu);
69 #endif
70 
71  };
72 }
73 
74 #ifndef __GCCXML__
75 namespace simb{
77  enum curr_type_{
78  kCC,
80  };
81 
83  enum int_type_{
85  kQE = 0,
86  kRes = 1,
87  kDIS = 2,
88  kCoh = 3,
95  kMEC = 10,
97  kEM = 12,
98  kWeakMix = 13,
99  kNuanceOffset = 1000,
147  };
148 }
149 
150 inline const simb::MCParticle& simb::MCNeutrino::Nu() const { return fNu; }
151 inline const simb::MCParticle& simb::MCNeutrino::Lepton() const { return fLepton; }
152 inline int simb::MCNeutrino::CCNC() const { return fCCNC; }
153 inline int simb::MCNeutrino::Mode() const { return fMode; }
155 inline int simb::MCNeutrino::Target() const { return fTarget; }
156 inline int simb::MCNeutrino::HitNuc() const { return fHitNuc; }
157 inline int simb::MCNeutrino::HitQuark() const { return fHitQuark; }
158 inline double simb::MCNeutrino::W() const { return fW; }
159 inline double simb::MCNeutrino::X() const { return fX; }
160 inline double simb::MCNeutrino::Y() const { return fY; }
161 inline double simb::MCNeutrino::QSqr() const { return fQSqr; }
162 #endif
163 
164 #endif //SIMB_MCNEUTRINO_H
165 
Float_t x
Definition: compare.C:6
neutral current quasi-elastic
Definition: MCNeutrino.h:101
resonant charged current, nubar p -> nubar n pi+
Definition: MCNeutrino.h:113
int CCNC() const
Definition: MCNeutrino.h:152
double QSqr() const
Definition: MCNeutrino.h:161
double Theta() const
angle between incoming and outgoing leptons, in radians
Definition: MCNeutrino.cxx:63
double fX
Bjorken x=Q^2/(2M*(E_neutrino-E_lepton)), unitless.
Definition: MCNeutrino.h:34
resonant neutral current, nu p -> nu p pi0
Definition: MCNeutrino.h:105
charged current deep inelastic scatter
Definition: MCNeutrino.h:138
resonant charged current, nubar p -> l+ p pi-
Definition: MCNeutrino.h:111
int fTarget
Nuclear target, as PDG code.
Definition: MCNeutrino.h:30
int HitQuark() const
Definition: MCNeutrino.h:157
neutrino electron elastic scatter
Definition: MCNeutrino.h:144
cout<< "-> Edep in the target
Definition: analysis.C:54
Float_t y
Definition: compare.C:6
const simb::MCParticle & Nu() const
Definition: MCNeutrino.h:150
resonant charged current, nubar p -> l+ n pi0
Definition: MCNeutrino.h:110
int fHitNuc
Hit nucleon (2212 (proton) or 2112 (neutron))
Definition: MCNeutrino.h:31
double Pt() const
transverse momentum of interaction, in GeV/c
Definition: MCNeutrino.cxx:74
int HitNuc() const
Definition: MCNeutrino.h:156
offset to account for adding in Nuance codes to this enum
Definition: MCNeutrino.h:99
charged current quasi-elastic
Definition: MCNeutrino.h:100
Particle class.
double fQSqr
Momentum transfer Q^2, in GeV^2.
Definition: MCNeutrino.h:36
resonant charged current, nubar p -> nubar p pi0
Definition: MCNeutrino.h:112
simb::MCParticle fNu
the incoming neutrino
Definition: MCNeutrino.h:25
resonant charged current, nu n -> l- p pi0
Definition: MCNeutrino.h:103
int InteractionType() const
Definition: MCNeutrino.h:154
simb::MCParticle fLepton
the outgoing lepton
Definition: MCNeutrino.h:26
const simb::MCParticle & Lepton() const
Definition: MCNeutrino.h:151
resonant neutral current, nu n -> nu n pi0
Definition: MCNeutrino.h:107
double W() const
Definition: MCNeutrino.h:158
resonant charged current, nu n -> l- n pi+
Definition: MCNeutrino.h:104
extension of nuance encoding for MEC / 2p2h
Definition: MCNeutrino.h:146
double Y() const
Definition: MCNeutrino.h:160
int fMode
Interaction mode (QE/1-pi/DIS...) see enum list.
Definition: MCNeutrino.h:27
resonant charged current, nubar n -> nubar p pi-
Definition: MCNeutrino.h:115
double fW
Hadronic invariant mass, in GeV.
Definition: MCNeutrino.h:33
int fCCNC
CC or NC interaction? see enum list.
Definition: MCNeutrino.h:29
double X() const
Definition: MCNeutrino.h:159
Framework includes.
charged current deep inelastic scatter
Definition: MCNeutrino.h:137
resonant charged current, nu p -> l- p pi+
Definition: MCNeutrino.h:102
friend std::ostream & operator<<(std::ostream &output, const simb::MCNeutrino &mcnu)
Definition: MCNeutrino.cxx:80
n.b.: this group is similar but not quite, entirely unlike GENIE ScatteringType convention ...
Definition: MCNeutrino.h:84
int Target() const
Definition: MCNeutrino.h:155
charged current coherent pion
Definition: MCNeutrino.h:143
resonant neutral current, nu n -> nu p pi-
Definition: MCNeutrino.h:108
double fY
Inelasticity y=1-(E_lepton/E_neutrino), unitless.
Definition: MCNeutrino.h:35
inverse muon decay
Definition: MCNeutrino.h:145
int_type_
Neutrino interaction categories.
Definition: MCNeutrino.h:83
resonant charged current, nubar n -> l+ n pi-
Definition: MCNeutrino.h:109
resonant charged current, nubar n -> nubar n pi0
Definition: MCNeutrino.h:114
int fHitQuark
For DIS events only, as PDG code.
Definition: MCNeutrino.h:32
curr_type_
Neutrino interaction categories.
Definition: MCNeutrino.h:77
Event generator information.
Definition: MCNeutrino.h:18
Float_t w
Definition: plot.C:23
int Mode() const
Definition: MCNeutrino.h:153
resonant neutral current, nu p -> nu p pi+
Definition: MCNeutrino.h:106
int fInteractionType
More detailed interaction type, see enum list below kNuanceOffset.
Definition: MCNeutrino.h:28