Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TEveVSDStructs.h
Go to the documentation of this file.
1// @(#)root/eve:$Id$
2// Authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 2007
3
4/*************************************************************************
5 * Copyright (C) 1995-2007, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12#ifndef ROOT_TEveVSDStructs
13#define ROOT_TEveVSDStructs
14
15#include "TObject.h"
16#include "TParticle.h"
17#include "TEveVector.h"
18
19/******************************************************************************/
20// VSD Structures
21/******************************************************************************/
22
23// Basic structures for Reve VSD concept. Design criteria:
24//
25// * provide basic cross-referencing functionality;
26//
27// * small memory/disk footprint (floats / count on compression in
28// split mode);
29//
30// * simple usage from tree selections;
31//
32// * placement in TClonesArray (composites are TObject derived);
33//
34// * minimal member-naming (impossible to make everybody happy).
35//
36
37
38/******************************************************************************/
39// TEveMCTrack
40/******************************************************************************/
41
42class TEveMCTrack : public TParticle // ?? Copy stuff over ??
43{
44public:
45 Int_t fLabel; // Label of the track
46 Int_t fIndex; // Index of the track (in some source array)
47 Int_t fEvaLabel; // Label of primary particle
48
49 Bool_t fDecayed; // True if decayed during tracking.
50 // ?? Perhaps end-of-tracking point/momentum would be better.
51 Float_t fTDecay; // Decay time
52 TEveVector fVDecay; // Decay vertex
53 TEveVector fPDecay; // Decay momentum
54
57 virtual ~TEveMCTrack() {}
58
60 { *((TParticle*)this) = p; return *this; }
61
62 void ResetPdgCode() { fPdgCode = 0; }
63
64 ClassDef(TEveMCTrack, 1); // Monte Carlo track (also used in VSD).
65};
66
67
68/******************************************************************************/
69// TEveHit
70/******************************************************************************/
71
72// Representation of a hit.
73
74// Members det_id (and fSubdetId) serve for cross-referencing into
75// geometry. Hits should be stored in fDetId (+some label ordering) in
76// order to maximize branch compression.
77
78
79class TEveHit : public TObject
80{
81public:
82 UShort_t fDetId; // Custom detector id.
83 UShort_t fSubdetId; // Custom sub-detector id.
84 Int_t fLabel; // Label of particle that produced the hit.
85 Int_t fEvaLabel; // Label of primary particle, ancestor of label.
86 TEveVector fV; // Hit position.
87
88 // Float_t charge; probably specific.
89
90 TEveHit() : fDetId(0), fSubdetId(0), fLabel(0), fEvaLabel(0), fV() {}
91 virtual ~TEveHit() {}
92
93 ClassDef(TEveHit, 1); // Monte Carlo hit (also used in VSD).
94};
95
96
97/******************************************************************************/
98// TEveCluster
99/******************************************************************************/
100
101// Base class for reconstructed clusters
102
103// ?? Should TEveHit and cluster have common base? No.
104
105class TEveCluster : public TObject
106{
107public:
108 UShort_t fDetId; // Custom detector id.
109 UShort_t fSubdetId; // Custom sub-detector id.
110 Int_t fLabel[3]; // Labels of particles that contributed hits.
111
112 // ?? Should include reconstructed track(s) using it? Rather not, separate.
113
114 TEveVector fV; // Vertex.
115 // TEveVector fW; // Cluster widths.
116 // Coord system? Errors and/or widths Wz, Wy?
117
118 TEveCluster() : fDetId(0), fSubdetId(0), fV() { fLabel[0] = fLabel[1] = fLabel[2] = 0; }
119 virtual ~TEveCluster() {}
120
121 ClassDef(TEveCluster, 1); // Reconstructed cluster (also used in VSD).
122};
123
124
125/******************************************************************************/
126// TEveRecTrack
127/******************************************************************************/
128template <typename TT>
129class TEveRecTrackT : public TObject
130{
131public:
132 Int_t fLabel; // Label of the track.
133 Int_t fIndex; // Index of the track (in some source array).
134 Int_t fStatus; // Status as exported from reconstruction.
135 Int_t fSign; // Charge of the track.
136 TEveVectorT<TT> fV; // Start vertex from reconstruction.
137 TEveVectorT<TT> fP; // Reconstructed momentum at start vertex.
138 TT fBeta; // Relativistic beta factor.
139 Double32_t fDcaXY; // dca xy to the primary vertex
140 Double32_t fDcaZ; // dca z to the primary vertex
144 // PID data missing
145
146 TEveRecTrackT() : fLabel(-1), fIndex(-1), fStatus(0), fSign(0), fV(), fP(), fBeta(0), fDcaXY(0), fDcaZ(0), fPVX(0), fPVY(0), fPVZ(0) {}
147 virtual ~TEveRecTrackT() {}
148
149 Float_t Pt() { return fP.Perp(); }
150
151 ClassDef(TEveRecTrackT, 2); // Template for reconstructed track (also used in VSD).
152};
153
157
158/******************************************************************************/
159// TEveRecKink
160/******************************************************************************/
161
162class TEveRecKink : public TObject
163{
164public:
165
166 TEveVector fVKink; // Kink vertex: reconstructed position of the kink
167 TEveVector fPMother; // Momentum of the mother track
168 TEveVector fVMother; // Vertex of the mother track
169 TEveVector fPDaughter; // Momentum of the daughter track
170 TEveVector fVDaughter; // Vertex of the daughter track
171 Double32_t fKinkAngle[3]; // three angles
172 Int_t fSign; // sign of the track
173 Int_t fStatus; // Status as exported from reconstruction
174
175 // Data from simulation
176 Int_t fKinkLabel[2]; // Labels of the mother and daughter tracks
177 Int_t fKinkIndex[2]; // Indices of the mother and daughter tracks
178 Int_t fKinkPdg[2]; // PDG code of mother and daughter.
179
181 {
182 fKinkAngle[0] = fKinkAngle[1] = fKinkAngle[2] = 0;
183 fKinkLabel[0] = fKinkLabel[1] = 0;
184 fKinkIndex[0] = fKinkIndex[1] = 0;
185 fKinkPdg[0] = fKinkPdg[1] = 0;
186 }
187 virtual ~TEveRecKink() {}
188
189 ClassDef(TEveRecKink, 1); // Reconstructed kink (also used in VSD).
190};
191
192
193/******************************************************************************/
194// TEveRecV0
195/******************************************************************************/
196
197class TEveRecV0 : public TObject
198{
199public:
201
202 TEveVector fVNeg; // Vertex of negative track.
203 TEveVector fPNeg; // Momentum of negative track.
204 TEveVector fVPos; // Vertex of positive track.
205 TEveVector fPPos; // Momentum of positive track.
206
207 TEveVector fVCa; // Point of closest approach.
208 TEveVector fV0Birth; // Reconstucted birth point of neutral particle.
209
210 // ? Data from simulation.
211 Int_t fLabel; // Neutral mother label read from kinematics.
212 Int_t fPdg; // PDG code of mother.
213 Int_t fDLabel[2]; // Daughter labels.
214
216 fVCa(), fV0Birth(), fLabel(0), fPdg(0)
217 { fDLabel[0] = fDLabel[1] = 0; }
218 virtual ~TEveRecV0() {}
219
220 ClassDef(TEveRecV0, 1); // Reconstructed V0 (also used in VSD).
221};
222
223
224/******************************************************************************/
225// TEveRecCascade
226/******************************************************************************/
227
229{
230public:
232
233 TEveVector fVBac; // Vertex of bachelor track.
234 TEveVector fPBac; // Momentum of bachelor track.
235
236 TEveVector fCascadeVCa; // Point of closest approach for Cascade.
237 TEveVector fCascadeBirth; // Reconstucted birth point of cascade particle.
238
239 // ? Data from simulation.
240 Int_t fLabel; // Cascade mother label read from kinematics.
241 Int_t fPdg; // PDG code of mother.
242 Int_t fDLabel; // Daughter label.
243
246 fLabel(0), fPdg(0), fDLabel(0) {}
247 virtual ~TEveRecCascade() {}
248
249 ClassDef(TEveRecCascade, 1); // Reconstructed Cascade (also used in VSD).
250};
251
252
253/******************************************************************************/
254// TEveMCRecCrossRef
255/******************************************************************************/
256
258{
259public:
260 Bool_t fIsRec; // Is reconstructed.
266
267 TEveMCRecCrossRef() : fIsRec(false), fHasV0(false), fHasKink(false),
268 fLabel(0), fNHits(0), fNClus(0) {}
270
271 ClassDef(TEveMCRecCrossRef, 1); // Cross-reference of sim/rec data per particle (also used in VSD).
272};
273
274
275/******************************************************************************/
276// Missing primary vertex class.
277/******************************************************************************/
278
279
280/******************************************************************************/
281/******************************************************************************/
282
283// This whole construction is somewhat doubtable. It requires
284// shameless copying of experiment data. What is good about this
285// scheme:
286//
287// 1) Filters can be applied at copy time so that only part of the
288// data is copied over.
289//
290// 2) Once the data is extracted it can be used without experiment
291// software. Thus, external service can provide this data and local
292// client can be really thin.
293//
294// 3) Some pretty advanced visualization schemes/selections can be
295// implemented in a general framework by providing data extractors
296// only. This is also good for PR or VIP displays.
297//
298// 4) These classes can be extended by particular implementations. The
299// container classes will use TClonesArray with user-specified element
300// class.
301
302// The common behaviour could be implemented entirely without usage of
303// a common base classes, by just specifying names of members that
304// retrieve specific data. This is fine as long as one only uses tree
305// selections but becomes painful for extraction of data into local
306// structures (could a) use interpreter but this is an overkill and
307// would cause serious trouble for multi-threaded environment; b) use
308// member offsets and data-types from the dictionary).
309
310#endif
unsigned short UShort_t
Definition RtypesCore.h:40
float Float_t
Definition RtypesCore.h:57
constexpr Bool_t kFALSE
Definition RtypesCore.h:101
#define ClassDef(name, id)
Definition Rtypes.h:337
TEveRecTrackT< Float_t > TEveRecTrack
TEveRecTrackT< Double_t > TEveRecTrackD
TEveRecTrackT< Float_t > TEveRecTrackF
winID h TVirtualViewer3D TVirtualGLPainter p
UShort_t fDetId
virtual ~TEveCluster()
Int_t fLabel[3]
TEveVector fV
UShort_t fSubdetId
UShort_t fSubdetId
virtual ~TEveHit()
TEveVector fV
Int_t fLabel
UShort_t fDetId
Int_t fEvaLabel
virtual ~TEveMCRecCrossRef()
TEveVector fVDecay
TEveMCTrack & operator=(const TParticle &p)
TEveVector fPDecay
void ResetPdgCode()
Float_t fTDecay
virtual ~TEveMCTrack()
TEveVector fCascadeVCa
virtual ~TEveRecCascade()
TEveVector fCascadeBirth
Double32_t fKinkAngle[3]
Int_t fKinkPdg[2]
TEveVector fVMother
Int_t fKinkLabel[2]
TEveVector fPMother
TEveVector fVKink
TEveVector fPDaughter
virtual ~TEveRecKink()
Int_t fKinkIndex[2]
TEveVector fVDaughter
virtual ~TEveRecTrackT()
TEveVectorT< TT > fP
Double32_t fDcaZ
Double32_t fPVY
Double32_t fPVX
TEveVectorT< TT > fV
Double32_t fDcaXY
Double32_t fPVZ
virtual ~TEveRecV0()
TEveVector fPPos
TEveVector fVPos
TEveVector fVNeg
TEveVector fVCa
TEveVector fV0Birth
Int_t fDLabel[2]
TEveVector fPNeg
Mother of all ROOT objects.
Definition TObject.h:41
Description of the dynamic properties of a particle.
Definition TParticle.h:26
Int_t fPdgCode
Definition TParticle.h:31