CMS 3D CMS Logo

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

Producer for PFRecTracks and PFSimParticles. More...

#include <PFSimParticleProducer.h>

Inheritance diagram for PFSimParticleProducer:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Types

typedef edm::Handle
< reco::PFRecTrackCollection
TrackHandle
 
- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
typedef WorkerT< EDProducerWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 

Public Member Functions

virtual void beginRun (const edm::Run &r, const edm::EventSetup &c) override
 
void getSimIDs (const TrackHandle &trackh, std::vector< unsigned > &recTrackSimID)
 
 PFSimParticleProducer (const edm::ParameterSet &)
 
virtual void produce (edm::Event &, const edm::EventSetup &) override
 
 ~PFSimParticleProducer ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndex indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Attributes

edm::InputTag inputTagEcalRecHitsEB_
 
edm::InputTag inputTagEcalRecHitsEE_
 
edm::InputTag inputTagFamosSimHits_
 
edm::InputTag inputTagRecTracks_
 
edm::InputTag inputTagSim_
 module label for retrieving input simtrack and simvertex More...
 
bool mctruthMatchingInfo_
 
FSimEventmySimEvent
 
edm::ParameterSet particleFilter_
 
bool processParticles_
 process particles on/off More...
 
bool verbose_
 verbose ? More...
 

Additional Inherited Members

- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDProducer
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Producer for PFRecTracks and PFSimParticles.

Author
Colin Bernet
Date
April 2007

Definition at line 38 of file PFSimParticleProducer.h.

Member Typedef Documentation

Definition at line 49 of file PFSimParticleProducer.h.

Constructor & Destructor Documentation

PFSimParticleProducer::PFSimParticleProducer ( const edm::ParameterSet iConfig)
explicit

Definition at line 67 of file PFSimParticleProducer.cc.

References edm::ParameterSet::getParameter(), and edm::ParameterSet::getUntrackedParameter().

68 {
69 
70 
72  iConfig.getUntrackedParameter<bool>("process_Particles",true);
73 
74 
76  = iConfig.getParameter<InputTag>("sim");
77 
78  //retrieving collections for MC Truth Matching
79 
80  //modif-beg
82  = iConfig.getUntrackedParameter<InputTag>("famosSimHits");
84  iConfig.getUntrackedParameter<bool>("MCTruthMatchingInfo",false);
85  //modif-end
86 
88  = iConfig.getParameter<InputTag>("RecTracks");
90  = iConfig.getParameter<InputTag>("ecalRecHitsEB");
92  = iConfig.getParameter<InputTag>("ecalRecHitsEE");
93 
94  verbose_ =
95  iConfig.getUntrackedParameter<bool>("verbose",false);
96 
97 
98  // register products
99  produces<reco::PFSimParticleCollection>();
100 
102  ( "ParticleFilter" );
103 
104  mySimEvent = new FSimEvent( particleFilter_ );
105 }
bool processParticles_
process particles on/off
edm::InputTag inputTagEcalRecHitsEE_
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::ParameterSet particleFilter_
edm::InputTag inputTagEcalRecHitsEB_
edm::InputTag inputTagSim_
module label for retrieving input simtrack and simvertex
edm::InputTag inputTagFamosSimHits_
PFSimParticleProducer::~PFSimParticleProducer ( )

Definition at line 109 of file PFSimParticleProducer.cc.

110 {
111  delete mySimEvent;
112 }

Member Function Documentation

void PFSimParticleProducer::beginRun ( const edm::Run r,
const edm::EventSetup c 
)
overridevirtual

Reimplemented from edm::EDProducer.

Definition at line 116 of file PFSimParticleProducer.cc.

References edm::EventSetup::getData(), and ParticleTable::instance().

118 {
119 
120  // init Particle data table (from Pythia)
122  // edm::ESHandle < DefaultConfig::ParticleDataTable > pdt;
123  es.getData(pdt);
125  mySimEvent->initializePdt(&(*pdt));
126 
127 }
static ParticleTable * instance()
Definition: ParticleTable.h:15
void initializePdt(const HepPDT::ParticleDataTable *aPdt)
Initialize the particle data table.
void PFSimParticleProducer::getSimIDs ( const TrackHandle trackh,
std::vector< unsigned > &  recTrackSimID 
)

Definition at line 536 of file PFSimParticleProducer.cc.

References i, edm::HandleBase::isValid(), SiTrackerGSMatchedRecHit2D::simtrackId(), and reco::PFRecTrack::trackRef().

538 {
539 
540  if( trackh.isValid() ) {
541 // cout << "Size=" << trackh->size() << endl;
542 
543  for(unsigned i=0;i<trackh->size(); i++) {
544 
545  reco::PFRecTrackRef ref( trackh,i );
546  const reco::PFRecTrack& PFT = *ref;
547  const reco::TrackRef trackref = PFT.trackRef();
548 
549 // double Pt = trackref->pt();
550 // double DPt = trackref->ptError();
551 // cout << " PFBlockProducer: PFrecTrack->Track Pt= "
552 // << Pt << " DPt = " << DPt << endl;
553 
554  trackingRecHit_iterator rhitbeg
555  = trackref->recHitsBegin();
556  trackingRecHit_iterator rhitend
557  = trackref->recHitsEnd();
558  for (trackingRecHit_iterator it = rhitbeg;
559  it != rhitend; it++){
560 
561  if( it->get()->isValid() ){
562 
563  const SiTrackerGSMatchedRecHit2D * rechit
564  = (const SiTrackerGSMatchedRecHit2D*) (it->get());
565 
566 // cout << "rechit"
567 // << " corresponding simId "
568 // << rechit->simtrackId()
569 // << endl;
570 
571  recTrackSimID.push_back( rechit->simtrackId() );
572  break;
573  }
574  }//loop track rechit
575  }//loop recTracks
576  }//track handle valid
577 }
reconstructed track used as an input to particle flow
Definition: PFRecTrack.h:22
int i
Definition: DBlmapReader.cc:9
const reco::TrackRef & trackRef() const
Definition: PFRecTrack.h:54
void PFSimParticleProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overridevirtual

Implements edm::EDProducer.

Definition at line 130 of file PFSimParticleProducer.cc.

References reco::PFTrajectoryPoint::BeamPipeOrEndVertex, FSimTrack::charge(), reco::PFTrajectoryPoint::ClosestApproach, FSimTrack::daughters(), cond::rpcobgas::detid, reco::PFTrajectoryPoint::ECALEntrance, FSimTrack::ecalEntrance(), FSimTrack::endVertex(), edm::EventID::event(), edm::hlt::Exception, newFWLiteAna::found, edm::Event::getByLabel(), ecalpyutils::hashedIndex(), reco::PFTrajectoryPoint::HCALEntrance, FSimTrack::hcalEntrance(), reco::PFTrajectoryPoint::HCALExit, FSimTrack::hcalExit(), FSimTrack::hoEntrance(), reco::PFTrajectoryPoint::HOLayer, i, edm::EventBase::id(), FSimTrack::id(), edm::HandleBase::isValid(), FSimTrack::layer1Entrance(), FSimTrack::layer2Entrance(), LogDebug, FSimTrack::momentum(), FSimTrack::mother(), FSimTrack::noEndVertex(), FSimTrack::noMother(), FSimTrack::onEcal(), FSimTrack::onHcal(), FSimTrack::onLayer1(), FSimTrack::onLayer2(), FSimVertex::position(), edm::Handle< T >::product(), reco::PFTrajectoryPoint::PS1, reco::PFTrajectoryPoint::PS2, edm::Event::put(), edm::EventID::run(), CoreSimTrack::type(), FSimTrack::vertex(), RawParticle::x(), RawParticle::y(), and RawParticle::z().

132 {
133 
134  LogDebug("PFSimParticleProducer")<<"START event: "<<iEvent.id().event()
135  <<" in run "<<iEvent.id().run()<<endl;
136 
137  //MC Truth Matching only with Famos and UnFoldedMode option to true!!
138 
139  //vector to store the trackIDs of rectracks corresponding
140  //to the simulated particle.
141  std::vector<unsigned> recTrackSimID;
142 
143  //In order to know which simparticule contribute to
144  //a given Ecal RecHit energy, we need to access
145  //the PCAloHit from FastSim.
146 
147  typedef std::pair<double, unsigned> hitSimID;
148  typedef std::list< std::pair<double, unsigned> >::iterator ITM;
149  std::vector< std::list <hitSimID> > caloHitsEBID(62000);
150  std::vector<double> caloHitsEBTotE(62000,0.0);
151 
153 
154  //getting the PCAloHit
156  // bool found_phit
157  // = iEvent.getByLabel("famosSimHits","EcalHitsEB",
158  // pcalohits);
159  //modif-beg
160  bool found_phit
161  = iEvent.getByLabel(inputTagFamosSimHits_,pcalohits);
162  //modif-end
163 
164  if(!found_phit) {
165  ostringstream err;
166  err<<"could not find pcaloHit "<<"famosSimHits:EcalHitsEB";
167  LogError("PFSimParticleProducer")<<err.str()<<endl;
168 
169  throw cms::Exception( "MissingProduct", err.str());
170  }
171  else {
172  assert( pcalohits.isValid() );
173 
174  // cout << "PFSimParticleProducer: number of pcalohits="
175  // << pcalohits->size() << endl;
176 
177  edm::PCaloHitContainer::const_iterator it
178  = pcalohits.product()->begin();
179  edm::PCaloHitContainer::const_iterator itend
180  = pcalohits.product()->end();
181 
182  //loop on the PCaloHit from FastSim Calorimetry
183  for(;it!=itend;++it)
184  {
185  EBDetId detid(it->id());
186 
187  // cout << detid << " " << detid.rawId()
188  // << " " << detid.hashedIndex()
189  // << " " << it->energy()
190  // << " " << it->id()
191  // << " trackId="
192  // << it->geantTrackId()
193  // << endl;
194 
195  if(it->energy() > 0.0) {
196  std::pair<double, unsigned> phitsimid
197  = make_pair(it->energy(),it->geantTrackId());
198  caloHitsEBID[detid.hashedIndex()].push_back(phitsimid);
199  caloHitsEBTotE[detid.hashedIndex()]
200  += it->energy(); //summing pcalhit energy
201  }//energy > 0
202 
203  }//loop PcaloHits
204  }//pcalohit handle access
205 
206  //Retrieving the PFRecTrack collection for
207  //Monte Carlo Truth Matching tool
209  try{
210  LogDebug("PFSimParticleProducer")<<"getting PFRecTracks"<<endl;
211  iEvent.getByLabel(inputTagRecTracks_, recTracks);
212 
213  } catch (cms::Exception& err) {
214  LogError("PFSimParticleProducer")<<err
215  <<" cannot get collection "
216  <<"particleFlowBlock"<<":"
217  <<""
218  <<endl;
219  }//pfrectrack handle access
220 
221  //getting the simID corresponding to
222  //each of the PFRecTracks
223  getSimIDs( recTracks, recTrackSimID );
224 
225  }//mctruthMatchingInfo_ //modif
226 
227  // deal with true particles
228  if( processParticles_) {
229 
230  auto_ptr< reco::PFSimParticleCollection >
231  pOutputPFSimParticleCollection(new reco::PFSimParticleCollection );
232 
233  Handle<vector<SimTrack> > simTracks;
234  bool found = iEvent.getByLabel(inputTagSim_,simTracks);
235  if(!found) {
236 
237  ostringstream err;
238  err<<"cannot find sim tracks: "<<inputTagSim_;
239  LogError("PFSimParticleProducer")<<err.str()<<endl;
240 
241  throw cms::Exception( "MissingProduct", err.str());
242  }
243 
244 
245 
246  Handle<vector<SimVertex> > simVertices;
247  found = iEvent.getByLabel(inputTagSim_,simVertices);
248  if(!found) {
249  LogError("PFSimParticleProducer")
250  <<"cannot find sim vertices: "<<inputTagSim_<<endl;
251  return;
252  }
253 
254 
255 
256  // for(unsigned it = 0; it<simTracks->size(); it++ ) {
257  // cout<<"\t track "<< (*simTracks)[it]<<" "
258  // <<(*simTracks)[it].momentum().vect().perp()<<" "
259  // <<(*simTracks)[it].momentum().e()<<endl;
260  // }
261 
262  mySimEvent->fill( *simTracks, *simVertices );
263 
264  if(verbose_)
265  mySimEvent->print();
266 
267  // const std::vector<FSimTrack>& fsimTracks = *(mySimEvent->tracks() );
268  for(unsigned i=0; i<mySimEvent->nTracks(); i++) {
269 
270  const FSimTrack& fst = mySimEvent->track(i);
271 
272  int motherId = -1;
273  if( ! fst.noMother() ) // a mother exist
274  motherId = fst.mother().id();
275 
276  //This is finding out the simID corresponding
277  //to the recTrack
278 // cout << "Particle " << i
279 // << " " << fst.genpartIndex()
280 // << " -------------------------------------" << endl;
281 
282  //GETTING THE TRACK ID
283  unsigned recTrackID = 99999;
284  vector<unsigned> recHitContrib; //modif
285  vector<double> recHitContribFrac; //modif
286 
287  if(mctruthMatchingInfo_){ //modif
288 
289  for(unsigned lo=0; lo<recTrackSimID.size();
290  lo++) {
291  if( i == recTrackSimID[lo] ) {
292 // cout << "Corresponding Rec Track "
293 // << lo << endl;
294  recTrackID = lo;
295  }//match track
296  }//loop rectrack
297 // if( recTrackID == 99999 )
298 // cout << "Sim Track not reconstructed pT=" <<
299 // fst.momentum().pt() << endl;
300 
301  // get the ecalBarrel rechits for MC truth matching tool
303  bool found = iEvent.getByLabel(inputTagEcalRecHitsEB_,
304  rhcHandle);
305  if(!found) {
306  ostringstream err;
307  err<<"could not find rechits "<< inputTagEcalRecHitsEB_;
308  LogError("PFSimParticleProducer")<<err.str()<<endl;
309 
310  throw cms::Exception( "MissingProduct", err.str());
311  }
312  else {
313  assert( rhcHandle.isValid() );
314 // cout << "PFSimParticleProducer: number of rechits="
315 // << rhcHandle->size() << endl;
316 
318  = rhcHandle.product()->begin();
320  = rhcHandle.product()->end();
321 
322  for(;it_rh!=itend_rh;++it_rh)
323  {
324  unsigned rhit_hi
325  = EBDetId(it_rh->id()).hashedIndex();
326  EBDetId detid(it_rh->id());
327 // cout << detid << " " << detid.rawId()
328 // << " " << detid.hashedIndex()
329 // << " " << it_rh->energy() << endl;
330 
331  ITM it_phit = caloHitsEBID[rhit_hi].begin();
332  ITM itend_phit = caloHitsEBID[rhit_hi].end();
333  for(;it_phit!=itend_phit;++it_phit)
334  {
335  if(i == it_phit->second)
336  {
337  //Alex (08/10/08) TO BE REMOVED, eliminating
338  //duplicated rechits
339  bool alreadyin = false;
340  for( unsigned ihit = 0; ihit < recHitContrib.size();
341  ++ihit )
342  if(detid.rawId() == recHitContrib[ihit])
343  alreadyin = true;
344 
345  if(!alreadyin){
346  double pcalofraction = 0.0;
347  if(caloHitsEBTotE[rhit_hi] != 0.0)
348  pcalofraction
349  = (it_phit->first/caloHitsEBTotE[rhit_hi])*100.0;
350 
351  //store info
352  recHitContrib.push_back(it_rh->id());
353  recHitContribFrac.push_back(pcalofraction);
354  }//selected rechits
355  }//matching
356  }//loop pcalohit
357 
358  }//loop rechits
359 
360  }//getting the rechits
361 
362  }//mctruthMatchingInfo_ //modif
363 
364 // cout << "This particule has " << recHitContrib.size()
365 // << " rechit contribution" << endl;
366 // for( unsigned ih = 0; ih < recHitContrib.size(); ++ih )
367 // cout << recHitContrib[ih]
368 // << " f=" << recHitContribFrac[ih] << " ";
369 // cout << endl;
370 
371  reco::PFSimParticle particle( fst.charge(),
372  fst.type(),
373  fst.id(),
374  motherId,
375  fst.daughters(),
376  recTrackID,
377  recHitContrib,
378  recHitContribFrac);
379 
380 
381  const FSimVertex& originVtx = fst.vertex();
382 
383  math::XYZPoint posOrig( originVtx.position().x(),
384  originVtx.position().y(),
385  originVtx.position().z() );
386 
387  math::XYZTLorentzVector momOrig( fst.momentum().px(),
388  fst.momentum().py(),
389  fst.momentum().pz(),
390  fst.momentum().e() );
392  pointOrig(-1,
394  posOrig, momOrig);
395 
396  // point 0 is origin vertex
397  particle.addPoint(pointOrig);
398 
399 
400  if( ! fst.noEndVertex() ) {
401  const FSimVertex& endVtx = fst.endVertex();
402 
403  math::XYZPoint posEnd( endVtx.position().x(),
404  endVtx.position().y(),
405  endVtx.position().z() );
406  // cout<<"end vertex : "
407  // <<endVtx.position().x()<<" "
408  // <<endVtx.position().y()<<endl;
409 
411 
413  pointEnd( -1,
415  posEnd, momEnd);
416 
417  particle.addPoint(pointEnd);
418  }
419  else { // add a dummy point
421  particle.addPoint(dummy);
422  }
423 
424 
425  if( fst.onLayer1() ) { // PS layer1
426  const RawParticle& rp = fst.layer1Entrance();
427 
428  math::XYZPoint posLayer1( rp.x(), rp.y(), rp.z() );
429  math::XYZTLorentzVector momLayer1( rp.px(), rp.py(), rp.pz(), rp.e() );
431  posLayer1, momLayer1);
432 
433  particle.addPoint( layer1Pt );
434 
435  // extrapolate to cluster depth
436  }
437  else { // add a dummy point
439  particle.addPoint(dummy);
440  }
441 
442  if( fst.onLayer2() ) { // PS layer2
443  const RawParticle& rp = fst.layer2Entrance();
444 
445  math::XYZPoint posLayer2( rp.x(), rp.y(), rp.z() );
446  math::XYZTLorentzVector momLayer2( rp.px(), rp.py(), rp.pz(), rp.e() );
448  posLayer2, momLayer2);
449 
450  particle.addPoint( layer2Pt );
451 
452  // extrapolate to cluster depth
453  }
454  else { // add a dummy point
456  particle.addPoint(dummy);
457  }
458 
459  if( fst.onEcal() ) {
460  const RawParticle& rp = fst.ecalEntrance();
461 
462  math::XYZPoint posECAL( rp.x(), rp.y(), rp.z() );
463  math::XYZTLorentzVector momECAL( rp.px(), rp.py(), rp.pz(), rp.e() );
464  reco::PFTrajectoryPoint ecalPt(-1,
466  posECAL, momECAL);
467 
468  particle.addPoint( ecalPt );
469 
470  // extrapolate to cluster depth
471  }
472  else { // add a dummy point
474  particle.addPoint(dummy);
475  }
476 
477  // add a dummy point for ECAL Shower max
479  particle.addPoint(dummy);
480 
481  if( fst.onHcal() ) {
482 
483  const RawParticle& rpin = fst.hcalEntrance();
484 
485  math::XYZPoint posHCALin( rpin.x(), rpin.y(), rpin.z() );
486  math::XYZTLorentzVector momHCALin( rpin.px(), rpin.py(), rpin.pz(),
487  rpin.e() );
488  reco::PFTrajectoryPoint hcalPtin(-1,
490  posHCALin, momHCALin);
491 
492  particle.addPoint( hcalPtin );
493 
494  const RawParticle& rpout = fst.hcalExit();
495 
496  math::XYZPoint posHCALout( rpout.x(), rpout.y(), rpout.z() );
497  math::XYZTLorentzVector momHCALout( rpout.px(), rpout.py(), rpout.pz(),
498  rpout.e() );
500  hcalPtout(-1, reco::PFTrajectoryPoint::HCALExit,
501  posHCALout, momHCALout);
502 
503  particle.addPoint( hcalPtout );
504 
505  const RawParticle& rpho = fst.hoEntrance();
506 
507  math::XYZPoint posHOEntrance( rpho.x(), rpho.y(), rpho.z() );
508  math::XYZTLorentzVector momHOEntrance( rpho.px(), rpho.py(), rpho.pz(),
509  rpho.e() );
512  posHOEntrance, momHOEntrance);
513 
514  particle.addPoint( hoPtin );
515 
516 
517 
518 
519  }
520  else { // add a dummy point
522  particle.addPoint(dummy);
523  }
524 
525  pOutputPFSimParticleCollection->push_back( particle );
526  }
527 
528  iEvent.put(pOutputPFSimParticleCollection);
529  }
530 
531 
532  LogDebug("PFSimParticleProducer")<<"STOP event: "<<iEvent.id().event()
533  <<" in run "<<iEvent.id().run()<<endl;
534 }
#define LogDebug(id)
RunNumber_t run() const
Definition: EventID.h:42
bool processParticles_
process particles on/off
EventNumber_t event() const
Definition: EventID.h:44
int i
Definition: DBlmapReader.cc:9
bool noEndVertex() const
no end vertex
float charge() const
charge
Definition: FSimTrack.h:47
const RawParticle & hoEntrance() const
The particle at HCAL exir.
Definition: FSimTrack.h:145
int onLayer2() const
Definition: FSimTrack.h:96
const FSimVertex & endVertex() const
end vertex
const RawParticle & hcalExit() const
The particle at HCAL exir.
Definition: FSimTrack.h:142
double y() const
y of vertex
Definition: RawParticle.h:270
double z() const
z of vertex
Definition: RawParticle.h:271
const XYZTLorentzVector & momentum() const
Temporary (until move of SimTrack to Mathcore) - No! Actually very useful.
Definition: FSimTrack.h:190
void fill(const HepMC::GenEvent &hev, edm::EventID &Id)
fill the FBaseSimEvent from the current HepMC::GenEvent
Definition: FSimEvent.cc:26
const RawParticle & layer1Entrance() const
The particle at Preshower Layer 1.
Definition: FSimTrack.h:127
std::vector< EcalRecHit >::const_iterator const_iterator
edm::InputTag inputTagEcalRecHitsEB_
std::vector< PFSimParticle > PFSimParticleCollection
collection of PFSimParticle objects
const std::vector< int > & daughters() const
Vector of daughter indices.
double x() const
x of vertex
Definition: RawParticle.h:269
int onEcal() const
Definition: FSimTrack.h:101
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
int hashedIndex(int ieta, int iphi)
Definition: EcalPyUtils.cc:42
edm::InputTag inputTagSim_
module label for retrieving input simtrack and simvertex
edm::InputTag inputTagFamosSimHits_
Point of closest approach from beam axis (initial point in the case of PFSimParticle) ...
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
const math::XYZTLorentzVector & position() const
Temporary (until CMSSW moves to Mathcore) - No ! Actually very useful.
Definition: FSimVertex.h:49
const RawParticle & ecalEntrance() const
The particle at ECAL entrance.
Definition: FSimTrack.h:133
std::list< std::pair< double, unsigned > >::iterator ITM
const FSimVertex & vertex() const
Origin vertex.
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
true particle for particle flow
Definition: PFSimParticle.h:19
unsigned int nTracks() const
Number of tracks.
Definition: FSimEvent.cc:48
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
T const * product() const
Definition: Handle.h:74
edm::EventID id() const
Definition: EventBase.h:56
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:40
int id() const
the index in FBaseSimEvent and other vectors
Definition: FSimTrack.h:86
const RawParticle & layer2Entrance() const
The particle at Preshower Layer 2.
Definition: FSimTrack.h:130
int onLayer1() const
Definition: FSimTrack.h:91
int onHcal() const
Definition: FSimTrack.h:106
bool noMother() const
no mother particle
const RawParticle & hcalEntrance() const
The particle at HCAL entrance.
Definition: FSimTrack.h:136
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
void getSimIDs(const TrackHandle &trackh, std::vector< unsigned > &recTrackSimID)
const FSimTrack & mother() const
mother
void print() const
print the FBaseSimEvent in an intelligible way
FSimTrack & track(int id) const
Return track with given Id.

Member Data Documentation

edm::InputTag PFSimParticleProducer::inputTagEcalRecHitsEB_
private

Definition at line 66 of file PFSimParticleProducer.h.

edm::InputTag PFSimParticleProducer::inputTagEcalRecHitsEE_
private

Definition at line 67 of file PFSimParticleProducer.h.

edm::InputTag PFSimParticleProducer::inputTagFamosSimHits_
private

Definition at line 62 of file PFSimParticleProducer.h.

edm::InputTag PFSimParticleProducer::inputTagRecTracks_
private

Definition at line 65 of file PFSimParticleProducer.h.

edm::InputTag PFSimParticleProducer::inputTagSim_
private

module label for retrieving input simtrack and simvertex

Definition at line 57 of file PFSimParticleProducer.h.

bool PFSimParticleProducer::mctruthMatchingInfo_
private

Definition at line 61 of file PFSimParticleProducer.h.

FSimEvent* PFSimParticleProducer::mySimEvent
private

Definition at line 72 of file PFSimParticleProducer.h.

edm::ParameterSet PFSimParticleProducer::particleFilter_
private

Definition at line 71 of file PFSimParticleProducer.h.

bool PFSimParticleProducer::processParticles_
private

process particles on/off

Definition at line 77 of file PFSimParticleProducer.h.

bool PFSimParticleProducer::verbose_
private

verbose ?

Definition at line 80 of file PFSimParticleProducer.h.