Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013 #include <algorithm>
00014 #include <memory>
00015
00016 #include "RecoTauTag/RecoTau/interface/RecoTauPiZeroPlugins.h"
00017
00018 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
00019 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00020 #include "DataFormats/VertexReco/interface/Vertex.h"
00021 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00022 #include "DataFormats/TauReco/interface/RecoTauPiZero.h"
00023 #include "DataFormats/JetReco/interface/PFJet.h"
00024 #include "CommonTools/CandUtils/interface/AddFourMomenta.h"
00025 #include "DataFormats/Math/interface/deltaPhi.h"
00026
00027 #include "RecoTauTag/RecoTau/interface/RecoTauCommonUtilities.h"
00028 #include "RecoTauTag/RecoTau/interface/RecoTauQualityCuts.h"
00029 #include "RecoTauTag/RecoTau/interface/RecoTauVertexAssociator.h"
00030 #include "RecoTauTag/RecoTau/interface/CombinatoricGenerator.h"
00031
00032 namespace reco { namespace tau {
00033
00034 namespace {
00035
00036 math::XYZTLorentzVector applyMassConstraint(
00037 const math::XYZTLorentzVector& vec,double mass) {
00038 double factor = sqrt(vec.energy()*vec.energy()-mass*mass)/vec.P();
00039 return math::XYZTLorentzVector(
00040 vec.px()*factor,vec.py()*factor,vec.pz()*factor,vec.energy());
00041 }
00042 }
00043
00044
00045 class RecoTauPiZeroStripPlugin : public RecoTauPiZeroBuilderPlugin {
00046 public:
00047 explicit RecoTauPiZeroStripPlugin(const edm::ParameterSet& pset);
00048 virtual ~RecoTauPiZeroStripPlugin() {}
00049
00050 return_type operator()(const reco::PFJet& jet) const;
00051
00052 virtual void beginEvent();
00053
00054 private:
00055 RecoTauQualityCuts qcuts_;
00056 RecoTauVertexAssociator vertexAssociator_;
00057
00058 std::vector<int> inputPdgIds_;
00059 double etaAssociationDistance_;
00060 double phiAssociationDistance_;
00061
00062
00063 bool combineStrips_;
00064 int maxStrips_;
00065 double combinatoricStripMassHypo_;
00066
00067 AddFourMomenta p4Builder_;
00068 };
00069
00070 RecoTauPiZeroStripPlugin::RecoTauPiZeroStripPlugin(
00071 const edm::ParameterSet& pset):RecoTauPiZeroBuilderPlugin(pset),
00072 qcuts_(pset.getParameterSet(
00073 "qualityCuts").getParameterSet("signalQualityCuts")),
00074 vertexAssociator_(pset.getParameter<edm::ParameterSet>("qualityCuts")) {
00075 inputPdgIds_ = pset.getParameter<std::vector<int> >(
00076 "stripCandidatesParticleIds");
00077 etaAssociationDistance_ = pset.getParameter<double>(
00078 "stripEtaAssociationDistance");
00079 phiAssociationDistance_ = pset.getParameter<double>(
00080 "stripPhiAssociationDistance");
00081 combineStrips_ = pset.getParameter<bool>("makeCombinatoricStrips");
00082 if (combineStrips_) {
00083 maxStrips_ = pset.getParameter<int>("maxInputStrips");
00084 combinatoricStripMassHypo_ =
00085 pset.getParameter<double>("stripMassWhenCombining");
00086 }
00087 }
00088
00089
00090 void RecoTauPiZeroStripPlugin::beginEvent() {
00091 vertexAssociator_.setEvent(*evt());
00092 }
00093
00094 RecoTauPiZeroStripPlugin::return_type RecoTauPiZeroStripPlugin::operator()(
00095 const reco::PFJet& jet) const {
00096
00097 typedef std::vector<reco::PFCandidatePtr> PFCandPtrs;
00098 typedef PFCandPtrs::iterator PFCandIter;
00099 PiZeroVector output;
00100
00101
00102 qcuts_.setPV(vertexAssociator_.associatedVertex(jet));
00103 PFCandPtrs candsVector = qcuts_.filterRefs(pfCandidates(jet, inputPdgIds_));
00104
00105
00106
00107 typedef std::list<reco::PFCandidatePtr> PFCandPtrList;
00108 typedef std::list<reco::PFCandidatePtr>::iterator PFCandPtrListIter;
00109 PFCandPtrList cands;
00110 cands.insert(cands.end(), candsVector.begin(), candsVector.end());
00111
00112 while (cands.size() > 0) {
00113
00114 PFCandidatePtr seed = cands.front();
00115 cands.pop_front();
00116
00117
00118 std::auto_ptr<RecoTauPiZero> strip(new RecoTauPiZero(
00119 *seed, RecoTauPiZero::kStrips));
00120 strip->addDaughter(seed);
00121
00122
00123 PFCandPtrListIter stripCand = cands.begin();
00124 while(stripCand != cands.end()) {
00125 if( fabs(strip->eta() - (*stripCand)->eta()) < etaAssociationDistance_
00126 && fabs(deltaPhi(*strip, **stripCand)) < phiAssociationDistance_ ) {
00127
00128 strip->addDaughter(*stripCand);
00129
00130 p4Builder_.set(*strip);
00131
00132
00133 stripCand = cands.erase(stripCand);
00134 } else {
00135
00136 ++stripCand;
00137 }
00138 }
00139
00140 if (strip->daughterPtr(0).isNonnull())
00141 strip->setVertex(strip->daughterPtr(0)->vertex());
00142 output.push_back(strip);
00143 }
00144
00145
00146 if (combineStrips_ && output.size() > 1) {
00147 PiZeroVector stripCombinations;
00148
00149 output.sort(output.begin(), output.end(),
00150 boost::bind(&RecoTauPiZero::pt, _1) >
00151 boost::bind(&RecoTauPiZero::pt, _2));
00152
00153 PiZeroVector::const_iterator end_iter = takeNElements(
00154 output.begin(), output.end(), maxStrips_);
00155
00156
00157 for (PiZeroVector::const_iterator first = output.begin();
00158 first != end_iter-1; ++first) {
00159 for (PiZeroVector::const_iterator second = first+1;
00160 second != end_iter; ++second) {
00161 Candidate::LorentzVector firstP4 = first->p4();
00162 Candidate::LorentzVector secondP4 = second->p4();
00163
00164 firstP4 = applyMassConstraint(firstP4, combinatoricStripMassHypo_);
00165 secondP4 = applyMassConstraint(secondP4, combinatoricStripMassHypo_);
00166 Candidate::LorentzVector totalP4 = firstP4 + secondP4;
00167
00168 std::auto_ptr<RecoTauPiZero> combinedStrips(
00169 new RecoTauPiZero(0, totalP4,
00170 Candidate::Point(0, 0, 0),
00171
00172 111, 10001, true, RecoTauPiZero::kUndefined));
00173
00174
00175 BOOST_FOREACH(const RecoTauPiZero::daughters::value_type& gamma,
00176 first->daughterPtrVector()) {
00177 combinedStrips->addDaughter(gamma);
00178 }
00179 BOOST_FOREACH(const RecoTauPiZero::daughters::value_type& gamma,
00180 second->daughterPtrVector()) {
00181 combinedStrips->addDaughter(gamma);
00182 }
00183
00184 if (combinedStrips->daughterPtr(0).isNonnull())
00185 combinedStrips->setVertex(combinedStrips->daughterPtr(0)->vertex());
00186
00187 stripCombinations.push_back(combinedStrips);
00188 }
00189 }
00190
00191
00192 output.transfer(output.end(), stripCombinations);
00193 }
00194
00195 return output.release();
00196 }
00197 }}
00198
00199 #include "FWCore/Framework/interface/MakerMacros.h"
00200 DEFINE_EDM_PLUGIN(RecoTauPiZeroBuilderPluginFactory,
00201 reco::tau::RecoTauPiZeroStripPlugin, "RecoTauPiZeroStripPlugin");