CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Types | Private Attributes
reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin Class Reference
Inheritance diagram for reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin:
reco::tau::PFRecoTauChargedHadronBuilderPlugin reco::tau::RecoTauEventHolderPlugin reco::tau::RecoTauNamedPlugin

Public Member Functions

virtual void beginEvent ()
 Hook called at the beginning of the event. More...
 
return_type operator() (const reco::PFJet &) const
 Build a collection of chargedHadrons from objects in the input jet. More...
 
 PFRecoTauChargedHadronFromPFCandidatePlugin (const edm::ParameterSet &, edm::ConsumesCollector &&iC)
 
virtual ~PFRecoTauChargedHadronFromPFCandidatePlugin ()
 
- Public Member Functions inherited from reco::tau::PFRecoTauChargedHadronBuilderPlugin
 PFRecoTauChargedHadronBuilderPlugin (const edm::ParameterSet &pset, edm::ConsumesCollector &&iC)
 
virtual ~PFRecoTauChargedHadronBuilderPlugin ()
 
- Public Member Functions inherited from reco::tau::RecoTauEventHolderPlugin
const edm::Eventevt () const
 
edm::Eventevt ()
 
const edm::EventSetupevtSetup () const
 
 RecoTauEventHolderPlugin (const edm::ParameterSet &pset)
 
void setup (edm::Event &, const edm::EventSetup &)
 
virtual ~RecoTauEventHolderPlugin ()
 
- Public Member Functions inherited from reco::tau::RecoTauNamedPlugin
const std::string & name () const
 
 RecoTauNamedPlugin (const edm::ParameterSet &pset)
 
virtual ~RecoTauNamedPlugin ()
 

Private Types

typedef std::vector
< reco::PFCandidatePtr
PFCandPtrs
 

Private Attributes

double dRmergeNeutralHadronWrtChargedHadron_
 
double dRmergeNeutralHadronWrtElectron_
 
double dRmergeNeutralHadronWrtNeutralHadron_
 
double dRmergeNeutralHadronWrtOther_
 
double dRmergePhotonWrtChargedHadron_
 
double dRmergePhotonWrtElectron_
 
double dRmergePhotonWrtNeutralHadron_
 
double dRmergePhotonWrtOther_
 
std::vector< int > inputPdgIds_
 
int maxUnmatchedBlockElementsNeutralHadron_
 
int maxUnmatchedBlockElementsPhoton_
 
int minBlockElementMatchesNeutralHadron_
 
int minBlockElementMatchesPhoton_
 
double minMergeChargedHadronPt_
 
double minMergeGammaEt_
 
double minMergeNeutralHadronEt_
 
RecoTauQualityCutsqcuts_
 
int verbosity_
 
RecoTauVertexAssociator vertexAssociator_
 

Additional Inherited Members

- Public Types inherited from reco::tau::PFRecoTauChargedHadronBuilderPlugin
typedef boost::ptr_vector
< PFRecoTauChargedHadron
ChargedHadronVector
 
typedef std::auto_ptr
< ChargedHadronVector
return_type
 

Detailed Description

Definition at line 41 of file PFRecoTauChargedHadronFromPFCandidatePlugin.cc.

Member Typedef Documentation

Constructor & Destructor Documentation

reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::PFRecoTauChargedHadronFromPFCandidatePlugin ( const edm::ParameterSet pset,
edm::ConsumesCollector &&  iC 
)
explicit

Definition at line 79 of file PFRecoTauChargedHadronFromPFCandidatePlugin.cc.

References dRmergeNeutralHadronWrtChargedHadron_, dRmergeNeutralHadronWrtElectron_, dRmergeNeutralHadronWrtNeutralHadron_, dRmergeNeutralHadronWrtOther_, dRmergePhotonWrtChargedHadron_, dRmergePhotonWrtElectron_, dRmergePhotonWrtNeutralHadron_, dRmergePhotonWrtOther_, edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), edm::ParameterSet::getParameterSet(), inputPdgIds_, maxUnmatchedBlockElementsNeutralHadron_, maxUnmatchedBlockElementsPhoton_, minBlockElementMatchesNeutralHadron_, minBlockElementMatchesPhoton_, minMergeChargedHadronPt_, minMergeGammaEt_, minMergeNeutralHadronEt_, qcuts_, and verbosity_.

80  : PFRecoTauChargedHadronBuilderPlugin(pset,std::move(iC)),
81  vertexAssociator_(pset.getParameter<edm::ParameterSet>("qualityCuts"),std::move(iC)),
82  qcuts_(0)
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 }
T getParameter(std::string const &) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
ParameterSet const & getParameterSet(std::string const &) const
PFRecoTauChargedHadronBuilderPlugin(const edm::ParameterSet &pset, edm::ConsumesCollector &&iC)
reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::~PFRecoTauChargedHadronFromPFCandidatePlugin ( )
virtual

Member Function Documentation

void reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::beginEvent ( )
virtual
PFRecoTauChargedHadronFromPFCandidatePlugin::return_type reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::operator() ( const reco::PFJet ) const
virtual

Build a collection of chargedHadrons from objects in the input jet.

Implements reco::tau::PFRecoTauChargedHadronBuilderPlugin.

Definition at line 166 of file PFRecoTauChargedHadronFromPFCandidatePlugin.cc.

References reco::tau::RecoTauVertexAssociator::associatedVertex(), gather_cfg::cout, reco::deltaR(), PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, dRmergeNeutralHadronWrtChargedHadron_, dRmergeNeutralHadronWrtElectron_, dRmergeNeutralHadronWrtNeutralHadron_, dRmergeNeutralHadronWrtOther_, dRmergePhotonWrtChargedHadron_, dRmergePhotonWrtElectron_, dRmergePhotonWrtNeutralHadron_, dRmergePhotonWrtOther_, reco::PFCandidate::e, reco::tau::RecoTauQualityCuts::filterCandRefs(), reco::PFCandidate::gamma, reco::PFJet::getPFConstituents(), reco::PFCandidate::h, reco::PFCandidate::h0, inputPdgIds_, reco::PFRecoTauChargedHadron::kChargedPFCandidate, reco::PFRecoTauChargedHadron::kPFNeutralHadron, reco::PFRecoTauChargedHadron::kUndefined, maxUnmatchedBlockElementsNeutralHadron_, maxUnmatchedBlockElementsPhoton_, minBlockElementMatchesNeutralHadron_, minBlockElementMatchesPhoton_, minMergeChargedHadronPt_, minMergeGammaEt_, minMergeNeutralHadronEt_, reco::tau::RecoTauNamedPlugin::name(), convertSQLitetoXML_cfg::output, reco::tau::pfCandidates(), qcuts_, edm::refToPtr(), reco::tau::setChargedHadronP4(), reco::tau::RecoTauQualityCuts::setPV(), verbosity_, and vertexAssociator_.

167 {
168  if ( verbosity_ ) {
169  std::cout << "<PFRecoTauChargedHadronFromPFCandidatePlugin::operator()>:" << std::endl;
170  std::cout << " pluginName = " << name() << std::endl;
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  std::cout << "processing PFCandidate: Pt = " << (*cand)->pt() << ", eta = " << (*cand)->eta() << ", phi = " << (*cand)->phi()
183  << " (type = " << getPFCandidateType((*cand)->particleId()) << ", charge = " << (*cand)->charge() << ")" << std::endl;
184  }
185 
187  if ( fabs((*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  chargedHadron->print(std::cout);
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 }
reco::VertexRef associatedVertex(const PFJet &tau) const
ParticleType
particle types
Definition: PFCandidate.h:43
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.
std::vector< reco::PFCandidatePtr > PFCandPtrs
boost::ptr_vector< PFRecoTauChargedHadron > ChargedHadronVector
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
void setPV(const reco::VertexRef &vtx) const
Update the primary vertex.
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
tuple cout
Definition: gather_cfg.py:121
const std::string & name() const
void setChargedHadronP4(reco::PFRecoTauChargedHadron &chargedHadron, double scaleFactor_neutralPFCands=1.0)

Member Data Documentation

double reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::dRmergeNeutralHadronWrtChargedHadron_
private
double reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::dRmergeNeutralHadronWrtElectron_
private
double reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::dRmergeNeutralHadronWrtNeutralHadron_
private
double reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::dRmergeNeutralHadronWrtOther_
private
double reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::dRmergePhotonWrtChargedHadron_
private
double reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::dRmergePhotonWrtElectron_
private
double reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::dRmergePhotonWrtNeutralHadron_
private
double reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::dRmergePhotonWrtOther_
private
std::vector<int> reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::inputPdgIds_
private
int reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::maxUnmatchedBlockElementsNeutralHadron_
private
int reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::maxUnmatchedBlockElementsPhoton_
private
int reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::minBlockElementMatchesNeutralHadron_
private
int reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::minBlockElementMatchesPhoton_
private
double reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::minMergeChargedHadronPt_
private
double reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::minMergeGammaEt_
private
double reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::minMergeNeutralHadronEt_
private
RecoTauQualityCuts* reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::qcuts_
private
int reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::verbosity_
private
RecoTauVertexAssociator reco::tau::PFRecoTauChargedHadronFromPFCandidatePlugin::vertexAssociator_
private

Definition at line 54 of file PFRecoTauChargedHadronFromPFCandidatePlugin.cc.

Referenced by beginEvent(), and operator()().