CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FSimEvent.cc
Go to the documentation of this file.
1 //FAMOS Headers
3 
4 //C++ Headers
5 
7  : FBaseSimEvent(kine), id_(edm::EventID(0,0,0)), weight_(0)
8 {}
9 
11 {}
12 
13 void
16  id_ = Id;
17 }
18 
19 void
20 FSimEvent::fill(const std::vector<SimTrack>& simTracks,
21  const std::vector<SimVertex>& simVertices) {
22  FBaseSimEvent::fill(simTracks,simVertices);
23  id_ = edm::EventID();
24 }
25 
27 FSimEvent::id() const {
28  return id_;
29 }
30 
31 float FSimEvent::weight() const {
32  return weight_;
33 }
34 
35 unsigned int
37  return FBaseSimEvent::nTracks();
38 }
39 
40 unsigned int
42  return FBaseSimEvent::nVertices();
43 }
44 
45 unsigned int
47  return FBaseSimEvent::nGenParts();
48 }
49 
50 void
52 {
53  for (unsigned int i=0; i<nTracks(); ++i) {
54  // SimTrack t = SimTrack(ip,p,iv,ig);
55  const SimTrack& t = embdTrack(i);
56  // Save all tracks
57  c.push_back(t);
58  // Save also some muons for later parameterization
59  if ( abs(t.type()) == 13 &&
60  t.momentum().perp2() > 1.0 &&
61  fabs(t.momentum().eta()) < 3.0 &&
62  track(i).noEndVertex() ) {
63  // Actually save the muon mother (and the attached muon) in case
64  if ( !track(i).noMother() && track(i).mother().closestDaughterId() == (int)i ) {
65  const SimTrack& T = embdTrack(track(i).mother().id());
66  m.push_back(T);
67  }
68  m.push_back(t);
69  }
70  }
71 }
72 
73 void
75 {
76  for (unsigned int i=0; i<nVertices(); ++i) {
77  // SimTrack t = SimTrack(ip,p,iv,ig);
78  c.push_back(embdVertex(i));
79  }
80 }
81 
82 
83 void
85 {
86 
87  for (unsigned int i=0; i<nVertices(); ++i) {
88  c.push_back(embdVertexType(i));
89  }
90 }
91 
92 
93 
94 
95 
96 
97 
int i
Definition: DBlmapReader.cc:9
bool noEndVertex() const
no end vertex
void load(edm::SimTrackContainer &c, edm::SimTrackContainer &m) const
Load containers of tracks (and muons) and vertices for the edm::Event.
Definition: FSimEvent.cc:51
void fill(const HepMC::GenEvent &hev, edm::EventID &Id)
fill the FBaseSimEvent from the current HepMC::GenEvent
Definition: FSimEvent.cc:14
FSimEvent(const edm::ParameterSet &kine)
Default constructor.
Definition: FSimEvent.cc:6
double weight_
Definition: FSimEvent.h:67
unsigned int nVertices() const
Number of vertices.
Definition: FBaseSimEvent.h:89
const FSimVertexType & embdVertexType(int i) const
return embedded vertex type with given id
unsigned int nTracks() const
Number of tracks.
Definition: FBaseSimEvent.h:84
int closestDaughterId() const
Get the index of the closest charged daughter.
Definition: FSimTrack.h:187
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< FSimVertexType > FSimVertexTypeCollection
collection of FSimVertexType objects
unsigned int nTracks() const
Number of tracks.
Definition: FSimEvent.cc:36
void fill(const HepMC::GenEvent &hev)
fill the FBaseSimEvent from the current HepMC::GenEvent
edm::EventID id() const
Method to return the EventId.
Definition: FSimEvent.cc:27
unsigned int nGenParts() const
Number of MC particles.
Definition: FSimEvent.cc:46
const SimVertex & embdVertex(int i) const
return embedded vertex with given id
virtual ~FSimEvent()
usual virtual destructor
Definition: FSimEvent.cc:10
unsigned int nVertices() const
Number of vertices.
Definition: FSimEvent.cc:41
std::vector< SimVertex > SimVertexContainer
const SimTrack & embdTrack(int i) const
return embedded track with given id
float weight() const
Method to return the event weight.
Definition: FSimEvent.cc:31
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:25
unsigned int nGenParts() const
Number of generator particles.
Definition: FBaseSimEvent.h:94
const math::XYZTLorentzVectorD & momentum() const
Definition: CoreSimTrack.h:22
bool noMother() const
no mother particle
long double T
std::vector< SimTrack > SimTrackContainer
const FSimTrack & mother() const
mother
FSimTrack & track(int id) const
Return track with given Id.
edm::EventID id_
Definition: FSimEvent.h:66