CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PhotonIDValueMapProducer.cc
Go to the documentation of this file.
3 
6 
8 
11 
14 
17 
19 
21 
22 #include <memory>
23 #include <vector>
24 
25 // Introduce these two types for brevity of casting
28 
29 // This template function finds whether theCandidate is in thefootprint
30 // collection. It is templated to be able to handle both reco and pat
31 // photons (from AOD and miniAOD, respectively).
32 template <class T, class U>
33 bool isInFootprint(const T& thefootprint, const U& theCandidate) {
34  for ( auto itr = thefootprint.begin(); itr != thefootprint.end(); ++itr ) {
35  if( itr->key() == theCandidate.key() ) return true;
36  }
37  return false;
38 }
39 
41 
42  public:
43 
46 
47  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
48 
49  private:
50 
51  virtual void produce(edm::Event&, const edm::EventSetup&) override;
52 
55  const std::vector<float> & values,
56  const std::string & label) const ;
57 
58  // This function computes charged hadron isolation
59  // with respect to multiple PVs and returns the worst
60  // of the found isolation values.
61  // The function implements the computation method taken directly
62  // from Run 1 code of H->gamma gamma, specifically from
63  // the class CiCPhotonID of the HiggsTo2photons anaysis code.
64  // Template is introduced to handle reco/pat photons and aod/miniAOD
65  // PF candidates collections
66  template <class T, class U>
67  float computeWorstPFChargedIsolation(const T& photon,
68  const U& pfCandidates,
70  bool isAOD,
71  float dRmax, float dxyMax, float dzMax);
72 
73  // Some helper functions that are needed to access info in
74  // AOD vs miniAOD
76  candidatePdgId(const edm::Ptr<reco::Candidate> candidate, bool isAOD);
77 
78  const reco::Track* getTrackPointer(const edm::Ptr<reco::Candidate> candidate, bool isAOD);
79 
80 
81  // The object that will compute 5x5 quantities
82  std::unique_ptr<noZS::EcalClusterLazyTools> lazyToolnoZS;
83 
84  // for AOD case
92 
93  // for miniAOD case
101 
102  // Cluster shapes
103  constexpr static char phoFull5x5SigmaIEtaIEta_[] = "phoFull5x5SigmaIEtaIEta";
104  constexpr static char phoFull5x5SigmaIEtaIPhi_[] = "phoFull5x5SigmaIEtaIPhi";
105  constexpr static char phoFull5x5E1x3_[] = "phoFull5x5E1x3";
106  constexpr static char phoFull5x5E2x2_[] = "phoFull5x5E2x2";
107  constexpr static char phoFull5x5E2x5Max_[] = "phoFull5x5E2x5Max";
108  constexpr static char phoFull5x5E5x5_[] = "phoFull5x5E5x5";
109  constexpr static char phoESEffSigmaRR_[] = "phoESEffSigmaRR";
110  // Isolations
111  constexpr static char phoChargedIsolation_[] = "phoChargedIsolation";
112  constexpr static char phoNeutralHadronIsolation_[] = "phoNeutralHadronIsolation";
113  constexpr static char phoPhotonIsolation_[] = "phoPhotonIsolation";
114  constexpr static char phoWorstChargedIsolation_[] = "phoWorstChargedIsolation";
115 
116 };
117 
118 // Cluster shapes
126 // Isolations
131 
133 
134  //
135  // Declare consummables, handle both AOD and miniAOD case
136  //
137  ebReducedRecHitCollection_ = mayConsume<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>
138  ("ebReducedRecHitCollection"));
139  ebReducedRecHitCollectionMiniAOD_ = mayConsume<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>
140  ("ebReducedRecHitCollectionMiniAOD"));
141 
142  eeReducedRecHitCollection_ = mayConsume<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>
143  ("eeReducedRecHitCollection"));
144  eeReducedRecHitCollectionMiniAOD_ = mayConsume<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>
145  ("eeReducedRecHitCollectionMiniAOD"));
146 
147  esReducedRecHitCollection_ = mayConsume<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>
148  ("esReducedRecHitCollection"));
149  esReducedRecHitCollectionMiniAOD_ = mayConsume<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>
150  ("esReducedRecHitCollectionMiniAOD"));
151 
152  vtxToken_ = mayConsume<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"));
153  vtxTokenMiniAOD_ = mayConsume<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("verticesMiniAOD"));
154 
155  // reco photons are castable into pat photons, so no need to handle reco/pat seprately
156  src_ = mayConsume<edm::View<reco::Photon> >(iConfig.getParameter<edm::InputTag>("src"));
157  srcMiniAOD_ = mayConsume<edm::View<reco::Photon> >(iConfig.getParameter<edm::InputTag>("srcMiniAOD"));
158 
159  // The particleBasedIsolation object is relevant only for AOD, RECO format
160  particleBasedIsolationToken_ = mayConsume<edm::ValueMap<std::vector<reco::PFCandidateRef > > >
161  (iConfig.getParameter<edm::InputTag>("particleBasedIsolation"));
162 
163  // AOD has reco::PFCandidate vector, and miniAOD has pat::PackedCandidate vector.
164  // Both inherit from reco::Candidate, so we go to the base class. We can cast into
165  // the full type later if needed. Since the collection names are different, we
166  // introduce both collections
167  pfCandidatesToken_ = mayConsume< edm::View<reco::Candidate> >(iConfig.getParameter<edm::InputTag>("pfCandidates"));
168  pfCandidatesTokenMiniAOD_ = mayConsume< edm::View<reco::Candidate> >(iConfig.getParameter<edm::InputTag>("pfCandidatesMiniAOD"));
169 
170  //
171  // Declare producibles
172  //
173  // Cluster shapes
174  produces<edm::ValueMap<float> >(phoFull5x5SigmaIEtaIEta_);
175  produces<edm::ValueMap<float> >(phoFull5x5SigmaIEtaIPhi_);
176  produces<edm::ValueMap<float> >(phoFull5x5E1x3_);
177  produces<edm::ValueMap<float> >(phoFull5x5E2x2_);
178  produces<edm::ValueMap<float> >(phoFull5x5E2x5Max_);
179  produces<edm::ValueMap<float> >(phoFull5x5E5x5_);
180  produces<edm::ValueMap<float> >(phoESEffSigmaRR_);
181  // Isolations
182  produces<edm::ValueMap<float> >(phoChargedIsolation_);
183  produces<edm::ValueMap<float> >(phoNeutralHadronIsolation_);
184  produces<edm::ValueMap<float> >(phoPhotonIsolation_);
185  produces<edm::ValueMap<float> >(phoWorstChargedIsolation_);
186 
187 }
188 
190 }
191 
193 
194  using namespace edm;
195 
196  // Constants
197  const float coneSizeDR = 0.3;
198  const float dxyMax = 0.1;
199  const float dzMax = 0.2;
200 
202 
203  bool isAOD = true;
204  iEvent.getByToken(src_, src);
205  if( !src.isValid() ){
206  isAOD = false;
207  iEvent.getByToken(srcMiniAOD_,src);
208  }
209  if( !src.isValid() ) {
210  throw cms::Exception("IllDefinedDataTier")
211  << "DataFormat does not contain a photon source!";
212  }
213 
214  // Configure Lazy Tools
215  if( isAOD )
216  lazyToolnoZS = std::unique_ptr<noZS::EcalClusterLazyTools>(new noZS::EcalClusterLazyTools(iEvent, iSetup,
220  else
221  lazyToolnoZS = std::unique_ptr<noZS::EcalClusterLazyTools>(new noZS::EcalClusterLazyTools(iEvent, iSetup,
225 
226  // Get PV
228  iEvent.getByToken(vtxToken_, vertices);
229  if( !vertices.isValid() )
231  if (vertices->empty()) return; // skip the event if no PV found
232  const reco::Vertex &pv = vertices->front();
233 
235  if( isAOD ){
236  // this exists only in AOD
237  iEvent.getByToken(particleBasedIsolationToken_, particleBasedIsolationMap);
238  }
239 
240  edm::Handle< edm::View<reco::Candidate> > pfCandidatesHandle;
241 
242  iEvent.getByToken(pfCandidatesToken_, pfCandidatesHandle);
243  if( !pfCandidatesHandle.isValid() )
244  iEvent.getByToken(pfCandidatesTokenMiniAOD_, pfCandidatesHandle);
245 
246  if( !isAOD && src->size() ) {
247  edm::Ptr<pat::Photon> test(src->ptrAt(0));
248  if( test.isNull() || !test.isAvailable() ) {
249  throw cms::Exception("InvalidConfiguration")
250  <<"DataFormat is detected as miniAOD but cannot cast to pat::Photon!";
251  }
252  }
253 
254  // size_t n = src->size();
255  // Cluster shapes
256  std::vector<float> phoFull5x5SigmaIEtaIEta;
257  std::vector<float> phoFull5x5SigmaIEtaIPhi;
258  std::vector<float> phoFull5x5E1x3;
259  std::vector<float> phoFull5x5E2x2;
260  std::vector<float> phoFull5x5E2x5Max;
261  std::vector<float> phoFull5x5E5x5;
262  std::vector<float> phoESEffSigmaRR;
263  // Isolations
264  std::vector<float> phoChargedIsolation;
265  std::vector<float> phoNeutralHadronIsolation;
266  std::vector<float> phoPhotonIsolation;
267  std::vector<float> phoWorstChargedIsolation;
268 
269  // reco::Photon::superCluster() is virtual so we can exploit polymorphism
270  for (unsigned idxpho = 0; idxpho < src->size(); ++idxpho) {
271  const auto& iPho = src->ptrAt(idxpho);
272 
273  //
274  // Compute full 5x5 quantities
275  //
276  const auto& theseed = *(iPho->superCluster()->seed());
277 
278  // For full5x5_sigmaIetaIeta, for 720 we use: lazy tools for AOD,
279  // and userFloats or lazy tools for miniAOD. From some point in 72X and on, one can
280  // retrieve the full5x5 directly from the object with ->full5x5_sigmaIetaIeta()
281  // for both formats.
282  float see = -999;
283  std::vector<float> vCov = lazyToolnoZS->localCovariances( theseed );
284  see = (isnan(vCov[0]) ? 0. : sqrt(vCov[0]));
285  float sep = vCov[1];
286  phoFull5x5SigmaIEtaIEta.push_back(see);
287  phoFull5x5SigmaIEtaIPhi.push_back(sep);
288 
289  phoFull5x5E1x3 .push_back(lazyToolnoZS-> e1x3 (theseed) );
290  phoFull5x5E2x2 .push_back(lazyToolnoZS-> e2x2 (theseed) );
291  phoFull5x5E2x5Max.push_back(lazyToolnoZS-> e2x5Max(theseed) );
292  phoFull5x5E5x5 .push_back(lazyToolnoZS-> e5x5 (theseed) );
293 
294  phoESEffSigmaRR .push_back(lazyToolnoZS->eseffsirir( *(iPho->superCluster()) ) );
295 
296  //
297  // Compute absolute uncorrected isolations with footprint removal
298  //
299 
300  // First, find photon direction with respect to the good PV
301  math::XYZVector photon_directionWrtVtx(iPho->superCluster()->x() - pv.x(),
302  iPho->superCluster()->y() - pv.y(),
303  iPho->superCluster()->z() - pv.z());
304 
305  // Zero the isolation sums
306  float chargedIsoSum = 0;
307  float neutralHadronIsoSum = 0;
308  float photonIsoSum = 0;
309 
310  // Loop over all PF candidates
311  for (unsigned idxcand = 0; idxcand < pfCandidatesHandle->size(); ++idxcand ){
312 
313  // Here, the type will be a simple reco::Candidate. We cast it
314  // for full PFCandidate or PackedCandidate below as necessary
315  const auto& iCand = pfCandidatesHandle->ptrAt(idxcand);
316 
317  // One would think that we should check that this iCand from
318  // the generic PF collection is not identical to the iPho photon
319  // for which we are computing the isolations. However, it turns out
320  // to be unnecessary. Below, in the function isInFootprint(), we drop
321  // this iCand if it is in the footprint, and this always removes
322  // the iCand if it matches the iPho.
323  // The explicit check at this point is not totally trivial because
324  // of non-triviality of implementation of this check for miniAOD (PackedCandidates
325  // of the PF collection do not contain the supercluser link, so can't use that).
326  // if( isAOD ){
327  // if( ((const recoCandPtr)iCand)->superClusterRef() == iPho->superCluster() ) continue;
328  // }
329 
330 
331  // Check if this candidate is within the isolation cone
332  float dR2 = deltaR2(photon_directionWrtVtx.Eta(),photon_directionWrtVtx.Phi(),
333  iCand->eta(), iCand->phi());
334  if( dR2 > coneSizeDR*coneSizeDR ) continue;
335 
336  // Check if this candidate is not in the footprint
337  bool inFootprint = false;
338  if(isAOD) {
339  inFootprint = isInFootprint( (*particleBasedIsolationMap)[iPho], iCand );
340  } else {
341  edm::Ptr<pat::Photon> patPhotonPtr(src->ptrAt(idxpho));
342  inFootprint = isInFootprint(patPhotonPtr->associatedPackedPFCandidates(), iCand);
343  }
344 
345  if( inFootprint ) continue;
346 
347  // Find candidate type
348  reco::PFCandidate::ParticleType thisCandidateType = candidatePdgId(iCand, isAOD);
349 
350  // Increment the appropriate isolation sum
351  if( thisCandidateType == reco::PFCandidate::h ){
352  // for charged hadrons, additionally check consistency
353  // with the PV
354  const reco::Track *theTrack = getTrackPointer( iCand, isAOD );
355 
356  float dxy = theTrack->dxy(pv.position());
357  if(fabs(dxy) > dxyMax) continue;
358 
359  float dz = theTrack->dz(pv.position());
360  if (fabs(dz) > dzMax) continue;
361 
362  // The candidate is eligible, increment the isolaiton
363  chargedIsoSum += iCand->pt();
364  }
365 
366  if( thisCandidateType == reco::PFCandidate::h0 )
367  neutralHadronIsoSum += iCand->pt();
368 
369  if( thisCandidateType == reco::PFCandidate::gamma )
370  photonIsoSum += iCand->pt();
371  }
372 
373  phoChargedIsolation .push_back( chargedIsoSum );
374  phoNeutralHadronIsolation.push_back( neutralHadronIsoSum );
375  phoPhotonIsolation .push_back( photonIsoSum );
376 
377  float worstChargedIso =
378  computeWorstPFChargedIsolation(iPho, pfCandidatesHandle, vertices,
379  isAOD, coneSizeDR, dxyMax, dzMax);
380  phoWorstChargedIsolation .push_back( worstChargedIso );
381 
382 
383  }
384 
385  // Cluster shapes
386  writeValueMap(iEvent, src, phoFull5x5SigmaIEtaIEta, phoFull5x5SigmaIEtaIEta_);
387  writeValueMap(iEvent, src, phoFull5x5SigmaIEtaIPhi, phoFull5x5SigmaIEtaIPhi_);
388  writeValueMap(iEvent, src, phoFull5x5E1x3 , phoFull5x5E1x3_);
389  writeValueMap(iEvent, src, phoFull5x5E2x2 , phoFull5x5E2x2_);
390  writeValueMap(iEvent, src, phoFull5x5E2x5Max, phoFull5x5E2x5Max_);
391  writeValueMap(iEvent, src, phoFull5x5E5x5 , phoFull5x5E5x5_);
392  writeValueMap(iEvent, src, phoESEffSigmaRR , phoESEffSigmaRR_);
393  // Isolations
394  writeValueMap(iEvent, src, phoChargedIsolation, phoChargedIsolation_);
395  writeValueMap(iEvent, src, phoNeutralHadronIsolation, phoNeutralHadronIsolation_);
396  writeValueMap(iEvent, src, phoPhotonIsolation, phoPhotonIsolation_);
397  writeValueMap(iEvent, src, phoWorstChargedIsolation, phoWorstChargedIsolation_);
398 }
399 
402  const std::vector<float> & values,
403  const std::string & label) const
404 {
405  using namespace edm;
406  using namespace std;
407  auto_ptr<ValueMap<float> > valMap(new ValueMap<float>());
408  edm::ValueMap<float>::Filler filler(*valMap);
409  filler.insert(handle, values.begin(), values.end());
410  filler.fill();
411  iEvent.put(valMap, label);
412 }
413 
415  //The following says we do not know what parameters are allowed so do no validation
416  // Please change this to state exactly what you do use, even if it is no parameters
418  desc.setUnknown();
419  descriptions.addDefault(desc);
420 }
421 
422 // Charged isolation with respect to the worst vertex. See more
423 // comments above at the function declaration.
424 template <class T, class U>
428  bool isAOD,
429  float dRmax, float dxyMax, float dzMax){
430 
431  float worstIsolation = 999;
432  std::vector<float> allIsolations;
433 
434  // Constants below: there are no vetos and no min pt requirement,
435  // just like in the original H->gamma gamma code.
436  const float dRvetoBarrel = 0.0;
437  const float dRvetoEndcap = 0.0;
438  const float ptMin = 0.0;
439 
440  float dRveto;
441  if (photon->isEB())
442  dRveto = dRvetoBarrel;
443  else
444  dRveto = dRvetoEndcap;
445 
446  //Calculate isolation sum separately for each vertex
447  for(unsigned int ivtx=0; ivtx<vertices->size(); ++ivtx) {
448 
449  // Shift the photon according to the vertex
450  reco::VertexRef vtx(vertices, ivtx);
451  math::XYZVector photon_directionWrtVtx(photon->superCluster()->x() - vtx->x(),
452  photon->superCluster()->y() - vtx->y(),
453  photon->superCluster()->z() - vtx->z());
454 
455  float sum = 0;
456  // Loop over the PFCandidates
457  for(unsigned i=0; i<pfCandidates->size(); i++) {
458 
459  const auto& iCand = pfCandidates->ptrAt(i);
460 
461  //require that PFCandidate is a charged hadron
462  reco::PFCandidate::ParticleType thisCandidateType = candidatePdgId(iCand, isAOD);
463  if (thisCandidateType != reco::PFCandidate::h)
464  continue;
465 
466  if (iCand->pt() < ptMin)
467  continue;
468 
469  const reco::Track *theTrack = getTrackPointer( iCand, isAOD );
470  float dxy = theTrack->dxy(vtx->position());
471  if( fabs(dxy) > dxyMax) continue;
472 
473  float dz = theTrack->dz(vtx->position());
474  if ( fabs(dz) > dzMax) continue;
475 
476  float dR2 = deltaR2(photon_directionWrtVtx.Eta(), photon_directionWrtVtx.Phi(),
477  iCand->eta(), iCand->phi());
478  if(dR2 > dRmax*dRmax || dR2 < dRveto*dRveto) continue;
479 
480  sum += iCand->pt();
481  }
482 
483  allIsolations.push_back(sum);
484  }
485 
486  if( allIsolations.size()>0 )
487  worstIsolation = * std::max_element( allIsolations.begin(), allIsolations.end() );
488 
489  return worstIsolation;
490 }
491 
494  bool isAOD){
495 
497  if( isAOD )
498  thisCandidateType = ( (const recoCandPtr)candidate)->particleId();
499  else {
500  // the neutral hadrons and charged hadrons can be of pdgId types
501  // only 130 (K0L) and +-211 (pi+-) in packed candidates
502  const int pdgId = ( (const patCandPtr)candidate)->pdgId();
503  if( pdgId == 22 )
504  thisCandidateType = reco::PFCandidate::gamma;
505  else if( abs(pdgId) == 130) // PDG ID for K0L
506  thisCandidateType = reco::PFCandidate::h0;
507  else if( abs(pdgId) == 211) // PDG ID for pi+-
508  thisCandidateType = reco::PFCandidate::h;
509  }
510  return thisCandidateType;
511 }
512 
513 const reco::Track*
515 
516  const reco::Track* theTrack = nullptr;
517  if( isAOD )
518  theTrack = &*( ((const recoCandPtr) candidate)->trackRef());
519  else
520  theTrack = &( ((const patCandPtr) candidate)->pseudoTrack());
521 
522  return theTrack;
523 }
524 
void writeValueMap(edm::Event &iEvent, const edm::Handle< edm::View< reco::Photon > > &handle, const std::vector< float > &values, const std::string &label) const
T getParameter(std::string const &) const
edm::EDGetTokenT< EcalRecHitCollection > esReducedRecHitCollectionMiniAOD_
const reco::Track * getTrackPointer(const edm::Ptr< reco::Candidate > candidate, bool isAOD)
int i
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< reco::VertexCollection > vtxToken_
ParticleType
particle types
Definition: PFCandidate.h:44
edm::EDGetTokenT< reco::VertexCollection > vtxTokenMiniAOD_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
std::unique_ptr< noZS::EcalClusterLazyTools > lazyToolnoZS
edm::EDGetTokenT< EcalRecHitCollection > eeReducedRecHitCollectionMiniAOD_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
reco::PFCandidate::ParticleType candidatePdgId(const edm::Ptr< reco::Candidate > candidate, bool isAOD)
edm::EDGetTokenT< EcalRecHitCollection > ebReducedRecHitCollection_
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
PhotonIDValueMapProducer(const edm::ParameterSet &)
edm::EDGetTokenT< edm::ValueMap< std::vector< reco::PFCandidateRef > > > particleBasedIsolationToken_
edm::EDGetTokenT< edm::ValueMap< std::vector< reco::PFCandidateRef > > > particleBasedIsolationTokenMiniAOD_
EcalClusterLazyToolsT< noZS::EcalClusterTools > EcalClusterLazyTools
edm::EDGetTokenT< EcalRecHitCollection > ebReducedRecHitCollectionMiniAOD_
#define constexpr
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< edm::View< reco::Candidate > > pfCandidatesTokenMiniAOD_
void addDefault(ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< EcalRecHitCollection > eeReducedRecHitCollection_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
bool isnan(float x)
Definition: math.h:13
T sqrt(T t)
Definition: SSEVec.h:48
bool isInFootprint(const T &thefootprint, const U &theCandidate)
edm::Ptr< pat::PackedCandidate > patCandPtr
tuple handle
Definition: patZpeak.py:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< edm::View< reco::Candidate > > pfCandidatesToken_
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
virtual void produce(edm::Event &, const edm::EventSetup &) override
float computeWorstPFChargedIsolation(const T &photon, const U &pfCandidates, const edm::Handle< reco::VertexCollection > vertices, bool isAOD, float dRmax, float dxyMax, float dzMax)
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:562
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
edm::EDGetTokenT< EcalRecHitCollection > esReducedRecHitCollection_
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:544
long double T
edm::Ptr< reco::PFCandidate > recoCandPtr