CMS 3D CMS Logo

PFRecoTauChargedHadronFromPFCandidatePlugin.cc
Go to the documentation of this file.
1 /*
2  * PFRecoTauChargedHadronFromPFCandidatePlugin
3  *
4  * Build PFRecoTauChargedHadron objects
5  * using charged PFCandidates and/or PFNeutralHadrons as input
6  *
7  * Author: Christian Veelken, LLR
8  *
9  */
10 
12 
27 
32 
33 #include <memory>
34 
35 namespace reco
36 {
37 
38 namespace tau
39 {
40 
42 {
43  public:
46  // Return type is auto_ptr<ChargedHadronVector>
47  return_type operator()(const reco::PFJet&) const override;
48  // Hook to update PV information
49  void beginEvent() override;
50 
51  private:
52  typedef std::vector<reco::PFCandidatePtr> PFCandPtrs;
53 
55 
57 
58  std::vector<int> inputPdgIds_; // type of candidates to clusterize
59 
75 
77 };
78 
81  vertexAssociator_(pset.getParameter<edm::ParameterSet>("qualityCuts"),std::move(iC)),
83 {
84  edm::ParameterSet qcuts_pset = pset.getParameterSet("qualityCuts").getParameterSet("signalQualityCuts");
85  qcuts_ = new RecoTauQualityCuts(qcuts_pset);
86 
87  inputPdgIds_ = pset.getParameter<std::vector<int> >("chargedHadronCandidatesParticleIds");
88 
89  dRmergeNeutralHadronWrtChargedHadron_ = pset.getParameter<double>("dRmergeNeutralHadronWrtChargedHadron");
90  dRmergeNeutralHadronWrtNeutralHadron_ = pset.getParameter<double>("dRmergeNeutralHadronWrtNeutralHadron");
91  dRmergeNeutralHadronWrtElectron_ = pset.getParameter<double>("dRmergeNeutralHadronWrtElectron");
92  dRmergeNeutralHadronWrtOther_ = pset.getParameter<double>("dRmergeNeutralHadronWrtOther");
93  minBlockElementMatchesNeutralHadron_ = pset.getParameter<int>("minBlockElementMatchesNeutralHadron");
94  maxUnmatchedBlockElementsNeutralHadron_ = pset.getParameter<int>("maxUnmatchedBlockElementsNeutralHadron");
95  dRmergePhotonWrtChargedHadron_ = pset.getParameter<double>("dRmergePhotonWrtChargedHadron");
96  dRmergePhotonWrtNeutralHadron_ = pset.getParameter<double>("dRmergePhotonWrtNeutralHadron");
97  dRmergePhotonWrtElectron_ = pset.getParameter<double>("dRmergePhotonWrtElectron");
98  dRmergePhotonWrtOther_ = pset.getParameter<double>("dRmergePhotonWrtOther");
99  minBlockElementMatchesPhoton_ = pset.getParameter<int>("minBlockElementMatchesPhoton");
100  maxUnmatchedBlockElementsPhoton_ = pset.getParameter<int>("maxUnmatchedBlockElementsPhoton");
101  minMergeNeutralHadronEt_ = pset.getParameter<double>("minMergeNeutralHadronEt");
102  minMergeGammaEt_ = pset.getParameter<double>("minMergeGammaEt");
103  minMergeChargedHadronPt_ = pset.getParameter<double>("minMergeChargedHadronPt");
104 
105  verbosity_ = ( pset.exists("verbosity") ) ?
106  pset.getParameter<int>("verbosity") : 0;
107 }
108 
110 {
111  delete qcuts_;
112 }
113 
114 // Update the primary vertex
116 {
118 }
119 
120 namespace
121 {
122  std::string getPFCandidateType(reco::PFCandidate::ParticleType pfCandidateType)
123  {
124  if ( pfCandidateType == reco::PFCandidate::X ) return "undefined";
125  else if ( pfCandidateType == reco::PFCandidate::h ) return "PFChargedHadron";
126  else if ( pfCandidateType == reco::PFCandidate::e ) return "PFElectron";
127  else if ( pfCandidateType == reco::PFCandidate::mu ) return "PFMuon";
128  else if ( pfCandidateType == reco::PFCandidate::gamma ) return "PFGamma";
129  else if ( pfCandidateType == reco::PFCandidate::h0 ) return "PFNeutralHadron";
130  else if ( pfCandidateType == reco::PFCandidate::h_HF ) return "HF_had";
131  else if ( pfCandidateType == reco::PFCandidate::egamma_HF ) return "HF_em";
132  else assert(0);
133  }
134 
135  bool isMatchedByBlockElement(const reco::PFCandidate& pfCandidate1, const reco::PFCandidate& pfCandidate2, int minMatches1, int minMatches2, int maxUnmatchedBlockElements1plus2)
136  {
137  reco::PFCandidate::ElementsInBlocks blockElements1 = pfCandidate1.elementsInBlocks();
138  int numBlocks1 = blockElements1.size();
139  reco::PFCandidate::ElementsInBlocks blockElements2 = pfCandidate2.elementsInBlocks();
140  int numBlocks2 = blockElements2.size();
141  int numBlocks_matched = 0;
142  for ( reco::PFCandidate::ElementsInBlocks::const_iterator blockElement1 = blockElements1.begin();
143  blockElement1 != blockElements1.end(); ++blockElement1 ) {
144  bool isMatched = false;
145  for ( reco::PFCandidate::ElementsInBlocks::const_iterator blockElement2 = blockElements2.begin();
146  blockElement2 != blockElements2.end(); ++blockElement2 ) {
147  if ( blockElement1->first.id() == blockElement2->first.id() &&
148  blockElement1->first.key() == blockElement2->first.key() &&
149  blockElement1->second == blockElement2->second ) {
150  isMatched = true;
151  }
152  }
153  if ( isMatched ) ++numBlocks_matched;
154  }
155  assert(numBlocks_matched <= numBlocks1);
156  assert(numBlocks_matched <= numBlocks2);
157  if ( numBlocks_matched >= minMatches1 && numBlocks_matched >= minMatches2 &&
158  ((numBlocks1 - numBlocks_matched) + (numBlocks2 - numBlocks_matched)) <= maxUnmatchedBlockElements1plus2 ) {
159  return true;
160  } else {
161  return false;
162  }
163  }
164 }
165 
167 {
168  if ( verbosity_ ) {
169  edm::LogPrint("TauChHadronFromPF") << "<PFRecoTauChargedHadronFromPFCandidatePlugin::operator()>:";
170  edm::LogPrint("TauChHadronFromPF") << " pluginName = " << name() ;
171  }
172 
174 
175  // Get the candidates passing our quality cuts
178 
179  for ( PFCandPtrs::iterator cand = candsVector.begin();
180  cand != candsVector.end(); ++cand ) {
181  if ( verbosity_ ) {
182  edm::LogPrint("TauChHadronFromPF") << "processing PFCandidate: Pt = " << (*cand)->pt() << ", eta = " << (*cand)->eta() << ", phi = " << (*cand)->phi()
183  << " (type = " << getPFCandidateType((*cand)->particleId()) << ", charge = " << (*cand)->charge() << ")" ;
184  }
185 
187  if ( std::abs((*cand)->charge()) > 0.5 ) algo = PFRecoTauChargedHadron::kChargedPFCandidate;
189  std::auto_ptr<PFRecoTauChargedHadron> chargedHadron(new PFRecoTauChargedHadron(**cand, algo));
190  if ( (*cand)->trackRef().isNonnull() ) chargedHadron->track_ = edm::refToPtr((*cand)->trackRef());
191  else if ( (*cand)->muonRef().isNonnull() && (*cand)->muonRef()->innerTrack().isNonnull() ) chargedHadron->track_ = edm::refToPtr((*cand)->muonRef()->innerTrack());
192  else if ( (*cand)->muonRef().isNonnull() && (*cand)->muonRef()->globalTrack().isNonnull() ) chargedHadron->track_ = edm::refToPtr((*cand)->muonRef()->globalTrack());
193  else if ( (*cand)->muonRef().isNonnull() && (*cand)->muonRef()->outerTrack().isNonnull() ) chargedHadron->track_ = edm::refToPtr((*cand)->muonRef()->outerTrack());
194  else if ( (*cand)->gsfTrackRef().isNonnull() ) chargedHadron->track_ = edm::refToPtr((*cand)->gsfTrackRef());
195  chargedHadron->chargedPFCandidate_ = (*cand);
196  chargedHadron->addDaughter(*cand);
197 
198  chargedHadron->positionAtECALEntrance_ = (*cand)->positionAtECALEntrance();
199 
200  reco::PFCandidate::ParticleType chargedPFCandidateType = chargedHadron->chargedPFCandidate_->particleId();
201 
202  if ( chargedHadron->pt() > minMergeChargedHadronPt_ ) {
203  std::vector<reco::PFCandidatePtr> jetConstituents = jet.getPFConstituents();
204  for ( std::vector<reco::PFCandidatePtr>::const_iterator jetConstituent = jetConstituents.begin();
205  jetConstituent != jetConstituents.end(); ++jetConstituent ) {
206  // CV: take care of not double-counting energy in case "charged" PFCandidate is in fact a PFNeutralHadron
207  if ( (*jetConstituent) == chargedHadron->chargedPFCandidate_ ) continue;
208 
209  reco::PFCandidate::ParticleType jetConstituentType = (*jetConstituent)->particleId();
210  if ( !(jetConstituentType == reco::PFCandidate::h0 || jetConstituentType == reco::PFCandidate::gamma) ) continue;
211 
212  double dR = deltaR((*jetConstituent)->positionAtECALEntrance(), chargedHadron->positionAtECALEntrance_);
213  double dRmerge = -1.;
214  int minBlockElementMatches = 1000;
215  int maxUnmatchedBlockElements = 0;
216  double minMergeEt = 1.e+6;
217  if ( jetConstituentType == reco::PFCandidate::h0 ) {
218  if ( chargedPFCandidateType == reco::PFCandidate::h ) dRmerge = dRmergeNeutralHadronWrtChargedHadron_;
219  else if ( chargedPFCandidateType == reco::PFCandidate::h0 ) dRmerge = dRmergeNeutralHadronWrtNeutralHadron_;
220  else if ( chargedPFCandidateType == reco::PFCandidate::e ) dRmerge = dRmergeNeutralHadronWrtElectron_;
221  else dRmerge = dRmergeNeutralHadronWrtOther_;
222  minBlockElementMatches = minBlockElementMatchesNeutralHadron_;
223  maxUnmatchedBlockElements = maxUnmatchedBlockElementsNeutralHadron_;
224  minMergeEt = minMergeNeutralHadronEt_;
225  } else if ( jetConstituentType == reco::PFCandidate::gamma ) {
226  if ( chargedPFCandidateType == reco::PFCandidate::h ) dRmerge = dRmergePhotonWrtChargedHadron_;
227  else if ( chargedPFCandidateType == reco::PFCandidate::h0 ) dRmerge = dRmergePhotonWrtNeutralHadron_;
228  else if ( chargedPFCandidateType == reco::PFCandidate::e ) dRmerge = dRmergePhotonWrtElectron_;
229  else dRmerge = dRmergePhotonWrtOther_;
230  minBlockElementMatches = minBlockElementMatchesPhoton_;
231  maxUnmatchedBlockElements = maxUnmatchedBlockElementsPhoton_;
232  minMergeEt = minMergeGammaEt_;
233  }
234  if ( (*jetConstituent)->et() > minMergeEt &&
235  (dR < dRmerge || isMatchedByBlockElement(**jetConstituent, *chargedHadron->chargedPFCandidate_, minBlockElementMatches, minBlockElementMatches, maxUnmatchedBlockElements)) ) {
236  chargedHadron->neutralPFCandidates_.push_back(*jetConstituent);
237  chargedHadron->addDaughter(*jetConstituent);
238  }
239  }
240  }
241 
242  setChargedHadronP4(*chargedHadron);
243 
244  if ( verbosity_ ) {
245  edm::LogPrint("TauChHadronFromPF") << *chargedHadron;
246  }
247  // Update the vertex
248  if ( chargedHadron->daughterPtr(0).isNonnull() ) chargedHadron->setVertex(chargedHadron->daughterPtr(0)->vertex());
249  output.push_back(chargedHadron);
250  }
251 
252  return output.release();
253 }
254 
255 }} // end namespace reco::tau
256 
258 
T getParameter(std::string const &) const
ParticleType
particle types
Definition: PFCandidate.h:45
Ptr< typename C::value_type > refToPtr(Ref< C, typename C::value_type, refhelper::FindUsingAdvance< C, typename C::value_type > > const &ref)
Definition: RefToPtr.h:18
Coll filterCandRefs(const Coll &refcoll, bool invert=false) const
Filter a ref vector of PFCandidates.
PFRecoTauChargedHadronFromPFCandidatePlugin(const edm::ParameterSet &, edm::ConsumesCollector &&iC)
reco::VertexRef associatedVertex(const PFJet &jet) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
void setEvent(const edm::Event &evt)
Load the vertices from the event.
#define nullptr
boost::ptr_vector< PFRecoTauChargedHadron > ChargedHadronVector
Jets made from PFObjects.
Definition: PFJet.h:21
std::vector< ElementInBlock > ElementsInBlocks
Definition: PFCandidate.h:387
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
void setPV(const reco::VertexRef &vtx) const
Update the primary vertex.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isMatched(TrackingRecHit const &hit)
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
return_type operator()(const reco::PFJet &) const override
Build a collection of chargedHadrons from objects in the input jet.
ParameterSet const & getParameterSet(std::string const &) const
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:40
fixed size matrix
HLT enums.
virtual std::vector< reco::PFCandidatePtr > getPFConstituents() const
get all constituents
Definition: PFJet.cc:52
#define DEFINE_EDM_PLUGIN(factory, type, name)
const std::string & name() const
void setChargedHadronP4(reco::PFRecoTauChargedHadron &chargedHadron, double scaleFactor_neutralPFCands=1.0)
const ElementsInBlocks & elementsInBlocks() const
Definition: PFCandidate.cc:691
void beginEvent() override
Hook called at the beginning of the event.
def move(src, dest)
Definition: eostools.py:511