CMS 3D CMS Logo

PFSimParticleProducer.cc
Go to the documentation of this file.
2 
5 
8 
11 
12 // include files used for reconstructed tracks
13 // #include "DataFormats/TrackReco/interface/Track.h"
14 // #include "DataFormats/TrackReco/interface/TrackFwd.h"
15 // #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
16 // #include "DataFormats/TrackingRecHit/interface/TrackingRecHitFwd.h"
17 // #include "DataFormats/TrackCandidate/interface/TrackCandidate.h"
18 // #include "DataFormats/TrackCandidate/interface/TrackCandidateCollection.h"
19 // #include "DataFormats/TrajectorySeed/interface/PropagationDirection.h"
20 
37 // #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
39 
43 
48 
58 
59 #include <set>
60 #include <sstream>
61 
62 using namespace std;
63 using namespace edm;
64 
66  processParticles_ = iConfig.getUntrackedParameter<bool>("process_Particles", true);
67 
68  inputTagSim_ = iConfig.getParameter<InputTag>("sim");
69  tokenSim_ = consumes<std::vector<SimTrack> >(inputTagSim_);
70  tokenSimVertices_ = consumes<std::vector<SimVertex> >(inputTagSim_);
71 
72  //retrieving collections for MC Truth Matching
73 
74  //modif-beg
75  inputTagFastSimProducer_ = iConfig.getUntrackedParameter<InputTag>("fastSimProducer");
76  tokenFastSimProducer_ = consumes<edm::PCaloHitContainer>(inputTagFastSimProducer_);
77  mctruthMatchingInfo_ = iConfig.getUntrackedParameter<bool>("MCTruthMatchingInfo", false);
78  //modif-end
79 
80  inputTagRecTracks_ = iConfig.getParameter<InputTag>("RecTracks");
81  tokenRecTracks_ = consumes<reco::PFRecTrackCollection>(inputTagRecTracks_);
82 
83  inputTagEcalRecHitsEB_ = iConfig.getParameter<InputTag>("ecalRecHitsEB");
84  tokenEcalRecHitsEB_ = consumes<EcalRecHitCollection>(inputTagEcalRecHitsEB_);
85  inputTagEcalRecHitsEE_ = iConfig.getParameter<InputTag>("ecalRecHitsEE");
86  tokenEcalRecHitsEE_ = consumes<EcalRecHitCollection>(inputTagEcalRecHitsEE_);
87 
88  verbose_ = iConfig.getUntrackedParameter<bool>("verbose", false);
89 
90  // register products
91  produces<reco::PFSimParticleCollection>();
92 
93  particleFilter_ = iConfig.getParameter<ParameterSet>("ParticleFilter");
94 
95  mySimEvent = new FSimEvent(particleFilter_);
96 }
97 
99 
101  // init Particle data table (from Pythia)
103  // edm::ESHandle < DefaultConfig::ParticleDataTable > pdt;
104  iSetup.getData(pdt);
105  mySimEvent->initializePdt(&(*pdt));
106 
107  LogDebug("PFSimParticleProducer") << "START event: " << iEvent.id().event() << " in run " << iEvent.id().run()
108  << endl;
109 
110  //MC Truth Matching only with Famos and UnFoldedMode option to true!!
111 
112  //vector to store the trackIDs of rectracks corresponding
113  //to the simulated particle.
114  std::vector<unsigned> recTrackSimID;
115 
116  //In order to know which simparticule contribute to
117  //a given Ecal RecHit energy, we need to access
118  //the PCAloHit from FastSim.
119 
120  typedef std::pair<double, unsigned> hitSimID;
121  typedef std::list<std::pair<double, unsigned> >::iterator ITM;
122  std::vector<std::list<hitSimID> > caloHitsEBID(62000);
123  std::vector<double> caloHitsEBTotE(62000, 0.0);
124 
125  if (mctruthMatchingInfo_) {
126  //getting the PCAloHit
128  // bool found_phit
129  // = iEvent.getByLabel("fastSimProducer","EcalHitsEB",
130  // pcalohits);
131  //modif-beg
132  bool found_phit = iEvent.getByToken(tokenFastSimProducer_, pcalohits);
133  //modif-end
134 
135  if (!found_phit) {
136  ostringstream err;
137  err << "could not find pcaloHit "
138  << "fastSimProducer:EcalHitsEB";
139  LogError("PFSimParticleProducer") << err.str() << endl;
140 
141  throw cms::Exception("MissingProduct", err.str());
142  } else {
143  assert(pcalohits.isValid());
144 
145  edm::PCaloHitContainer::const_iterator it = pcalohits.product()->begin();
146  edm::PCaloHitContainer::const_iterator itend = pcalohits.product()->end();
147 
148  //loop on the PCaloHit from FastSim Calorimetry
149  for (; it != itend; ++it) {
150  EBDetId detid(it->id());
151 
152  if (it->energy() > 0.0) {
153  std::pair<double, unsigned> phitsimid = make_pair(it->energy(), it->geantTrackId());
154  caloHitsEBID[detid.hashedIndex()].push_back(phitsimid);
155  caloHitsEBTotE[detid.hashedIndex()] += it->energy(); //summing pcalhit energy
156  } //energy > 0
157 
158  } //loop PcaloHits
159  } //pcalohit handle access
160 
161  //Retrieving the PFRecTrack collection for
162  //Monte Carlo Truth Matching tool
164  try {
165  LogDebug("PFSimParticleProducer") << "getting PFRecTracks" << endl;
166  iEvent.getByToken(tokenRecTracks_, recTracks);
167 
168  } catch (cms::Exception& err) {
169  LogError("PFSimParticleProducer") << err << " cannot get collection "
170  << "particleFlowBlock"
171  << ":"
172  << "" << endl;
173  } //pfrectrack handle access
174 
175  //getting the simID corresponding to
176  //each of the PFRecTracks
177  getSimIDs(recTracks, recTrackSimID);
178 
179  } //mctruthMatchingInfo_ //modif
180 
181  // deal with true particles
182  if (processParticles_) {
183  auto pOutputPFSimParticleCollection = std::make_unique<reco::PFSimParticleCollection>();
184 
186  bool found = iEvent.getByToken(tokenSim_, simTracks);
187  if (!found) {
188  ostringstream err;
189  err << "cannot find sim tracks: " << inputTagSim_;
190  LogError("PFSimParticleProducer") << err.str() << endl;
191 
192  throw cms::Exception("MissingProduct", err.str());
193  }
194 
196  found = iEvent.getByToken(tokenSimVertices_, simVertices);
197  if (!found) {
198  LogError("PFSimParticleProducer") << "cannot find sim vertices: " << inputTagSim_ << endl;
199  return;
200  }
201 
202  mySimEvent->fill(*simTracks, *simVertices);
203 
204  if (verbose_)
205  mySimEvent->print();
206 
207  for (unsigned i = 0; i < mySimEvent->nTracks(); i++) {
208  const FSimTrack& fst = mySimEvent->track(i);
209 
210  int motherId = -1;
211  if (!fst.noMother()) // a mother exist
212  motherId = fst.mother().id();
213 
214  //This is finding out the simID corresponding
215  //to the recTrack
216 
217  //GETTING THE TRACK ID
218  unsigned recTrackID = 99999;
219  vector<unsigned> recHitContrib; //modif
220  vector<double> recHitContribFrac; //modif
221 
222  if (mctruthMatchingInfo_) { //modif
223 
224  for (unsigned lo = 0; lo < recTrackSimID.size(); lo++) {
225  if (i == recTrackSimID[lo]) {
226  recTrackID = lo;
227  } //match track
228  } //loop rectrack
229 
230  // get the ecalBarrel rechits for MC truth matching tool
232  bool found = iEvent.getByToken(tokenEcalRecHitsEB_, rhcHandle);
233  if (!found) {
234  ostringstream err;
235  err << "could not find rechits " << inputTagEcalRecHitsEB_;
236  LogError("PFSimParticleProducer") << err.str() << endl;
237 
238  throw cms::Exception("MissingProduct", err.str());
239  } else {
240  assert(rhcHandle.isValid());
241 
242  EBRecHitCollection::const_iterator it_rh = rhcHandle.product()->begin();
243  EBRecHitCollection::const_iterator itend_rh = rhcHandle.product()->end();
244 
245  for (; it_rh != itend_rh; ++it_rh) {
246  unsigned rhit_hi = EBDetId(it_rh->id()).hashedIndex();
247  EBDetId detid(it_rh->id());
248 
249  ITM it_phit = caloHitsEBID[rhit_hi].begin();
250  ITM itend_phit = caloHitsEBID[rhit_hi].end();
251  for (; it_phit != itend_phit; ++it_phit) {
252  if (i == it_phit->second) {
253  //Alex (08/10/08) TO BE REMOVED, eliminating
254  //duplicated rechits
255  bool alreadyin = false;
256  for (unsigned ihit = 0; ihit < recHitContrib.size(); ++ihit)
257  if (detid.rawId() == recHitContrib[ihit])
258  alreadyin = true;
259 
260  if (!alreadyin) {
261  double pcalofraction = 0.0;
262  if (caloHitsEBTotE[rhit_hi] != 0.0)
263  pcalofraction = (it_phit->first / caloHitsEBTotE[rhit_hi]) * 100.0;
264 
265  //store info
266  recHitContrib.push_back(it_rh->id());
267  recHitContribFrac.push_back(pcalofraction);
268  } //selected rechits
269  } //matching
270  } //loop pcalohit
271 
272  } //loop rechits
273 
274  } //getting the rechits
275 
276  } //mctruthMatchingInfo_ //modif
277 
278  reco::PFSimParticle particle(
279  fst.charge(), fst.type(), fst.id(), motherId, fst.daughters(), recTrackID, recHitContrib, recHitContribFrac);
280 
281  const FSimVertex& originVtx = fst.vertex();
282 
283  math::XYZPoint posOrig(originVtx.position().x(), originVtx.position().y(), originVtx.position().z());
284 
285  math::XYZTLorentzVector momOrig(
286  fst.momentum().px(), fst.momentum().py(), fst.momentum().pz(), fst.momentum().e());
287  reco::PFTrajectoryPoint pointOrig(-1, reco::PFTrajectoryPoint::ClosestApproach, posOrig, momOrig);
288 
289  // point 0 is origin vertex
290  particle.addPoint(pointOrig);
291 
292  if (!fst.noEndVertex()) {
293  const FSimVertex& endVtx = fst.endVertex();
294 
295  math::XYZPoint posEnd(endVtx.position().x(), endVtx.position().y(), endVtx.position().z());
296 
298 
300 
301  particle.addPoint(pointEnd);
302  } else { // add a dummy point
304  particle.addPoint(dummy);
305  }
306 
307  if (fst.onLayer1()) { // PS layer1
308  const RawParticle& rp = fst.layer1Entrance();
309 
310  math::XYZPoint posLayer1(rp.x(), rp.y(), rp.z());
311  math::XYZTLorentzVector momLayer1(rp.px(), rp.py(), rp.pz(), rp.e());
312  reco::PFTrajectoryPoint layer1Pt(-1, reco::PFTrajectoryPoint::PS1, posLayer1, momLayer1);
313 
314  particle.addPoint(layer1Pt);
315 
316  // extrapolate to cluster depth
317  } else { // add a dummy point
319  particle.addPoint(dummy);
320  }
321 
322  if (fst.onLayer2()) { // PS layer2
323  const RawParticle& rp = fst.layer2Entrance();
324 
325  math::XYZPoint posLayer2(rp.x(), rp.y(), rp.z());
326  math::XYZTLorentzVector momLayer2(rp.px(), rp.py(), rp.pz(), rp.e());
327  reco::PFTrajectoryPoint layer2Pt(-1, reco::PFTrajectoryPoint::PS2, posLayer2, momLayer2);
328 
329  particle.addPoint(layer2Pt);
330 
331  // extrapolate to cluster depth
332  } else { // add a dummy point
334  particle.addPoint(dummy);
335  }
336 
337  if (fst.onEcal()) {
338  const RawParticle& rp = fst.ecalEntrance();
339 
340  math::XYZPoint posECAL(rp.x(), rp.y(), rp.z());
341  math::XYZTLorentzVector momECAL(rp.px(), rp.py(), rp.pz(), rp.e());
343 
344  particle.addPoint(ecalPt);
345 
346  // extrapolate to cluster depth
347  } else { // add a dummy point
349  particle.addPoint(dummy);
350  }
351 
352  // add a dummy point for ECAL Shower max
354  particle.addPoint(dummy);
355 
356  if (fst.onHcal()) {
357  const RawParticle& rpin = fst.hcalEntrance();
358 
359  math::XYZPoint posHCALin(rpin.x(), rpin.y(), rpin.z());
360  math::XYZTLorentzVector momHCALin(rpin.px(), rpin.py(), rpin.pz(), rpin.e());
361  reco::PFTrajectoryPoint hcalPtin(-1, reco::PFTrajectoryPoint::HCALEntrance, posHCALin, momHCALin);
362 
363  particle.addPoint(hcalPtin);
364 
365  const RawParticle& rpout = fst.hcalExit();
366 
367  math::XYZPoint posHCALout(rpout.x(), rpout.y(), rpout.z());
368  math::XYZTLorentzVector momHCALout(rpout.px(), rpout.py(), rpout.pz(), rpout.e());
369  reco::PFTrajectoryPoint hcalPtout(-1, reco::PFTrajectoryPoint::HCALExit, posHCALout, momHCALout);
370 
371  particle.addPoint(hcalPtout);
372 
373  const RawParticle& rpho = fst.hoEntrance();
374 
375  math::XYZPoint posHOEntrance(rpho.x(), rpho.y(), rpho.z());
376  math::XYZTLorentzVector momHOEntrance(rpho.px(), rpho.py(), rpho.pz(), rpho.e());
377  reco::PFTrajectoryPoint hoPtin(-1, reco::PFTrajectoryPoint::HOLayer, posHOEntrance, momHOEntrance);
378 
379  particle.addPoint(hoPtin);
380 
381  } else { // add a dummy point
383  particle.addPoint(dummy);
384  }
385 
386  pOutputPFSimParticleCollection->push_back(particle);
387  }
388 
389  iEvent.put(std::move(pOutputPFSimParticleCollection));
390  }
391 
392  LogDebug("PFSimParticleProducer") << "STOP event: " << iEvent.id().event() << " in run " << iEvent.id().run() << endl;
393 }
394 
395 void PFSimParticleProducer::getSimIDs(const TrackHandle& trackh, std::vector<unsigned>& recTrackSimID) {
396  if (trackh.isValid()) {
397  for (unsigned i = 0; i < trackh->size(); i++) {
398  const reco::PFRecTrackRef ref(trackh, i);
399 
400  for (auto const& hit : ref->trackRef()->recHits()) {
401  if (hit->isValid()) {
402  auto rechit = dynamic_cast<const FastTrackerRecHit*>(hit);
403 
404  for (unsigned int st_index = 0; st_index < rechit->nSimTrackIds(); ++st_index) {
405  recTrackSimID.push_back(rechit->simTrackId(st_index));
406  }
407  break;
408  }
409  } //loop track rechit
410  } //loop recTracks
411  } //track handle valid
412 }
PFSimParticleProducer(const edm::ParameterSet &)
#define LogDebug(id)
RunNumber_t run() const
Definition: EventID.h:38
double pz() const
z of the momentum
Definition: RawParticle.h:302
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:40
T getUntrackedParameter(std::string const &, T const &) const
bool noEndVertex() const
no end vertex
float charge() const
charge
Definition: FSimTrack.h:56
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
const RawParticle & hoEntrance() const
The particle at HCAL exir.
Definition: FSimTrack.h:155
int onLayer2() const
Definition: FSimTrack.h:106
const FSimVertex & endVertex() const
end vertex
const RawParticle & hcalExit() const
The particle at HCAL exir.
Definition: FSimTrack.h:152
double y() const
y of vertex
Definition: RawParticle.h:283
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
double z() const
z of vertex
Definition: RawParticle.h:284
const XYZTLorentzVector & momentum() const
Temporary (until move of SimTrack to Mathcore) - No! Actually very useful.
Definition: FSimTrack.h:209
const RawParticle & layer1Entrance() const
The particle at Preshower Layer 1.
Definition: FSimTrack.h:137
std::vector< EcalRecHit >::const_iterator const_iterator
const std::vector< int > & daughters() const
Vector of daughter indices.
double x() const
x of vertex
Definition: RawParticle.h:282
int onEcal() const
Definition: FSimTrack.h:111
void produce(edm::Event &, const edm::EventSetup &) override
bool getData(T &iHolder) const
Definition: EventSetup.h:113
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int hashedIndex(int ieta, int iphi)
Definition: EcalPyUtils.cc:38
int iEvent
Definition: GenABIO.cc:224
double e() const
energy of the momentum
Definition: RawParticle.h:305
Point of closest approach from beam axis (initial point in the case of PFSimParticle) ...
const math::XYZTLorentzVector & position() const
Temporary (until CMSSW moves to Mathcore) - No ! Actually very useful.
Definition: FSimVertex.h:48
const RawParticle & ecalEntrance() const
The particle at ECAL entrance.
Definition: FSimTrack.h:143
bool isValid() const
Definition: HandleBase.h:70
true particle for particle flow
Definition: PFSimParticle.h:19
const_iterator end() const
T const * product() const
Definition: Handle.h:69
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
edm::EventID id() const
Definition: EventBase.h:59
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:22
HLT enums.
int id() const
the index in FBaseSimEvent and other vectors
Definition: FSimTrack.h:96
const RawParticle & layer2Entrance() const
The particle at Preshower Layer 2.
Definition: FSimTrack.h:140
int onLayer1() const
Definition: FSimTrack.h:101
int onHcal() const
Definition: FSimTrack.h:116
bool noMother() const
no mother particle
const RawParticle & hcalEntrance() const
The particle at HCAL entrance.
Definition: FSimTrack.h:146
double px() const
x of the momentum
Definition: RawParticle.h:296
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
const FSimVertex vertex() const
Origin vertex.
double py() const
y of the momentum
Definition: RawParticle.h:299
void getSimIDs(const TrackHandle &trackh, std::vector< unsigned > &recTrackSimID)
const FSimTrack & mother() const
mother
def move(src, dest)
Definition: eostools.py:511
const_iterator begin() const