CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoTauTag/RecoTau/plugins/RecoTauPiZeroStripPlugin.cc

Go to the documentation of this file.
00001 /*
00002  * RecoTauPiZeroStripPlugin
00003  *
00004  * Merges PFGammas in a PFJet into Candidate piZeros defined as
00005  * strips in eta-phi.
00006  *
00007  * Author: Michail Bachtis (University of Wisconsin)
00008  *
00009  * Code modifications: Evan Friis (UC Davis)
00010  *
00011  * $Id $
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/CombinatoricGenerator.h"
00030 
00031 namespace reco { namespace tau {
00032 
00033 namespace {
00034 // Apply a hypothesis on the mass of the strips.
00035 math::XYZTLorentzVector applyMassConstraint(
00036     const math::XYZTLorentzVector& vec,double mass) {
00037   double factor = sqrt(vec.energy()*vec.energy()-mass*mass)/vec.P();
00038   return math::XYZTLorentzVector(
00039       vec.px()*factor,vec.py()*factor,vec.pz()*factor,vec.energy());
00040 }
00041 }
00042 
00043 
00044 class RecoTauPiZeroStripPlugin : public RecoTauPiZeroBuilderPlugin {
00045   public:
00046     explicit RecoTauPiZeroStripPlugin(const edm::ParameterSet& pset);
00047     virtual ~RecoTauPiZeroStripPlugin() {}
00048     // Return type is auto_ptr<PiZeroVector>
00049     return_type operator()(const reco::PFJet& jet) const;
00050     // Hook to update PV information
00051     virtual void beginEvent();
00052 
00053   private:
00054     // PV needed for quality cuts
00055     edm::InputTag pvSrc_;
00056     RecoTauQualityCuts qcuts_;
00057 
00058     std::vector<int> inputPdgIds_; //type of candidates to clusterize
00059     double etaAssociationDistance_;//eta Clustering Association Distance
00060     double phiAssociationDistance_;//phi Clustering Association Distance
00061 
00062     // Parameters for build strip combinations
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.getParameter<edm::ParameterSet>("qualityCuts"))
00073 {
00074   pvSrc_ = pset.getParameter<edm::InputTag>("primaryVertexSrc");
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 // Update the primary vertex
00090 void RecoTauPiZeroStripPlugin::beginEvent() {
00091   edm::Handle<reco::VertexCollection> pvHandle;
00092   evt()->getByLabel(pvSrc_, pvHandle);
00093   if (pvHandle->size()) {
00094     qcuts_.setPV(reco::VertexRef(pvHandle, 0));
00095   }
00096 }
00097 
00098 RecoTauPiZeroStripPlugin::return_type RecoTauPiZeroStripPlugin::operator()(
00099     const reco::PFJet& jet) const {
00100   // Get list of gamma candidates
00101   typedef std::vector<reco::PFCandidatePtr> PFCandPtrs;
00102   typedef PFCandPtrs::iterator PFCandIter;
00103   PiZeroVector output;
00104 
00105   // Get the candidates passing our quality cuts
00106   PFCandPtrs candsVector = qcuts_.filterRefs(pfCandidates(jet, inputPdgIds_));
00107   //PFCandPtrs candsVector = qcuts_.filterRefs(pfGammas(jet));
00108 
00109   // Convert to stl::list to allow fast deletions
00110   typedef std::list<reco::PFCandidatePtr> PFCandPtrList;
00111   typedef std::list<reco::PFCandidatePtr>::iterator PFCandPtrListIter;
00112   PFCandPtrList cands;
00113   cands.insert(cands.end(), candsVector.begin(), candsVector.end());
00114 
00115   while (cands.size() > 0) {
00116     // Seed this new strip, and delete it from future strips
00117     PFCandidatePtr seed = cands.front();
00118     cands.pop_front();
00119 
00120     // Add a new candidate to our collection using this seed
00121     std::auto_ptr<RecoTauPiZero> strip(new RecoTauPiZero(
00122             *seed, RecoTauPiZero::kStrips));
00123     strip->addDaughter(seed);
00124 
00125     // Find all other objects in the strip
00126     PFCandPtrListIter stripCand = cands.begin();
00127     while(stripCand != cands.end()) {
00128       if( fabs(strip->eta() - (*stripCand)->eta()) < etaAssociationDistance_
00129           && fabs(deltaPhi(*strip, **stripCand)) < phiAssociationDistance_ ) {
00130         // Add candidate to strip
00131         strip->addDaughter(*stripCand);
00132         // Update the strips four momenta
00133         p4Builder_.set(*strip);
00134         // Delete this candidate from future strips and move on to
00135         // the next potential candidate
00136         stripCand = cands.erase(stripCand);
00137       } else {
00138         // This candidate isn't compatabile - just move to the next candidate
00139         ++stripCand;
00140       }
00141     }
00142     // Update the vertex
00143     if (strip->daughterPtr(0).isNonnull())
00144       strip->setVertex(strip->daughterPtr(0)->vertex());
00145     output.push_back(strip);
00146   }
00147 
00148   // Check if we want to combine our strips
00149   if (combineStrips_ && output.size() > 1) {
00150     PiZeroVector stripCombinations;
00151     // Sort the output by descending pt
00152     output.sort(output.begin(), output.end(),
00153         boost::bind(&RecoTauPiZero::pt, _1) >
00154         boost::bind(&RecoTauPiZero::pt, _2));
00155     // Get the end of interesting set of strips to try and combine
00156     PiZeroVector::const_iterator end_iter = takeNElements(
00157         output.begin(), output.end(), maxStrips_);
00158 
00159     // Look at all the combinations
00160     for (PiZeroVector::const_iterator first = output.begin();
00161         first != end_iter-1; ++first) {
00162       for (PiZeroVector::const_iterator second = first+1;
00163           second != end_iter; ++second) {
00164         Candidate::LorentzVector firstP4 = first->p4();
00165         Candidate::LorentzVector secondP4 = second->p4();
00166         // If we assume a certain mass for each strip apply it here.
00167         firstP4 = applyMassConstraint(firstP4, combinatoricStripMassHypo_);
00168         secondP4 = applyMassConstraint(secondP4, combinatoricStripMassHypo_);
00169         Candidate::LorentzVector totalP4 = firstP4 + secondP4;
00170         // Make our new combined strip
00171         std::auto_ptr<RecoTauPiZero> combinedStrips(
00172             new RecoTauPiZero(0, totalP4,
00173               Candidate::Point(0, 0, 0),
00174               //111, 10001, true, RecoTauPiZero::kCombinatoricStrips));
00175               111, 10001, true, RecoTauPiZero::kUndefined));
00176 
00177         // Now loop over the strip members
00178         BOOST_FOREACH(const RecoTauPiZero::daughters::value_type& gamma,
00179             first->daughterPtrVector()) {
00180           combinedStrips->addDaughter(gamma);
00181         }
00182         BOOST_FOREACH(const RecoTauPiZero::daughters::value_type& gamma,
00183             second->daughterPtrVector()) {
00184           combinedStrips->addDaughter(gamma);
00185         }
00186         // Update the vertex
00187         if (combinedStrips->daughterPtr(0).isNonnull())
00188           combinedStrips->setVertex(combinedStrips->daughterPtr(0)->vertex());
00189         // Add to our collection of combined strips
00190         stripCombinations.push_back(combinedStrips);
00191       }
00192     }
00193     // When done doing all the combinations, add the combined strips to the
00194     // output.
00195     output.transfer(output.end(), stripCombinations);
00196   }
00197 
00198   return output.release();
00199 }
00200 }} // end namespace reco::tau
00201 
00202 #include "FWCore/Framework/interface/MakerMacros.h"
00203 DEFINE_EDM_PLUGIN(RecoTauPiZeroBuilderPluginFactory,
00204     reco::tau::RecoTauPiZeroStripPlugin, "RecoTauPiZeroStripPlugin");