CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/RecoTauTag/RecoTau/plugins/PFRecoTauDiscriminationByIsolation.cc

Go to the documentation of this file.
00001 #include <functional>
00002 #include <boost/foreach.hpp>
00003 #include "RecoTauTag/RecoTau/interface/TauDiscriminationProducerBase.h"
00004 #include "DataFormats/Candidate/interface/LeafCandidate.h"
00005 #include "RecoTauTag/RecoTau/interface/RecoTauQualityCuts.h"
00006 #include "RecoTauTag/RecoTau/interface/RecoTauVertexAssociator.h"
00007 #include "RecoTauTag/RecoTau/interface/ConeTools.h"
00008 #include "DataFormats/VertexReco/interface/Vertex.h"
00009 
00010 #include "TMath.h"
00011 #include "TFormula.h"
00012 
00013 /* class PFRecoTauDiscriminationByIsolation
00014  * created : Jul 23 2007,
00015  * revised : Thu Aug 13 14:44:40 PDT 2009
00016  * contributors : Ludovic Houchu (IPHC, Strasbourg),
00017  *               Christian Veelken (UC Davis),
00018 *                Evan K. Friis (UC Davis)
00019  *               Michalis Bachtis (UW Madison)
00020  */
00021 
00022 using namespace reco;
00023 using namespace std;
00024 
00025 class PFRecoTauDiscriminationByIsolation :
00026   public PFTauDiscriminationProducerBase  {
00027   public:
00028     explicit PFRecoTauDiscriminationByIsolation(const edm::ParameterSet& pset):
00029       PFTauDiscriminationProducerBase(pset),
00030       qualityCutsPSet_(pset.getParameter<edm::ParameterSet>("qualityCuts")) {
00031 
00032         includeTracks_ = pset.getParameter<bool>(
00033             "ApplyDiscriminationByTrackerIsolation");
00034         includeGammas_ = pset.getParameter<bool>(
00035             "ApplyDiscriminationByECALIsolation");
00036 
00037         applyOccupancyCut_ = pset.getParameter<bool>("applyOccupancyCut");
00038         maximumOccupancy_ = pset.getParameter<uint32_t>("maximumOccupancy");
00039 
00040         applySumPtCut_ = pset.getParameter<bool>("applySumPtCut");
00041         maximumSumPt_ = pset.getParameter<double>("maximumSumPtCut");
00042 
00043         applyRelativeSumPtCut_ = pset.getParameter<bool>(
00044             "applyRelativeSumPtCut");
00045         maximumRelativeSumPt_ = pset.getParameter<double>(
00046             "relativeSumPtCut");
00047 
00048         storeRawOccupancy_ = pset.exists("storeRawOccupancy") ?
00049           pset.getParameter<bool>("storeRawOccupancy") : false;
00050         storeRawSumPt_ = pset.exists("storeRawSumPt") ?
00051           pset.getParameter<bool>("storeRawSumPt") : false;
00052 
00053         // Sanity check on requested options.  We can't apply cuts and store the
00054         // raw output at the same time
00055         if (applySumPtCut_ || applyOccupancyCut_ || applyRelativeSumPtCut_) {
00056           if (storeRawSumPt_ || storeRawOccupancy_) {
00057             throw cms::Exception("BadIsoConfig") <<
00058               "A 'store raw' and a 'apply cut' option have been set to true "
00059               << "simultaneously.  These options are mutually exclusive.";
00060           }
00061         }
00062 
00063         // Can only store one type
00064         if (storeRawSumPt_ && storeRawOccupancy_) {
00065             throw cms::Exception("BadIsoConfig") <<
00066               "Both 'store sum pt' and 'store occupancy' options are set."
00067               << " These options are mutually exclusive.";
00068         }
00069 
00070         if (pset.exists("customOuterCone")) {
00071           customIsoCone_ = pset.getParameter<double>("customOuterCone");
00072         } else {
00073           customIsoCone_ = -1;
00074         }
00075 
00076         // Get the quality cuts specific to the isolation region
00077         edm::ParameterSet isolationQCuts = qualityCutsPSet_.getParameterSet(
00078             "isolationQualityCuts");
00079 
00080         qcuts_.reset(new tau::RecoTauQualityCuts(isolationQCuts));
00081 
00082         vertexAssociator_.reset(
00083             new tau::RecoTauVertexAssociator(qualityCutsPSet_));
00084 
00085         applyDeltaBeta_ = pset.exists("applyDeltaBetaCorrection") ?
00086           pset.getParameter<bool>("applyDeltaBetaCorrection") : false;
00087 
00088         if (applyDeltaBeta_) {
00089           // Factorize the isolation QCuts into those that are used to
00090           // select PU and those that are not.
00091           std::pair<edm::ParameterSet, edm::ParameterSet> puFactorizedIsoQCuts =
00092             reco::tau::factorizePUQCuts(isolationQCuts);
00093 
00094           // Determine the pt threshold for the PU tracks
00095           // First check if the user specifies explicitly the cut.
00096           if (pset.exists("deltaBetaPUTrackPtCutOverride")) {
00097             puFactorizedIsoQCuts.second.addParameter<double>(
00098                 "minTrackPt",
00099                 pset.getParameter<double>("deltaBetaPUTrackPtCutOverride"));
00100           } else {
00101             // Secondly take it from the minGammaEt
00102             puFactorizedIsoQCuts.second.addParameter<double>(
00103                 "minTrackPt",
00104                 isolationQCuts.getParameter<double>("minGammaEt"));
00105           }
00106 
00107           pileupQcutsPUTrackSelection_.reset(new tau::RecoTauQualityCuts(
00108                 puFactorizedIsoQCuts.first));
00109 
00110           pileupQcutsGeneralQCuts_.reset(new tau::RecoTauQualityCuts(
00111                 puFactorizedIsoQCuts.second));
00112 
00113           pfCandSrc_ = pset.getParameter<edm::InputTag>("particleFlowSrc");
00114           vertexSrc_ = pset.getParameter<edm::InputTag>("vertexSrc");
00115           deltaBetaCollectionCone_ = pset.getParameter<double>(
00116               "isoConeSizeForDeltaBeta");
00117           std::string deltaBetaFactorFormula =
00118             pset.getParameter<string>("deltaBetaFactor");
00119           deltaBetaFormula_.reset(
00120               new TFormula("DB_corr", deltaBetaFactorFormula.c_str()));
00121         }
00122 
00123         applyRhoCorrection_ = pset.exists("applyRhoCorrection") ?
00124           pset.getParameter<bool>("applyRhoCorrection") : false;
00125         if (applyRhoCorrection_) {
00126           rhoProducer_ = pset.getParameter<edm::InputTag>("rhoProducer");
00127           rhoConeSize_ = pset.getParameter<double>("rhoConeSize");
00128           rhoUEOffsetCorrection_ =
00129             pset.getParameter<double>("rhoUEOffsetCorrection");
00130         }
00131       }
00132 
00133     ~PFRecoTauDiscriminationByIsolation(){}
00134 
00135     void beginEvent(const edm::Event& evt, const edm::EventSetup& evtSetup);
00136     double discriminate(const PFTauRef& pfTau);
00137 
00138   private:
00139     edm::ParameterSet qualityCutsPSet_;
00140     std::auto_ptr<tau::RecoTauQualityCuts> qcuts_;
00141 
00142     // Inverted QCut which selects tracks with bad DZ/trackWeight
00143     std::auto_ptr<tau::RecoTauQualityCuts> pileupQcutsPUTrackSelection_;
00144     std::auto_ptr<tau::RecoTauQualityCuts> pileupQcutsGeneralQCuts_;
00145 
00146     std::auto_ptr<tau::RecoTauVertexAssociator> vertexAssociator_;
00147 
00148     bool includeTracks_;
00149     bool includeGammas_;
00150     bool applyOccupancyCut_;
00151     uint32_t maximumOccupancy_;
00152     bool applySumPtCut_;
00153     double maximumSumPt_;
00154     bool applyRelativeSumPtCut_;
00155     double maximumRelativeSumPt_;
00156     double customIsoCone_;
00157 
00158     // Options to store the raw value in the discriminator instead of
00159     // boolean float
00160     bool storeRawOccupancy_;
00161     bool storeRawSumPt_;
00162 
00163     /* **********************************************************************
00164        **** Pileup Subtraction Parameters ***********************************
00165        **********************************************************************/
00166 
00167     // Delta Beta correction
00168     bool applyDeltaBeta_;
00169     edm::InputTag pfCandSrc_;
00170     // Keep track of how many vertices are in the event
00171     edm::InputTag vertexSrc_;
00172     std::vector<reco::PFCandidateRef> chargedPFCandidatesInEvent_;
00173     // Size of cone used to collect PU tracks
00174     double deltaBetaCollectionCone_;
00175     std::auto_ptr<TFormula> deltaBetaFormula_;
00176     double deltaBetaFactorThisEvent_;
00177 
00178     // Rho correction
00179     bool applyRhoCorrection_;
00180     edm::InputTag rhoProducer_;
00181     double rhoConeSize_;
00182     double rhoUEOffsetCorrection_;
00183     double rhoCorrectionThisEvent_;
00184     double rhoThisEvent_;
00185   };
00186 
00187 void PFRecoTauDiscriminationByIsolation::beginEvent(const edm::Event& event,
00188     const edm::EventSetup& eventSetup) {
00189 
00190   // NB: The use of the PV in this context is necessitated by its use in
00191   // applying quality cuts to the different objects in the isolation cone
00192   // The vertex associator contains the logic to select the appropriate vertex
00193   // We need to pass it the event so it can load the vertices.
00194   vertexAssociator_->setEvent(event);
00195 
00196   // If we are applying the delta beta correction, we need to get the PF
00197   // candidates from the event so we can find the PU tracks.
00198   chargedPFCandidatesInEvent_.clear();
00199   if (applyDeltaBeta_) {
00200     // Collect all the PF pile up tracks
00201     edm::Handle<reco::PFCandidateCollection> pfCandHandle_;
00202     event.getByLabel(pfCandSrc_, pfCandHandle_);
00203     chargedPFCandidatesInEvent_.reserve(pfCandHandle_->size());
00204     for (size_t i = 0; i < pfCandHandle_->size(); ++i) {
00205       reco::PFCandidateRef pfCand(pfCandHandle_, i);
00206       if (pfCand->charge() != 0)
00207         chargedPFCandidatesInEvent_.push_back(pfCand);
00208     }
00209     // Count all the vertices in the event, to parameterize the DB
00210     // correction factor
00211     edm::Handle<reco::VertexCollection> vertices;
00212     event.getByLabel(vertexSrc_, vertices);
00213     size_t nVtxThisEvent = vertices->size();
00214     deltaBetaFactorThisEvent_ = deltaBetaFormula_->Eval(nVtxThisEvent);
00215   }
00216 
00217   if (applyRhoCorrection_) {
00218     edm::Handle<double> rhoHandle_;
00219     event.getByLabel(rhoProducer_, rhoHandle_);
00220     rhoThisEvent_ = (*rhoHandle_ - rhoUEOffsetCorrection_)*
00221       (3.14159)*rhoConeSize_*rhoConeSize_;
00222   }
00223 }
00224 
00225 double
00226 PFRecoTauDiscriminationByIsolation::discriminate(const PFTauRef& pfTau) {
00227   // collect the objects we are working with (ie tracks, tracks+gammas, etc)
00228   std::vector<PFCandidateRef> isoCharged;
00229   std::vector<PFCandidateRef> isoNeutral;
00230   std::vector<PFCandidateRef> isoPU;
00231 
00232   // Get the primary vertex associated to this tau
00233   reco::VertexRef pv = vertexAssociator_->associatedVertex(*pfTau);
00234   // Let the quality cuts know which the vertex to use when applying selections
00235   // on dz, etc.
00236   qcuts_->setPV(pv);
00237   qcuts_->setLeadTrack(pfTau->leadPFChargedHadrCand());
00238   if (applyDeltaBeta_) {
00239     pileupQcutsGeneralQCuts_->setPV(pv);
00240     pileupQcutsGeneralQCuts_->setLeadTrack(pfTau->leadPFChargedHadrCand());
00241     pileupQcutsPUTrackSelection_->setPV(pv);
00242     pileupQcutsPUTrackSelection_->setLeadTrack(pfTau->leadPFChargedHadrCand());
00243   }
00244   // Load the tracks if they are being used.
00245   if (includeTracks_) {
00246     BOOST_FOREACH(const reco::PFCandidateRef& cand,
00247         pfTau->isolationPFChargedHadrCands()) {
00248       if (qcuts_->filterRef(cand))
00249         isoCharged.push_back(cand);
00250     }
00251   }
00252 
00253   if (includeGammas_) {
00254     BOOST_FOREACH(const reco::PFCandidateRef& cand,
00255         pfTau->isolationPFGammaCands()) {
00256       if (qcuts_->filterRef(cand))
00257         isoNeutral.push_back(cand);
00258     }
00259   }
00260   typedef reco::tau::cone::DeltaRPtrFilter<PFCandidateRef> DRFilter;
00261 
00262   // If desired, get PU tracks.
00263   if (applyDeltaBeta_) {
00264     // First select by inverted the DZ/track weight cuts. True = invert
00265     //std::cout << "Initial PFCands: " << chargedPFCandidatesInEvent_.size()
00266     //  << std::endl;
00267 
00268     std::vector<PFCandidateRef> allPU =
00269       pileupQcutsPUTrackSelection_->filterRefs(
00270           chargedPFCandidatesInEvent_, true);
00271 
00272     //std::cout << "After track cuts: " << allPU.size() << std::endl;
00273 
00274     // Now apply the rest of the cuts, like pt, and TIP, tracker hits, etc
00275     std::vector<PFCandidateRef> cleanPU =
00276       pileupQcutsGeneralQCuts_->filterRefs(allPU);
00277 
00278     //std::cout << "After cleaning cuts: " << cleanPU.size() << std::endl;
00279 
00280     // Only select PU tracks inside the isolation cone.
00281     DRFilter deltaBetaFilter(pfTau->p4(), 0, deltaBetaCollectionCone_);
00282     BOOST_FOREACH(const reco::PFCandidateRef& cand, cleanPU) {
00283       if (deltaBetaFilter(cand)) {
00284         isoPU.push_back(cand);
00285       }
00286     }
00287     //std::cout << "After cone cuts: " << isoPU.size() << std::endl;
00288   }
00289 
00290   // Check if we want a custom iso cone
00291   if (customIsoCone_ >= 0.) {
00292     DRFilter filter(pfTau->p4(), 0, customIsoCone_);
00293     std::vector<PFCandidateRef> isoCharged_filter;
00294     std::vector<PFCandidateRef> isoNeutral_filter;
00295     // Remove all the objects not in our iso cone
00296      BOOST_FOREACH(const PFCandidateRef& isoObject, isoCharged) {
00297       if(filter(isoObject)) isoCharged_filter.push_back(isoObject);
00298     }
00299     BOOST_FOREACH(const PFCandidateRef& isoObject, isoNeutral) {
00300       if(filter(isoObject)) isoNeutral_filter.push_back(isoObject);
00301     }
00302     isoCharged.clear();
00303     isoCharged=isoCharged_filter;
00304     isoNeutral.clear();
00305     isoNeutral=isoNeutral_filter;
00306 
00307   }
00308 
00309   bool failsOccupancyCut     = false;
00310   bool failsSumPtCut         = false;
00311   bool failsRelativeSumPtCut = false;
00312 
00313   //--- nObjects requirement
00314   int neutrals = isoNeutral.size();
00315 
00316   if (applyDeltaBeta_) {
00317     neutrals -= TMath::Nint(deltaBetaFactorThisEvent_*isoPU.size());
00318   }
00319   if(neutrals < 0) {
00320     neutrals=0;
00321   }
00322 
00323   size_t nOccupants = isoCharged.size() + neutrals;
00324 
00325   failsOccupancyCut = ( nOccupants > maximumOccupancy_ );
00326 
00327   double totalPt=0.0;
00328   //--- Sum PT requirement
00329   if( applySumPtCut_ || applyRelativeSumPtCut_ || storeRawSumPt_) {
00330     double chargedPt=0.0;
00331     double puPt=0.0;
00332     double neutralPt=0.0;
00333     BOOST_FOREACH(const PFCandidateRef& isoObject, isoCharged) {
00334       chargedPt += isoObject->pt();
00335     }
00336     BOOST_FOREACH(const PFCandidateRef& isoObject, isoNeutral) {
00337       neutralPt += isoObject->pt();
00338     }
00339     BOOST_FOREACH(const PFCandidateRef& isoObject, isoPU) {
00340       puPt += isoObject->pt();
00341     }
00342 
00343     if (applyDeltaBeta_) {
00344       neutralPt -= deltaBetaFactorThisEvent_*puPt;
00345     }
00346 
00347     if (applyRhoCorrection_) {
00348       neutralPt -= rhoThisEvent_;
00349     }
00350 
00351     if (neutralPt < 0.0) {
00352       neutralPt = 0.0;
00353     }
00354 
00355     totalPt = chargedPt+neutralPt;
00356 
00357     failsSumPtCut = (totalPt > maximumSumPt_);
00358 
00359     //--- Relative Sum PT requirement
00360     failsRelativeSumPtCut = (
00361         (pfTau->pt() > 0 ? totalPt/pfTau->pt() : 0 ) > maximumRelativeSumPt_ );
00362   }
00363 
00364   bool fails = (applyOccupancyCut_ && failsOccupancyCut) ||
00365     (applySumPtCut_ && failsSumPtCut) ||
00366     (applyRelativeSumPtCut_ && failsRelativeSumPtCut);
00367 
00368   // We did error checking in the constructor, so this is safe.
00369   if (storeRawSumPt_)
00370     return totalPt;
00371   else if (storeRawOccupancy_)
00372     return nOccupants;
00373   else
00374     return (fails ? 0. : 1.);
00375 }
00376 
00377 DEFINE_FWK_MODULE(PFRecoTauDiscriminationByIsolation);