CMS 3D CMS Logo

PuppiPhoton.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 
4 // user include files
26 //Main File
28 
29 // ------------------------------------------------------------------------------------------
31  tokenPFCandidates_ = consumes<CandidateView>(iConfig.getParameter<edm::InputTag>("candName"));
32  tokenPuppiCandidates_ = consumes<CandidateView>(iConfig.getParameter<edm::InputTag>("puppiCandName"));
33  tokenPhotonCandidates_ = consumes<CandidateView>(iConfig.getParameter<edm::InputTag>("photonName"));
34  usePhotonId_ = (iConfig.getParameter<edm::InputTag>("photonId")).label().size() != 0;
35  if(usePhotonId_)
36  tokenPhotonId_ = consumes<edm::ValueMap<bool> >(iConfig.getParameter<edm::InputTag>("photonId"));
37  runOnMiniAOD_ = iConfig.getParameter<bool>("runOnMiniAOD");
38  if(!runOnMiniAOD_)
39  reco2pf_ = consumes<edm::ValueMap<std::vector<reco::PFCandidateRef> > >(iConfig.getParameter<edm::InputTag>("recoToPFMap"));
40  useValueMap_ = iConfig.getParameter<bool>("useValueMap");
41  if(useValueMap_)
42  tokenWeights_ = consumes<edm::ValueMap<float> >(iConfig.getParameter<edm::InputTag>("weightsName"));
43 
44  pt_ = iConfig.getParameter<double>("pt");
45  eta_ = iConfig.getParameter<double>("eta");
46  dRMatch_ = iConfig.getParameter<std::vector<double> > ("dRMatch");
47  pdgIds_ = iConfig.getParameter<std::vector<int32_t> >("pdgids");
48  usePFRef_ = iConfig.getParameter<bool>("useRefs");
49  weight_ = iConfig.getParameter<double>("weight");
50  produces<PFOutputCollection>();
51  produces< edm::ValueMap<reco::CandidatePtr> >();
52 }
53 // ------------------------------------------------------------------------------------------
55 // ------------------------------------------------------------------------------------------
57 
58  edm::Handle<CandidateView> hPhoProduct;
59  iEvent.getByToken(tokenPhotonCandidates_,hPhoProduct);
60  const CandidateView *phoCol = hPhoProduct.product();
61 
63  if(usePhotonId_) iEvent.getByToken(tokenPhotonId_,photonId);
64  int iC = -1;
65  std::vector<const reco::Candidate*> phoCands;
66  std::vector<uint16_t> phoIndx;
67 
69  if(!runOnMiniAOD_) iEvent.getByToken(reco2pf_, reco2pf);
70 
71  // Get PFCandidate Collection
72  edm::Handle<CandidateView> hPFProduct;
73  iEvent.getByToken(tokenPFCandidates_,hPFProduct);
74  const CandidateView *pfCol = hPFProduct.product();
75 
76  edm::Handle<CandidateView> hPuppiProduct;
77  iEvent.getByToken(tokenPuppiCandidates_,hPuppiProduct);
78  const CandidateView *pupCol = hPuppiProduct.product();
79  for(CandidateView::const_iterator itPho = phoCol->begin(); itPho!=phoCol->end(); itPho++) {
80  iC++;
81  bool passObject = false;
82  if(itPho->isPhoton() && usePhotonId_) passObject = (*photonId) [phoCol->ptrAt(iC)];
83  if(itPho->pt() < pt_) continue;
84  if(!passObject && usePhotonId_) continue;
85  if(runOnMiniAOD_) {
86  const pat::Photon *pPho = dynamic_cast<const pat::Photon*>(&(*itPho));
87  if(pPho != 0) {
89  if(fabs(ref->eta()) < eta_ ) {
90  phoIndx.push_back(ref.key());
91  phoCands.push_back(&(*(pfCol->ptrAt(ref.key()))));
92  }
93  }
94  continue;
95  }
96  const pat::Electron *pElectron = dynamic_cast<const pat::Electron*>(&(*itPho));
97  if(pElectron != 0) {
99  if(fabs(ref->eta()) < eta_ ) {
100  phoIndx.push_back(ref.key());
101  phoCands.push_back(&(*(pfCol->ptrAt(ref.key()))));
102  }
103  }
104  } else {
105  for( const edm::Ref<std::vector<reco::PFCandidate> > & ref : (*reco2pf)[phoCol->ptrAt(iC)] ) {
106  if(fabs(ref->eta()) < eta_ ) {
107  phoIndx.push_back(ref.key());
108  phoCands.push_back(&(*(pfCol->ptrAt(ref.key()))));
109  }
110  }
111  }
112  }
113  //Get Weights
114  edm::Handle<edm::ValueMap<float> > pupWeights;
115  if(useValueMap_)
116  iEvent.getByToken(tokenWeights_,pupWeights);
117  std::unique_ptr<edm::ValueMap<LorentzVector> > p4PupOut(new edm::ValueMap<LorentzVector>());
118  LorentzVectorCollection puppiP4s;
119  std::vector<reco::CandidatePtr> values(hPFProduct->size());
120  int iPF = 0;
121  std::vector<float> lWeights;
122  static const reco::PFCandidate dummySinceTranslateIsNotStatic;
123  corrCandidates_.reset( new PFOutputCollection );
124  std::set<int> foundPhoIndex;
125  for(CandidateView::const_iterator itPF = pupCol->begin(); itPF!=pupCol->end(); itPF++) {
126  auto id = dummySinceTranslateIsNotStatic.translatePdgIdToType(itPF->pdgId());
127  const reco::PFCandidate *pPF = dynamic_cast<const reco::PFCandidate*>(&(*itPF));
128  reco::PFCandidate pCand( pPF ? *pPF : reco::PFCandidate(itPF->charge(), itPF->p4(), id) );
129  LorentzVector pVec = itPF->p4();
130  float pWeight = 1.;
131  if(useValueMap_) pWeight = (*pupWeights)[pupCol->ptrAt(iPF)];
132  if(!usePFRef_) {
133  int iPho = -1;
134  for(std::vector<const reco::Candidate*>::iterator itPho = phoCands.begin(); itPho!=phoCands.end(); itPho++) {
135  iPho++;
136  if((!matchPFCandidate(&(*itPF),*itPho))||(foundPhoIndex.count(iPho)!=0)) continue;
137  pWeight = weight_;
138  if(!useValueMap_ && itPF->pt() != 0) pWeight = pWeight*(phoCands[iPho]->pt()/itPF->pt());
139  if(!useValueMap_ && itPF->pt() == 0) pVec.SetPxPyPzE(phoCands[iPho]->px()*pWeight,phoCands[iPho]->py()*pWeight,phoCands[iPho]->pz()*pWeight,phoCands[iPho]->energy()*pWeight);
140  foundPhoIndex.insert(iPho);
141  }
142  } else {
143  int iPho = -1;
144  for(std::vector<uint16_t>::const_iterator itPho = phoIndx.begin(); itPho!=phoIndx.end(); itPho++) {
145  iPho++;
146  if(pupCol->refAt(iPF).key() != *itPho) continue;
147  pWeight = weight_;
148  if(!useValueMap_ && itPF->pt() != 0) pWeight = pWeight*(phoCands[iPho]->pt()/itPF->pt());
149  if(!useValueMap_ && itPF->pt() == 0) pVec.SetPxPyPzE(phoCands[iPho]->px()*pWeight,phoCands[iPho]->py()*pWeight,phoCands[iPho]->pz()*pWeight,phoCands[iPho]->energy()*pWeight);
150  foundPhoIndex.insert(iPho);
151  }
152  }
153  if(itPF->pt() != 0) pVec.SetPxPyPzE(itPF->px()*pWeight,itPF->py()*pWeight,itPF->pz()*pWeight,itPF->energy()*pWeight);
154 
155  lWeights.push_back(pWeight);
156  pCand.setP4(pVec);
157  puppiP4s.push_back( pVec );
158  pCand.setSourceCandidatePtr( itPF->sourceCandidatePtr(0) );
159  corrCandidates_->push_back(pCand);
160  iPF++;
161  }
162  //Add the missing pfcandidates
163  for(unsigned int iPho = 0; iPho < phoCands.size(); iPho++) {
164  if(foundPhoIndex.count(iPho)!=0) continue;
165  auto id = dummySinceTranslateIsNotStatic.translatePdgIdToType(phoCands[iPho]->pdgId());
166  reco::PFCandidate pCand(reco::PFCandidate(phoCands[iPho]->charge(), phoCands[iPho]->p4(),id) );
167  pCand.setSourceCandidatePtr( phoCands[iPho]->sourceCandidatePtr(0) );
168  LorentzVector pVec = phoCands[iPho]->p4();
169  pVec.SetPxPyPzE(phoCands[iPho]->px()*weight_,phoCands[iPho]->py()*weight_,phoCands[iPho]->pz()*weight_,phoCands[iPho]->energy()*weight_);
170  pCand.setP4(pVec);
171  lWeights.push_back(weight_);
172  puppiP4s.push_back( pVec );
173  corrCandidates_->push_back(pCand);
174  }
175  //Fill it into the event
177  for(unsigned int ic=0, nc = pupCol->size(); ic < nc; ++ic) {
178  reco::CandidatePtr pkref( oh, ic );
179  values[ic] = pkref;
180  }
181  std::unique_ptr<edm::ValueMap<reco::CandidatePtr> > pfMap_p(new edm::ValueMap<reco::CandidatePtr>());
183  filler.insert(hPFProduct, values.begin(), values.end());
184  filler.fill();
185  iEvent.put(std::move(pfMap_p));
186 }
187 // ------------------------------------------------------------------------------------------
189  double lDR = deltaR(iPF->eta(),iPF->phi(),iPho->eta(),iPho->phi());
190  for(unsigned int i0 = 0; i0 < pdgIds_.size(); i0++) {
191  if(std::abs(iPF->pdgId()) == pdgIds_[i0] && lDR < dRMatch_[i0]) return true;
192  }
193  return false;
194 }
195 // ------------------------------------------------------------------------------------------
197  //The following says we do not know what parameters are allowed so do no validation
198  // Please change this to state exactly what you do use, even if it is no parameters
200  desc.setUnknown();
201  descriptions.addDefault(desc);
202 }
203 //define this as a plug-in
std::vector< int32_t > pdgIds_
Definition: PuppiPhoton.h:46
double pt_
Definition: PuppiPhoton.h:40
T getParameter(std::string const &) const
Analysis-level Photon class.
Definition: Photon.h:47
edm::EDGetTokenT< CandidateView > tokenPhotonCandidates_
Definition: PuppiPhoton.h:36
PuppiPhoton(const edm::ParameterSet &)
Definition: PuppiPhoton.cc:30
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
bool runOnMiniAOD_
Definition: PuppiPhoton.h:43
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
Ptr< value_type > ptrAt(size_type i) const
std::vector< reco::PFCandidate > PFOutputCollection
Definition: PuppiPhoton.h:28
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
edm::RefVector< pat::PackedCandidateCollection > associatedPackedPFCandidates() const
References to PFCandidates linked to this object (e.g. for isolation vetos or masking before jet recl...
size_type size() const
void setSourceCandidatePtr(const PFCandidatePtr &ptr)
Definition: PFCandidate.h:125
edm::EDGetTokenT< CandidateView > tokenPuppiCandidates_
Definition: PuppiPhoton.h:35
bool usePFRef_
Definition: PuppiPhoton.h:42
std::unique_ptr< PFOutputCollection > corrCandidates_
Definition: PuppiPhoton.h:47
edm::EDGetTokenT< edm::ValueMap< std::vector< reco::PFCandidateRef > > > reco2pf_
Definition: PuppiPhoton.h:37
RefToBase< value_type > refAt(size_type i) const
int iEvent
Definition: GenABIO.cc:230
void addDefault(ParameterSetDescription const &psetDescription)
const_iterator begin() const
virtual int pdgId() const =0
PDG identifier.
double p4[4]
Definition: TauolaWrapper.h:92
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual void produce(edm::Event &, const edm::EventSetup &)
Definition: PuppiPhoton.cc:56
bool matchPFCandidate(const reco::Candidate *iPF, const reco::Candidate *iPho)
Definition: PuppiPhoton.cc:188
bool usePhotonId_
Definition: PuppiPhoton.h:44
double weight_
Definition: PuppiPhoton.h:48
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
double eta_
Definition: PuppiPhoton.h:41
virtual double eta() const =0
momentum pseudorapidity
T const * product() const
Definition: Handle.h:81
edm::EDGetTokenT< edm::ValueMap< float > > tokenWeights_
Definition: PuppiPhoton.h:38
Analysis-level electron class.
Definition: Electron.h:52
virtual void setP4(const LorentzVector &p4) final
set 4-momentum
ParticleType translatePdgIdToType(int pdgid) const
Definition: PFCandidate.cc:229
std::vector< double > dRMatch_
Definition: PuppiPhoton.h:45
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
edm::EDGetTokenT< edm::ValueMap< bool > > tokenPhotonId_
Definition: PuppiPhoton.h:39
const_iterator end() const
math::XYZTLorentzVector LorentzVector
Definition: PuppiPhoton.h:25
edm::EDGetTokenT< CandidateView > tokenPFCandidates_
Definition: PuppiPhoton.h:34
virtual double phi() const =0
momentum azimuthal angle
edm::RefVector< pat::PackedCandidateCollection > associatedPackedPFCandidates() const
References to PFCandidates linked to this object (e.g. for isolation vetos or masking before jet recl...
std::vector< LorentzVector > LorentzVectorCollection
Definition: PuppiPhoton.h:26
def move(src, dest)
Definition: eostools.py:510
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: PuppiPhoton.cc:196
bool useValueMap_
Definition: PuppiPhoton.h:49