CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
FSimEvent Class Reference

#include <FSimEvent.h>

Inheritance diagram for FSimEvent:
FBaseSimEvent

Public Member Functions

void fill (const HepMC::GenEvent &hev, edm::EventID &Id)
 fill the FBaseSimEvent from the current HepMC::GenEvent More...
 
void fill (const reco::GenParticleCollection &parts, edm::EventID &Id)
 fill the FBaseSimEvent from the current reco::GenParticleCollection More...
 
void fill (const std::vector< SimTrack > &simTracks, const std::vector< SimVertex > &simVertices)
 fill the FBaseSimEvent from the SimTrack's and SimVert'ices More...
 
 FSimEvent (const edm::ParameterSet &kine)
 Default constructor. More...
 
 FSimEvent (const edm::ParameterSet &vtx, const edm::ParameterSet &kine, const RandomEngine *engine)
 
edm::EventID id () const
 Method to return the EventId. More...
 
void load (edm::SimTrackContainer &c, edm::SimTrackContainer &m) const
 Load containers of tracks (and muons) and vertices for the edm::Event. More...
 
void load (edm::SimVertexContainer &c) const
 
void load (FSimVertexTypeCollection &c) const
 
unsigned int nGenParts () const
 Number of MC particles. More...
 
unsigned int nTracks () const
 Number of tracks. More...
 
unsigned int nVertices () const
 Number of vertices. More...
 
float weight () const
 Method to return the event weight. More...
 
virtual ~FSimEvent ()
 usual virtual destructor More...
 
- Public Member Functions inherited from FBaseSimEvent
void addChargedTrack (int id)
 Add an id in the vector of charged tracks id's. More...
 
void addParticles (const HepMC::GenEvent &hev)
 Add the particles and their vertices to the list. More...
 
void addParticles (const reco::GenParticleCollection &myGenParticles)
 
int addSimTrack (const RawParticle *p, int iv, int ig=-1, const HepMC::GenVertex *ev=0)
 Add a new track to the Event and to the various lists. More...
 
int addSimVertex (const XYZTLorentzVector &decayVertex, int im=-1, FSimVertexType::VertexType type=FSimVertexType::ANY)
 Add a new vertex to the Event and to the various lists. More...
 
int chargedTrack (int id) const
 return "reconstructed" charged tracks index. More...
 
void clear ()
 clear the FBaseSimEvent content before the next event More...
 
const HepMC::GenParticle * embdGenpart (int i) const
 return MC track with a given id More...
 
const SimTrackembdTrack (int i) const
 return embedded track with given id More...
 
const SimVertexembdVertex (int i) const
 return embedded vertex with given id More...
 
const FSimVertexTypeembdVertexType (int i) const
 return embedded vertex type with given id More...
 
 FBaseSimEvent (const edm::ParameterSet &kine)
 Default constructor. More...
 
 FBaseSimEvent (const edm::ParameterSet &vtx, const edm::ParameterSet &kine, const RandomEngine *engine)
 
void fill (const HepMC::GenEvent &hev)
 fill the FBaseSimEvent from the current HepMC::GenEvent More...
 
void fill (const reco::GenParticleCollection &hev)
 fill the FBaseSimEvent from the current reco::GenParticleCollection More...
 
void fill (const std::vector< SimTrack > &, const std::vector< SimVertex > &)
 fill the FBaseSimEvent from SimTrack's and SimVert'ices More...
 
const KineParticleFilterfilter () const
 
void initializePdt (const HepPDT::ParticleDataTable *aPdt)
 Initialize the particle data table. More...
 
unsigned int nChargedTracks () const
 Number of "reconstructed" charged tracks. More...
 
unsigned int nGenParts () const
 Number of generator particles. More...
 
unsigned int nTracks () const
 Number of tracks. More...
 
unsigned int nVertices () const
 Number of vertices. More...
 
void print () const
 print the FBaseSimEvent in an intelligible way More...
 
void printMCTruth (const HepMC::GenEvent &hev)
 print the original MCTruth event More...
 
void setBeamSpot (const math::XYZPoint &aBeamSpot)
 Set the beam spot position. More...
 
PrimaryVertexGeneratorthePrimaryVertexGenerator () const
 
const HepPDT::ParticleDataTabletheTable () const
 Get the pointer to the particle data table. More...
 
FSimTracktrack (int id) const
 Return track with given Id. More...
 
FSimVertexvertex (int id) const
 Return vertex with given Id. More...
 
FSimVertexTypevertexType (int id) const
 Return vertex with given Id. More...
 
 ~FBaseSimEvent ()
 usual virtual destructor More...
 

Private Attributes

edm::EventID id_
 
double weight_
 

Additional Inherited Members

- Protected Member Functions inherited from FBaseSimEvent
std::vector
< HepMC::GenParticle * > * 
genparts () const
 The pointer to the vector of GenParticle's. More...
 
std::vector< FSimTrack > * tracks () const
 The pointer to the vector of FSimTrack's. More...
 
std::vector< FSimVertex > * vertices () const
 The pointer to the vector of FSimVertex's. More...
 

Detailed Description

Definition at line 31 of file FSimEvent.h.

Constructor & Destructor Documentation

FSimEvent::FSimEvent ( const edm::ParameterSet kine)

Default constructor.

Definition at line 6 of file FSimEvent.cc.

7  : FBaseSimEvent(kine), id_(edm::EventID(0,0,0)), weight_(0)
8 {}
FBaseSimEvent(const edm::ParameterSet &kine)
Default constructor.
double weight_
Definition: FSimEvent.h:76
edm::EventID id_
Definition: FSimEvent.h:75
FSimEvent::FSimEvent ( const edm::ParameterSet vtx,
const edm::ParameterSet kine,
const RandomEngine engine 
)

Definition at line 10 of file FSimEvent.cc.

13  : FBaseSimEvent(vtx,kine,engine), id_(edm::EventID(0,0,0)), weight_(0)
14 {}
FBaseSimEvent(const edm::ParameterSet &kine)
Default constructor.
double weight_
Definition: FSimEvent.h:76
edm::EventID id_
Definition: FSimEvent.h:75
FSimEvent::~FSimEvent ( )
virtual

usual virtual destructor

Definition at line 16 of file FSimEvent.cc.

17 {}

Member Function Documentation

void FSimEvent::fill ( const HepMC::GenEvent &  hev,
edm::EventID Id 
)

fill the FBaseSimEvent from the current HepMC::GenEvent

Definition at line 26 of file FSimEvent.cc.

References FBaseSimEvent::fill(), and id_.

Referenced by FamosManager::reconstruct().

26  {
27  FBaseSimEvent::fill(hev);
28  id_ = Id;
29 }
void fill(const HepMC::GenEvent &hev)
fill the FBaseSimEvent from the current HepMC::GenEvent
edm::EventID id_
Definition: FSimEvent.h:75
void FSimEvent::fill ( const reco::GenParticleCollection parts,
edm::EventID Id 
)

fill the FBaseSimEvent from the current reco::GenParticleCollection

Definition at line 20 of file FSimEvent.cc.

References FBaseSimEvent::fill(), and id_.

20  {
22  id_ = Id;
23 }
void fill(const HepMC::GenEvent &hev)
fill the FBaseSimEvent from the current HepMC::GenEvent
edm::EventID id_
Definition: FSimEvent.h:75
void FSimEvent::fill ( const std::vector< SimTrack > &  simTracks,
const std::vector< SimVertex > &  simVertices 
)

fill the FBaseSimEvent from the SimTrack's and SimVert'ices

Definition at line 32 of file FSimEvent.cc.

References FBaseSimEvent::fill(), and id_.

33  {
34  FBaseSimEvent::fill(simTracks,simVertices);
35  id_ = edm::EventID();
36 }
void fill(const HepMC::GenEvent &hev)
fill the FBaseSimEvent from the current HepMC::GenEvent
edm::EventID id_
Definition: FSimEvent.h:75
edm::EventID FSimEvent::id ( void  ) const

Method to return the EventId.

Definition at line 39 of file FSimEvent.cc.

References id_.

Referenced by CalorimetryManager::reconstruct().

39  {
40  return id_;
41 }
edm::EventID id_
Definition: FSimEvent.h:75
void FSimEvent::load ( edm::SimTrackContainer c,
edm::SimTrackContainer m 
) const

Load containers of tracks (and muons) and vertices for the edm::Event.

Definition at line 63 of file FSimEvent.cc.

References abs, FSimTrack::closestDaughterId(), FBaseSimEvent::embdTrack(), i, CoreSimTrack::momentum(), FSimTrack::mother(), FSimTrack::noEndVertex(), nTracks(), lumiQTWidget::t, FBaseSimEvent::track(), and CoreSimTrack::type().

64 {
65  for (unsigned int i=0; i<nTracks(); ++i) {
66  // SimTrack t = SimTrack(ip,p,iv,ig);
67  const SimTrack& t = embdTrack(i);
68  // Save all tracks
69  c.push_back(t);
70  // Save also some muons for later parameterization
71  if ( abs(t.type()) == 13 &&
72  t.momentum().perp2() > 1.0 &&
73  fabs(t.momentum().eta()) < 3.0 &&
74  track(i).noEndVertex() ) {
75  // Actually save the muon mother (and the attached muon) in case
76  if ( track(i).mother().closestDaughterId() == (int)i ) {
77  const SimTrack& T = embdTrack(track(i).mother().id());
78  m.push_back(T);
79  }
80  m.push_back(t);
81  }
82  }
83 }
int i
Definition: DBlmapReader.cc:9
bool noEndVertex() const
no end vertex
#define abs(x)
Definition: mlp_lapack.h:159
int closestDaughterId() const
Get the index of the closest charged daughter.
Definition: FSimTrack.h:187
unsigned int nTracks() const
Number of tracks.
Definition: FSimEvent.cc:48
const SimTrack & embdTrack(int i) const
return embedded track with given id
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:40
const math::XYZTLorentzVectorD & momentum() const
particle info...
Definition: CoreSimTrack.h:36
long double T
const FSimTrack & mother() const
mother
FSimTrack & track(int id) const
Return track with given Id.
void FSimEvent::load ( edm::SimVertexContainer c) const

Definition at line 86 of file FSimEvent.cc.

References FBaseSimEvent::embdVertex(), i, and nVertices().

87 {
88  for (unsigned int i=0; i<nVertices(); ++i) {
89  // SimTrack t = SimTrack(ip,p,iv,ig);
90  c.push_back(embdVertex(i));
91  }
92 }
int i
Definition: DBlmapReader.cc:9
const SimVertex & embdVertex(int i) const
return embedded vertex with given id
unsigned int nVertices() const
Number of vertices.
Definition: FSimEvent.cc:53
void FSimEvent::load ( FSimVertexTypeCollection c) const

Definition at line 96 of file FSimEvent.cc.

References FBaseSimEvent::embdVertexType(), i, and nVertices().

97 {
98 
99  for (unsigned int i=0; i<nVertices(); ++i) {
100  c.push_back(embdVertexType(i));
101  }
102 }
int i
Definition: DBlmapReader.cc:9
const FSimVertexType & embdVertexType(int i) const
return embedded vertex type with given id
unsigned int nVertices() const
Number of vertices.
Definition: FSimEvent.cc:53
unsigned int FSimEvent::nGenParts ( ) const

Number of MC particles.

Definition at line 58 of file FSimEvent.cc.

References FBaseSimEvent::nGenParts().

Referenced by FamosManager::reconstruct().

58  {
59  return FBaseSimEvent::nGenParts();
60 }
unsigned int nGenParts() const
Number of generator particles.
unsigned int FSimEvent::nTracks ( ) const

Number of tracks.

Definition at line 48 of file FSimEvent.cc.

References FBaseSimEvent::nTracks().

Referenced by load(), FamosManager::reconstruct(), CalorimetryManager::reconstruct(), and TrajectoryManager::reconstruct().

48  {
49  return FBaseSimEvent::nTracks();
50 }
unsigned int nTracks() const
Number of tracks.
Definition: FBaseSimEvent.h:95
unsigned int FSimEvent::nVertices ( ) const

Number of vertices.

Definition at line 53 of file FSimEvent.cc.

References FBaseSimEvent::nVertices().

Referenced by load(), and FamosManager::reconstruct().

53  {
54  return FBaseSimEvent::nVertices();
55 }
unsigned int nVertices() const
Number of vertices.
float FSimEvent::weight ( void  ) const

Method to return the event weight.

Definition at line 43 of file FSimEvent.cc.

References weight_.

Referenced by FamosManager::reconstruct().

43  {
44  return weight_;
45 }
double weight_
Definition: FSimEvent.h:76

Member Data Documentation

edm::EventID FSimEvent::id_
private

Definition at line 75 of file FSimEvent.h.

Referenced by fill(), and id().

double FSimEvent::weight_
private

Definition at line 76 of file FSimEvent.h.

Referenced by weight().