CMS 3D CMS Logo

PFTau3ProngReco.cc
Go to the documentation of this file.
1 /* class PFTau3ProngReco
2  * EDProducer of the
3  * author: Ian M. Nugent
4  * The idea of the fully reconstructing the tau using a kinematic fit comes from
5  * Lars Perchalla and Philip Sauerland Theses under Achim Stahl supervision. This
6  * code is a result of the continuation of this work by Ian M. Nugent and Vladimir Cherepanov.
7  */
8 
18 
23 
39 
45 
48 #include <TFormula.h>
49 
50 #include <memory>
51 
62 #include "TLorentzVector.h"
63 
65 
66 using namespace reco;
67 using namespace edm;
68 using namespace std;
69 using namespace tauImpactParameter;
70 
71 class PFTau3ProngReco : public EDProducer {
72 public:
73  enum Alg { useKalmanFit = 0, useTrackHelix };
74 
75  struct DiscCutPair {
76  DiscCutPair() : cutFormula_(nullptr) {}
77  ~DiscCutPair() { delete cutFormula_; }
80  double cut_;
81  TFormula* cutFormula_;
82  };
83  typedef std::vector<DiscCutPair*> DiscCutPairVec;
84 
85  explicit PFTau3ProngReco(const edm::ParameterSet& iConfig);
86  ~PFTau3ProngReco() override;
87  void produce(edm::Event&, const edm::EventSetup&) override;
88 
89 private:
94  std::unique_ptr<StringCutObjectSelector<reco::PFTau>> cut_;
95  int ndfPVT_;
97 };
98 
100  : PFTauTag_(iConfig.getParameter<edm::InputTag>("PFTauTag")),
101  PFTauTIPTag_(iConfig.getParameter<edm::InputTag>("PFTauTIPTag")),
102  Algorithm_(iConfig.getParameter<int>("Algorithm")),
103  ndfPVT_(iConfig.getUntrackedParameter("ndfPVT", (int)5)) {
105  std::vector<edm::ParameterSet> discriminators =
106  iConfig.getParameter<std::vector<edm::ParameterSet>>("discriminators");
107  // Build each of our cuts
108  for (auto const& pset : discriminators) {
109  DiscCutPair* newCut = new DiscCutPair();
110  newCut->inputTag_ = pset.getParameter<edm::InputTag>("discriminator");
111  if (pset.existsAs<std::string>("selectionCut"))
112  newCut->cutFormula_ = new TFormula("selectionCut", pset.getParameter<std::string>("selectionCut").data());
113  else
114  newCut->cut_ = pset.getParameter<double>("selectionCut");
115  discriminators_.push_back(newCut);
116  }
117  // Build a string cut if desired
118  if (iConfig.exists("cut"))
119  cut_ = std::make_unique<StringCutObjectSelector<reco::PFTau>>(iConfig.getParameter<std::string>("cut"));
121  produces<edm::AssociationVector<PFTauRefProd, std::vector<reco::PFTau3ProngSummaryRef>>>();
122  produces<PFTau3ProngSummaryCollection>("PFTau3ProngSummary");
123 }
124 
126 
127 namespace {
128  inline const reco::Track* getTrack(const reco::Candidate& cand) {
129  const reco::PFCandidate* pfCandPtr = dynamic_cast<const reco::PFCandidate*>(&cand);
130  if (pfCandPtr) {
131  if (pfCandPtr->trackRef().isNonnull())
132  return pfCandPtr->trackRef().get();
133  else if (pfCandPtr->gsfTrackRef().isNonnull())
134  return pfCandPtr->gsfTrackRef().get();
135  else
136  return nullptr;
137  }
138  const pat::PackedCandidate* packedCand = dynamic_cast<const pat::PackedCandidate*>(&cand);
139  if (packedCand != nullptr)
140  return packedCand->bestTrack();
141 
142  return nullptr;
143  }
144 } // namespace
145 
147  // Obtain Collections
148  edm::ESHandle<TransientTrackBuilder> transTrackBuilder;
149  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", transTrackBuilder);
150 
152  iEvent.getByLabel(PFTauTag_, Tau);
153 
155  iEvent.getByLabel(PFTauTIPTag_, TIPAV);
156 
157  auto AVPFTau3PS = std::make_unique<edm::AssociationVector<PFTauRefProd, std::vector<reco::PFTau3ProngSummaryRef>>>(
158  PFTauRefProd(Tau));
159  auto PFTau3PSCollection_out = std::make_unique<PFTau3ProngSummaryCollection>();
160  reco::PFTau3ProngSummaryRefProd PFTau3RefProd_out =
161  iEvent.getRefBeforePut<reco::PFTau3ProngSummaryCollection>("PFTau3ProngSummary");
162 
163  // Load each discriminator
164  for (auto& disc : discriminators_) {
165  iEvent.getByLabel(disc->inputTag_, disc->handle_);
166  }
167 
168  // For each Tau Run Algorithim
169  if (Tau.isValid()) {
170  for (reco::PFTauCollection::size_type iPFTau = 0; iPFTau < Tau->size(); iPFTau++) {
171  reco::PFTauRef tau(Tau, iPFTau);
172  reco::PFTau3ProngSummary PFTau3PS;
174  // Check if it passed all the discrimiantors
175  bool passed(true);
176  for (auto const& disc : discriminators_) {
177  // Check this discriminator passes
178  bool passedDisc = true;
179  if (disc->cutFormula_)
180  passedDisc = (disc->cutFormula_->Eval((*disc->handle_)[tau]) > 0.5);
181  else
182  passedDisc = ((*disc->handle_)[tau] > disc->cut_);
183  if (!passedDisc) {
184  passed = false;
185  break;
186  }
187  }
188  if (passed && cut_.get()) {
189  passed = (*cut_)(*tau);
190  }
191  if (passed) {
192  PDGInfo pdgInfo;
193  const reco::PFTauTransverseImpactParameterRef theTIP = TIPAV->value(tau.key());
194  const reco::VertexRef primaryVertex = theTIP->primaryVertex();
196  // Now compute the 3 prong Tau
197  bool SecondaryVtxOK(false);
199  if (theTIP->hasSecondaryVertex() && primaryVertex->ndof() > ndfPVT_) {
200  const VertexRef secVtx = theTIP->secondaryVertex();
201  GlobalPoint sv(secVtx->position().x(), secVtx->position().y(), secVtx->position().z());
202  double vtxchi2(0), vtxndf(1);
203  if (useKalmanFit == Algorithm_) {
204  vtxchi2 = secVtx->chi2();
205  vtxndf = secVtx->ndof();
206  const std::vector<reco::Track>& selectedTracks = secVtx->refittedTracks();
207  std::vector<reco::TransientTrack> transTrkVect;
208  for (unsigned int i = 0; i != selectedTracks.size(); i++)
209  transTrkVect.push_back(transTrackBuilder->build(selectedTracks[i]));
211  float piMassSigma(1.e-6), piChi(0.0), piNdf(0.0);
212  std::vector<RefCountedKinematicParticle> pions;
213  for (unsigned int i = 0; i < transTrkVect.size(); i++)
214  pions.push_back(kinFactory.particle(transTrkVect[i], pdgInfo.pi_mass(), piChi, piNdf, sv, piMassSigma));
215  RefCountedKinematicTree jpTree = kpvFitter_.fit(pions);
216  jpTree->movePointerToTheTop();
217  const KinematicParameters parameters = jpTree->currentParticle()->currentState().kinematicParameters();
218  AlgebraicSymMatrix77 cov = jpTree->currentParticle()->currentState().kinematicParametersError().matrix();
219  // get pions
220  double c(0);
221  std::vector<reco::Track> Tracks;
222  std::vector<LorentzVectorParticle> ReFitPions;
223  for (unsigned int i = 0; i < transTrkVect.size(); i++) {
224  c += transTrkVect[i].charge();
225  ReFitPions.push_back(ParticleBuilder::createLorentzVectorParticle(transTrkVect[i], *secVtx, true, true));
226  }
227  // now covert a1 into LorentzVectorParticle
228  TVectorT<double> a1_par(LorentzVectorParticle::NLorentzandVertexPar);
229  TMatrixTSym<double> a1_cov(LorentzVectorParticle::NLorentzandVertexPar);
230  for (int i = 0; i < LorentzVectorParticle::NLorentzandVertexPar; i++) {
231  a1_par(i) = parameters(i);
232  for (int j = 0; j < LorentzVectorParticle::NLorentzandVertexPar; j++) {
233  a1_cov(i, j) = cov(i, j);
234  }
235  }
237  a1_par, a1_cov, abs(PdtPdgMini::a_1_plus) * c, c, transTrackBuilder->field()->inInverseGeV(sv).z());
238  SecondaryVtxOK = true;
239  PFTau3PS = reco::PFTau3ProngSummary(theTIP, a1.p4(), vtxchi2, vtxndf);
240  } else if (useTrackHelix == Algorithm_) {
241  // use Track Helix
242  std::vector<TrackParticle> pions;
243  GlobalPoint pvpoint(
244  primaryVertex->position().x(), primaryVertex->position().y(), primaryVertex->position().z());
245  const std::vector<edm::Ptr<reco::Candidate>> cands = tau->signalChargedHadrCands();
246  for (std::vector<edm::Ptr<reco::Candidate>>::const_iterator iter = cands.begin(); iter != cands.end();
247  ++iter) {
248  const reco::Track* track = getTrack(**iter);
249  if (track != nullptr) {
250  reco::TransientTrack transTrk = transTrackBuilder->build(*track);
251  pions.push_back(ParticleBuilder::createTrackParticle(transTrk, pvpoint, true, true));
252  }
253  }
254  TVector3 pv(secVtx->position().x(), secVtx->position().y(), secVtx->position().z());
255  Chi2VertexFitter chi2v(pions, pv);
256  SecondaryVtxOK = chi2v.fit();
257  double c(0);
258  for (unsigned int i = 0; i < pions.size(); i++) {
259  c += pions[i].charge();
260  }
261  int pdgid = abs(PdtPdgMini::a_1_plus) * c;
262  a1 = chi2v.getMother(pdgid);
263  PFTau3PS = reco::PFTau3ProngSummary(theTIP, a1.p4(), vtxchi2, vtxndf);
264  }
265  }
266  if (SecondaryVtxOK) {
267  // Tau Solver
268  TVector3 pv(primaryVertex->position().x(), primaryVertex->position().y(), primaryVertex->position().z());
269  TMatrixTSym<double> pvcov(LorentzVectorParticle::NVertex);
271  primaryVertex->fill(pvCov);
272  for (int i = 0; i < LorentzVectorParticle::NVertex; i++) {
273  for (int j = 0; j < LorentzVectorParticle::NVertex; j++) {
274  pvcov(i, j) = pvCov(i, j);
275  }
276  }
277  for (unsigned int i = 0; i < PFTau3ProngSummary::nsolutions; i++) {
278  TauA1NuConstrainedFitter TauA1NU(i, a1, pv, pvcov);
279  bool isFitOK = TauA1NU.fit();
280  if (isFitOK) {
281  LorentzVectorParticle theTau = TauA1NU.getMother();
282  std::vector<LorentzVectorParticle> daughter = TauA1NU.getRefitDaughters();
283  std::vector<TLorentzVector> daughter_p4;
284  std::vector<int> daughter_charge, daughter_PDGID;
285  for (unsigned int d = 0; d < daughter.size(); d++) {
286  daughter_p4.push_back(daughter[d].p4());
287  daughter_charge.push_back((int)daughter[d].charge());
288  daughter_PDGID.push_back(daughter[d].pdgId());
289  }
290  PFTau3PS.AddSolution(i, theTau.p4(), daughter_p4, daughter_charge, daughter_PDGID, (isFitOK), 0.0, -999);
291  }
292  }
293  }
294  }
295  reco::PFTau3ProngSummaryRef PFTau3PSRef =
296  reco::PFTau3ProngSummaryRef(PFTau3RefProd_out, PFTau3PSCollection_out->size());
297  PFTau3PSCollection_out->push_back(PFTau3PS);
298  AVPFTau3PS->setValue(iPFTau, PFTau3PSRef);
299  }
300  }
301  iEvent.put(std::move(PFTau3PSCollection_out), "PFTau3ProngSummary");
302  iEvent.put(std::move(AVPFTau3PS));
303 }
304 
PdtPdgMini.h
reco::PFCandidate::trackRef
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:408
edm::RefProd
Definition: EDProductfwd.h:25
PFTau3ProngReco::Algorithm_
int Algorithm_
Definition: PFTau3ProngReco.cc:92
RefProd.h
mps_fire.i
i
Definition: mps_fire.py:355
PFTauFwd.h
Muon.h
MessageLogger.h
tauImpactParameter::LorentzVectorParticle
Definition: LorentzVectorParticle.h:17
EDProducer.h
ESHandle.h
PFTauDiscriminator.h
metsig::tau
Definition: SignAlgoResolutions.h:49
KinematicParticleFactoryFromTransientTrack
Definition: KinematicParticleFactoryFromTransientTrack.h:16
PFTau3ProngReco::discriminators_
DiscCutPairVec discriminators_
Definition: PFTau3ProngReco.cc:93
reco::PFTau3ProngSummaryCollection
std::vector< reco::PFTau3ProngSummary > PFTau3ProngSummaryCollection
collection of PFTau3ProngSummary objects
Definition: PFTau3ProngSummaryFwd.h:9
L1TRate_Offline_cfi.Tau
Tau
Definition: L1TRate_Offline_cfi.py:43
edm
HLT enums.
Definition: AlignableModifier.h:19
TrackParticle.h
PDWG_TauSkim_cff.discriminators
discriminators
Definition: PDWG_TauSkim_cff.py:6
TransientTrackBuilder::field
const MagneticField * field() const
Definition: TransientTrackBuilder.h:70
reco::PFTauRefProd
edm::RefProd< PFTauCollection > PFTauRefProd
references to PFTau collection
Definition: PFTauFwd.h:15
PFTau3ProngReco::DiscCutPair::handle_
edm::Handle< reco::PFTauDiscriminator > handle_
Definition: PFTau3ProngReco.cc:78
PFTau3ProngReco::DiscCutPair
Definition: PFTau3ProngReco.cc:75
tauImpactParameter
Definition: Chi2VertexFitter.h:14
edm::Ref::get
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:232
AlgebraicSymMatrix77
ROOT::Math::SMatrix< double, 7, 7, ROOT::Math::MatRepSym< double, 7 > > AlgebraicSymMatrix77
Definition: Matrices.h:9
ReferenceCountingPointer< KinematicTree >
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
Tau
Definition: Tau.py:1
Association.h
Particle.h
PdtPdgMini::a_1_plus
Definition: PdtPdgMini.h:89
edm::Handle< reco::PFTauDiscriminator >
KinematicParticleVertexFitter.h
TtSemiLepEvtBuilder_cfi.disc
disc
Definition: TtSemiLepEvtBuilder_cfi.py:60
reco::PFTau3ProngSummary
Definition: PFTau3ProngSummary.h:26
edm::Ref< PFTauCollection >
parameters
parameters
Definition: BeamSpot_PayloadInspector.cc:14
GenParticle.h
PFTau3ProngSummaryFwd.h
PV3DBase::z
T z() const
Definition: PV3DBase.h:61
BeamMonitor_cff.primaryVertex
primaryVertex
hltOfflineBeamSpot for HLTMON
Definition: BeamMonitor_cff.py:7
MakerMacros.h
trigger::size_type
uint16_t size_type
Definition: TriggerTypeDefs.h:18
Track.h
edm::EventSetup::get
T get() const
Definition: EventSetup.h:73
TrackFwd.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
BeamSpot.h
pat::PackedCandidate::bestTrack
const reco::Track * bestTrack() const override
return a pointer to the track if present. otherwise, return a null pointer
Definition: PackedCandidate.h:779
PFTau3ProngReco::DiscCutPair::inputTag_
edm::InputTag inputTag_
Definition: PFTau3ProngReco.cc:79
tauImpactParameter::TauA1NuConstrainedFitter
Definition: TauA1NuConstrainedFitter.h:21
MuonFwd.h
tauImpactParameter::PDGInfo
Definition: PDGInfo.h:11
pfDeepBoostedJetPreprocessParams_cfi.sv
sv
Definition: pfDeepBoostedJetPreprocessParams_cfi.py:226
TransientTrackRecord
Definition: TransientTrackRecord.h:11
reco::Track
Definition: Track.h:27
edm::ESHandle< TransientTrackBuilder >
PFTau3ProngReco::useTrackHelix
Definition: PFTau3ProngReco.cc:73
Point3DBase< float, GlobalTag >
tauImpactParameter::LorentzVectorParticle::p4
TLorentzVector p4() const
Definition: LorentzVectorParticle.h:55
PFTau3ProngReco::PFTau3ProngReco
PFTau3ProngReco(const edm::ParameterSet &iConfig)
Definition: PFTau3ProngReco.cc:99
LorentzVectorParticle.h
KinematicParticleFactoryFromTransientTrack::particle
RefCountedKinematicParticle particle(const reco::TransientTrack &initialTrack, const ParticleMass &massGuess, float chiSquared, float degreesOfFr, float &m_sigma) const
Definition: KinematicParticleFactoryFromTransientTrack.cc:15
RefToBase.h
KinematicParticleVertexFitter::fit
RefCountedKinematicTree fit(const std::vector< RefCountedKinematicParticle > &particles) const
Definition: KinematicParticleVertexFitter.cc:47
MagneticField::inInverseGeV
GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
Definition: MagneticField.h:36
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
Vertex.h
PFTau3ProngReco::PFTauTIPTag_
edm::InputTag PFTauTIPTag_
Definition: PFTau3ProngReco.cc:91
ALCARECOTkAlJpsiMuMu_cff.charge
charge
Definition: ALCARECOTkAlJpsiMuMu_cff.py:47
PFTau3ProngReco::DiscCutPairVec
std::vector< DiscCutPair * > DiscCutPairVec
Definition: PFTau3ProngReco.cc:83
edm::ParameterSet::exists
bool exists(std::string const &parameterName) const
checks if a parameter exists
Definition: ParameterSet.cc:674
PFTau3ProngReco::DiscCutPair::cut_
double cut_
Definition: PFTau3ProngReco.cc:80
ElectronFwd.h
TransientTrackBuilder.h
PFTau3ProngReco::PFTauTag_
edm::InputTag PFTauTag_
Definition: PFTau3ProngReco.cc:90
edm::ParameterSet
Definition: ParameterSet.h:36
reco::PFCandidate::gsfTrackRef
reco::GsfTrackRef gsfTrackRef() const
Definition: PFCandidate.cc:440
Event.h
pat::PackedCandidate
Definition: PackedCandidate.h:22
PackedCandidate.h
edm::Ref::isNonnull
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
tauImpactParameter::PDGInfo::pi_mass
static double pi_mass()
Definition: PDGInfo.h:13
PFTau3ProngReco_cfi.PFTau3ProngReco
PFTau3ProngReco
Definition: PFTau3ProngReco_cfi.py:3
cand
Definition: decayParser.h:34
KinematicParameters
Definition: KinematicParameters.h:15
createfilelist.int
int
Definition: createfilelist.py:10
MetAnalyzer.pv
def pv(vc)
Definition: MetAnalyzer.py:7
iEvent
int iEvent
Definition: GenABIO.cc:224
gather_cfg.Tracks
Tracks
Definition: gather_cfg.py:251
p4
double p4[4]
Definition: TauolaWrapper.h:92
GsfTrack.h
PFTau3ProngReco::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: PFTau3ProngReco.cc:146
PFTau3ProngReco::DiscCutPair::cutFormula_
TFormula * cutFormula_
Definition: PFTau3ProngReco.cc:81
EgammaValidation_cff.pdgId
pdgId
Definition: EgammaValidation_cff.py:118
edm::EventSetup
Definition: EventSetup.h:57
HltBtagPostValidation_cff.c
c
Definition: HltBtagPostValidation_cff.py:31
TransientTrackRecord.h
get
#define get
tauImpactParameter::TauA1NuConstrainedFitter::fit
bool fit()
Definition: TauA1NuConstrainedFitter.cc:222
edm::Ptr< reco::Candidate >
reco::Candidate
Definition: Candidate.h:27
VertexFwd.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
reco::PFTau3ProngSummaryRef
edm::Ref< PFTau3ProngSummaryCollection > PFTau3ProngSummaryRef
presistent reference to a PFTau3ProngSummary
Definition: PFTau3ProngSummaryFwd.h:13
getTrack
static const TrackGhostTrackState * getTrack(const BasicGhostTrackState *basic)
Definition: GhostTrackState.cc:49
PFTau3ProngReco::cut_
std::unique_ptr< StringCutObjectSelector< reco::PFTau > > cut_
Definition: PFTau3ProngReco.cc:94
TransientVertex.h
Electron.h
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
PFTau3ProngReco::DiscCutPair::~DiscCutPair
~DiscCutPair()
Definition: PFTau3ProngReco.cc:77
StringCutObjectSelector.h
GsfTrackFwd.h
reco::TransientTrack
Definition: TransientTrack.h:19
TauA1NuConstrainedFitter.h
Frameworkfwd.h
KinematicParticleVertexFitter
Definition: KinematicParticleVertexFitter.h:25
StringCutObjectSelector< reco::PFTau >
HLTMuonOfflineAnalyzer_cfi.NVertex
NVertex
parameters for fixed-width plots
Definition: HLTMuonOfflineAnalyzer_cfi.py:93
math::Error::type
ErrorD< N >::type type
Definition: Error.h:32
Chi2VertexFitter.h
PFTau.h
TriggerAnalyzer.passed
passed
Definition: TriggerAnalyzer.py:62
HLT_2018_cff.cands
cands
Definition: HLT_2018_cff.py:13762
PFTau3ProngReco::Alg
Alg
Definition: PFTau3ProngReco.cc:73
KinematicParticleFactoryFromTransientTrack.h
EventSetup.h
edm::EDProducer
Definition: EDProducer.h:36
Exception.h
reco::PFCandidate
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:40
TransientTrackBuilder::build
reco::TransientTrack build(const reco::Track *p) const
Definition: TransientTrackBuilder.cc:20
HLT_2018_cff.track
track
Definition: HLT_2018_cff.py:10352
ztail.d
d
Definition: ztail.py:151
ParticleBuilder.h
AssociationVector.h
tauImpactParameter::TauA1NuConstrainedFitter::getRefitDaughters
std::vector< LorentzVectorParticle > getRefitDaughters()
Definition: TauA1NuConstrainedFitter.cc:176
PDGInfo.h
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
tauImpactParameter::Chi2VertexFitter
Definition: Chi2VertexFitter.h:16
ParameterSet.h
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
edm::Event
Definition: Event.h:73
EgammaValidation_cff.pdgid
pdgid
Definition: EgammaValidation_cff.py:30
AdaptiveVertexFitter.h
PFTau3ProngReco
Definition: PFTau3ProngReco.cc:71
tauImpactParameter::TauA1NuConstrainedFitter::getMother
LorentzVectorParticle getMother()
Definition: TauA1NuConstrainedFitter.cc:195
PFTau3ProngSummary.h
edm::InputTag
Definition: InputTag.h:15
PFTau3ProngReco::kpvFitter_
KinematicParticleVertexFitter kpvFitter_
Definition: PFTau3ProngReco.cc:96
reco::PFTau3ProngSummary::AddSolution
virtual bool AddSolution(unsigned int solution, const TLorentzVector &tau, const std::vector< TLorentzVector > &daughter_p4, const std::vector< int > &daughter_charge, const std::vector< int > &daughter_PDGID, bool has3ProngSolution, double solutionChi2, double thetaGJsig)
Definition: PFTau3ProngSummary.cc:65
PFTau3ProngReco::ndfPVT_
int ndfPVT_
Definition: PFTau3ProngReco.cc:95
PFTau3ProngReco::DiscCutPair::DiscCutPair
DiscCutPair()
Definition: PFTau3ProngReco.cc:76
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27
TrackCollections2monitor_cff.selectedTracks
selectedTracks
Definition: TrackCollections2monitor_cff.py:32
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37
PFTau3ProngReco::~PFTau3ProngReco
~PFTau3ProngReco() override
Definition: PFTau3ProngReco.cc:125
PFTau3ProngReco::useKalmanFit
Definition: PFTau3ProngReco.cc:73