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 
114  inline bool propagated() const { return prop; }
115 
117  inline const RawParticle& layer1Entrance() const { return Layer1_Entrance; }
118 
120  inline const RawParticle& layer2Entrance() const { return Layer2_Entrance; }
121 
123  inline const RawParticle& ecalEntrance() const { return ECAL_Entrance; }
124 
126  inline const RawParticle& hcalEntrance() const { return HCAL_Entrance; }
127 
129  inline const RawParticle& vfcalEntrance() const { return VFCAL_Entrance; }
130 
132  inline void setEndVertex(int endv) { endv_ = endv; }
133 
135  void setPropagate();
136 
138  void setLayer1(const RawParticle& pp, int success);
139 
141  void setLayer2(const RawParticle& pp, int success);
142 
144  void setEcal(const RawParticle& pp,int success);
145 
147  void setHcal(const RawParticle& pp, int success);
148 
150  void setVFcal(const RawParticle& pp, int success);
151 
153  // void addRecHit(const FamosBasicRecHit* hit, unsigned layer);
154 
156  // void addSimHit(const RawParticle& pp, unsigned layer);
157 
159  inline void addDaughter(int i) { daugh_.push_back(i); }
160 
162  inline void setClosestDaughterId(int id) { closestDaughterId_ = id; }
163 
165  inline int closestDaughterId() const { return closestDaughterId_; }
166 
168  const XYZTLorentzVector& momentum() const { return momentum_; }
169 
171  inline void setMomentum(const math::XYZTLorentzVector& newMomentum) {momentum_ = newMomentum; }
172 
174  inline const SimTrack& simTrack() const { return *this; }
175 
177  inline double decayTime() const { return properDecayTime; }
178 
179  private:
180 
181  // HepMC::GenParticle* me_;
182 
184  // int embd_; // The index in the SimTrack vector
185  int id_; // The index in the FSimTrackVector
186 
187  int endv_; // The index of the end vertex in FSimVertex
188 
189  int layer1;// 1 if the particle was propagated to preshower layer1
190  int layer2;// 1 if the particle was propagated to preshower layer2
191  int ecal; // 1 if the particle was propagated to ECAL/HCAL barrel
192  int hcal; // 2 if the particle was propagated to ECAL/HCAL endcap
193  int vfcal; // 1 if the particle was propagated to VFCAL
194 
195  bool prop; // true if the propagation to the calorimeters was done
196 
197  RawParticle Layer1_Entrance; // the particle at preshower Layer1
198  RawParticle Layer2_Entrance; // the particle at preshower Layer2
199  RawParticle ECAL_Entrance; // the particle at ECAL entrance
200  RawParticle HCAL_Entrance; // the particle at HCAL entrance
201  RawParticle VFCAL_Entrance; // the particle at VFCAL entrance
202 
203  std::vector<int> daugh_; // The indices of the daughters in FSimTrack
204  int closestDaughterId_; // The index of the closest daughter id
205 
206  const HepPDT::ParticleData* info_; // The PDG info
207 
209 
210  double properDecayTime; // The proper decay time (default is -1)
211 
212 };
213 
214 #include<iosfwd>
215 std::ostream& operator <<(std::ostream& o , const FSimTrack& t);
216 
217 #include "FastSimulation/Event/interface/FSimTrack.icc"
218 
219 
220 
221 #endif // FSimTrack_H
222 
223 
224 
225 
226 
float dt
Definition: AMPTWrapper.h:126
int i
Definition: DBlmapReader.cc:9
bool noEndVertex() const
no end vertex
int endv_
Definition: FSimTrack.h:187
float charge() const
charge
Definition: FSimTrack.h:47
const SimTrack & simTrack() const
Simply returns the SimTrack.
Definition: FSimTrack.h:174
const HepPDT::ParticleData * particleInfo() const
particle info...
Definition: FSimTrack.h:42
const RawParticle & vfcalEntrance() const
The particle at VFCAL entrance.
Definition: FSimTrack.h:129
int onLayer2() const
Definition: FSimTrack.h:96
int layer2
Definition: FSimTrack.h:190
const FSimVertex & endVertex() const
end vertex
int ecal
Definition: FSimTrack.h:191
tuple pp
Definition: createTree.py:15
RawParticle Layer1_Entrance
Definition: FSimTrack.h:197
const HepPDT::ParticleData * info_
Definition: FSimTrack.h:206
RawParticle VFCAL_Entrance
Definition: FSimTrack.h:201
const XYZTLorentzVector & momentum() const
Temporary (until move of SimTrack to Mathcore) - No! Actually very useful.
Definition: FSimTrack.h:168
virtual ~FSimTrack()
Destructor.
Definition: FSimTrack.cc:30
const RawParticle & layer1Entrance() const
The particle at Preshower Layer 1.
Definition: FSimTrack.h:117
const FSimTrack & daughter(int i) const
Ith daughter.
void setEcal(const RawParticle &pp, int success)
Set the ecal variables.
Definition: FSimTrack.cc:69
void setPropagate()
The particle has been propgated through the tracker.
Definition: FSimTrack.cc:49
int nDaughters() const
Number of daughters.
void setClosestDaughterId(int id)
Set the index of the closest charged daughter.
Definition: FSimTrack.h:162
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.
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:30
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:114
int closestDaughterId_
Definition: FSimTrack.h:204
int closestDaughterId() const
Get the index of the closest charged daughter.
Definition: FSimTrack.h:165
const RawParticle & ecalEntrance() const
The particle at ECAL entrance.
Definition: FSimTrack.h:123
int onVFcal() const
Definition: FSimTrack.h:111
const FSimVertex & vertex() const
Origin vertex.
int vfcal
Definition: FSimTrack.h:193
FSimTrack()
Default constructor.
Definition: FSimTrack.cc:9
HepPDT::ParticleData ParticleData
int layer1
Definition: FSimTrack.h:189
FBaseSimEvent * mom_
Definition: FSimTrack.h:183
const HepMC::GenParticle * genParticle() const
The original GenParticle.
RawParticle HCAL_Entrance
Definition: FSimTrack.h:200
void setEndVertex(int endv)
Set the end vertex.
Definition: FSimTrack.h:132
XYZTLorentzVector momentum_
Definition: FSimTrack.h:208
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:159
int hcal
Definition: FSimTrack.h:192
const RawParticle & layer2Entrance() const
The particle at Preshower Layer 2.
Definition: FSimTrack.h:120
int onLayer1() const
Definition: FSimTrack.h:91
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:126
bool noDaughter() const
no daughters
RawParticle Layer2_Entrance
Definition: FSimTrack.h:198
void setHcal(const RawParticle &pp, int success)
Set the hcal variables.
Definition: FSimTrack.cc:76
double properDecayTime
Definition: FSimTrack.h:210
bool prop
Definition: FSimTrack.h:195
RawParticle ECAL_Entrance
Definition: FSimTrack.h:199
void setMomentum(const math::XYZTLorentzVector &newMomentum)
Reset the momentum (to be used with care)
Definition: FSimTrack.h:171
double decayTime() const
Return the pre-defined decay time.
Definition: FSimTrack.h:177
const FSimTrack & mother() const
mother
std::vector< int > daugh_
Definition: FSimTrack.h:203
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:15