MCNeutrino.cxx
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 /// \file MCNeutrino.cxx
3 /// \brief Simple MC truth class, holds a vector of TParticles
4 ///
5 /// \version $Id: MCNeutrino.cxx,v 1.5 2012-10-15 20:36:27 brebel Exp $
6 /// \author jpaley@indiana.edu
7 ////////////////////////////////////////////////////////////////////////
10 #include "TVector3.h"
11 #include <iostream>
12 #include <climits>
13 
14 namespace simb{
15 
16  //......................................................................
18  : fNu ()
19  , fLepton ()
20  , fMode (std::numeric_limits<int>::min())
21  , fInteractionType(std::numeric_limits<int>::min())
22  , fCCNC (std::numeric_limits<int>::min())
23  , fTarget (std::numeric_limits<int>::min())
24  , fHitNuc (std::numeric_limits<int>::min())
25  , fHitQuark (std::numeric_limits<int>::min())
26  , fW (std::numeric_limits<double>::min())
27  , fX (std::numeric_limits<double>::min())
28  , fY (std::numeric_limits<double>::min())
29  , fQSqr (std::numeric_limits<double>::min())
30  {
31  }
32 
33  //......................................................................
34  ///nu is the incoming neutrino and lep is the outgoing lepton
36  MCParticle &lep,
37  int CCNC,
38  int mode,
39  int interactionType,
40  int target,
41  int nucleon,
42  int quark,
43  double w,
44  double x,
45  double y,
46  double qsqr)
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  }
61 
62  //......................................................................
63  double MCNeutrino::Theta() const
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  }
72 
73  //......................................................................
74  double MCNeutrino::Pt() const
75  {
76  return fNu.Pt();
77  }
78 
79  //......................................................................
80  std::ostream& operator<< (std::ostream& output, const simb::MCNeutrino &mcnu)
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  }
97 
98 }
99 ////////////////////////////////////////////////////////////////////////
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
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
double Py(const int i=0) const
Definition: MCParticle.h:231
int fTarget
Nuclear target, as PDG code.
Definition: MCNeutrino.h:30
int HitQuark() const
Definition: MCNeutrino.h:153
const simb::MCParticle & Nu() const
Definition: MCNeutrino.h:146
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
double Px(const int i=0) const
Definition: MCParticle.h:230
int HitNuc() const
Definition: MCNeutrino.h:152
STL namespace.
Particle class.
double fQSqr
Momentum transfer Q^2, in GeV^2.
Definition: MCNeutrino.h:36
simb::MCParticle fNu
the incoming neutrino
Definition: MCNeutrino.h:25
int InteractionType() const
Definition: MCNeutrino.h:150
double Pt(const int i=0) const
Definition: MCParticle.h:236
simb::MCParticle fLepton
the outgoing lepton
Definition: MCNeutrino.h:26
double W() const
Definition: MCNeutrino.h:154
double Y() const
Definition: MCNeutrino.h:156
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 X() const
Definition: MCNeutrino.h:155
Base utilities and modules for event generation and detector simulation.
friend std::ostream & operator<<(std::ostream &output, const simb::MCNeutrino &mcnu)
Definition: MCNeutrino.cxx:80
int Target() const
Definition: MCNeutrino.h:151
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
double Pz(const int i=0) const
Definition: MCParticle.h:232
list x
Definition: train.py:276
int fHitQuark
For DIS events only, as PDG code.
Definition: MCNeutrino.h:32
Event generator information.
Definition: MCNeutrino.h:18
int Mode() const
Definition: MCNeutrino.h:149
int fInteractionType
More detailed interaction type, see enum list below kNuanceOffset.
Definition: MCNeutrino.h:28
QTextStream & endl(QTextStream &s)