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  eta_ = iConfig.getParameter<double>("eta");
37  dRMatch_ = iConfig.getParameter<std::vector<double> > ("dRMatch");
38  pdgIds_ = iConfig.getParameter<std::vector<int32_t> >("pdgids");
39  usePFRef_ = iConfig.getParameter<bool>("useRefs");
40  weight_ = iConfig.getParameter<double>("weight");
41  useValueMap_ = iConfig.getParameter<bool>("useValueMap");
42  tokenWeights_ = consumes<edm::ValueMap<float> >(iConfig.getParameter<edm::InputTag>("weightsName"));
43 
44  usePhotonId_ = (iConfig.getParameter<edm::InputTag>("photonId")).label().size() == 0;
45  produces<PFOutputCollection>();
46  produces< edm::ValueMap<reco::CandidatePtr> >();
47 }
48 // ------------------------------------------------------------------------------------------
50 // ------------------------------------------------------------------------------------------
52 
53  edm::Handle<CandidateView> hPhoProduct;
54  iEvent.getByToken(tokenPhotonCandidates_,hPhoProduct);
55  const CandidateView *phoCol = hPhoProduct.product();
56 
58  if(usePhotonId_) iEvent.getByToken(tokenPhotonId_,photonId);
59  int iC = -1;
60  std::vector<const reco::Candidate*> phoCands;
61  std::vector<uint16_t> phoIndx;
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  for(CandidateView::const_iterator itPho = phoCol->begin(); itPho!=phoCol->end(); itPho++) {
72  iC++;
73  bool passObject = false;
74  if(itPho->isPhoton() && usePhotonId_) passObject = (*photonId) [phoCol->ptrAt(iC)];
75  if(itPho->pt() < pt_) continue;
76  if(!passObject && usePhotonId_) continue;
77  //if(!usePFRef_ && fabs(itPho->eta()) < eta_) phoCands.push_back(&(*itPho)); ===> should add a flag useAODRef_ in place of usePFRef_
78  //if(!usePFRef_) continue;
79  const pat::Photon *pPho = dynamic_cast<const pat::Photon*>(&(*itPho));
80  if(pPho != 0) {
82  if(fabs(pfCol->ptrAt(ref.key())->eta()) < eta_ ) {
83  phoIndx.push_back(ref.key());
84  phoCands.push_back(&(*(pfCol->ptrAt(ref.key()))));
85  }
86  }
87  continue;
88  }
89  const pat::Electron *pElectron = dynamic_cast<const pat::Electron*>(&(*itPho));
90  if(pElectron != 0) {
92  if(fabs(pfCol->ptrAt(ref.key())->eta()) < eta_ ) {
93  phoIndx.push_back(ref.key());
94  phoCands.push_back(&(*(pfCol->ptrAt(ref.key()))));
95  }
96  }
97  }
98  //Get Weights
100  iEvent.getByToken(tokenWeights_,pupWeights);
101  std::auto_ptr<edm::ValueMap<LorentzVector> > p4PupOut(new edm::ValueMap<LorentzVector>());
102  LorentzVectorCollection puppiP4s;
103  std::vector<reco::CandidatePtr> values(hPFProduct->size());
104  int iPF = 0;
105  std::vector<float> lWeights;
106  static const reco::PFCandidate dummySinceTranslateIsNotStatic;
107  corrCandidates_.reset( new PFOutputCollection );
108  std::vector<int> foundPhoIndex;
109  for(CandidateView::const_iterator itPF = pupCol->begin(); itPF!=pupCol->end(); itPF++) {
110  auto id = dummySinceTranslateIsNotStatic.translatePdgIdToType(itPF->pdgId());
111  const reco::PFCandidate *pPF = dynamic_cast<const reco::PFCandidate*>(&(*itPF));
112  reco::PFCandidate pCand( pPF ? *pPF : reco::PFCandidate(itPF->charge(), itPF->p4(), id) );
113  LorentzVector pVec = itPF->p4();
114  float pWeight = 1.;
115  if(useValueMap_) pWeight = (*pupWeights)[pupCol->ptrAt(iPF)];
116  if(!usePFRef_) {
117  int iPho = -1;
118  for(std::vector<const reco::Candidate*>::iterator itPho = phoCands.begin(); itPho!=phoCands.end(); itPho++) {
119  iPho++;
120  if(!matchPFCandidate(&(*itPF),*itPho)) continue;
121  pWeight = weight_;
122  if(!useValueMap_ && itPF->pt() != 0) pWeight = pWeight*(phoCands[iPho]->pt()/itPF->pt());
123  if(!useValueMap_ && itPF->pt() == 0) pVec.SetPxPyPzE(phoCands[iPho]->px()*pWeight,phoCands[iPho]->py()*pWeight,phoCands[iPho]->pz()*pWeight,phoCands[iPho]->energy()*pWeight);
124  foundPhoIndex.push_back(iPho);
125  }
126  } else {
127  int iPho = -1;
128  for(std::vector<uint16_t>::const_iterator itPho = phoIndx.begin(); itPho!=phoIndx.end(); itPho++) {
129  iPho++;
130  if(pupCol->refAt(iPF).key() != *itPho) continue;
131  pWeight = weight_;
132  if(!useValueMap_ && itPF->pt() != 0) pWeight = pWeight*(phoCands[iPho]->pt()/itPF->pt());
133  if(!useValueMap_ && itPF->pt() == 0) pVec.SetPxPyPzE(phoCands[iPho]->px()*pWeight,phoCands[iPho]->py()*pWeight,phoCands[iPho]->pz()*pWeight,phoCands[iPho]->energy()*pWeight);
134  foundPhoIndex.push_back(iPho);
135  }
136  }
137  if(itPF->pt() != 0) pVec.SetPxPyPzE(itPF->px()*pWeight,itPF->py()*pWeight,itPF->pz()*pWeight,itPF->energy()*pWeight);
138 
139  lWeights.push_back(pWeight);
140  pCand.setP4(pVec);
141  puppiP4s.push_back( pVec );
142  pCand.setSourceCandidatePtr( itPF->sourceCandidatePtr(0) );
143  corrCandidates_->push_back(pCand);
144  iPF++;
145  }
146  //Add the missing pfcandidates
147  for(unsigned int iPho = 0; iPho < phoCands.size(); iPho++) {
148  bool pFound = false;
149  for(unsigned int jPho = 0; jPho < foundPhoIndex.size(); jPho++) {
150  if(foundPhoIndex[jPho] == int(iPho)) {
151  pFound = true;
152  break;
153  }
154  }
155  if(pFound) continue;
156  auto id = dummySinceTranslateIsNotStatic.translatePdgIdToType(phoCands[iPho]->pdgId());
157  reco::PFCandidate pCand(reco::PFCandidate(phoCands[iPho]->charge(), phoCands[iPho]->p4(),id) );
158  pCand.setSourceCandidatePtr( phoCands[iPho]->sourceCandidatePtr(0) );
159  LorentzVector pVec = phoCands[iPho]->p4();
160  pVec.SetPxPyPzE(phoCands[iPho]->px()*weight_,phoCands[iPho]->py()*weight_,phoCands[iPho]->pz()*weight_,phoCands[iPho]->energy()*weight_);
161  pCand.setP4(pVec);
162  lWeights.push_back(weight_);
163  puppiP4s.push_back( pVec );
164  corrCandidates_->push_back(pCand);
165  }
166  //Fill it into the event
168  for(unsigned int ic=0, nc = pupCol->size(); ic < nc; ++ic) {
169  reco::CandidatePtr pkref( oh, ic );
170  values[ic] = pkref;
171  }
172  std::auto_ptr<edm::ValueMap<reco::CandidatePtr> > pfMap_p(new edm::ValueMap<reco::CandidatePtr>());
174  filler.insert(hPFProduct, values.begin(), values.end());
175  filler.fill();
176  iEvent.put(pfMap_p);
177 }
178 // ------------------------------------------------------------------------------------------
180  double lDR = deltaR(iPF->eta(),iPF->phi(),iPho->eta(),iPho->phi());
181  for(unsigned int i0 = 0; i0 < pdgIds_.size(); i0++) {
182  if(std::abs(iPF->pdgId()) == pdgIds_[i0] && lDR < dRMatch_[i0]) return true;
183  }
184  return false;
185 }
186 // ------------------------------------------------------------------------------------------
188  //The following says we do not know what parameters are allowed so do no validation
189  // Please change this to state exactly what you do use, even if it is no parameters
191  desc.setUnknown();
192  descriptions.addDefault(desc);
193 }
194 //define this as a plug-in
std::vector< int32_t > pdgIds_
Definition: PuppiPhoton.h:44
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:45
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:41
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:51
bool matchPFCandidate(const reco::Candidate *iPF, const reco::Candidate *iPho)
Definition: PuppiPhoton.cc:179
bool usePhotonId_
Definition: PuppiPhoton.h:42
double weight_
Definition: PuppiPhoton.h:46
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
virtual int pdgId() const =0
PDG identifier.
double eta_
Definition: PuppiPhoton.h:40
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:43
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:187
bool useValueMap_
Definition: PuppiPhoton.h:47
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity