CMS 3D CMS Logo

FSimTrack.cc
Go to the documentation of this file.
2 
3 //C++ Headers
4 #include <iomanip>
5 #include <string>
6 
7 //using namespace HepPDT;
8 
10  SimTrack(), mom_(0), id_(-1), endv_(-1),
11  layer1(0), layer2(0), ecal(0), hcal(0), vfcal(0), hcalexit(0), hoentr(0),
12  prop(false), closestDaughterId_(-1), info_(0),
13  properDecayTime(1E99) {;}
14 
16  int iv, int ig, int id,
17  FBaseSimEvent* mom,
18  double dt) :
19  // SimTrack(p->pid(),*p,iv,ig), // to uncomment once Mathcore is installed
20  SimTrack(p->pid(),p->momentum(),iv,ig),
21  mom_(mom), id_(id), endv_(-1),
22  layer1(0), layer2(0), ecal(0), hcal(0), vfcal(0), hcalexit(0), hoentr(0), prop(false),
24  properDecayTime(dt)
25 {
26  setTrackId(id);
27  info_ = mom_->theTable()->particle(HepPDT::ParticleID(type()));
28 }
29 
31 
32 bool
34  // If there is no end vertex, nothing to compare to
35  if ( noEndVertex() ) return true;
36  // If the particle immediately decays, no need to propagate
37  if ( (endVertex().position()-vertex().position()).Vect().Mag2() < 1e-4 )
38  return false;
39  // If the end vertex has a larger radius, not yet there
40  if ( endVertex().position().Perp2() > pos.Perp2()+1e-10 ) return true;
41  // If the end vertex has a larger z, not yet there
42  if ( fabs(endVertex().position().Z()) > fabs(pos.Z())+1e-5 ) return true;
43  // Otherwise, the end vertex is overtaken already
44  return false;
45 }
46 
48 void
50  prop=true;
51 }
52 
54 void
57  layer1=success;
58 }
59 
61 void
64  layer2=success;
65 }
66 
68 void
71  ecal=success;
72 }
73 
75 void
78  hcal=success;
79 }
80 
82 void
85  vfcal=success;
86 }
87 
89 void
91  HCAL_Exit=pp;
93 }
94 
95 void
97  HO_Entrance=pp;
98  hoentr=success;
99 }
100 
101 
102 
103 
104 
105 
106 std::ostream& operator <<(std::ostream& o , const FSimTrack& t) {
107 
108  std::string name = t.particleInfo() ? t.particleInfo()->name() : "Unknown";
109  XYZTLorentzVector momentum1 = t.momentum();
110  XYZVector vertex1 = t.vertex().position().Vect();
111  int vertexId1 = t.vertex().id();
112 
113  o.setf(std::ios::fixed, std::ios::floatfield);
114  o.setf(std::ios::right, std::ios::adjustfield);
115 
116  o << std::setw(4) << t.id() << " "
117  << std::setw(4) << t.genpartIndex() << " "
118  << name;
119 
120  for(unsigned int k=0;k<11-name.length() && k<12; k++) o << " ";
121 
122  o << std::setw(6) << std::setprecision(2) << momentum1.eta() << " "
123  << std::setw(6) << std::setprecision(2) << momentum1.phi() << " "
124  << std::setw(6) << std::setprecision(2) << momentum1.pt() << " "
125  << std::setw(6) << std::setprecision(2) << momentum1.e() << " "
126  << std::setw(4) << vertexId1 << " "
127  << std::setw(6) << std::setprecision(1) << vertex1.x() << " "
128  << std::setw(6) << std::setprecision(1) << vertex1.y() << " "
129  << std::setw(6) << std::setprecision(1) << vertex1.z() << " "
130  << std::setw(4) << (t.noMother() ? -1 :t.mother().id()) << " ";
131 
132  if ( !t.noEndVertex() ) {
133  XYZTLorentzVector vertex2 = t.endVertex().position();
134  int vertexId2 = t.endVertex().id();
135 
136  o << std::setw(4) << vertexId2 << " "
137  << std::setw(6) << std::setprecision(2) << vertex2.eta() << " "
138  << std::setw(6) << std::setprecision(2) << vertex2.phi() << " "
139  << std::setw(5) << std::setprecision(1) << vertex2.pt() << " "
140  << std::setw(6) << std::setprecision(1) << vertex2.z() << " ";
141  for (int i=0; i<t.nDaughters(); ++i)
142  o << std::setw(4) << t.daughter(i).id() << " ";
143 
144  } else {
145 
146  if ( t.onLayer1() ) {
147 
148  XYZTLorentzVector vertex2 = t.layer1Entrance().vertex();
149 
150  o << std::setw(4) << -t.onLayer1() << " "
151  << std::setw(6) << std::setprecision(2) << vertex2.eta() << " "
152  << std::setw(6) << std::setprecision(2) << vertex2.phi() << " "
153  << std::setw(5) << std::setprecision(1) << vertex2.pt() << " "
154  << std::setw(6) << std::setprecision(1) << vertex2.z() << " "
155  << std::setw(6) << std::setprecision(2) << t.layer1Entrance().pt() << " "
156  << std::setw(6) << std::setprecision(2) << t.layer1Entrance().e() << " ";
157 
158  } else if ( t.onEcal() ) {
159 
160  XYZTLorentzVector vertex2 = t.ecalEntrance().vertex();
161 
162  o << std::setw(4) << -t.onEcal() << " "
163  << std::setw(6) << std::setprecision(2) << vertex2.eta() << " "
164  << std::setw(6) << std::setprecision(2) << vertex2.phi() << " "
165  << std::setw(5) << std::setprecision(1) << vertex2.pt() << " "
166  << std::setw(6) << std::setprecision(1) << vertex2.z() << " "
167  << std::setw(6) << std::setprecision(2) << t.ecalEntrance().pt() << " "
168  << std::setw(6) << std::setprecision(2) << t.ecalEntrance().e() << " ";
169  }
170  }
171  return o;
172 }
const double Z[kNumberCalorimeter]
int id() const
the index in FBaseSimEvent
Definition: FSimVertex.h:44
float dt
Definition: AMPTWrapper.h:126
int i
Definition: DBlmapReader.cc:9
bool noEndVertex() const
no end vertex
int endv_
Definition: FSimTrack.h:209
const HepPDT::ParticleData * particleInfo() const
particle info...
Definition: FSimTrack.h:42
int layer2
Definition: FSimTrack.h:212
const FSimVertex & endVertex() const
end vertex
int ecal
Definition: FSimTrack.h:213
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
const HepPDT::ParticleDataTable * theTable() const
Get the pointer to the particle data table.
Definition: FBaseSimEvent.h:57
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 setLayer2(const RawParticle &pp, int success)
Set the preshower layer2 variables.
Definition: FSimTrack.cc:62
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
void setLayer1(const RawParticle &pp, int success)
Set the preshower layer1 variables.
Definition: FSimTrack.cc:55
math::XYZVector XYZVector
const math::XYZTLorentzVector & position() const
Temporary (until CMSSW moves to Mathcore) - No ! Actually very useful.
Definition: FSimVertex.h:49
int genpartIndex() const
index of the corresponding Generator particle in the Event container (-1 if no Genpart) ...
Definition: SimTrack.h:33
int closestDaughterId_
Definition: FSimTrack.h:232
const RawParticle & ecalEntrance() const
The particle at ECAL entrance.
Definition: FSimTrack.h:133
const FSimVertex & vertex() const
Origin vertex.
int vfcal
Definition: FSimTrack.h:215
FSimTrack()
Default constructor.
Definition: FSimTrack.cc:9
int layer1
Definition: FSimTrack.h:211
int k[5][pyjets_maxn]
FBaseSimEvent * mom_
Definition: FSimTrack.h:205
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:285
RawParticle HCAL_Entrance
Definition: FSimTrack.h:225
void setTrackId(unsigned int t)
Definition: CoreSimTrack.h:33
void setHcalExit(const RawParticle &pp, int success)
Set the hcal exit variables.
Definition: FSimTrack.cc:90
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:25
XYZTLorentzVector momentum_
Definition: FSimTrack.h:236
static int position[264][3]
Definition: ReadPGInfo.cc:509
int id() const
the index in FBaseSimEvent and other vectors
Definition: FSimTrack.h:86
int onLayer1() const
Definition: FSimTrack.h:91
void setHO(const RawParticle &pp, int success)
Set the ho variables.
Definition: FSimTrack.cc:96
void setVFcal(const RawParticle &pp, int success)
Set the hcal variables.
Definition: FSimTrack.cc:83
bool noMother() const
no mother particle
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
const FSimTrack & mother() const
mother
std::ostream & operator<<(std::ostream &o, const FSimTrack &t)
Definition: FSimTrack.cc:106
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:15