CMS 3D CMS Logo

Classes | Namespaces | Functions
FSimTrack.h File Reference
#include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
#include "SimDataFormats/Track/interface/SimTrack.h"
#include "SimDataFormats/Vertex/interface/SimVertex.h"
#include "CommonTools/BaseParticlePropagator/interface/RawParticle.h"
#include <vector>
#include <iosfwd>
#include "FastSimulation/Event/interface/FSimTrack.icc"

Go to the source code of this file.

Classes

class  FSimTrack
 

Namespaces

 HepMC
 

Functions

std::ostream & operator<< (std::ostream &o, const FSimTrack &t)
 

Function Documentation

std::ostream& operator<< ( std::ostream &  o,
const FSimTrack t 
)

Definition at line 118 of file FSimTrack.cc.

References FSimTrack::daughter(), RawParticle::e(), FSimTrack::ecalEntrance(), FSimTrack::endVertex(), alignBH_cfg::fixed, SimTrack::genpartIndex(), FSimTrack::hcalEntrance(), mps_fire::i, FSimVertex::id(), FSimTrack::id(), FSimTrack::layer1Entrance(), FSimTrack::layer2Entrance(), FSimTrack::momentum(), FSimTrack::mother(), dataset::name, FSimTrack::nDaughters(), FSimTrack::noEndVertex(), FSimTrack::noMother(), connectstrParser::o, FSimTrack::onEcal(), FSimTrack::onHcal(), FSimTrack::onLayer1(), FSimTrack::onLayer2(), FSimTrack::onVFcal(), FSimTrack::particleInfo(), FSimVertex::position(), RawParticle::pt(), AlCaHLTBitMon_QueryRunRegistry::string, FSimTrack::vertex(), RawParticle::vertex(), and FSimTrack::vfcalEntrance().

118  {
119  std::string name = t.particleInfo() ? t.particleInfo()->name() : "Unknown";
120  const XYZTLorentzVector& momentum1 = t.momentum();
121  XYZVector vertex1 = t.vertex().position().Vect();
122  int vertexId1 = t.vertex().id();
123 
124  o.setf(std::ios::fixed, std::ios::floatfield);
125  o.setf(std::ios::right, std::ios::adjustfield);
126 
127  o << std::setw(4) << t.id() << "; "
128  << std::setw(4) << t.genpartIndex() << "; "
129  << name << "\n\t\t";
130 
131  //for(unsigned int k=0;k<11-name.length() && k<12; k++) o << "; ";
132 
133  o << "Position: "
134  << std::setw(6) << std::setprecision(2) << momentum1.eta() << "; "
135  << std::setw(6) << std::setprecision(2) << momentum1.phi() << "; "
136  << std::setw(6) << std::setprecision(2) << momentum1.pt() << "; "
137  << std::setw(6) << std::setprecision(2) << momentum1.e() << "; "
138  << std::setw(4) << vertexId1 << "; "
139  << std::setw(6) << std::setprecision(1) << vertex1.x() << "; "
140  << std::setw(6) << std::setprecision(1) << vertex1.y() << "; "
141  << std::setw(6) << std::setprecision(1) << vertex1.z() << "; "
142  << std::setw(4) << (t.particleInfo() ? (t.noMother() ? -1 :t.mother().id()) :-1) << "; " << "\n\t\t";
143 
144  if ( !t.noEndVertex() ) {
145  XYZTLorentzVector vertex2 = t.endVertex().position();
146  int vertexId2 = t.endVertex().id();
147 
148  o << "Decayed particle: "
149  << std::setw(4) << vertexId2 << "; "
150  << std::setw(6) << std::setprecision(2) << vertex2.eta() << "; "
151  << std::setw(6) << std::setprecision(2) << vertex2.phi() << "; "
152  << std::setw(5) << std::setprecision(1) << vertex2.pt() << "; "
153  << std::setw(6) << std::setprecision(1) << vertex2.z() << "; ";
154  for (int i=0; i<t.nDaughters(); ++i)
155  o << std::setw(4) << t.daughter(i).id() << "; " << "\n\t\t";
156 
157  } else {
158 
159  if ( t.onLayer1() ) {
160 
161  XYZTLorentzVector vertex2 = t.layer1Entrance().vertex();
162 
163  o << "Layer 1: "
164  << std::setw(4) << -t.onLayer1() << "; "
165  << std::setw(6) << std::setprecision(2) << vertex2.eta() << "; "
166  << std::setw(6) << std::setprecision(2) << vertex2.phi() << "; "
167  << std::setw(5) << std::setprecision(1) << vertex2.pt() << "; "
168  << std::setw(6) << std::setprecision(1) << vertex2.z() << "; "
169  << std::setw(6) << std::setprecision(2) << t.layer1Entrance().pt() << "; "
170  << std::setw(6) << std::setprecision(2) << t.layer1Entrance().e() << "; " << "\n\t\t";
171 
172  } if ( t.onLayer2() ) {
173 
174  XYZTLorentzVector vertex2 = t.layer2Entrance().vertex();
175 
176  o << "Layer 2: "
177  << std::setw(4) << -t.onLayer2() << "; "
178  << std::setw(6) << std::setprecision(2) << vertex2.eta() << "; "
179  << std::setw(6) << std::setprecision(2) << vertex2.phi() << "; "
180  << std::setw(5) << std::setprecision(1) << vertex2.pt() << "; "
181  << std::setw(6) << std::setprecision(1) << vertex2.z() << "; "
182  << std::setw(6) << std::setprecision(2) << t.layer2Entrance().pt() << "; "
183  << std::setw(6) << std::setprecision(2) << t.layer2Entrance().e() << "; " << "\n\t\t";
184 
185  }
186  //if ( t.onEcal() ) {
187 
188  XYZTLorentzVector vertex2 = t.ecalEntrance().vertex();
189 
190  o << "ECAL: "
191  << std::setw(4) << -t.onEcal() << "; "
192  << std::setw(6) << std::setprecision(2) << vertex2.eta() << "; "
193  << std::setw(6) << std::setprecision(2) << vertex2.phi() << "; "
194  << std::setw(5) << std::setprecision(1) << vertex2.pt() << "; "
195  << std::setw(6) << std::setprecision(1) << vertex2.z() << "; "
196  << std::setw(6) << std::setprecision(2) << t.ecalEntrance().pt() << "; "
197  << std::setw(6) << std::setprecision(2) << t.ecalEntrance().e() << "; " << "\n\t\t";
198  //}
199  //if ( t.onHcal() ) {
200 
201  vertex2 = t.hcalEntrance().vertex();
202 
203  o << "HCAL: "
204  << std::setw(4) << -t.onHcal() << "; "
205  << std::setw(6) << std::setprecision(2) << vertex2.eta() << "; "
206  << std::setw(6) << std::setprecision(2) << vertex2.phi() << "; "
207  << std::setw(5) << std::setprecision(1) << vertex2.pt() << "; "
208  << std::setw(6) << std::setprecision(1) << vertex2.z() << "; "
209  << std::setw(6) << std::setprecision(2) << t.hcalEntrance().pt() << "; "
210  << std::setw(6) << std::setprecision(2) << t.hcalEntrance().e() << "; " << "\n\t\t";
211  //}
212  //if ( t.onVFcal() ) {
213 
214  vertex2 = t.vfcalEntrance().vertex();
215 
216  o << "VFCAL: "
217  << std::setw(4) << -t.onVFcal() << "; "
218  << std::setw(6) << std::setprecision(2) << vertex2.eta() << "; "
219  << std::setw(6) << std::setprecision(2) << vertex2.phi() << "; "
220  << std::setw(5) << std::setprecision(1) << vertex2.pt() << "; "
221  << std::setw(6) << std::setprecision(1) << vertex2.z() << "; "
222  << std::setw(6) << std::setprecision(2) << t.vfcalEntrance().pt() << "; "
223  << std::setw(6) << std::setprecision(2) << t.vfcalEntrance().e() << "; " << "\n\t\t";
224  //}
225  }
226  return o;
227 }
int id() const
the index in FBaseSimEvent
Definition: FSimVertex.h:44
bool noEndVertex() const
no end vertex
const HepPDT::ParticleData * particleInfo() const
particle info...
Definition: FSimTrack.h:46
double pt() const
transverse momentum
Definition: RawParticle.h:328
const RawParticle & vfcalEntrance() const
The particle at VFCAL entrance.
Definition: FSimTrack.h:144
int onLayer2() const
Definition: FSimTrack.h:101
const FSimVertex & endVertex() const
end vertex
const XYZTLorentzVector & momentum() const
Temporary (until move of SimTrack to Mathcore) - No! Actually very useful.
Definition: FSimTrack.h:204
const RawParticle & layer1Entrance() const
The particle at Preshower Layer 1.
Definition: FSimTrack.h:132
const FSimTrack & daughter(int i) const
Ith daughter.
int nDaughters() const
Number of daughters.
int onEcal() const
Definition: FSimTrack.h:106
double e() const
energy of the momentum
Definition: RawParticle.h:324
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:34
const RawParticle & ecalEntrance() const
The particle at ECAL entrance.
Definition: FSimTrack.h:138
int onVFcal() const
Definition: FSimTrack.h:116
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:339
int id() const
the index in FBaseSimEvent and other vectors
Definition: FSimTrack.h:91
const RawParticle & layer2Entrance() const
The particle at Preshower Layer 2.
Definition: FSimTrack.h:135
int onLayer1() const
Definition: FSimTrack.h:96
int onHcal() const
Definition: FSimTrack.h:111
bool noMother() const
no mother particle
const RawParticle & hcalEntrance() const
The particle at HCAL entrance.
Definition: FSimTrack.h:141
const FSimVertex vertex() const
Origin vertex.
math::XYZVector XYZVector
Definition: RawParticle.h:28
const FSimTrack & mother() const
mother
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:27