CMS 3D CMS Logo

PFLinker.cc
Go to the documentation of this file.
1 
20 
21 namespace edm {
22  class EventSetup;
23 } // namespace edm
24 
26 public:
27  explicit PFLinker(const edm::ParameterSet&);
28 
29  ~PFLinker() override;
30 
31  void produce(edm::Event&, const edm::EventSetup&) override;
32  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
33 
34 private:
35  template <typename TYPE>
39  edm::Handle<TYPE>& inputObjCollection,
40  const std::map<edm::Ref<TYPE>, reco::PFCandidatePtr>& mapToTheCandidate,
41  const edm::OrphanHandle<reco::PFCandidateCollection>& newPFCandColl) const;
42 
43 private:
45  std::vector<edm::EDGetTokenT<reco::PFCandidateCollection>> inputTagPFCandidates_;
46 
49 
52 
59 
62 
65 
68 
71 
74 
77 };
78 
80 
83  desc.add<std::vector<edm::InputTag>>("PFCandidate", {edm::InputTag("particleFlow")});
84  desc.add<edm::InputTag>("GsfElectrons", {"gedGsfElectrons"});
85  desc.add<edm::InputTag>("Photons", {"gedPhotons"});
86  desc.add<edm::InputTag>("Muons", {"muons", "muons1stStep2muonsMap"});
87  desc.add<bool>("ProducePFCandidates", true);
88  desc.add<bool>("FillMuonRefs", true);
89  desc.add<std::string>("OutputPF", "");
90  desc.add<std::string>("ValueMapElectrons", "electrons");
91  desc.add<std::string>("ValueMapPhotons", "photons");
92  desc.add<std::string>("ValueMapMerged", "all");
93  desc.add<bool>("forceElectronsInHGCAL", false);
94  descriptions.addWithDefaultLabel(desc);
95 }
96 
98  // vector of InputTag; more than 1 is not for RECO, it is for analysis
99 
100  std::vector<edm::InputTag> tags = iConfig.getParameter<std::vector<edm::InputTag>>("PFCandidate");
101  for (unsigned int i = 0; i < tags.size(); ++i)
102  inputTagPFCandidates_.push_back(consumes<reco::PFCandidateCollection>(tags[i]));
103 
104  inputTagGsfElectrons_ = consumes<reco::GsfElectronCollection>(iConfig.getParameter<edm::InputTag>("GsfElectrons"));
105 
106  inputTagPhotons_ = consumes<reco::PhotonCollection>(iConfig.getParameter<edm::InputTag>("Photons"));
107 
108  muonTag_ = iConfig.getParameter<edm::InputTag>("Muons");
109  inputTagMuons_ = consumes<reco::MuonCollection>(edm::InputTag(muonTag_.label()));
110  inputTagMuonMap_ = consumes<reco::MuonToMuonMap>(muonTag_);
111 
112  nameOutputPF_ = iConfig.getParameter<std::string>("OutputPF");
113 
114  nameOutputElectronsPF_ = iConfig.getParameter<std::string>("ValueMapElectrons");
115 
116  nameOutputPhotonsPF_ = iConfig.getParameter<std::string>("ValueMapPhotons");
117 
118  producePFCandidates_ = iConfig.getParameter<bool>("ProducePFCandidates");
119 
120  nameOutputMergedPF_ = iConfig.getParameter<std::string>("ValueMapMerged");
121 
122  fillMuonRefs_ = iConfig.getParameter<bool>("FillMuonRefs");
123 
124  forceElectronsInHGCAL_ = iConfig.getParameter<bool>("forceElectronsInHGCAL");
125 
126  // should not produce PFCandidates and read seve
127  if (producePFCandidates_ && inputTagPFCandidates_.size() > 1) {
128  edm::LogError("PFLinker")
129  << " cannot read several collections of PFCandidates and produce a new collection at the same time. "
130  << std::endl;
131  assert(false);
132  }
133  if (producePFCandidates_) {
134  produces<reco::PFCandidateCollection>(nameOutputPF_);
135  }
136  produces<edm::ValueMap<reco::PFCandidatePtr>>(nameOutputElectronsPF_);
137  produces<edm::ValueMap<reco::PFCandidatePtr>>(nameOutputPhotonsPF_);
138  produces<edm::ValueMap<reco::PFCandidatePtr>>(nameOutputMergedPF_);
139  if (fillMuonRefs_)
140  produces<edm::ValueMap<reco::PFCandidatePtr>>(muonTag_.label());
141 }
142 
144 
146  auto pfCandidates_p = std::make_unique<reco::PFCandidateCollection>();
147 
148  auto gsfElectrons = iEvent.getHandle(inputTagGsfElectrons_);
149 
150  std::map<reco::GsfElectronRef, reco::PFCandidatePtr> electronCandidateMap;
151 
152  auto photons = iEvent.getHandle(inputTagPhotons_);
153  std::map<reco::PhotonRef, reco::PFCandidatePtr> photonCandidateMap;
154 
156  if (fillMuonRefs_)
157  muonMap = iEvent.getHandle(inputTagMuonMap_);
158  std::map<reco::MuonRef, reco::PFCandidatePtr> muonCandidateMap;
159 
160  unsigned nColPF = inputTagPFCandidates_.size();
161 
162  for (unsigned icol = 0; icol < nColPF; ++icol) {
163  auto pfCandidates = iEvent.getHandle(inputTagPFCandidates_[icol]);
164  unsigned ncand = pfCandidates->size();
165 
166  for (unsigned i = 0; i < ncand; ++i) {
168  reco::PFCandidate cand(candPtr);
169 
170  bool isphoton = cand.particleId() == reco::PFCandidate::gamma && cand.mva_nothing_gamma() > 0.;
171  bool iselectron = cand.particleId() == reco::PFCandidate::e;
172  // PFCandidates may have a valid MuonRef though they are not muons.
173  bool hasNonNullMuonRef = cand.muonRef().isNonnull() && fillMuonRefs_;
174 
175  // if not an electron or a photon or a muon just fill the PFCandidate collection
176  if (!(isphoton || iselectron || hasNonNullMuonRef)) {
177  pfCandidates_p->push_back(cand);
178  continue;
179  }
180 
181  if (hasNonNullMuonRef) {
182  reco::MuonRef muRef = (*muonMap)[cand.muonRef()];
183  cand.setMuonRef(muRef);
184  muonCandidateMap[muRef] = candPtr;
185  }
186 
187  // if it is an electron. Find the GsfElectron with the same GsfTrack
188  if (iselectron) {
189  const reco::GsfTrackRef& gsfTrackRef(cand.gsfTrackRef());
190  auto itcheck = find_if(gsfElectrons->begin(), gsfElectrons->end(), [&gsfTrackRef](const auto& ele) {
191  return (ele.gsfTrack() == gsfTrackRef);
192  });
193  if (itcheck == gsfElectrons->end()) {
194  if (!forceElectronsInHGCAL_) {
195  std::ostringstream err;
196  err << " Problem in PFLinker: no GsfElectron " << std::endl;
197  edm::LogError("PFLinker") << err.str();
198  } else {
199  LogDebug("PFLinker") << "Forcing an electron pfCandidate at: " << cand.eta() << " in HGCAL" << std::endl;
200  pfCandidates_p->push_back(cand);
201  }
202  continue; // Watch out ! Continue
203  }
204  reco::GsfElectronRef electronRef(gsfElectrons, itcheck - gsfElectrons->begin());
205  cand.setGsfElectronRef(electronRef);
206  cand.setSuperClusterRef(electronRef->superCluster());
207  // update energy information since now it is done post-particleFlowTmp
208  cand.setEcalEnergy(electronRef->superCluster()->rawEnergy(), electronRef->ecalEnergy());
209  cand.setDeltaP(electronRef->p4Error(reco::GsfElectron::P4_COMBINATION));
210  cand.setP4(electronRef->p4(reco::GsfElectron::P4_COMBINATION));
211  electronCandidateMap[electronRef] = candPtr;
212  }
213 
214  // if it is a photon, find the one with the same PF super-cluster
215  if (isphoton) {
216  const reco::SuperClusterRef& scRef(cand.superClusterRef());
217  auto itcheck = find_if(
218  photons->begin(), photons->end(), [&scRef](const auto& photon) { return photon.superCluster() == scRef; });
219  if (itcheck == photons->end()) {
220  std::ostringstream err;
221  err << " Problem in PFLinker: no Photon " << std::endl;
222  edm::LogError("PFLinker") << err.str();
223  continue; // Watch out ! Continue
224  }
225  reco::PhotonRef photonRef(photons, itcheck - photons->begin());
226  cand.setPhotonRef(photonRef);
227  cand.setSuperClusterRef(photonRef->superCluster());
228  // update energy information since now it is done post-particleFlowTmp
229  cand.setEcalEnergy(photonRef->superCluster()->rawEnergy(),
230  photonRef->getCorrectedEnergy(reco::Photon::regression2));
231  cand.setDeltaP(photonRef->getCorrectedEnergyError(reco::Photon::regression2));
232  cand.setP4(photonRef->p4(reco::Photon::regression2));
233  photonCandidateMap[photonRef] = candPtr;
234  }
235 
236  pfCandidates_p->push_back(cand);
237  }
238  // save the PFCandidates and get a valid handle
239  }
240  const edm::OrphanHandle<reco::PFCandidateCollection> pfCandidateRefProd =
241  (producePFCandidates_) ? iEvent.put(std::move(pfCandidates_p), nameOutputPF_)
243 
244  // now make the valuemaps
245 
246  edm::ValueMap<reco::PFCandidatePtr> pfMapGsfElectrons = fillValueMap<reco::GsfElectronCollection>(
247  iEvent, nameOutputElectronsPF_, gsfElectrons, electronCandidateMap, pfCandidateRefProd);
248 
249  edm::ValueMap<reco::PFCandidatePtr> pfMapPhotons = fillValueMap<reco::PhotonCollection>(
250  iEvent, nameOutputPhotonsPF_, photons, photonCandidateMap, pfCandidateRefProd);
251 
253 
254  if (fillMuonRefs_) {
255  auto muons = iEvent.getHandle(inputTagMuons_);
256 
257  pfMapMuons =
258  fillValueMap<reco::MuonCollection>(iEvent, muonTag_.label(), muons, muonCandidateMap, pfCandidateRefProd);
259  }
260 
261  auto pfMapMerged = std::make_unique<edm::ValueMap<reco::PFCandidatePtr>>();
262  edm::ValueMap<reco::PFCandidatePtr>::Filler pfMapMergedFiller(*pfMapMerged);
263 
264  *pfMapMerged += pfMapGsfElectrons;
265  *pfMapMerged += pfMapPhotons;
266  if (fillMuonRefs_)
267  *pfMapMerged += pfMapMuons;
268 
269  iEvent.put(std::move(pfMapMerged), nameOutputMergedPF_);
270 }
271 
272 template <typename TYPE>
274  edm::Event& event,
276  edm::Handle<TYPE>& inputObjCollection,
277  const std::map<edm::Ref<TYPE>, reco::PFCandidatePtr>& mapToTheCandidate,
278  const edm::OrphanHandle<reco::PFCandidateCollection>& newPFCandColl) const {
279  auto pfMap_p = std::make_unique<edm::ValueMap<reco::PFCandidatePtr>>();
281 
282  typedef typename std::map<edm::Ref<TYPE>, reco::PFCandidatePtr>::const_iterator MapTYPE_it;
283 
284  unsigned nObj = inputObjCollection->size();
285  std::vector<reco::PFCandidatePtr> values(nObj);
286 
287  for (unsigned iobj = 0; iobj < nObj; ++iobj) {
288  edm::Ref<TYPE> objRef(inputObjCollection, iobj);
289  MapTYPE_it itcheck = mapToTheCandidate.find(objRef);
290 
291  reco::PFCandidatePtr candPtr;
292 
293  if (itcheck != mapToTheCandidate.end())
294  candPtr = producePFCandidates_ ? reco::PFCandidatePtr(newPFCandColl, itcheck->second.key()) : itcheck->second;
295 
296  values[iobj] = candPtr;
297  }
298 
299  filler.insert(inputObjCollection, values.begin(), values.end());
300  filler.fill();
301  edm::ValueMap<reco::PFCandidatePtr> returnValue = *pfMap_p;
302  event.put(std::move(pfMap_p), label);
303  return returnValue;
304 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
edm::InputTag muonTag_
Input Muons.
Definition: PFLinker.cc:54
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
edm::EDGetTokenT< reco::MuonCollection > inputTagMuons_
Definition: PFLinker.cc:55
void produce(edm::Event &, const edm::EventSetup &) override
Definition: PFLinker.cc:145
bool producePFCandidates_
Flags - if true: References will be towards new collection ; if false to the original one...
Definition: PFLinker.cc:70
std::string const & label() const
Definition: InputTag.h:36
Log< level::Error, false > LogError
edm::EDGetTokenT< reco::PhotonCollection > inputTagPhotons_
Input Photons.
Definition: PFLinker.cc:51
muons
the two sets of parameters below are mutually exclusive, depending if RECO or ALCARECO is used the us...
Definition: DiMuonV_cfg.py:214
assert(be >=bs)
PFLinker(const edm::ParameterSet &)
Definition: PFLinker.cc:97
edm::EDGetTokenT< reco::MuonToMuonMap > inputTagMuonMap_
Definition: PFLinker.cc:56
char const * label
int iEvent
Definition: GenABIO.cc:224
std::string nameOutputElectronsPF_
name of output ValueMap electrons
Definition: PFLinker.cc:61
std::string nameOutputPhotonsPF_
name of output ValueMap photons
Definition: PFLinker.cc:64
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::string nameOutputPF_
name of output collection of PFCandidate
Definition: PFLinker.cc:58
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: PFLinker.cc:81
bool fillMuonRefs_
Set muon refs and produce the value map?
Definition: PFLinker.cc:73
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
HLT enums.
bool forceElectronsInHGCAL_
Put Electrons within HGCAL coming from SimPFProducer.
Definition: PFLinker.cc:76
~PFLinker() override
Definition: PFLinker.cc:143
edm::EDGetTokenT< reco::GsfElectronCollection > inputTagGsfElectrons_
Input GsfElectrons.
Definition: PFLinker.cc:48
def move(src, dest)
Definition: eostools.py:511
std::vector< edm::EDGetTokenT< reco::PFCandidateCollection > > inputTagPFCandidates_
Input PFCandidates.
Definition: PFLinker.cc:45
Definition: event.py:1
std::string nameOutputMergedPF_
name of output merged ValueMap
Definition: PFLinker.cc:67
#define LogDebug(id)
edm::ValueMap< reco::PFCandidatePtr > fillValueMap(edm::Event &event, std::string label, edm::Handle< TYPE > &inputObjCollection, const std::map< edm::Ref< TYPE >, reco::PFCandidatePtr > &mapToTheCandidate, const edm::OrphanHandle< reco::PFCandidateCollection > &newPFCandColl) const
Definition: PFLinker.cc:273