CMS 3D CMS Logo

PFRecoTauDecayModeDeterminator Class Reference

#include <RecoTauTag/RecoTau/interface/PFRecoTauDecayModeDeterminator.h>

Inheritance diagram for PFRecoTauDecayModeDeterminator:

edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Types

typedef list< CompositeCandidatecompCandList
typedef list
< CompositeCandidate >
::reverse_iterator 
compCandRevIter

Public Member Functions

void mergePiZeroes (compCandList &, compCandRevIter)
 PFRecoTauDecayModeDeterminator (const ParameterSet &iConfig)
virtual void produce (Event &, const EventSetup &)
 ~PFRecoTauDecayModeDeterminator ()

Protected Attributes

const double chargedPionMass
const double neutralPionMass

Private Attributes

AddFourMomenta addP4
TauTagTools::sortByAscendingPt
< CompositeCandidate
candAscendingSorter
TauTagTools::sortByDescendingPt
< CompositeCandidate
candDescendingSorter
bool filterPhotons_
bool filterTwoProngs_
uint32_t maxPhotonsToMerge_
double maxPiZeroMass_
bool mergeLowPtPhotonsFirst_
double minPtFractionForGammas_
double minPtFractionForSecondProng_
InputTag PFTauProducer_
bool refitTracks_
PFCandCommonVertexFitterBasevertexFitter_


Detailed Description

Definition at line 47 of file PFRecoTauDecayModeDeterminator.h.


Member Typedef Documentation

typedef list<CompositeCandidate> PFRecoTauDecayModeDeterminator::compCandList

Definition at line 50 of file PFRecoTauDecayModeDeterminator.h.

typedef list<CompositeCandidate>::reverse_iterator PFRecoTauDecayModeDeterminator::compCandRevIter

Definition at line 51 of file PFRecoTauDecayModeDeterminator.h.


Constructor & Destructor Documentation

PFRecoTauDecayModeDeterminator::PFRecoTauDecayModeDeterminator ( const ParameterSet iConfig  )  [explicit]

Definition at line 3 of file PFRecoTauDecayModeDeterminator.cc.

References filterPhotons_, filterTwoProngs_, edm::ParameterSet::getParameter(), maxPhotonsToMerge_, maxPiZeroMass_, mergeLowPtPhotonsFirst_, minPtFractionForGammas_, minPtFractionForSecondProng_, PFTauProducer_, refitTracks_, and vertexFitter_.

00003                                                                                          :chargedPionMass(0.13957),neutralPionMass(0.13497){
00004   PFTauProducer_                = iConfig.getParameter<InputTag>("PFTauProducer");
00005   maxPhotonsToMerge_            = iConfig.getParameter<uint32_t>("maxPhotonsToMerge");
00006   maxPiZeroMass_                = iConfig.getParameter<double>("maxPiZeroMass");             
00007   mergeLowPtPhotonsFirst_       = iConfig.getParameter<bool>("mergeLowPtPhotonsFirst");
00008   refitTracks_                  = iConfig.getParameter<bool>("refitTracks");
00009   filterTwoProngs_              = iConfig.getParameter<bool>("filterTwoProngs");
00010   filterPhotons_                = iConfig.getParameter<bool>("filterPhotons");
00011   minPtFractionForSecondProng_  = iConfig.getParameter<double>("minPtFractionForSecondProng");
00012   minPtFractionForGammas_       = iConfig.getParameter<double>("minPtFractionForThirdGamma");
00013   //setup vertex fitter
00014   vertexFitter_ = new PFCandCommonVertexFitter<KalmanVertexFitter>(iConfig);
00015   produces<PFTauDecayModeAssociation>();      
00016 }

PFRecoTauDecayModeDeterminator::~PFRecoTauDecayModeDeterminator (  ) 

Definition at line 18 of file PFRecoTauDecayModeDeterminator.cc.

00019 {
00020 //   delete vertexFitter_;  //now a very small memory leak, fix me later
00021 }


Member Function Documentation

void PFRecoTauDecayModeDeterminator::mergePiZeroes ( compCandList input,
compCandRevIter  seed 
)

Definition at line 27 of file PFRecoTauDecayModeDeterminator.cc.

References funct::abs(), addP4, pyDBSguiBaseClass::base, maxPhotonsToMerge_, maxPiZeroMass_, neutralPionMass, and AddFourMomenta::set().

Referenced by produce().

00028 {
00029    //uses std::list instead of vector, so that iterators can be deleted in situ
00030    //we go backwards for historical reasons ;)
00031    if(seed == input.rend())
00032       return;
00033    compCandRevIter bestMatchSoFar;
00034    LorentzVector combinationCandidate;
00035    float closestInvariantMassDifference = maxPiZeroMass_ + 1;
00036    bool foundACompatibleMatch = false;
00037    //find the best match to make a pi0
00038    compCandRevIter potentialMatch = seed;
00039    ++potentialMatch;
00040    for(; potentialMatch != input.rend(); ++potentialMatch)
00041    {
00042       // see how close this combination comes to the pion mass
00043       LorentzVector seedFourVector              = seed->p4();
00044       LorentzVector toAddFourVect               = potentialMatch->p4();
00045       combinationCandidate                      = seedFourVector + toAddFourVect;
00046       float combinationCandidateMass            = combinationCandidate.M();
00047       float differenceToTruePiZeroMass          = abs(combinationCandidateMass - neutralPionMass);
00048       if(combinationCandidateMass < maxPiZeroMass_ && differenceToTruePiZeroMass < closestInvariantMassDifference)
00049       {
00050          closestInvariantMassDifference = differenceToTruePiZeroMass;
00051          bestMatchSoFar = potentialMatch;
00052          foundACompatibleMatch = true;
00053       }
00054    }
00055    //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
00056    if(foundACompatibleMatch && seed->numberOfDaughters() < maxPhotonsToMerge_)
00057    {
00058       //combine match into Seed and update four vector
00059       if(bestMatchSoFar->numberOfDaughters() > 0)
00060       {
00061          const Candidate* photonToAdd = (*bestMatchSoFar).daughter(0);
00062          seed->addDaughter(*photonToAdd);
00063       }
00064       addP4.set(*seed);
00065       //remove match as it is now contained in the seed 
00066       input.erase( (++bestMatchSoFar).base() );  //convert to normal iterator, after correct for offset
00067       mergePiZeroes(input, seed);
00068    } else
00069    {
00070       // otherwise move to next highest object and recurse
00071       addP4.set(*seed);
00072       ++seed;
00073       mergePiZeroes(input, seed);
00074    }
00075 }

void PFRecoTauDecayModeDeterminator::produce ( Event iEvent,
const EventSetup iSetup 
) [virtual]

Implements edm::EDProducer.

Definition at line 77 of file PFRecoTauDecayModeDeterminator.cc.

References reco::CompositeCandidate::addDaughter(), addP4, edm::RefVector< C, T, F >::begin(), candAscendingSorter, candDescendingSorter, chargedPionMass, reco::CompositeCandidate::daughter(), edm::RefVector< C, T, F >::end(), edm::EventSetup::get(), edm::Event::getByLabel(), mergeLowPtPhotonsFirst_, mergePiZeroes(), neutralPionMass, reco::CompositeCandidate::numberOfDaughters(), reco::Particle::p4(), PFTauProducer_, edm::ESHandle< T >::product(), edm::Event::put(), refitTracks_, HLT_VtxMuL3::result, AddFourMomenta::set(), PFCandCommonVertexFitterBase::set(), reco::Particle::setMass(), reco::PFTauDecayMode::setPFTauRef(), reco::PFTau::signalPFChargedHadrCands(), reco::PFTau::signalPFGammaCands(), and vertexFitter_.

00077                                                                                   {
00078 
00079   ESHandle<TransientTrackBuilder> myTransientTrackBuilder;
00080   ESHandle<MagneticField> myMF;
00081 
00082   if (refitTracks_)
00083   {
00084      iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",myTransientTrackBuilder);
00085      iSetup.get<IdealMagneticFieldRecord>().get(myMF);
00086   }
00087 
00088   Handle<PFTauCollection> thePFTauCollection;
00089   iEvent.getByLabel(PFTauProducer_,thePFTauCollection);
00090 
00091   auto_ptr<PFTauDecayModeAssociation> result(new PFTauDecayModeAssociation(PFTauRefProd(thePFTauCollection)));
00092 
00093   size_t numberOfPFTaus = thePFTauCollection->size();
00094   for(size_t iPFTau = 0; iPFTau < numberOfPFTaus; ++iPFTau)
00095   {
00096      //get the reference to the PFTau
00097      PFTauRef           pfTauRef(thePFTauCollection, iPFTau);
00098      PFTau              myPFTau = *pfTauRef;
00099      
00100      //get the charged & neutral collections corresponding to this PFTau
00101      const PFCandidateRefVector& theChargedHadronCandidates = myPFTau.signalPFChargedHadrCands();
00102      const PFCandidateRefVector& theGammaCandidates         = myPFTau.signalPFGammaCands();
00103 
00104      LorentzVector totalFourVector;
00105 
00106      //shallow clone everything
00107      vector<ShallowCloneCandidate>    chargedCandidates;
00108      list<CompositeCandidate>         gammaCandidates;
00109      VertexCompositeCandidate chargedCandsToAdd;  //move me if wanted for filter
00110 
00111      for( PFCandidateRefVector::const_iterator iCharged  = theChargedHadronCandidates.begin();
00112                                                iCharged != theChargedHadronCandidates.end();
00113                                              ++iCharged)
00114      {
00115         // copy as shallow clone, and asssume mass of pi+
00116         chargedCandsToAdd.addDaughter(ShallowCloneCandidate(CandidateBaseRef(*iCharged)));
00117         Candidate* justAdded = chargedCandsToAdd.daughter(chargedCandsToAdd.numberOfDaughters()-1);
00118         totalFourVector += justAdded->p4();
00119         justAdded->setMass(chargedPionMass);
00120      }
00121      for( PFCandidateRefVector::const_iterator iGamma  = theGammaCandidates.begin();
00122                                                iGamma != theGammaCandidates.end();
00123                                              ++iGamma)
00124      {
00125         CompositeCandidate potentialPiZero;
00126         potentialPiZero.addDaughter(ShallowCloneCandidate(CandidateBaseRef(*iGamma)));
00127         addP4.set(potentialPiZero);
00128         totalFourVector += potentialPiZero.p4();
00129         gammaCandidates.push_back(potentialPiZero);
00130      }
00131 
00132      //sort the photons by pt before passing to merger
00133      if (mergeLowPtPhotonsFirst_)
00134         gammaCandidates.sort(candDescendingSorter);
00135      else
00136         gammaCandidates.sort(candAscendingSorter);
00137 
00138      mergePiZeroes(gammaCandidates, gammaCandidates.rbegin());
00139 
00140      CompositeCandidate mergedPiZerosToAdd;
00141      for( list<CompositeCandidate>::iterator iGamma  = gammaCandidates.begin();
00142                                              iGamma != gammaCandidates.end();
00143                                            ++iGamma)
00144      {
00145         iGamma->setMass(neutralPionMass);
00146         mergedPiZerosToAdd.addDaughter(*iGamma);
00147      }
00148      addP4.set(mergedPiZerosToAdd);
00149 
00150      CompositeCandidate filteredStuff; //empty for now.
00151 
00152      // apply vertex fitting.
00153      if (refitTracks_ && chargedCandsToAdd.numberOfDaughters() > 1)
00154      {
00155         vertexFitter_->set(myMF.product());
00156         vertexFitter_->set(chargedCandsToAdd);  //refits tracks, adds vertexing info
00157      }
00158 
00159      addP4.set(chargedCandsToAdd);
00160 
00161      /*
00162      LorentzVector refitFourVector = chargedCandsToAdd.p4() + mergedPiZerosToAdd.p4();
00163 
00164      edm::LogInfo("PFTauDecayModeDeterminator") << "Found nCharged: " << chargedCandsToAdd.numberOfDaughters()
00165                                   << " and nNeutral: " << mergedPiZerosToAdd.numberOfDaughters()
00166                                   << " Former mass: " << totalFourVector.mass() 
00167                                   << " New mass: " << refitFourVector.mass();
00168      */
00169 
00170      PFTauDecayMode myDecayModeTau(chargedCandsToAdd, mergedPiZerosToAdd, filteredStuff);
00171      myDecayModeTau.setPFTauRef(pfTauRef);
00172      result->setValue(iPFTau, myDecayModeTau);
00173   }
00174   iEvent.put(result);
00175 }


Member Data Documentation

AddFourMomenta PFRecoTauDecayModeDeterminator::addP4 [private]

Definition at line 66 of file PFRecoTauDecayModeDeterminator.h.

Referenced by mergePiZeroes(), and produce().

TauTagTools::sortByAscendingPt<CompositeCandidate> PFRecoTauDecayModeDeterminator::candAscendingSorter [private]

Definition at line 76 of file PFRecoTauDecayModeDeterminator.h.

Referenced by produce().

TauTagTools::sortByDescendingPt<CompositeCandidate> PFRecoTauDecayModeDeterminator::candDescendingSorter [private]

Definition at line 75 of file PFRecoTauDecayModeDeterminator.h.

Referenced by produce().

const double PFRecoTauDecayModeDeterminator::chargedPionMass [protected]

Definition at line 60 of file PFRecoTauDecayModeDeterminator.h.

Referenced by produce().

bool PFRecoTauDecayModeDeterminator::filterPhotons_ [private]

Definition at line 72 of file PFRecoTauDecayModeDeterminator.h.

Referenced by PFRecoTauDecayModeDeterminator().

bool PFRecoTauDecayModeDeterminator::filterTwoProngs_ [private]

Definition at line 71 of file PFRecoTauDecayModeDeterminator.h.

Referenced by PFRecoTauDecayModeDeterminator().

uint32_t PFRecoTauDecayModeDeterminator::maxPhotonsToMerge_ [private]

Definition at line 67 of file PFRecoTauDecayModeDeterminator.h.

Referenced by mergePiZeroes(), and PFRecoTauDecayModeDeterminator().

double PFRecoTauDecayModeDeterminator::maxPiZeroMass_ [private]

Definition at line 68 of file PFRecoTauDecayModeDeterminator.h.

Referenced by mergePiZeroes(), and PFRecoTauDecayModeDeterminator().

bool PFRecoTauDecayModeDeterminator::mergeLowPtPhotonsFirst_ [private]

Definition at line 69 of file PFRecoTauDecayModeDeterminator.h.

Referenced by PFRecoTauDecayModeDeterminator(), and produce().

double PFRecoTauDecayModeDeterminator::minPtFractionForGammas_ [private]

Definition at line 74 of file PFRecoTauDecayModeDeterminator.h.

Referenced by PFRecoTauDecayModeDeterminator().

double PFRecoTauDecayModeDeterminator::minPtFractionForSecondProng_ [private]

Definition at line 73 of file PFRecoTauDecayModeDeterminator.h.

Referenced by PFRecoTauDecayModeDeterminator().

const double PFRecoTauDecayModeDeterminator::neutralPionMass [protected]

Definition at line 61 of file PFRecoTauDecayModeDeterminator.h.

Referenced by mergePiZeroes(), and produce().

InputTag PFRecoTauDecayModeDeterminator::PFTauProducer_ [private]

Definition at line 65 of file PFRecoTauDecayModeDeterminator.h.

Referenced by PFRecoTauDecayModeDeterminator(), and produce().

bool PFRecoTauDecayModeDeterminator::refitTracks_ [private]

Definition at line 70 of file PFRecoTauDecayModeDeterminator.h.

Referenced by PFRecoTauDecayModeDeterminator(), and produce().

PFCandCommonVertexFitterBase* PFRecoTauDecayModeDeterminator::vertexFitter_ [private]

Definition at line 64 of file PFRecoTauDecayModeDeterminator.h.

Referenced by PFRecoTauDecayModeDeterminator(), and produce().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:29:45 2009 for CMSSW by  doxygen 1.5.4