CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FSimTrack.h
Go to the documentation of this file.
1 #ifndef FastSimulation_Event_FSimTrack_H
2 #define FastSimulation_Event_FSimTrack_H
3 
4 // HepPDT Headers
6 
7 // CMSSW Headers
9 
10 // FAMOS headers
12 
13 #include <vector>
14 
15 class FSimVertex;
16 class FBaseSimEvent;
17 
18 namespace HepMC {
19  class GenParticle;
20  class GenVertex;
21 }
22 
29 class FSimTrack : public SimTrack {
30 
31  public:
33  FSimTrack();
34 
36  FSimTrack(const RawParticle* p, int iv, int ig, int id, FBaseSimEvent* mom, double dt=-1.);
37 
39  virtual ~FSimTrack();
40 
42  inline const HepPDT::ParticleData* particleInfo() const {
43  return info_;
44  }
45 
47  inline float charge() const {
48  return particleInfo()->charge();
49  }
50 
51 
53  inline const FSimVertex& vertex() const;
54 
56  inline const FSimVertex& endVertex() const;
57 
59  inline const FSimTrack& mother() const;
60 
62  inline const FSimTrack& daughter(int i) const;
63 
65  inline int nDaughters() const;
66 
68  inline const std::vector<int>& daughters() const;
69 
71  inline bool noEndVertex() const;
72 
74  bool notYetToEndVertex(const XYZTLorentzVector& pos) const;
75 
77  inline bool noMother() const;
78 
80  inline bool noDaughter() const;
81 
83  inline const HepMC::GenParticle* genParticle() const;
84 
86  inline int id() const { return id_; }
87 
91  inline int onLayer1() const { return layer1; }
92 
96  inline int onLayer2() const { return layer2; }
97 
101  inline int onEcal() const { return ecal; }
102 
106  inline int onHcal() const { return hcal; }
107 
111  inline int onVFcal() const { return vfcal; }
112 
116  inline int outHcal() const { return hcalexit; }
117 
118  //The particle was propagated to the HO front face
121  inline int onHO() const { return hoentr; }
122 
124  inline bool propagated() const { return prop; }
125 
127  inline const RawParticle& layer1Entrance() const { return Layer1_Entrance; }
128 
130  inline const RawParticle& layer2Entrance() const { return Layer2_Entrance; }
131 
133  inline const RawParticle& ecalEntrance() const { return ECAL_Entrance; }
134 
136  inline const RawParticle& hcalEntrance() const { return HCAL_Entrance; }
137 
139  inline const RawParticle& vfcalEntrance() const { return VFCAL_Entrance; }
140 
142  inline const RawParticle& hcalExit() const { return HCAL_Exit; }
143 
145  inline const RawParticle& hoEntrance() const { return HO_Entrance; }
146 
148  inline void setEndVertex(int endv) { endv_ = endv; }
149 
151  void setPropagate();
152 
154  void setLayer1(const RawParticle& pp, int success);
155 
157  void setLayer2(const RawParticle& pp, int success);
158 
160  void setEcal(const RawParticle& pp,int success);
161 
163  void setHcal(const RawParticle& pp, int success);
164 
166  void setVFcal(const RawParticle& pp, int success);
167 
169  void setHcalExit(const RawParticle& pp, int success);
170 
172  void setHO(const RawParticle& pp, int success);
173 
175  // void addRecHit(const FamosBasicRecHit* hit, unsigned layer);
176 
178  // void addSimHit(const RawParticle& pp, unsigned layer);
179 
181  inline void addDaughter(int i) { daugh_.push_back(i); }
182 
184  inline void setClosestDaughterId(int id) { closestDaughterId_ = id; }
185 
187  inline int closestDaughterId() const { return closestDaughterId_; }
188 
190  const XYZTLorentzVector& momentum() const { return momentum_; }
191 
193  inline void setMomentum(const math::XYZTLorentzVector& newMomentum) {momentum_ = newMomentum; }
194 
196  inline const SimTrack& simTrack() const { return *this; }
197 
199  inline double decayTime() const { return properDecayTime; }
200 
201  private:
202 
203  // HepMC::GenParticle* me_;
204 
206  // int embd_; // The index in the SimTrack vector
207  int id_; // The index in the FSimTrackVector
208 
209  int endv_; // The index of the end vertex in FSimVertex
210 
211  int layer1;// 1 if the particle was propagated to preshower layer1
212  int layer2;// 1 if the particle was propagated to preshower layer2
213  int ecal; // 1 if the particle was propagated to ECAL/HCAL barrel
214  int hcal; // 2 if the particle was propagated to ECAL/HCAL endcap
215  int vfcal; // 1 if the particle was propagated to VFCAL
216  int hcalexit; // 2 if the particle was propagated to HCAL Exit point
217  int hoentr; // 1 if the particle was propagated to HO
218 
219 
220  bool prop; // true if the propagation to the calorimeters was done
221 
222  RawParticle Layer1_Entrance; // the particle at preshower Layer1
223  RawParticle Layer2_Entrance; // the particle at preshower Layer2
224  RawParticle ECAL_Entrance; // the particle at ECAL entrance
225  RawParticle HCAL_Entrance; // the particle at HCAL entrance
226  RawParticle VFCAL_Entrance; // the particle at VFCAL entrance
227  RawParticle HCAL_Exit; // the particle at HCAL ezit point
228  RawParticle HO_Entrance; // the particle at HO entrance
229 
230 
231  std::vector<int> daugh_; // The indices of the daughters in FSimTrack
232  int closestDaughterId_; // The index of the closest daughter id
233 
234  const HepPDT::ParticleData* info_; // The PDG info
235 
237 
238  double properDecayTime; // The proper decay time (default is -1)
239 
240 };
241 
242 #include<iosfwd>
243 std::ostream& operator <<(std::ostream& o , const FSimTrack& t);
244 
245 #include "FastSimulation/Event/interface/FSimTrack.icc"
246 
247 
248 
249 #endif // FSimTrack_H
250 
251 
252 
253 
254 
tuple t
Definition: tree.py:139
float dt
Definition: AMPTWrapper.h:126
int i
Definition: DBlmapReader.cc:9
bool noEndVertex() const
no end vertex
int endv_
Definition: FSimTrack.h:209
float charge() const
charge
Definition: FSimTrack.h:47
const SimTrack & simTrack() const
Simply returns the SimTrack.
Definition: FSimTrack.h:196
const HepPDT::ParticleData * particleInfo() const
particle info...
Definition: FSimTrack.h:42
const RawParticle & vfcalEntrance() const
The particle at VFCAL entrance.
Definition: FSimTrack.h:139
const RawParticle & hoEntrance() const
The particle at HCAL exir.
Definition: FSimTrack.h:145
int onLayer2() const
Definition: FSimTrack.h:96
int layer2
Definition: FSimTrack.h:212
const FSimVertex & endVertex() const
end vertex
const RawParticle & hcalExit() const
The particle at HCAL exir.
Definition: FSimTrack.h:142
int ecal
Definition: FSimTrack.h:213
tuple pp
Definition: createTree.py:15
RawParticle Layer1_Entrance
Definition: FSimTrack.h:222
const HepPDT::ParticleData * info_
Definition: FSimTrack.h:234
RawParticle VFCAL_Entrance
Definition: FSimTrack.h:226
const XYZTLorentzVector & momentum() const
Temporary (until move of SimTrack to Mathcore) - No! Actually very useful.
Definition: FSimTrack.h:190
virtual ~FSimTrack()
Destructor.
Definition: FSimTrack.cc:30
int hoentr
Definition: FSimTrack.h:217
const RawParticle & layer1Entrance() const
The particle at Preshower Layer 1.
Definition: FSimTrack.h:127
const FSimTrack & daughter(int i) const
Ith daughter.
void setEcal(const RawParticle &pp, int success)
Set the ecal variables.
Definition: FSimTrack.cc:69
RawParticle HO_Entrance
Definition: FSimTrack.h:228
void setPropagate()
The particle has been propgated through the tracker.
Definition: FSimTrack.cc:49
RawParticle HCAL_Exit
Definition: FSimTrack.h:227
int nDaughters() const
Number of daughters.
void setClosestDaughterId(int id)
Set the index of the closest charged daughter.
Definition: FSimTrack.h:184
void setLayer2(const RawParticle &pp, int success)
Set the preshower layer2 variables.
Definition: FSimTrack.cc:62
const std::vector< int > & daughters() const
Vector of daughter indices.
int outHcal() const
Definition: FSimTrack.h:116
std::ostream & operator<<(std::ostream &out, const ALILine &li)
Definition: ALILine.cc:187
bool notYetToEndVertex(const XYZTLorentzVector &pos) const
Compare the end vertex position with another position.
Definition: FSimTrack.cc:33
int onEcal() const
Definition: FSimTrack.h:101
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
void setLayer1(const RawParticle &pp, int success)
Set the preshower layer1 variables.
Definition: FSimTrack.cc:55
bool propagated() const
The particle was tentatively propagated to calorimeters.
Definition: FSimTrack.h:124
int closestDaughterId_
Definition: FSimTrack.h:232
int closestDaughterId() const
Get the index of the closest charged daughter.
Definition: FSimTrack.h:187
const RawParticle & ecalEntrance() const
The particle at ECAL entrance.
Definition: FSimTrack.h:133
int onHO() const
Definition: FSimTrack.h:121
int onVFcal() const
Definition: FSimTrack.h:111
const FSimVertex & vertex() const
Origin vertex.
int vfcal
Definition: FSimTrack.h:215
FSimTrack()
Default constructor.
Definition: FSimTrack.cc:9
HepPDT::ParticleData ParticleData
int layer1
Definition: FSimTrack.h:211
FBaseSimEvent * mom_
Definition: FSimTrack.h:205
const HepMC::GenParticle * genParticle() const
The original GenParticle.
RawParticle HCAL_Entrance
Definition: FSimTrack.h:225
void setHcalExit(const RawParticle &pp, int success)
Set the hcal exit variables.
Definition: FSimTrack.cc:90
void setEndVertex(int endv)
Set the end vertex.
Definition: FSimTrack.h:148
XYZTLorentzVector momentum_
Definition: FSimTrack.h:236
int id() const
the index in FBaseSimEvent and other vectors
Definition: FSimTrack.h:86
void addDaughter(int i)
Add a RecHit for a track on a layer.
Definition: FSimTrack.h:181
int hcal
Definition: FSimTrack.h:214
const RawParticle & layer2Entrance() const
The particle at Preshower Layer 2.
Definition: FSimTrack.h:130
int onLayer1() const
Definition: FSimTrack.h:91
void setHO(const RawParticle &pp, int success)
Set the ho variables.
Definition: FSimTrack.cc:96
int onHcal() const
Definition: FSimTrack.h:106
void setVFcal(const RawParticle &pp, int success)
Set the hcal variables.
Definition: FSimTrack.cc:83
bool noMother() const
no mother particle
const RawParticle & hcalEntrance() const
The particle at HCAL entrance.
Definition: FSimTrack.h:136
bool noDaughter() const
no daughters
RawParticle Layer2_Entrance
Definition: FSimTrack.h:223
void setHcal(const RawParticle &pp, int success)
Set the hcal variables.
Definition: FSimTrack.cc:76
double properDecayTime
Definition: FSimTrack.h:238
int hcalexit
Definition: FSimTrack.h:216
bool prop
Definition: FSimTrack.h:220
RawParticle ECAL_Entrance
Definition: FSimTrack.h:224
void setMomentum(const math::XYZTLorentzVector &newMomentum)
Reset the momentum (to be used with care)
Definition: FSimTrack.h:193
double decayTime() const
Return the pre-defined decay time.
Definition: FSimTrack.h:199
const FSimTrack & mother() const
mother
std::vector< int > daugh_
Definition: FSimTrack.h:231
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:15