Public Member Functions | Private Attributes | Friends | List of all members
simb::MCNeutrino Class Reference

Event generator information. More...

#include <MCNeutrino.h>

Public Member Functions

 MCNeutrino ()
 
 MCNeutrino (simb::MCParticle &nu, simb::MCParticle &lep, int CCNC, int mode, int interactionType, int target, int nucleon, int quark, double w, double x, double y, double qsqr)
 nu is the incoming neutrino and lep is the outgoing lepton More...
 
const simb::MCParticleNu () const
 
const simb::MCParticleLepton () const
 
int CCNC () const
 
int Mode () const
 
int InteractionType () const
 
int Target () const
 
int HitNuc () const
 
int HitQuark () const
 
double W () const
 
double X () const
 
double Y () const
 
double QSqr () const
 
double Pt () const
 transverse momentum of interaction, in GeV/c More...
 
double Theta () const
 angle between incoming and outgoing leptons, in radians More...
 

Private Attributes

simb::MCParticle fNu
 the incoming neutrino More...
 
simb::MCParticle fLepton
 the outgoing lepton More...
 
int fMode
 Interaction mode (QE/1-pi/DIS...) see enum list. More...
 
int fInteractionType
 More detailed interaction type, see enum list below kNuanceOffset. More...
 
int fCCNC
 CC or NC interaction? see enum list. More...
 
int fTarget
 Nuclear target, as PDG code. More...
 
int fHitNuc
 Hit nucleon (2212 (proton) or 2112 (neutron)) More...
 
int fHitQuark
 For DIS events only, as PDG code. More...
 
double fW
 Hadronic invariant mass, in GeV. More...
 
double fX
 Bjorken x=Q^2/(2M*(E_neutrino-E_lepton)), unitless. More...
 
double fY
 Inelasticity y=1-(E_lepton/E_neutrino), unitless. More...
 
double fQSqr
 Momentum transfer Q^2, in GeV^2. More...
 

Friends

std::ostream & operator<< (std::ostream &output, const simb::MCNeutrino &mcnu)
 

Detailed Description

Event generator information.

Definition at line 18 of file MCNeutrino.h.

Constructor & Destructor Documentation

simb::MCNeutrino::MCNeutrino ( )

Definition at line 17 of file MCNeutrino.cxx.

18  : fNu ()
19  , fLepton ()
30  {
31  }
double fX
Bjorken x=Q^2/(2M*(E_neutrino-E_lepton)), unitless.
Definition: MCNeutrino.h:34
int fTarget
Nuclear target, as PDG code.
Definition: MCNeutrino.h:30
int fHitNuc
Hit nucleon (2212 (proton) or 2112 (neutron))
Definition: MCNeutrino.h:31
double fQSqr
Momentum transfer Q^2, in GeV^2.
Definition: MCNeutrino.h:36
simb::MCParticle fNu
the incoming neutrino
Definition: MCNeutrino.h:25
simb::MCParticle fLepton
the outgoing lepton
Definition: MCNeutrino.h:26
int fMode
Interaction mode (QE/1-pi/DIS...) see enum list.
Definition: MCNeutrino.h:27
double fW
Hadronic invariant mass, in GeV.
Definition: MCNeutrino.h:33
int fCCNC
CC or NC interaction? see enum list.
Definition: MCNeutrino.h:29
T min(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
Definition: statistics.h:55
double fY
Inelasticity y=1-(E_lepton/E_neutrino), unitless.
Definition: MCNeutrino.h:35
int fHitQuark
For DIS events only, as PDG code.
Definition: MCNeutrino.h:32
int fInteractionType
More detailed interaction type, see enum list below kNuanceOffset.
Definition: MCNeutrino.h:28
simb::MCNeutrino::MCNeutrino ( simb::MCParticle nu,
simb::MCParticle lep,
int  CCNC,
int  mode,
int  interactionType,
int  target,
int  nucleon,
int  quark,
double  w,
double  x,
double  y,
double  qsqr 
)

nu is the incoming neutrino and lep is the outgoing lepton

Definition at line 35 of file MCNeutrino.cxx.

47  : fNu(nu)
48  , fLepton(lep)
49  , fMode(mode)
50  , fInteractionType(interactionType)
51  , fCCNC(CCNC)
52  , fTarget(target)
53  , fHitNuc(nucleon)
54  , fHitQuark(quark)
55  , fW(w)
56  , fX(x)
57  , fY(y)
58  , fQSqr(qsqr)
59  {
60  }
int CCNC() const
Definition: MCNeutrino.h:148
double fX
Bjorken x=Q^2/(2M*(E_neutrino-E_lepton)), unitless.
Definition: MCNeutrino.h:34
int fTarget
Nuclear target, as PDG code.
Definition: MCNeutrino.h:30
int fHitNuc
Hit nucleon (2212 (proton) or 2112 (neutron))
Definition: MCNeutrino.h:31
double fQSqr
Momentum transfer Q^2, in GeV^2.
Definition: MCNeutrino.h:36
simb::MCParticle fNu
the incoming neutrino
Definition: MCNeutrino.h:25
simb::MCParticle fLepton
the outgoing lepton
Definition: MCNeutrino.h:26
int fMode
Interaction mode (QE/1-pi/DIS...) see enum list.
Definition: MCNeutrino.h:27
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 fY
Inelasticity y=1-(E_lepton/E_neutrino), unitless.
Definition: MCNeutrino.h:35
list x
Definition: train.py:276
int fHitQuark
For DIS events only, as PDG code.
Definition: MCNeutrino.h:32
int fInteractionType
More detailed interaction type, see enum list below kNuanceOffset.
Definition: MCNeutrino.h:28

Member Function Documentation

int simb::MCNeutrino::CCNC ( ) const
inline

Definition at line 148 of file MCNeutrino.h.

148 { return fCCNC; }
int fCCNC
CC or NC interaction? see enum list.
Definition: MCNeutrino.h:29
int simb::MCNeutrino::HitNuc ( ) const
inline

Definition at line 152 of file MCNeutrino.h.

152 { return fHitNuc; }
int fHitNuc
Hit nucleon (2212 (proton) or 2112 (neutron))
Definition: MCNeutrino.h:31
int simb::MCNeutrino::HitQuark ( ) const
inline

Definition at line 153 of file MCNeutrino.h.

153 { return fHitQuark; }
int fHitQuark
For DIS events only, as PDG code.
Definition: MCNeutrino.h:32
int simb::MCNeutrino::InteractionType ( ) const
inline

Definition at line 150 of file MCNeutrino.h.

150 { return fInteractionType; }
int fInteractionType
More detailed interaction type, see enum list below kNuanceOffset.
Definition: MCNeutrino.h:28
const simb::MCParticle & simb::MCNeutrino::Lepton ( ) const
inline

Definition at line 147 of file MCNeutrino.h.

147 { return fLepton; }
simb::MCParticle fLepton
the outgoing lepton
Definition: MCNeutrino.h:26
int simb::MCNeutrino::Mode ( ) const
inline

Definition at line 149 of file MCNeutrino.h.

149 { return fMode; }
int fMode
Interaction mode (QE/1-pi/DIS...) see enum list.
Definition: MCNeutrino.h:27
const simb::MCParticle & simb::MCNeutrino::Nu ( ) const
inline

Definition at line 146 of file MCNeutrino.h.

146 { return fNu; }
simb::MCParticle fNu
the incoming neutrino
Definition: MCNeutrino.h:25
double simb::MCNeutrino::Pt ( ) const

transverse momentum of interaction, in GeV/c

Definition at line 74 of file MCNeutrino.cxx.

75  {
76  return fNu.Pt();
77  }
simb::MCParticle fNu
the incoming neutrino
Definition: MCNeutrino.h:25
double Pt(const int i=0) const
Definition: MCParticle.h:236
double simb::MCNeutrino::QSqr ( ) const
inline

Definition at line 157 of file MCNeutrino.h.

157 { return fQSqr; }
double fQSqr
Momentum transfer Q^2, in GeV^2.
Definition: MCNeutrino.h:36
int simb::MCNeutrino::Target ( ) const
inline

Definition at line 151 of file MCNeutrino.h.

151 { return fTarget; }
int fTarget
Nuclear target, as PDG code.
Definition: MCNeutrino.h:30
double simb::MCNeutrino::Theta ( ) const

angle between incoming and outgoing leptons, in radians

make TVector3 objects for the momenta of the incoming neutrino and outgoing lepton

Definition at line 63 of file MCNeutrino.cxx.

64  {
65  ///make TVector3 objects for the momenta of the incoming neutrino
66  ///and outgoing lepton
67  TVector3 in(fNu.Px(), fNu.Py(), fNu.Pz());
68  TVector3 out(fLepton.Px(), fLepton.Py(), fLepton.Pz());
69 
70  return in.Angle(out);
71  }
double Py(const int i=0) const
Definition: MCParticle.h:231
double Px(const int i=0) const
Definition: MCParticle.h:230
simb::MCParticle fNu
the incoming neutrino
Definition: MCNeutrino.h:25
simb::MCParticle fLepton
the outgoing lepton
Definition: MCNeutrino.h:26
double Pz(const int i=0) const
Definition: MCParticle.h:232
double simb::MCNeutrino::W ( ) const
inline

Definition at line 154 of file MCNeutrino.h.

154 { return fW; }
double fW
Hadronic invariant mass, in GeV.
Definition: MCNeutrino.h:33
double simb::MCNeutrino::X ( ) const
inline

Definition at line 155 of file MCNeutrino.h.

155 { return fX; }
double fX
Bjorken x=Q^2/(2M*(E_neutrino-E_lepton)), unitless.
Definition: MCNeutrino.h:34
double simb::MCNeutrino::Y ( ) const
inline

Definition at line 156 of file MCNeutrino.h.

156 { return fY; }
double fY
Inelasticity y=1-(E_lepton/E_neutrino), unitless.
Definition: MCNeutrino.h:35

Friends And Related Function Documentation

std::ostream& operator<< ( std::ostream &  output,
const simb::MCNeutrino mcnu 
)
friend

Definition at line 80 of file MCNeutrino.cxx.

81  {
82  output << " neutrino = " << mcnu.Nu().PdgCode() << std::endl
83  << " neutrino energy = " << mcnu.Nu().E() << std::endl
84  << " CCNC = " << mcnu.CCNC() << std::endl
85  << " mode = " << mcnu.Mode() << std::endl
86  << " interaction type = " << mcnu.InteractionType() << std::endl
87  << " target = " << mcnu.Target() << std::endl
88  << " nucleon = " << mcnu.HitNuc() << std::endl
89  << " quark = " << mcnu.HitQuark() << std::endl
90  << " W = " << mcnu.W() << std::endl
91  << " X = " << mcnu.X() << std::endl
92  << " Y = " << mcnu.Y() << std::endl
93  << " Q^2 = " << mcnu.QSqr() << std::endl;
94 
95  return output;
96  }
double E(const int i=0) const
Definition: MCParticle.h:233
int PdgCode() const
Definition: MCParticle.h:212
int CCNC() const
Definition: MCNeutrino.h:148
double QSqr() const
Definition: MCNeutrino.h:157
int HitQuark() const
Definition: MCNeutrino.h:153
const simb::MCParticle & Nu() const
Definition: MCNeutrino.h:146
int HitNuc() const
Definition: MCNeutrino.h:152
int InteractionType() const
Definition: MCNeutrino.h:150
double W() const
Definition: MCNeutrino.h:154
double Y() const
Definition: MCNeutrino.h:156
double X() const
Definition: MCNeutrino.h:155
int Target() const
Definition: MCNeutrino.h:151
int Mode() const
Definition: MCNeutrino.h:149
QTextStream & endl(QTextStream &s)

Member Data Documentation

int simb::MCNeutrino::fCCNC
private

CC or NC interaction? see enum list.

Definition at line 29 of file MCNeutrino.h.

int simb::MCNeutrino::fHitNuc
private

Hit nucleon (2212 (proton) or 2112 (neutron))

Definition at line 31 of file MCNeutrino.h.

int simb::MCNeutrino::fHitQuark
private

For DIS events only, as PDG code.

Definition at line 32 of file MCNeutrino.h.

int simb::MCNeutrino::fInteractionType
private

More detailed interaction type, see enum list below kNuanceOffset.

Definition at line 28 of file MCNeutrino.h.

simb::MCParticle simb::MCNeutrino::fLepton
private

the outgoing lepton

Definition at line 26 of file MCNeutrino.h.

int simb::MCNeutrino::fMode
private

Interaction mode (QE/1-pi/DIS...) see enum list.

Definition at line 27 of file MCNeutrino.h.

simb::MCParticle simb::MCNeutrino::fNu
private

the incoming neutrino

Definition at line 25 of file MCNeutrino.h.

double simb::MCNeutrino::fQSqr
private

Momentum transfer Q^2, in GeV^2.

Definition at line 36 of file MCNeutrino.h.

int simb::MCNeutrino::fTarget
private

Nuclear target, as PDG code.

Definition at line 30 of file MCNeutrino.h.

double simb::MCNeutrino::fW
private

Hadronic invariant mass, in GeV.

Definition at line 33 of file MCNeutrino.h.

double simb::MCNeutrino::fX
private

Bjorken x=Q^2/(2M*(E_neutrino-E_lepton)), unitless.

Definition at line 34 of file MCNeutrino.h.

double simb::MCNeutrino::fY
private

Inelasticity y=1-(E_lepton/E_neutrino), unitless.

Definition at line 35 of file MCNeutrino.h.


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