CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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  tokenPhotonId_ = consumes<edm::ValueMap<bool> >(iConfig.getParameter<edm::InputTag>("photonId"));
35  pt_ = iConfig.getParameter<double>("pt");
36  dRMatch_ = iConfig.getParameter<std::vector<double> > ("dRMatch");
37  pdgIds_ = iConfig.getParameter<std::vector<int32_t> >("pdgids");
38  usePFRef_ = iConfig.getParameter<bool>("useRefs");
39  weight_ = iConfig.getParameter<double>("weight");
40  useValueMap_ = iConfig.getParameter<bool>("useValueMap");
41  tokenWeights_ = consumes<edm::ValueMap<float> >(iConfig.getParameter<edm::InputTag>("weightsName"));
42 
43  usePhotonId_ = (iConfig.getParameter<edm::InputTag>("photonId")).label().size() == 0;
44  produces<PFOutputCollection>();
45  produces< edm::ValueMap<reco::CandidatePtr> >();
46 }
47 // ------------------------------------------------------------------------------------------
49 // ------------------------------------------------------------------------------------------
51 
52  edm::Handle<CandidateView> hPhoProduct;
53  iEvent.getByToken(tokenPhotonCandidates_,hPhoProduct);
54  const CandidateView *phoCol = hPhoProduct.product();
55 
57  if(usePhotonId_) iEvent.getByToken(tokenPhotonId_,photonId);
58  int iC = -1;
59  std::vector<const reco::Candidate*> phoCands;
60  std::vector<uint16_t> phoIndx;
61 
62  // Get PFCandidate Collection
63  edm::Handle<CandidateView> hPFProduct;
64  iEvent.getByToken(tokenPFCandidates_,hPFProduct);
65  const CandidateView *pfCol = hPFProduct.product();
66 
67  edm::Handle<CandidateView> hPuppiProduct;
68  iEvent.getByToken(tokenPuppiCandidates_,hPuppiProduct);
69  const CandidateView *pupCol = hPuppiProduct.product();
70  for(CandidateView::const_iterator itPho = phoCol->begin(); itPho!=phoCol->end(); itPho++) {
71  iC++;
72  bool passObject = false;
73  if(itPho->isPhoton() && usePhotonId_) passObject = (*photonId) [phoCol->ptrAt(iC)];
74  if(itPho->pt() < pt_) continue;
75  if(!passObject && usePhotonId_) continue;
76  if(!usePFRef_) phoCands.push_back(&(*itPho));
77  if(!usePFRef_) continue;
78  const pat::Photon *pPho = dynamic_cast<const pat::Photon*>(&(*itPho));
79  if(pPho != 0) {
81  if(matchPFCandidate(&(*(pfCol->ptrAt(ref.key()))),&(*itPho))) {
82  phoIndx.push_back(ref.key());
83  phoCands.push_back(&(*(pfCol->ptrAt(ref.key()))));
84  }
85  }
86  continue;
87  }
88  const pat::Electron *pElectron = dynamic_cast<const pat::Electron*>(&(*itPho));
89  if(pElectron != 0) {
91  if(matchPFCandidate(&(*(pfCol->ptrAt(ref.key()))),&(*itPho))) {
92  phoIndx.push_back(ref.key());
93  phoCands.push_back(&(*(pfCol->ptrAt(ref.key()))));
94  }
95  }
96  }
97  //Get Weights
99  iEvent.getByToken(tokenWeights_,pupWeights);
100  std::auto_ptr<edm::ValueMap<LorentzVector> > p4PupOut(new edm::ValueMap<LorentzVector>());
101  LorentzVectorCollection puppiP4s;
102  std::vector<reco::CandidatePtr> values(hPFProduct->size());
103  int iPF = 0;
104  std::vector<float> lWeights;
105  static const reco::PFCandidate dummySinceTranslateIsNotStatic;
106  corrCandidates_.reset( new PFOutputCollection );
107  std::vector<int> foundPhoIndex;
108  for(CandidateView::const_iterator itPF = pupCol->begin(); itPF!=pupCol->end(); itPF++) {
109  auto id = dummySinceTranslateIsNotStatic.translatePdgIdToType(itPF->pdgId());
110  const reco::PFCandidate *pPF = dynamic_cast<const reco::PFCandidate*>(&(*itPF));
111  reco::PFCandidate pCand( pPF ? *pPF : reco::PFCandidate(itPF->charge(), itPF->p4(), id) );
112  LorentzVector pVec = itPF->p4();
113  float pWeight = 1.;
114  if(useValueMap_) pWeight = (*pupWeights)[pupCol->ptrAt(iPF)];
115  if(!usePFRef_) {
116  for(std::vector<const reco::Candidate*>::iterator itPho = phoCands.begin(); itPho!=phoCands.end(); itPho++) {
117  if(matchPFCandidate(&(*itPF),*itPho)) pWeight = weight_;
118  }
119  } else {
120  int iPho = -1;
121  for(std::vector<uint16_t>::const_iterator itPho = phoIndx.begin(); itPho!=phoIndx.end(); itPho++) {
122  iPho++;
123  if(pupCol->refAt(iPF).key() != *itPho) continue;
124  pWeight = weight_;
125  if(!useValueMap_) {
126  double pCorr = phoCands[iPho]->pt()/itPF->pt();
127  pWeight = pWeight*pCorr;
128  }
129  foundPhoIndex.push_back(iPho);
130  }
131  }
132  pVec.SetPxPyPzE(itPF->px()*pWeight,itPF->py()*pWeight,itPF->pz()*pWeight,itPF->energy()*pWeight);
133  lWeights.push_back(pWeight);
134  pCand.setP4(pVec);
135  puppiP4s.push_back( pVec );
136  pCand.setSourceCandidatePtr( itPF->sourceCandidatePtr(0) );
137  corrCandidates_->push_back(pCand);
138  iPF++;
139  }
140  //Add the missing pfcandidates
141  for(unsigned int iPho = 0; iPho < phoCands.size(); iPho++) {
142  bool pFound = false;
143  for(unsigned int jPho = 0; jPho < foundPhoIndex.size(); jPho++) {
144  if(foundPhoIndex[jPho] == int(iPho)) {
145  pFound = true;
146  break;
147  }
148  }
149  if(pFound) continue;
150  auto id = dummySinceTranslateIsNotStatic.translatePdgIdToType(phoCands[iPho]->pdgId());
151  reco::PFCandidate pCand(reco::PFCandidate(phoCands[iPho]->charge(), phoCands[iPho]->p4(),id) );
152  pCand.setSourceCandidatePtr( phoCands[iPho]->sourceCandidatePtr(0) );
153  LorentzVector pVec = phoCands[iPho]->p4();
154  pVec.SetPxPyPzE(phoCands[iPho]->px()*weight_,phoCands[iPho]->py()*weight_,phoCands[iPho]->pz()*weight_,phoCands[iPho]->energy()*weight_);
155  pCand.setP4(pVec);
156  lWeights.push_back(weight_);
157  puppiP4s.push_back( pVec );
158  corrCandidates_->push_back(pCand);
159  }
160  //Fill it into the event
162  for(unsigned int ic=0, nc = pupCol->size(); ic < nc; ++ic) {
163  reco::CandidatePtr pkref( oh, ic );
164  values[ic] = pkref;
165  }
166  std::auto_ptr<edm::ValueMap<reco::CandidatePtr> > pfMap_p(new edm::ValueMap<reco::CandidatePtr>());
168  filler.insert(hPFProduct, values.begin(), values.end());
169  filler.fill();
170  iEvent.put(pfMap_p);
171 }
172 // ------------------------------------------------------------------------------------------
174  double lDR = deltaR(iPF->eta(),iPF->phi(),iPho->eta(),iPho->phi());
175  for(unsigned int i0 = 0; i0 < pdgIds_.size(); i0++) {
176  if(std::abs(iPF->pdgId()) == pdgIds_[i0] && lDR < dRMatch_[i0]) return true;
177  }
178  return false;
179 }
180 // ------------------------------------------------------------------------------------------
182  //The following says we do not know what parameters are allowed so do no validation
183  // Please change this to state exactly what you do use, even if it is no parameters
185  desc.setUnknown();
186  descriptions.addDefault(desc);
187 }
188 //define this as a plug-in
std::vector< int32_t > pdgIds_
Definition: PuppiPhoton.h:43
double pt_
Definition: PuppiPhoton.h:39
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
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
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
std::auto_ptr< PFOutputCollection > corrCandidates_
Definition: PuppiPhoton.h:44
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:40
RefToBase< value_type > refAt(size_type i) const
int iEvent
Definition: GenABIO.cc:230
void addDefault(ParameterSetDescription const &psetDescription)
const_iterator begin() const
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
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:50
bool matchPFCandidate(const reco::Candidate *iPF, const reco::Candidate *iPho)
Definition: PuppiPhoton.cc:173
bool usePhotonId_
Definition: PuppiPhoton.h:41
double weight_
Definition: PuppiPhoton.h:45
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
virtual int pdgId() const =0
PDG identifier.
T const * product() const
Definition: Handle.h:81
edm::EDGetTokenT< edm::ValueMap< float > > tokenWeights_
Definition: PuppiPhoton.h:37
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:224
std::vector< double > dRMatch_
Definition: PuppiPhoton.h:42
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
edm::EDGetTokenT< edm::ValueMap< bool > > tokenPhotonId_
Definition: PuppiPhoton.h:38
const_iterator end() const
math::XYZTLorentzVector LorentzVector
Definition: PuppiPhoton.h:25
edm::EDGetTokenT< CandidateView > tokenPFCandidates_
Definition: PuppiPhoton.h:34
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
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: PuppiPhoton.cc:181
bool useValueMap_
Definition: PuppiPhoton.h:46
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity