eventAnalysis  7.0-49-g0ac7482
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
TTruthVerticesModule.hxx
Go to the documentation of this file.
1 ///
2 /// For questions or suggestions about this module please contact the
3 /// current responsible and CC in the eventAnalysis package manager.
4 ///
5 /// 26-Nov-2012: Current responsible for this module is,
6 /// Steve Dennis (stephen.dennis [*a*t*] warwick.ac.uk)
7 ///
8 
9 #ifndef TTruthVerticesModule_hxx_seen
10 #define TTruthVerticesModule_hxx_seen
11 
12 // Root
13 #include <TObject.h>
14 
15 // oaEvent
16 #include <TND280Event.hxx>
17 
18 // Geant4
19 #include "TG4PrimaryVertex.hxx"
20 
21 #include "eventAnalysisEnums.hxx"
22 
23 // eventAnalysis
24 #include "EeventAnalysis.hxx"
26 
27 namespace ND {
28 class TTruthVerticesModule;
29 OA_EXCEPTION(ETruthVerticesModule, EeventAnalysis);
30 OA_EXCEPTION(EUndefinedParticleCategoryVertex, ETruthVerticesModule);
31 }
32 
33 ///\brief eventAnalysis module for storing the truth information for primary
34 /// vertices in events
36  public:
37  class TTruthVertex;
38 
39  TTruthVerticesModule(const char* name = "Vertices",
40  const char* title = "True Primary Vertex information");
41  virtual ~TTruthVerticesModule();
42 
43  virtual Bool_t IsEnabledByDefault() const { return kTRUE; }
44  virtual void InitializeBranches();
45  virtual bool FillTree(ND::TND280Event&);
46  virtual Bool_t ProcessFirstEvent(ND::TND280Event&);
47 
48  private:
49  Bool_t IsGeantinoVtx(ND::TG4PrimaryVertex vtx);
50 
51  ///\brief The maximum number of vertices that the module will be saved in
52  /// a single event.
53  const unsigned int fMaxNVertices;
54 
55  public:
56  //--------------------------------------------------------------
57  // Tree Entries
58  //--------------------------------------------------------------
59 
60  Int_t fNVtx; ///< The total number of vertices recorded
61 
62 
63  Int_t fNVtxFGD1; ///< The number of vertices recorded in FGD 1
64 
65 
66  Int_t fNVtxFGD2; ///< The number of vertices recorded in FGD 2
67 
68 
69  Int_t fNVtxSFG; ///< The number of vertices recorded in SFG
70 
71 
72  Int_t fNVtxP0D; ///< The number of vertices recorded in the P0D
73 
74 
75  Int_t fNVtxDsECal; ///< The number of vertices recorded in the Downstream ECal
76 
77 
78  Int_t fNVtxBrlECal; ///< The number of vertices recorded in the Barrel ECal
79 
80 
81  Int_t fNVtxP0DECal; ///< The number of vertices recorded in the P0D ECal
82 
83 
84  Int_t fNVtxTPC1; ///< The number of vertices recorded in TPC 1
85 
86 
87  Int_t fNVtxTPC2; ///< The number of vertices recorded in TPC 2
88 
89 
90  Int_t fNVtxTPC3; ///< The number of vertices recorded in TPC 3
91 
92 
93  Int_t fNVtxSMRD; ///< The number of vertices recorded in the SMRD
94 
95 
96  Int_t fNVtxRestOfOffAxis; ///< The number of vertices recorded in the rest of the off-axis detector
97 
98 
99  Int_t fNVtxINGRID; ///< The number of vertices recorded in INGRID
100 
101 
102  Int_t fNVtxOther;///< The number of vertices recorded anywhere which does not fall into the other available categories
103 
104  TClonesArray* fVertices;///< The TClonesArray storing the TTruthVertex objects holding the information relating to each vertex.
105  ///<\details The array is sorted by the vertices ID number to enable more
106  ///< efficient retrieval by analysis tools.
107 
108  private:
109 };
110 
111 ///\brief Class used by the Truth Vertices Module to store information relating
112 /// to an individual true primary vertex.
114  public:
115  TTruthVertex();
116  virtual ~TTruthVertex(){};
117 
118  ///\brief Compare the values of the vertices' IDs so that a
119  /// TClonesArray can be sorted in order of increasing ID.
120  Int_t Compare(const TObject* obj) const;
121 
122  ///\brief Make the object sortable so that a TClonesArray can be
123  /// sorted in ID order.
124  Bool_t IsSortable() const { return kTRUE; }
125 
126  //--------------------------------------------------------------
127  // Vertex information
128  //--------------------------------------------------------------
129 
130  Int_t ID; ///< A number which uniquely identifies this vertex within the event. This ID is the interface between the Truth Vertices module and other eventAnalysis modules. Other modules should use this number to reference trajectories.
131 
132  TLorentzVector Position; ///< Position and time of the vertex
133 
134  std::string Generator; ///< The generator that created the vertex. eg: "genie:mean@free-spill"
135 
136  std::string ReactionCode; ///< The Reaction code according to the generator For Genie this will be of the form: "nu:14;tgt:1000260560;N:2112;proc:Weak[CC],QES;" For Neut it will be an integer, see definitions here: http://www.t2k.org/asg/xsec/niwgdocs/neut_xsecs/neutmodesC.h/view
137 
138  Int_t Subdetector; ///< Subdetector which the vertex occurs in.
139 
140 
141  Double_t Fiducial; ///< The distance inside the local fiducial volume [mm]. Not currently set for any detector other than the P0D.
142 
143  Int_t NPrimaryTraj; ///< The number of primary trajectories (ie: the number of primary particles exiting the interaction vertex).
144 
145  std::vector<Int_t> PrimaryTrajIDs; ///<\brief ID numbers which uniquely identify the trajectories of the primary particles of the vertex within the event.
146 
147  Int_t NeutrinoPDG; ///< The PDG number of the incoming neutrino. Set to 0 in the absence of a neutrino.
148 
149  TLorentzVector NeutrinoMomentum; ///< The four-momentum of the incoming neutrino. Set to (-999999.9, -999999.9, -999999.9, -999999.9) in the absence of a neutrino.
150 
151  Int_t TargetPDG; ///< The (extended for nuclei) PDG number of the target. Set to 0 in the absence of a target.
152 
153  TLorentzVector TargetMomentum; ///< The four-momentum of the target. Set to (-999999.9, -999999.9, -999999.9, -999999.9) in the absence of a target.
154 
155  private:
157 };
158 
159 #endif
std::string Generator
The generator that created the vertex. eg: &quot;genie:mean@free-spill&quot;.
const unsigned int fMaxNVertices
The maximum number of vertices that the module will be saved in a single event.
Int_t NeutrinoPDG
The PDG number of the incoming neutrino. Set to 0 in the absence of a neutrino.
Int_t fNVtxTPC1
The number of vertices recorded in TPC 1.
OA_EXCEPTION(EeventAnalysis, EoaCore)
Generate a base exception for the Analysis library.
eventAnalysis module for storing the truth information for primary vertices in events ...
Int_t fNVtxDsECal
The number of vertices recorded in the Downstream ECal.
Int_t fNVtxTPC3
The number of vertices recorded in TPC 3.
TLorentzVector TargetMomentum
The four-momentum of the target. Set to (-999999.9, -999999.9, -999999.9, -999999.9) in the absence of a target.
TClonesArray * fVertices
The TClonesArray storing the TTruthVertex objects holding the information relating to each vertex...
Bool_t IsGeantinoVtx(ND::TG4PrimaryVertex vtx)
std::string ReactionCode
The Reaction code according to the generator For Genie this will be of the form: &quot;nu:14;tgt:100026056...
Int_t fNVtxFGD2
The number of vertices recorded in FGD 2.
Int_t fNVtxINGRID
The number of vertices recorded in INGRID.
Int_t TargetPDG
The (extended for nuclei) PDG number of the target. Set to 0 in the absence of a target.
Int_t fNVtxBrlECal
The number of vertices recorded in the Barrel ECal.
Class used by the Truth Vertices Module to store information relating to an individual true primary v...
Bool_t IsSortable() const
Make the object sortable so that a TClonesArray can be sorted in ID order.
virtual Bool_t IsEnabledByDefault() const
Is the module is enabled by default.
Int_t fNVtxTPC2
The number of vertices recorded in TPC 2.
Int_t fNVtxSMRD
The number of vertices recorded in the SMRD.
Int_t fNVtxOther
The number of vertices recorded anywhere which does not fall into the other available categories...
Int_t fNVtxP0DECal
The number of vertices recorded in the P0D ECal.
Int_t fNVtxFGD1
The number of vertices recorded in FGD 1.
A base class for analysis output modules which contain event truth information.
Int_t fNVtxP0D
The number of vertices recorded in the P0D.
Int_t Compare(const TObject *obj) const
Compare the values of the vertices&#39; IDs so that a TClonesArray can be sorted in order of increasing I...
Int_t fNVtxRestOfOffAxis
The number of vertices recorded in the rest of the off-axis detector.
TLorentzVector NeutrinoMomentum
The four-momentum of the incoming neutrino. Set to (-999999.9, -999999.9, -999999.9, -999999.9) in the absence of a neutrino.
virtual void InitializeBranches()
Initialize Branches. Don&#39;t do anything else in this function.
Int_t Subdetector
Subdetector which the vertex occurs in.
std::vector< Int_t > PrimaryTrajIDs
ID numbers which uniquely identify the trajectories of the primary particles of the vertex within the...
Int_t fNVtx
The total number of vertices recorded.
virtual Bool_t ProcessFirstEvent(ND::TND280Event &)
Is called after the first event is loaded in.
ClassDef(ND::TTruthVerticesModule::TTruthVertex, 1)
Int_t ID
A number which uniquely identifies this vertex within the event. This ID is the interface between the...
Double_t Fiducial
The distance inside the local fiducial volume [mm]. Not currently set for any detector other than the...
Int_t fNVtxSFG
The number of vertices recorded in SFG.
TLorentzVector Position
Position and time of the vertex.
TTruthVerticesModule(const char *name="Vertices", const char *title="True Primary Vertex information")
Int_t NPrimaryTraj
The number of primary trajectories (ie: the number of primary particles exiting the interaction verte...
virtual bool FillTree(ND::TND280Event &)
Fill all the stuff that goes in the output tree.

Package Summary
Package Name: eventAnalysis
Package Version: 7.0-49-g0ac7482
Package Manager:

Generated on Mon Mar 25 2024 14:43:59 for eventAnalysis by doxygen 1.8.5