File indexing completed on 2024-04-06 12:27:52
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013 #include "RecoTauTag/RecoTau/interface/RecoTauBuilderPlugins.h"
0014 #include "RecoTauTag/RecoTau/interface/RecoTauCommonUtilities.h"
0015
0016 namespace reco {
0017 namespace tau {
0018
0019 namespace {
0020
0021 std::vector<CandidatePtr> deleteFrom(const CandidatePtr& ptr, const std::vector<CandidatePtr>& collection) {
0022 std::vector<CandidatePtr> output;
0023 for (std::vector<CandidatePtr>::const_iterator cand = collection.begin(); cand != collection.end(); ++cand) {
0024 if ((*cand) != ptr)
0025 output.push_back(*cand);
0026 }
0027 return output;
0028 }
0029 }
0030
0031 class RecoTauTwoProngFilter : public RecoTauModifierPlugin {
0032 public:
0033 explicit RecoTauTwoProngFilter(const edm::ParameterSet& pset, edm::ConsumesCollector&& iC);
0034 ~RecoTauTwoProngFilter() override {}
0035 void operator()(PFTau&) const override;
0036
0037 private:
0038 double minPtFractionForSecondProng_;
0039 };
0040
0041 RecoTauTwoProngFilter::RecoTauTwoProngFilter(const edm::ParameterSet& pset, edm::ConsumesCollector&& iC)
0042 : RecoTauModifierPlugin(pset, std::move(iC)) {
0043 minPtFractionForSecondProng_ = pset.getParameter<double>("minPtFractionForSecondProng");
0044 }
0045
0046 void RecoTauTwoProngFilter::operator()(PFTau& tau) const {
0047 if (tau.signalChargedHadrCands().size() == 2) {
0048 const std::vector<CandidatePtr>& signalCharged = tau.signalChargedHadrCands();
0049 size_t indexOfHighestPt = (signalCharged[0]->pt() > signalCharged[1]->pt()) ? 0 : 1;
0050 size_t indexOfLowerPt = (indexOfHighestPt) ? 0 : 1;
0051 double ratio = signalCharged[indexOfLowerPt]->pt() / signalCharged[indexOfHighestPt]->pt();
0052
0053 if (ratio < minPtFractionForSecondProng_) {
0054 CandidatePtr keep = signalCharged[indexOfHighestPt];
0055 CandidatePtr filter = signalCharged[indexOfLowerPt];
0056
0057 std::vector<CandidatePtr> newSignalCharged;
0058 newSignalCharged.push_back(keep);
0059 std::vector<CandidatePtr> newSignal = deleteFrom(filter, tau.signalCands());
0060
0061
0062 std::vector<CandidatePtr> newIsolationCharged = tau.isolationChargedHadrCands();
0063 newIsolationCharged.push_back(filter);
0064 std::vector<CandidatePtr> newIsolation = tau.isolationCands();
0065 newIsolation.push_back(filter);
0066
0067
0068 tau.setP4(tau.p4() - filter->p4());
0069 tau.setisolationPFChargedHadrCandsPtSum(tau.isolationPFChargedHadrCandsPtSum() - filter->pt());
0070 tau.setCharge(tau.charge() - filter->charge());
0071
0072 tau.setsignalChargedHadrCands(newSignalCharged);
0073 tau.setsignalCands(newSignal);
0074 tau.setisolationChargedHadrCands(newIsolationCharged);
0075 tau.setisolationCands(newIsolation);
0076 }
0077 }
0078 }
0079 }
0080 }
0081 #include "FWCore/Framework/interface/MakerMacros.h"
0082 DEFINE_EDM_PLUGIN(RecoTauModifierPluginFactory, reco::tau::RecoTauTwoProngFilter, "RecoTauTwoProngFilter");