CMS 3D CMS Logo

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

Go to the documentation of this file.
00001 /* class PFRecoTauDecayModeDeterminator
00002  *
00003  * Takes PFCandidates from PFTau and reconstructs tau decay mode.
00004  * Notably, merges photons (PFGammas) into pi zeros.
00005  * PFChargedHadrons are assumed to be charged pions.
00006  * Output candidate collections are owned (as shallow clones) by this object.
00007  * 
00008  * author: Evan K. Friis, UC Davis (evan.klose.friis@cern.ch) 
00009  */
00010 
00011 #include "DataFormats/TauReco/interface/PFTauTagInfo.h"
00012 #include "DataFormats/TauReco/interface/PFTauDecayMode.h"
00013 #include "DataFormats/TauReco/interface/PFTauDecayModeAssociation.h"
00014 #include "DataFormats/TauReco/interface/PFTau.h"
00015 #include "DataFormats/Candidate/interface/Candidate.h"
00016 #include "DataFormats/Candidate/interface/ShallowCloneCandidate.h"
00017 #include "RecoTauTag/TauTagTools/interface/TauTagTools.h"
00018 
00019 #include "FWCore/Framework/interface/EventSetup.h"
00020 #include "FWCore/Framework/interface/ESHandle.h"
00021 #include "FWCore/Framework/interface/Event.h"
00022 #include "FWCore/Framework/interface/MakerMacros.h"
00023 #include "FWCore/Framework/interface/Frameworkfwd.h"
00024 #include "FWCore/Framework/interface/EDProducer.h"
00025 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00026 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00027 
00028 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00029 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00030 #include "RecoTauTag/TauTagTools/interface/PFCandCommonVertexFitter.h"
00031 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
00032 #include "CommonTools/CandUtils/interface/AddFourMomenta.h"
00033 
00034 #include "CLHEP/Random/RandGauss.h"
00035 
00036 #include <memory>
00037 #include <algorithm>
00038 
00039 using namespace reco;
00040 using namespace edm;
00041 using namespace std;
00042 typedef reco::Particle::LorentzVector LorentzVector;
00043 
00044 class PFRecoTauDecayModeDeterminator : public EDProducer {
00045  public:
00046 
00047   typedef std::list<CompositeCandidate>                    compCandList;
00048   typedef std::list<CompositeCandidate>::reverse_iterator  compCandRevIter;
00049 
00050   void mergePiZeroes(compCandList&, compCandRevIter);
00051   void mergePiZeroesByBestMatch(compCandList&);
00052 
00053   explicit PFRecoTauDecayModeDeterminator(const edm::ParameterSet& iConfig);
00054   ~PFRecoTauDecayModeDeterminator();
00055   virtual void produce(edm::Event&,const edm::EventSetup&);
00056 
00057  protected:
00058   const double chargedPionMass;
00059   const double neutralPionMass;
00060 
00061   struct  gammaMatchContainer {
00062      double matchQuality;
00063      size_t firstIndex;
00064      size_t secondIndex;
00065   };
00066 
00067   static bool gammaMatchSorter (const gammaMatchContainer& first, const gammaMatchContainer& second);
00068 
00069  private:
00070   PFCandCommonVertexFitterBase* vertexFitter_;
00071   edm::InputTag              PFTauProducer_;
00072   AddFourMomenta        addP4;
00073   uint32_t              maxPhotonsToMerge_;             //number of photons allowed in a merged pi0
00074   double                maxPiZeroMass_;             
00075   bool                  mergeLowPtPhotonsFirst_;
00076   bool                  mergeByBestMatch_;
00077   bool                  setChargedPionMass_;
00078   bool                  setPi0Mass_;
00079   bool                  setMergedPi0Mass_;
00080   bool                  refitTracks_;
00081   bool                  filterTwoProngs_;
00082   bool                  filterPhotons_;  
00083   double                minPtFractionForSecondProng_;   //2 prongs whose second prong falls under 
00084   double                minPtFractionSinglePhotons_; 
00085   double                minPtFractionPiZeroes_; 
00086   TauTagTools::sortByDescendingPt<CompositeCandidate>   candDescendingSorter;
00087   TauTagTools::sortByAscendingPt<CompositeCandidate>    candAscendingSorter;
00088 };
00089 
00090 PFRecoTauDecayModeDeterminator::PFRecoTauDecayModeDeterminator(const edm::ParameterSet& iConfig):chargedPionMass(0.13957),neutralPionMass(0.13497){
00091   PFTauProducer_                = iConfig.getParameter<edm::InputTag>("PFTauProducer");
00092   maxPhotonsToMerge_            = iConfig.getParameter<uint32_t>("maxPhotonsToMerge");
00093   maxPiZeroMass_                = iConfig.getParameter<double>("maxPiZeroMass");             
00094   mergeLowPtPhotonsFirst_       = iConfig.getParameter<bool>("mergeLowPtPhotonsFirst");
00095   mergeByBestMatch_             = iConfig.getParameter<bool>("mergeByBestMatch");
00096   setChargedPionMass_           = iConfig.getParameter<bool>("setChargedPionMass");
00097   setPi0Mass_                   = iConfig.getParameter<bool>("setPi0Mass");
00098   setMergedPi0Mass_             = iConfig.getParameter<bool>("setMergedPi0Mass");
00099   refitTracks_                  = iConfig.getParameter<bool>("refitTracks");
00100   filterTwoProngs_              = iConfig.getParameter<bool>("filterTwoProngs");
00101   filterPhotons_                = iConfig.getParameter<bool>("filterPhotons");
00102   minPtFractionForSecondProng_  = iConfig.getParameter<double>("minPtFractionForSecondProng");
00103   minPtFractionSinglePhotons_   = iConfig.getParameter<double>("minPtFractionSinglePhotons");
00104   minPtFractionPiZeroes_        = iConfig.getParameter<double>("minPtFractionPiZeroes");
00105   //setup vertex fitter
00106   vertexFitter_ = new PFCandCommonVertexFitter<KalmanVertexFitter>(iConfig);
00107   produces<PFTauDecayModeAssociation>();      
00108 }
00109 
00110 PFRecoTauDecayModeDeterminator::~PFRecoTauDecayModeDeterminator()
00111 {
00112 //   delete vertexFitter_;  //now a very small memory leak, fix me later
00113 }
00114 
00115 
00116 /* 
00117  * ******************************************************************
00118    **     Merges a list of photons in to Pi0 candidates            **
00119    ******************************************************************
00120  */
00121 void PFRecoTauDecayModeDeterminator::mergePiZeroes(compCandList& input, compCandRevIter seed)
00122 {
00123    //uses std::list instead of vector, so that iterators can be deleted in situ
00124    //we go backwards for historical reasons ;)
00125    if(seed == input.rend())
00126       return;
00127    compCandRevIter bestMatchSoFar;
00128    LorentzVector combinationCandidate;
00129    float closestInvariantMassDifference = maxPiZeroMass_ + 1;
00130    bool foundACompatibleMatch = false;
00131    //find the best match to make a pi0
00132    compCandRevIter potentialMatch = seed;
00133    ++potentialMatch;
00134    for(; potentialMatch != input.rend(); ++potentialMatch)
00135    {
00136       // see how close this combination comes to the pion mass
00137       LorentzVector seedFourVector              = seed->p4();
00138       LorentzVector toAddFourVect               = potentialMatch->p4();
00139       combinationCandidate                      = seedFourVector + toAddFourVect;
00140       float combinationCandidateMass            = combinationCandidate.M();
00141       float differenceToTruePiZeroMass          = std::abs(combinationCandidateMass - neutralPionMass);
00142       if(combinationCandidateMass < maxPiZeroMass_ && differenceToTruePiZeroMass < closestInvariantMassDifference)
00143       {
00144          closestInvariantMassDifference = differenceToTruePiZeroMass;
00145          bestMatchSoFar = potentialMatch;
00146          foundACompatibleMatch = true;
00147       }
00148    }
00149    //if we found a combination that might make a pi0, combine it into the seed gamma, erase it, then see if we can add anymore
00150    if(foundACompatibleMatch && seed->numberOfDaughters() < maxPhotonsToMerge_)
00151    {
00152       //combine match into Seed and update four vector
00153       if(bestMatchSoFar->numberOfDaughters() > 0)
00154       {
00155          const Candidate* photonToAdd = (*bestMatchSoFar).daughter(0);
00156          seed->addDaughter(*photonToAdd);
00157       }
00158       addP4.set(*seed);
00159       //remove match as it is now contained in the seed 
00160       input.erase( (++bestMatchSoFar).base() );  //convert to normal iterator, after correct for offset
00161       mergePiZeroes(input, seed);
00162    } else
00163    {
00164       // otherwise move to next highest object and recurse
00165       addP4.set(*seed);
00166       ++seed;
00167       mergePiZeroes(input, seed);
00168    }
00169 }
00170 
00171 bool 
00172 PFRecoTauDecayModeDeterminator::gammaMatchSorter(const gammaMatchContainer& first,
00173                                                  const gammaMatchContainer& second)
00174 {
00175    return (first.matchQuality < second.matchQuality);
00176 }
00177 
00178 void PFRecoTauDecayModeDeterminator::mergePiZeroesByBestMatch(compCandList& input)
00179 {
00180    if(!input.size()) //nothing to merge... (NOTE: this line is necessary, as for size_t x, x in [0, +inf), x < -1 = true)
00181       return;
00182 
00183    std::vector<compCandList::iterator> gammas;       // iterators to all the gammas.  needed as we are using a list for compatability
00184                                                 // with the original merging algorithm, and this implementation requires random access
00185    std::vector<gammaMatchContainer> matches;
00186 
00187    // populate the list of gammas
00188    for(compCandList::iterator iGamma = input.begin(); iGamma != input.end(); ++iGamma)
00189       gammas.push_back(iGamma);
00190 
00191 
00192    for(size_t gammaA = 0; gammaA < gammas.size()-1; ++gammaA)
00193    {
00194       for(size_t gammaB = gammaA+1; gammaB < gammas.size(); ++gammaB)
00195       {
00196          //construct invariant mass of this pair
00197          LorentzVector piZeroAB = gammas[gammaA]->p4() + gammas[gammaB]->p4();
00198          //different to true pizero mass
00199          double piZeroABMass               = piZeroAB.M();
00200          double differenceToTruePiZeroMass = std::abs(piZeroABMass - neutralPionMass);
00201 
00202          if(piZeroABMass < maxPiZeroMass_)
00203          {
00204             gammaMatchContainer   aMatch;
00205             aMatch.matchQuality = differenceToTruePiZeroMass;
00206             aMatch.firstIndex   = gammaA;
00207             aMatch.secondIndex  = gammaB;
00208             matches.push_back(aMatch);
00209          }
00210       }
00211    }
00212 
00213    sort(matches.begin(), matches.end(), gammaMatchSorter);
00214    //the pairs whose mass is closest to the true pi0 mass are now at the beginning
00215    //of this vector
00216 
00217    for(std::vector<gammaMatchContainer>::iterator iMatch  = matches.begin(); 
00218                                              iMatch != matches.end();
00219                                            ++iMatch)
00220    {
00221       size_t gammaA = iMatch->firstIndex;
00222       size_t gammaB = iMatch->secondIndex;
00223       //check to see that both gammas in this match have not been used (ie their iterators set to input.end())
00224       if( gammas[gammaA] != input.end() && gammas[gammaB] != input.end() )
00225       {
00226          //merge the second gamma into the first; loop occurs in case of multiple gamma merging option
00227          for(size_t bDaughter = 0; bDaughter < gammas[gammaB]->numberOfDaughters(); ++bDaughter)
00228             gammas[gammaA]->addDaughter( *(gammas[gammaB]->daughter(bDaughter)) );
00229          //update the four vector information
00230          addP4.set(*gammas[gammaA]);
00231          //delete gammaB from the list of photons/pi zeroes, as it has been merged into gammaA
00232          input.erase(gammas[gammaB]);
00233          //mark both as "merged"
00234          gammas[gammaA] = input.end();
00235          gammas[gammaB] = input.end();
00236       } // else this match contains a photon that has already been merged
00237    }
00238 
00239 }
00240 
00241 void PFRecoTauDecayModeDeterminator::produce(edm::Event& iEvent,const edm::EventSetup& iSetup){
00242 
00243   edm::ESHandle<TransientTrackBuilder> myTransientTrackBuilder;
00244   edm::ESHandle<MagneticField> myMF;
00245 
00246   if (refitTracks_)
00247   {
00248      iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",myTransientTrackBuilder);
00249      iSetup.get<IdealMagneticFieldRecord>().get(myMF);
00250   }
00251 
00252   edm::Handle<PFTauCollection> thePFTauCollection;
00253   iEvent.getByLabel(PFTauProducer_,thePFTauCollection);
00254 
00255   auto_ptr<PFTauDecayModeAssociation> result(new PFTauDecayModeAssociation(PFTauRefProd(thePFTauCollection)));
00256 
00257   size_t numberOfPFTaus = thePFTauCollection->size();
00258   for(size_t iPFTau = 0; iPFTau < numberOfPFTaus; ++iPFTau)
00259   {
00260      //get the reference to the PFTau
00261      PFTauRef           pfTauRef(thePFTauCollection, iPFTau);
00262      PFTau              myPFTau = *pfTauRef;
00263      
00264      //get the charged & neutral collections corresponding to this PFTau
00265      const PFCandidateRefVector& theChargedHadronCandidates = myPFTau.signalPFChargedHadrCands();
00266      const PFCandidateRefVector& theGammaCandidates         = myPFTau.signalPFGammaCands();
00267 
00268      LorentzVector totalFourVector;                       //contains non-filtered stuff only.
00269 
00270      //shallow clone everything
00271      std::vector<ShallowCloneCandidate>    chargedCandidates;
00272      std::list<CompositeCandidate>         gammaCandidates;
00273      VertexCompositeCandidate         chargedCandsToAdd;  
00274      CompositeCandidate               filteredStuff;      //empty for now.
00275 
00276      bool needToProcessTracks = true;
00277      if (filterTwoProngs_ && theChargedHadronCandidates.size() == 2)
00278      {
00279         size_t indexOfHighestPt = (theChargedHadronCandidates[0]->pt() > theChargedHadronCandidates[1]->pt()) ? 0 : 1;
00280         size_t indexOfLowerPt   = ( indexOfHighestPt ) ? 0 : 1; 
00281         //maybe include a like signed requirement?? (future)
00282         double highPt = theChargedHadronCandidates[indexOfHighestPt]->pt();
00283         double lowPt  = theChargedHadronCandidates[indexOfLowerPt]->pt();
00284         if (lowPt/highPt < minPtFractionForSecondProng_)  //if it is super low, filter it!
00285         {
00286            needToProcessTracks = false;  //we are doing it here instead
00287            chargedCandsToAdd.addDaughter(ShallowCloneCandidate(CandidateBaseRef(theChargedHadronCandidates[indexOfHighestPt])));
00288            Candidate* justAdded = chargedCandsToAdd.daughter(chargedCandsToAdd.numberOfDaughters()-1);
00289            totalFourVector += justAdded->p4();
00290            if(setChargedPionMass_)
00291               justAdded->setMass(chargedPionMass);
00292            //add the two prong to the list of filtered stuff (to be added to the isolation collection later)
00293            filteredStuff.addDaughter(ShallowCloneCandidate(CandidateBaseRef(theChargedHadronCandidates[indexOfLowerPt])));
00294         }
00295      }
00296 
00297      if(needToProcessTracks) //not a two prong, filter is turned off, or 2nd prong passes cuts
00298      {
00299         for( PFCandidateRefVector::const_iterator iCharged  = theChargedHadronCandidates.begin();
00300               iCharged != theChargedHadronCandidates.end();
00301               ++iCharged)
00302         {
00303            // copy as shallow clone, and asssume mass of pi+
00304            chargedCandsToAdd.addDaughter(ShallowCloneCandidate(CandidateBaseRef(*iCharged)));
00305            Candidate* justAdded = chargedCandsToAdd.daughter(chargedCandsToAdd.numberOfDaughters()-1);
00306            totalFourVector += justAdded->p4();
00307            if(setChargedPionMass_)
00308               justAdded->setMass(chargedPionMass);
00309         }
00310      }
00311 
00312      for( PFCandidateRefVector::const_iterator iGamma  = theGammaCandidates.begin();
00313                                                iGamma != theGammaCandidates.end();
00314                                              ++iGamma)
00315      {
00316         CompositeCandidate potentialPiZero;
00317         potentialPiZero.addDaughter(ShallowCloneCandidate(CandidateBaseRef(*iGamma)));
00318         addP4.set(potentialPiZero);
00319         totalFourVector += potentialPiZero.p4();
00320         gammaCandidates.push_back(potentialPiZero);
00321      }
00322 
00323      //sort the photons by pt before passing to merger
00324      if (mergeLowPtPhotonsFirst_)
00325         gammaCandidates.sort(candDescendingSorter);
00326      else
00327         gammaCandidates.sort(candAscendingSorter);
00328 
00329      if (mergeByBestMatch_)
00330         mergePiZeroesByBestMatch(gammaCandidates);
00331      else
00332         mergePiZeroes(gammaCandidates, gammaCandidates.rbegin());
00333 
00334      if (filterPhotons_)
00335      {
00336         //sort by pt, from high to low.
00337         gammaCandidates.sort(candAscendingSorter);
00338 
00339         compCandRevIter wimp = gammaCandidates.rbegin();
00340 
00341         bool doneFiltering = false;
00342         while(!doneFiltering && wimp != gammaCandidates.rend())
00343         {
00344            double ptFraction          = wimp->pt()/totalFourVector.pt();
00345            size_t numberOfPhotons     = wimp->numberOfDaughters();
00346 
00347            //check if it is a single photon or has been merged
00348            if ( (numberOfPhotons == 1 && ptFraction < minPtFractionSinglePhotons_) ||
00349                 (numberOfPhotons  > 1 && ptFraction < minPtFractionPiZeroes_     )    )
00350            {
00351               //remove
00352               totalFourVector -= wimp->p4();
00353               for(size_t iDaughter = 0; iDaughter < numberOfPhotons; ++iDaughter)
00354               {
00355                  filteredStuff.addDaughter(ShallowCloneCandidate(CandidateBaseRef( wimp->daughter(iDaughter)->masterClone() )));
00356               }
00357 
00358               //move to the next photon to filter
00359               ++wimp;
00360            } else
00361            {
00362               //if this pizero passes the filter, we are done looking
00363               doneFiltering = true;
00364            }
00365         }
00366         //delete the filtered objects
00367         gammaCandidates.erase(wimp.base(), gammaCandidates.end());
00368      }
00369 
00370 
00371      CompositeCandidate mergedPiZerosToAdd;
00372      for( std::list<CompositeCandidate>::iterator iGamma  = gammaCandidates.begin();
00373                                              iGamma != gammaCandidates.end();
00374                                            ++iGamma)
00375      {
00376         if (setPi0Mass_) // set mass as pi 0
00377         {
00378            if (iGamma->numberOfDaughters() == 1) // for merged gamma pairs, check if user wants to keep ECAL mass
00379               iGamma->setMass(neutralPionMass);
00380            else if (setMergedPi0Mass_)
00381               iGamma->setMass(neutralPionMass);
00382         }
00383         mergedPiZerosToAdd.addDaughter(*iGamma);
00384      }
00385 
00386      // apply vertex fitting.
00387      if (refitTracks_ && chargedCandsToAdd.numberOfDaughters() > 1)
00388      {
00389         vertexFitter_->set(myMF.product());
00390         vertexFitter_->set(chargedCandsToAdd);  //refits tracks, adds vertexing info
00391      }
00392 
00393      // correctly set the four vectors of the composite candidates
00394      addP4.set(chargedCandsToAdd);
00395      addP4.set(mergedPiZerosToAdd);
00396      addP4.set(filteredStuff);
00397 
00398      /*
00399      LorentzVector refitFourVector = chargedCandsToAdd.p4() + mergedPiZerosToAdd.p4();
00400 
00401      edm::LogInfo("PFTauDecayModeDeterminator") << "Found nCharged: " << chargedCandsToAdd.numberOfDaughters()
00402                                   << " and nNeutral: " << mergedPiZerosToAdd.numberOfDaughters()
00403                                   << " Former mass: " << totalFourVector.mass() 
00404                                   << " New mass: " << refitFourVector.mass();
00405      */
00406 
00407      PFTauDecayMode myDecayModeTau(chargedCandsToAdd, mergedPiZerosToAdd, filteredStuff);
00408      myDecayModeTau.setPFTauRef(pfTauRef);
00409      result->setValue(iPFTau, myDecayModeTau);
00410   }
00411   iEvent.put(result);
00412 }
00413 DEFINE_FWK_MODULE(PFRecoTauDecayModeDeterminator);