CMS 3D CMS Logo

SimPFProducer.cc
Go to the documentation of this file.
1 // This producer assigns vertex times (with a specified resolution) to tracks.
2 // The times are produced as valuemaps associated to tracks, so the track dataformat doesn't
3 // need to be modified.
4 
5 #include "CLHEP/Units/SystemOfUnits.h"
41 
42 #include <unordered_map>
43 #include <memory>
44 
46 public:
48 
49  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
50 
51  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
52 
53 private:
54  // parameters
56  const bool useTiming_;
57  const bool useTimingQuality_;
59 
60  // inputs
71  // tracking particle associators by order of preference
72  const std::vector<edm::EDGetTokenT<reco::TrackToTrackingParticleAssociator>> associators_;
73 };
74 
77 
79  // simPFProducer
81  desc.add<edm::InputTag>("simClustersSrc", {"particleFlowClusterHGCalFromSimCl"});
82  desc.add<edm::InputTag>("trackSrc", {"generalTracks"});
83  desc.add<std::vector<edm::InputTag>>("associators",
84  {
85  {"quickTrackAssociatorByHits"},
86  });
87  desc.add<edm::InputTag>("pfRecTrackSrc", {"hgcalTrackCollection", "TracksInHGCal"});
88  desc.add<edm::InputTag>("trackingParticleSrc", {"mix", "MergedTrackTruth"});
89  desc.add<double>("neutralEMThreshold", 0.25);
90  desc.add<edm::InputTag>("caloParticlesSrc", {"mix", "MergedCaloTruth"});
91  desc.add<double>("superClusterThreshold", 4.0);
92  desc.add<edm::InputTag>("simClusterTruthSrc", {"mix", "MergedCaloTruth"});
93  desc.add<edm::InputTag>("muonSrc", {"muons1stStep"});
94  desc.add<double>("neutralHADThreshold", 0.25);
95  desc.add<edm::InputTag>("gsfTrackSrc", {"electronGsfTracks"});
96 
97  // if useTiming_
98  desc.add<bool>("useTiming", false);
99  desc.add<bool>("useTimingQuality", false);
100  desc.add<edm::InputTag>("trackTimeValueMap", {});
101  desc.add<edm::InputTag>("trackTimeErrorMap", {});
102  desc.add<edm::InputTag>("trackTimeQualityMap", {});
103  desc.add<double>("timingQualityThreshold", 0);
104  desc.add<edm::InputTag>("gsfTrackTimeValueMap", {});
105  desc.add<edm::InputTag>("gsfTrackTimeErrorMap", {});
106  desc.add<edm::InputTag>("gsfTrackTimeQualityMap", {});
107 
108  descriptions.addWithDefaultLabel(desc);
109 }
110 
111 namespace {
112 
113  template <typename T>
114  uint64_t hashSimInfo(const T& simTruth, size_t i = 0) {
115  uint64_t evtid = simTruth.eventId().rawId();
116  uint64_t trackid = simTruth.g4Tracks()[i].trackId();
117  return ((evtid << 3) + 23401923) ^ trackid;
118  };
119 } // namespace
120 
122  : superClusterThreshold_(conf.getParameter<double>("superClusterThreshold")),
123  neutralEMThreshold_(conf.getParameter<double>("neutralEMThreshold")),
124  neutralHADThreshold_(conf.getParameter<double>("neutralHADThreshold")),
125  useTiming_(conf.getParameter<bool>("useTiming")),
126  useTimingQuality_(conf.getParameter<bool>("useTimingQuality")),
127  timingQualityThreshold_(useTimingQuality_ ? conf.getParameter<double>("timingQualityThreshold") : -99.),
128  pfRecTracks_(consumes<edm::View<reco::PFRecTrack>>(conf.getParameter<edm::InputTag>("pfRecTrackSrc"))),
129  tracks_(consumes<edm::View<reco::Track>>(conf.getParameter<edm::InputTag>("trackSrc"))),
130  gsfTracks_(consumes<edm::View<reco::Track>>(conf.getParameter<edm::InputTag>("gsfTrackSrc"))),
131  muons_(consumes<reco::MuonCollection>(conf.getParameter<edm::InputTag>("muonSrc"))),
132  srcTrackTime_(useTiming_ ? consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("trackTimeValueMap"))
133  : edm::EDGetTokenT<edm::ValueMap<float>>()),
134  srcTrackTimeError_(useTiming_
135  ? consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("trackTimeErrorMap"))
136  : edm::EDGetTokenT<edm::ValueMap<float>>()),
137  srcTrackTimeQuality_(useTimingQuality_
138  ? consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("trackTimeQualityMap"))
139  : edm::EDGetTokenT<edm::ValueMap<float>>()),
140  srcGsfTrackTime_(useTiming_
141  ? consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("gsfTrackTimeValueMap"))
142  : edm::EDGetTokenT<edm::ValueMap<float>>()),
143  srcGsfTrackTimeError_(
144  useTiming_ ? consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("gsfTrackTimeErrorMap"))
145  : edm::EDGetTokenT<edm::ValueMap<float>>()),
146  srcGsfTrackTimeQuality_(
147  useTimingQuality_ ? consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("gsfTrackTimeQualityMap"))
148  : edm::EDGetTokenT<edm::ValueMap<float>>()),
149  trackingParticles_(consumes<TrackingParticleCollection>(conf.getParameter<edm::InputTag>("trackingParticleSrc"))),
150  simClustersTruth_(consumes<SimClusterCollection>(conf.getParameter<edm::InputTag>("simClusterTruthSrc"))),
151  caloParticles_(consumes<CaloParticleCollection>(conf.getParameter<edm::InputTag>("caloParticlesSrc"))),
152  simClusters_(consumes<std::vector<reco::PFCluster>>(conf.getParameter<edm::InputTag>("simClustersSrc"))),
153  associators_(edm::vector_transform(
154  conf.getParameter<std::vector<edm::InputTag>>("associators"),
155  [this](const edm::InputTag& tag) { return this->consumes<reco::TrackToTrackingParticleAssociator>(tag); })) {
156  produces<reco::PFBlockCollection>();
157  produces<reco::SuperClusterCollection>("perfect");
158  produces<reco::PFCandidateCollection>();
159 }
160 
162  //get associators
163  std::vector<edm::Handle<reco::TrackToTrackingParticleAssociator>> associators;
164  for (const auto& token : associators_) {
165  associators.emplace_back();
166  auto& back = associators.back();
167  evt.getByToken(token, back);
168  }
169 
170  //get PFRecTrack
171  edm::Handle<edm::View<reco::PFRecTrack>> PFTrackCollectionH;
172  evt.getByToken(pfRecTracks_, PFTrackCollectionH);
173  const edm::View<reco::PFRecTrack> PFTrackCollection = *PFTrackCollectionH;
174  std::unordered_set<unsigned> PFTrackToGeneralTrack;
175  for (unsigned i = 0; i < PFTrackCollection.size(); ++i) {
176  const auto ptr = PFTrackCollection.ptrAt(i);
177  PFTrackToGeneralTrack.insert(ptr->trackRef().key());
178  }
179 
180  //get track collections
181  edm::Handle<edm::View<reco::Track>> TrackCollectionH;
182  evt.getByToken(tracks_, TrackCollectionH);
183  const edm::View<reco::Track>& TrackCollection = *TrackCollectionH;
184 
186  evt.getByToken(muons_, muons);
187  std::unordered_set<unsigned> MuonTrackToGeneralTrack;
188  for (auto const& mu : *muons.product()) {
189  reco::TrackRef muTrkRef = mu.track();
190  if (muTrkRef.isNonnull())
191  MuonTrackToGeneralTrack.insert(muTrkRef.key());
192  }
193 
194  // get timing, if enabled
195  edm::Handle<edm::ValueMap<float>> trackTimeH, trackTimeErrH, trackTimeQualH, gsfTrackTimeH, gsfTrackTimeErrH,
196  gsfTrackTimeQualH;
197  if (useTiming_) {
198  evt.getByToken(srcTrackTime_, trackTimeH);
199  evt.getByToken(srcTrackTimeError_, trackTimeErrH);
200  evt.getByToken(srcGsfTrackTime_, gsfTrackTimeH);
201  evt.getByToken(srcGsfTrackTimeError_, gsfTrackTimeErrH);
202 
203  if (useTimingQuality_) {
204  evt.getByToken(srcTrackTimeQuality_, trackTimeQualH);
205  evt.getByToken(srcGsfTrackTimeQuality_, gsfTrackTimeQualH);
206  }
207  }
208 
209  //get tracking particle collections
211  evt.getByToken(trackingParticles_, TPCollectionH);
212  //const TrackingParticleCollection& TPCollection = *TPCollectionH;
213 
214  // grab phony clustering information
215  edm::Handle<SimClusterCollection> SimClustersTruthH;
216  evt.getByToken(simClustersTruth_, SimClustersTruthH);
217  const SimClusterCollection& SimClustersTruth = *SimClustersTruthH;
218 
220  evt.getByToken(caloParticles_, CaloParticlesH);
221  const CaloParticleCollection& CaloParticles = *CaloParticlesH;
222 
224  evt.getByToken(simClusters_, SimClustersH);
225  const std::vector<reco::PFCluster>& SimClusters = *SimClustersH;
226 
227  std::unordered_map<uint64_t, size_t> hashToSimCluster;
228 
229  for (unsigned i = 0; i < SimClustersTruth.size(); ++i) {
230  const auto& simTruth = SimClustersTruth[i];
231  hashToSimCluster[hashSimInfo(simTruth)] = i;
232  }
233 
234  // associate the reco tracks / gsf Tracks
235  std::vector<reco::RecoToSimCollection> associatedTracks, associatedTracksGsf;
236  for (const auto& associator : associators) {
237  associatedTracks.emplace_back(associator->associateRecoToSim(TrackCollectionH, TPCollectionH));
238  //associatedTracksGsf.emplace_back(associator->associateRecoToSim(GsfTrackCollectionH, TPCollectionH));
239  }
240 
241  // make blocks out of calo particles so we can have cluster references
242  // likewise fill out superclusters
243  auto superclusters = std::make_unique<reco::SuperClusterCollection>();
244  auto blocks = std::make_unique<reco::PFBlockCollection>();
245  std::unordered_map<size_t, size_t> simCluster2Block;
246  std::unordered_map<size_t, size_t> simCluster2BlockIndex;
247  std::unordered_multimap<size_t, size_t> caloParticle2SimCluster;
248  std::vector<int> caloParticle2SuperCluster;
249  for (unsigned icp = 0; icp < CaloParticles.size(); ++icp) {
250  blocks->emplace_back();
251  auto& block = blocks->back();
252  const auto& simclusters = CaloParticles[icp].simClusters();
253  double pttot = 0.0;
254  double etot = 0.0;
255  std::vector<size_t> good_simclusters;
256  for (unsigned isc = 0; isc < simclusters.size(); ++isc) {
257  auto simc = simclusters[isc];
258  auto pdgId = std::abs(simc->pdgId());
259  edm::Ref<std::vector<reco::PFCluster>> clusterRef(SimClustersH, simc.key());
260  if (((pdgId == 22 || pdgId == 11) && clusterRef->energy() > neutralEMThreshold_) ||
261  clusterRef->energy() > neutralHADThreshold_) {
262  good_simclusters.push_back(isc);
263  etot += clusterRef->energy();
264  pttot += clusterRef->pt();
265  auto bec = std::make_unique<reco::PFBlockElementCluster>(clusterRef, reco::PFBlockElement::HGCAL);
266  block.addElement(bec.get());
267  simCluster2Block[simc.key()] = icp;
268  simCluster2BlockIndex[simc.key()] = bec->index();
269  caloParticle2SimCluster.emplace(CaloParticles[icp].g4Tracks()[0].trackId(), simc.key());
270  }
271  }
272 
273  auto pdgId = std::abs(CaloParticles[icp].pdgId());
274 
275  caloParticle2SuperCluster.push_back(-1);
276  if ((pdgId == 22 || pdgId == 11) && pttot > superClusterThreshold_) {
277  caloParticle2SuperCluster[icp] = superclusters->size();
278 
279  math::XYZPoint seedpos; // take seed pos as supercluster point
282  for (auto idx : good_simclusters) {
283  edm::Ptr<reco::PFCluster> ptr(SimClustersH, simclusters[idx].key());
284  clusters.push_back(ptr);
285  if (seed.isNull() || seed->energy() < ptr->energy()) {
286  seed = ptr;
287  seedpos = ptr->position();
288  }
289  }
290  superclusters->emplace_back(etot, seedpos, seed, clusters);
291  }
292  }
293 
294  auto blocksHandle = evt.put(std::move(blocks));
295  auto superClustersHandle = evt.put(std::move(superclusters), "perfect");
296 
297  // list tracks so we can mark them as used and/or fight over them
298  std::vector<bool> usedTrack(TrackCollection.size(), false),
299  //usedGsfTrack(GsfTrackCollection.size(),false),
300  usedSimCluster(SimClusters.size(), false);
301 
302  auto candidates = std::make_unique<reco::PFCandidateCollection>();
303  // in good particle flow fashion, start from the tracks and go out
304  for (unsigned itk = 0; itk < TrackCollection.size(); ++itk) {
305  auto tkRef = TrackCollection.refAt(itk);
306  // skip tracks not selected by PF
307  if (PFTrackToGeneralTrack.count(itk) == 0)
308  continue;
309  reco::RecoToSimCollection::const_iterator assoc_tps = associatedTracks.back().end();
310  for (const auto& association : associatedTracks) {
311  assoc_tps = association.find(tkRef);
312  if (assoc_tps != association.end())
313  break;
314  }
315  if (assoc_tps == associatedTracks.back().end())
316  continue;
317  // assured now that we are matched to a set of tracks
318  const auto& matches = assoc_tps->val;
319 
320  const auto absPdgId = std::abs(matches[0].first->pdgId());
321  const auto charge = tkRef->charge();
322  const auto three_mom = tkRef->momentum();
323  constexpr double mpion2 = 0.13957 * 0.13957;
324  double energy = std::sqrt(three_mom.mag2() + mpion2);
325  math::XYZTLorentzVector trk_p4(three_mom.x(), three_mom.y(), three_mom.z(), energy);
326 
328 
329  switch (absPdgId) {
330  case 11:
331  part_type = reco::PFCandidate::e;
332  break;
333  case 13:
334  part_type = reco::PFCandidate::mu;
335  break;
336  default:
337  part_type = reco::PFCandidate::h;
338  }
339 
340  candidates->emplace_back(charge, trk_p4, part_type);
341  auto& candidate = candidates->back();
342 
343  candidate.setTrackRef(tkRef.castTo<reco::TrackRef>());
344 
345  if (useTiming_) {
346  // check if track-mtd match is of sufficient quality
347  const bool assocQuality = useTimingQuality_ ? (*trackTimeQualH)[tkRef] > timingQualityThreshold_ : true;
348  if (assocQuality) {
349  candidate.setTime((*trackTimeH)[tkRef], (*trackTimeErrH)[tkRef]);
350  } else {
351  candidate.setTime(0., -1.);
352  }
353  }
354 
355  // bind to cluster if there is one and try to gather conversions, etc
356  for (const auto& match : matches) {
357  uint64_t hash = hashSimInfo(*(match.first));
358  if (hashToSimCluster.count(hash)) {
359  auto simcHash = hashToSimCluster[hash];
360 
361  if (!usedSimCluster[simcHash]) {
362  if (simCluster2Block.count(simcHash) && simCluster2BlockIndex.count(simcHash)) {
363  size_t block = simCluster2Block.find(simcHash)->second;
364  size_t blockIdx = simCluster2BlockIndex.find(simcHash)->second;
365  edm::Ref<reco::PFBlockCollection> blockRef(blocksHandle, block);
366  candidate.addElementInBlock(blockRef, blockIdx);
367  usedSimCluster[simcHash] = true;
368  }
369  }
370  if (absPdgId == 11) { // collect brems/conv. brems
371  if (simCluster2Block.count(simcHash)) {
372  auto block_index = simCluster2Block.find(simcHash)->second;
373  auto supercluster_index = caloParticle2SuperCluster[block_index];
374  if (supercluster_index != -1) {
375  edm::Ref<reco::PFBlockCollection> blockRef(blocksHandle, block_index);
376  for (const auto& elem : blockRef->elements()) {
377  const auto& ref = elem.clusterRef();
378  if (!usedSimCluster[ref.key()]) {
379  candidate.addElementInBlock(blockRef, elem.index());
380  usedSimCluster[ref.key()] = true;
381  }
382  }
383 
384  //*TODO* cluster time is not reliable at the moment, so just keep time from the track if available
385  }
386  }
387  }
388  }
389  // Now try to include also electrons that have been reconstructed using
390  // the GraphCaloParticles. In particular, recover the cases in which the
391  // tracking particle associated to the CaloParticle has not left any hits
392  // in the calorimeters or, if it had, the cluster has been skipped due to
393  // threshold requirements.
394  if (caloParticle2SimCluster.count(match.first->g4Tracks()[0].trackId())) {
395  auto range = caloParticle2SimCluster.equal_range(match.first->g4Tracks()[0].trackId());
396  for (auto it = range.first; it != range.second; ++it) {
397  if (!usedSimCluster[it->second]) {
398  usedSimCluster[it->second] = true;
399  if (simCluster2Block.find(it->second) != simCluster2Block.end()) {
400  size_t block = simCluster2Block.find(it->second)->second;
401  size_t blockIdx = simCluster2BlockIndex.find(it->second)->second;
402  edm::Ref<reco::PFBlockCollection> blockRef(blocksHandle, block);
403  candidate.addElementInBlock(blockRef, blockIdx);
404  }
405  }
406  }
407  }
408  }
409  usedTrack[tkRef.key()] = true;
410  // remove tracks already used by muons
411  if (MuonTrackToGeneralTrack.count(itk) || absPdgId == 13)
412  candidates->pop_back();
413  }
414 
415  // now loop over the non-collected clusters in blocks
416  // and turn them into neutral hadrons or photons
417  const auto& theblocks = *blocksHandle;
418  for (unsigned ibl = 0; ibl < theblocks.size(); ++ibl) {
419  reco::PFBlockRef blref(blocksHandle, ibl);
420  const auto& elements = theblocks[ibl].elements();
421  for (const auto& elem : elements) {
422  const auto& ref = elem.clusterRef();
423  const auto& simtruth = SimClustersTruth[ref.key()];
425  if (!usedSimCluster[ref.key()]) {
426  auto absPdgId = std::abs(simtruth.pdgId());
427  switch (absPdgId) {
428  case 11:
429  case 22:
430  part_type = reco::PFCandidate::gamma;
431  break;
432  default:
433  part_type = reco::PFCandidate::h0;
434  }
435  const auto three_mom = (ref->position() - math::XYZPoint(0, 0, 0)).unit() * ref->correctedEnergy();
436  math::XYZTLorentzVector clu_p4(three_mom.x(), three_mom.y(), three_mom.z(), ref->correctedEnergy());
437  candidates->emplace_back(0, clu_p4, part_type);
438  auto& candidate = candidates->back();
439  candidate.addElementInBlock(blref, elem.index());
440  }
441  }
442  }
443 
444  evt.put(std::move(candidates));
445 }
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:154
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
const edm::EDGetTokenT< edm::ValueMap< float > > srcTrackTimeError_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
ParticleType
particle types
Definition: PFCandidate.h:44
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:232
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:526
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
muons
the two sets of parameters below are mutually exclusive, depending if RECO or ALCARECO is used the us...
Definition: DiMuonV_cfg.py:214
key_type key() const
Accessor for product key.
Definition: Ref.h:244
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
std::tuple< layerClusterToCaloParticle, caloParticleToLayerCluster > association
const double neutralHADThreshold_
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
const edm::EDGetTokenT< edm::ValueMap< float > > srcTrackTimeQuality_
const double timingQualityThreshold_
const edm::EDGetTokenT< CaloParticleCollection > caloParticles_
const std::vector< edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > > associators_
T sqrt(T t)
Definition: SSEVec.h:23
double energy() const
cluster energy
Definition: PFCluster.h:74
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
key
prepare the HTCondor submission files and eventually submit them
const edm::EDGetTokenT< edm::View< reco::Track > > tracks_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::EDGetTokenT< edm::ValueMap< float > > srcGsfTrackTime_
const dim3 blockIdx
Definition: cudaCompat.h:32
Basic3DVector unit() const
const edm::EDGetTokenT< edm::View< reco::PFRecTrack > > pfRecTracks_
const edm::EDGetTokenT< SimClusterCollection > simClustersTruth_
const edm::EDGetTokenT< TrackingParticleCollection > trackingParticles_
unsigned long long uint64_t
Definition: Time.h:13
const bool useTimingQuality_
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
key_type index() const
Definition: Ref.h:247
SimPFProducer(const edm::ParameterSet &)
const double superClusterThreshold_
const edm::EDGetTokenT< reco::MuonCollection > muons_
const edm::EDGetTokenT< edm::ValueMap< float > > srcGsfTrackTimeQuality_
std::vector< l1t::PFTrack > PFTrackCollection
Definition: PFTrack.h:89
const edm::EDGetTokenT< edm::ValueMap< float > > srcTrackTime_
fixed size matrix
HLT enums.
constexpr float mpion2
Definition: Common.h:41
std::vector< CaloParticle > CaloParticleCollection
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const edm::EDGetTokenT< edm::ValueMap< float > > srcGsfTrackTimeError_
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
const bool useTiming_
const double neutralEMThreshold_
std::vector< TrackingParticle > TrackingParticleCollection
std::vector< SimCluster > SimClusterCollection
const edm::EDGetTokenT< std::vector< reco::PFCluster > > simClusters_
long double T
def move(src, dest)
Definition: eostools.py:511
const edm::EDGetTokenT< edm::View< reco::Track > > gsfTracks_