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