CMS 3D CMS Logo

PFTauDiscriminantManager.cc

Go to the documentation of this file.
00001 #include "RecoTauTag/TauTagTools/interface/PFTauDiscriminantManager.h"
00002 
00003 namespace PFTauDiscriminants
00004 {
00005 using namespace std;
00006 
00007 typedef std::vector<const reco::Candidate*> candPtrVector;
00008 
00009 PFTauDiscriminantManager::PFTauDiscriminantManager()
00010 {
00011    iAmSignal_           = false;
00012    iAmNull_             = false;
00013    eventWeight_         = 1.0;
00014    currentTauDecayMode_ = NULL;
00015    eventData_           = NULL;
00016    mainTrack_           = NULL;
00017 }
00018 
00019 void
00020 PFTauDiscriminantManager::addDiscriminant(Discriminant* const discriminant)
00021 {
00022    if (!discriminant)
00023    {
00024       edm::LogError("PFTauDiscriminantManager") << "Error adding a discriminant, null pointer!";
00025       return;
00026    }
00027    string            discriminantName = discriminant->name();
00028    myDiscriminants_.insert(make_pair(discriminantName, discriminant));
00029 }
00030 
00031 void 
00032 PFTauDiscriminantManager::clearCache()
00033 {
00034    mainTrack_ = NULL;
00035    signalObjectsSortedByPt_.clear();
00036    signalObjectsSortedByDR_.clear();
00037    outlierObjectsSortedByPt_.clear();
00038    outlierObjectsSortedByDR_.clear();
00039 }
00040 
00041 bool
00042 PFTauDiscriminantManager::setEventData(const reco::PFTauDecayMode& theTau, const edm::Event& iEvent, const double& eventWeight, bool prePass, bool preFail)
00043 {
00044    currentTauDecayMode_ = &theTau;
00045    eventData_           = &iEvent;
00046    eventWeight_         = eventWeight;
00047    iAmNull_             = false;
00048    prePass_             = prePass;
00049    preFail_             = preFail;
00050    //reset cached collections
00051    clearCache();
00052    
00053    for(discriminantHolder::iterator aDiscriminant  = myDiscriminants_.begin();
00054                                     aDiscriminant != myDiscriminants_.end();
00055                                   ++aDiscriminant)
00056    {
00057       Discriminant* const theAlgo  = aDiscriminant->second;
00058       if (!theAlgo)
00059       {
00060          string theName  = aDiscriminant->first;
00061          edm::LogError("PFTauDiscriminantManager") << "Error filling discriminant " << theName <<", null pointer!";
00062          return false;
00063       }
00064       theAlgo->compute(this);
00065    }
00066    return true;
00067 }
00068 
00069 bool
00070 PFTauDiscriminantManager::setNullResult(const edm::Event& iEvent, const double& eventWeight)
00071 {
00072    currentTauDecayMode_ = NULL;
00073    eventData_           = &iEvent;
00074    eventWeight_         = eventWeight;
00075    iAmNull_             = true;
00076    prePass_             = false;
00077    preFail_             = false;
00078    //reset cached collections
00079    clearCache();
00080 
00081    for(discriminantHolder::iterator aDiscriminant  = myDiscriminants_.begin();
00082                                     aDiscriminant != myDiscriminants_.end();
00083                                   ++aDiscriminant)
00084    {
00085       Discriminant* const theAlgo  = aDiscriminant->second;
00086       if (!theAlgo)
00087       {
00088          string theName  = aDiscriminant->first;
00089          edm::LogError("PFTauDiscriminantManager") << "Error filling discriminant " << theName <<", null pointer!";
00090          return false;
00091       }
00092       theAlgo->setNullResult(this);
00093    }
00094    return true;
00095 }
00096 
00097 
00098 
00099 void PFTauDiscriminantManager::fillSignalObjects(candPtrVector& toFill)
00100 {
00101    toFill.clear();
00102    if (currentTauDecayMode_ == NULL)
00103    {
00104       edm::LogError("PFTauDiscriminantManager") << "Trying to get signal objects from null PFTauDecayMode object!  Returning empty vector...";
00105       return;
00106    }
00107    candPtrVector tempChargedVector = currentTauDecayMode_->chargedPionCandidates();
00108    candPtrVector tempNeutralVector = currentTauDecayMode_->neutralPionCandidates();
00109    toFill.insert(toFill.end(), tempChargedVector.begin(), tempChargedVector.end());
00110    toFill.insert(toFill.end(), tempNeutralVector.begin(), tempNeutralVector.end());
00111 }
00112 
00113 void PFTauDiscriminantManager::fillOutlierObjects(candPtrVector& toFill)
00114 {
00115    toFill.clear();
00116    if (currentTauDecayMode_ == NULL)
00117    {
00118       edm::LogError("PFTauDiscriminantManager") << "Trying to get QCD objects from null PFTauDecayMode object!  Returning empty vector...";
00119       return;
00120    }
00121    // get associated PFTau from PFTauDecayMode
00122    const PFTau* originalTau = currentTauDecayMode_->pfTauRef().get();
00123 
00124    if(originalTau) //this may be null by design if there is no associated PFTau (e.g. if DecayMode is constructed from MC truth)
00125    {
00126       const PFCandidateRefVector& theOutliers = originalTau->isolationPFCands();
00127       for(PFCandidateRefVector::const_iterator iIsoCand  = theOutliers.begin();
00128                                                iIsoCand != theOutliers.end();
00129                                              ++iIsoCand)
00130       {
00131          const PFCandidate* pfCand = iIsoCand->get();
00132          const Candidate* castedCand = static_cast<const Candidate*>(pfCand);
00133          if (castedCand)
00134             toFill.push_back(castedCand);
00135       }
00136    }
00137 }
00138 
00139 const reco::Candidate*
00140 PFTauDiscriminantManager::mainTrack() 
00141 {
00142    if (mainTrack_ == NULL) //otherwise already cached or d.n.e
00143    {
00144       if (!this->getDecayMode())
00145       {
00146          edm::LogError("PFTauDiscriminantManager") << "In ::mainTrack(), trying to access a null PFTauDecayMode - returning null pointer for main track";
00147          return NULL;
00148       }
00149 
00150       vector<const reco::Candidate*> myChargedCandidates = getDecayMode()->chargedPionCandidates();
00151       size_t nTracks = myChargedCandidates.size();
00152       if (!nTracks) 
00153       {
00154          // ...removing this warning for now, not sure what to do about this case (as it shoudl be permissible to pass a jet->pftau->pfTauDecayMode of all gammas??)
00155          //edm::LogError("PFTauDiscriminantManager") << "In ::mainTrack(), associated PFTauDecayMode has no associated tracks, returning null pointer.";
00156          return NULL;
00157       }
00158 
00159       //if there are more than three tracks, only take the top three, by Pt
00160       TauTagTools::sortByDescendingPt<reco::Candidate> ptSorter;
00161       sort(myChargedCandidates.begin(), myChargedCandidates.end(), ptSorter);
00162       size_t maxTracks = (nTracks > 3) ? 3 : nTracks;
00163       int    charge    = 0;
00164 
00165       if (maxTracks < 3) //two or one track, returning higher Pt track
00166          mainTrack_ = myChargedCandidates[0];
00167       else
00168       {
00169          for(size_t iTrack = 0; iTrack < maxTracks; ++iTrack)
00170             charge += myChargedCandidates[iTrack]->charge();
00171 
00172          for(size_t iTrack = 0; iTrack < maxTracks; ++iTrack)
00173          {
00174             int currentCharge = myChargedCandidates[iTrack]->charge();
00175             if (currentCharge != charge)
00176             {
00177                mainTrack_ = myChargedCandidates[iTrack];
00178                break;
00179             }
00180          }
00181       }
00182    }
00183    return mainTrack_;
00184 }
00185    
00186 
00187 
00188 candPtrVector
00189 PFTauDiscriminantManager::filterByCharge(const candPtrVector& input, bool isCharged) const
00190 {
00191    candPtrVector output;
00192    for(candPtrVector::const_iterator iCandidate  = input.begin();
00193                                      iCandidate != input.end();
00194                                    ++iCandidate)
00195    {
00196       bool chargeType = (*iCandidate)->charge();
00197       if( chargeType == isCharged ) 
00198          output.push_back(*iCandidate);
00199    }
00200    return output;
00201 }
00202 
00203 const std::vector<const reco::Candidate*>&
00204 PFTauDiscriminantManager::signalObjectsSortedByPt()
00205 {
00206    // return already computed vector if has already been computed or is empty (due to null tau)
00207    if(!signalObjectsSortedByPt_.empty() || iAmNull_)  
00208    {
00209       return signalObjectsSortedByPt_;
00210    }
00211    else
00212    {  
00213       TauTagTools::sortByAscendingPt<reco::Candidate> mySorter;
00214       fillSignalObjects(signalObjectsSortedByPt_);
00215       sort(signalObjectsSortedByPt_.begin(), signalObjectsSortedByPt_.end(), mySorter);
00216    }
00217    return signalObjectsSortedByPt_;
00218 }
00219 
00220 const std::vector<const reco::Candidate*>&
00221 PFTauDiscriminantManager::signalObjectsSortedByDR()
00222 {
00223    // return already computed vector if has already been computed or is empty (due to null tau)
00224    if(!signalObjectsSortedByDR_.empty() || iAmNull_)  
00225    {
00226       return signalObjectsSortedByDR_;
00227    }
00228    else
00229    {  
00230       if (currentTauDecayMode_ == NULL)
00231       {
00232          edm::LogError("PFTauDiscriminantManager") << "Trying to get signal objects from null PFTauDecayMode object!  Returning empty vector...";
00233          return signalObjectsSortedByDR_;
00234       }
00235       math::XYZVector signalAxisVector = currentTauDecayMode_->momentum();
00236       TauTagTools::sortByOpeningAngleAscending<reco::Candidate> mySorter(signalAxisVector, TauTagTools::computeDeltaR);
00237       fillSignalObjects(signalObjectsSortedByDR_);
00238       sort(signalObjectsSortedByDR_.begin(), signalObjectsSortedByDR_.end(), mySorter);
00239    }
00240    return signalObjectsSortedByDR_;
00241 }
00242 
00243 const std::vector<const reco::Candidate*>&
00244 PFTauDiscriminantManager::outlierObjectsSortedByPt()
00245 {
00246    if(!outlierObjectsSortedByPt_.empty() || iAmNull_)
00247    {
00248       return outlierObjectsSortedByPt_;
00249    }
00250    else
00251    {
00252       fillOutlierObjects(outlierObjectsSortedByPt_);
00253       TauTagTools::sortByAscendingPt<reco::Candidate> mySorter;
00254       sort(outlierObjectsSortedByPt_.begin(), outlierObjectsSortedByPt_.end(), mySorter);
00255    }
00256    return outlierObjectsSortedByPt_;
00257 }
00258 
00259 const std::vector<const reco::Candidate*>&
00260 PFTauDiscriminantManager::outlierObjectsSortedByDR()
00261 {
00262    if(!outlierObjectsSortedByDR_.empty() || iAmNull_)
00263    {
00264       return outlierObjectsSortedByDR_;
00265    }
00266    else
00267    {
00268       if (currentTauDecayMode_ == NULL)
00269       {
00270          edm::LogError("PFTauDiscriminantManager") << "Trying to get outlier objects from null PFTauDecayMode object!  Returning empty vector...";
00271          return outlierObjectsSortedByDR_;
00272       }
00273       math::XYZVector signalAxisVector = currentTauDecayMode_->momentum();
00274       fillOutlierObjects(outlierObjectsSortedByDR_);
00275       TauTagTools::sortByOpeningAngleAscending<reco::Candidate> mySorter(signalAxisVector, TauTagTools::computeDeltaR);
00276       sort(outlierObjectsSortedByDR_.begin(), outlierObjectsSortedByDR_.end(), mySorter);
00277    }
00278    return outlierObjectsSortedByDR_;
00279 }
00280 
00281 
00282 bool
00283 PFTauDiscriminantManager::branchTree(TTree* treeToBranch)
00284 {
00285    if(!treeToBranch)
00286    {
00287       edm::LogError("PFTauDiscriminantManager") << "Error: trying to branch ttree - TTree pointer is null!";
00288       return false;
00289    }
00290    //add magic variables _TARGET_ (for sig/bkg) and _WEIGHT_, and ISNULL for non-existence
00291    //treeToBranch->Branch("__TARGET__", &iAmSignal_,  "__TARGET__/O");  //needs bugfix in MVA framework code..
00292    treeToBranch->Branch("__ISNULL__", &iAmNull_,    "__ISNULL__/O");
00293    treeToBranch->Branch("__WEIGHT__", &eventWeight_,"__WEIGHT__/D");
00294    treeToBranch->Branch("__PREPASS__", &prePass_,"__PREPASS__/O");
00295    treeToBranch->Branch("__PREFAIL__", &preFail_,"__PREFAIL__/O");
00296 
00297    //loop over all the variables and make a branch for each one
00298    for(discriminantHolder::iterator iVariable  = myDiscriminants_.begin();
00299                                     iVariable != myDiscriminants_.end();
00300                                   ++iVariable)
00301    {
00302       Discriminant * theDiscriminant = iVariable->second;
00303       edm::LogInfo("PFTauDiscriminantManager") << "Branching for discriminant w/ name: " << theDiscriminant->name();
00304       theDiscriminant->branchTree(treeToBranch);
00305    }
00306    return true;
00307 }
00308 
00309 void
00310 PFTauDiscriminantManager::buildMVAComputerLink(std::vector<PhysicsTools::Variable::Value>& toFill)
00311 {
00312    for(discriminantHolder::iterator iVariable  = myDiscriminants_.begin();
00313                                     iVariable != myDiscriminants_.end();
00314                                   ++iVariable)
00315    {
00316       Discriminant * theDiscriminant = iVariable->second;
00317       theDiscriminant->fillMVA(toFill);
00318    }
00319 }
00320 
00321 
00322 
00323 PFTauDiscriminantManager::~PFTauDiscriminantManager()
00324 {
00325 }
00326 
00327 
00328 
00329 } //end namespace
00330 
00331 
00332 

Generated on Tue Jun 9 17:45:10 2009 for CMSSW by  doxygen 1.5.4