CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/DataFormats/TauReco/src/PFTau.cc

Go to the documentation of this file.
00001 #include "DataFormats/TauReco/interface/PFTau.h"
00002 #include "DataFormats/Common/interface/RefToPtr.h"
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 
00005 //using namespace std;
00006 namespace reco {
00007 
00008 PFTau::PFTau() {
00009     leadPFChargedHadrCandsignedSipt_=NAN;
00010     isolationPFChargedHadrCandsPtSum_=NAN;
00011     isolationPFGammaCandsEtSum_=NAN;
00012     maximumHCALPFClusterEt_=NAN;
00013     emFraction_ = NAN;
00014     hcalTotOverPLead_ = NAN;
00015     hcalMaxOverPLead_ = NAN;
00016     hcal3x3OverPLead_ = NAN;
00017     ecalStripSumEOverPLead_= NAN;
00018     bremsRecoveryEOverPLead_ = NAN;
00019     electronPreIDOutput_ = NAN;
00020     electronPreIDDecision_= NAN;
00021     caloComp_ = NAN;
00022     segComp_ = NAN;
00023     muonDecision_ = NAN;
00024 }
00025 
00026 PFTau::PFTau(Charge q,const LorentzVector& p4,const Point& vtx) : BaseTau(q,p4,vtx)
00027 {
00028    leadPFChargedHadrCandsignedSipt_=NAN;
00029    isolationPFChargedHadrCandsPtSum_=NAN;
00030    isolationPFGammaCandsEtSum_=NAN;
00031    maximumHCALPFClusterEt_=NAN;
00032 
00033    emFraction_ = NAN;
00034    hcalTotOverPLead_ = NAN;
00035    hcalMaxOverPLead_ = NAN;
00036    hcal3x3OverPLead_ = NAN;
00037    ecalStripSumEOverPLead_= NAN;
00038    bremsRecoveryEOverPLead_ = NAN;
00039    electronPreIDOutput_ = NAN;
00040    electronPreIDDecision_= NAN;
00041 
00042    caloComp_ = NAN;
00043    segComp_ = NAN;
00044    muonDecision_ = NAN;
00045 }
00046 
00047 PFTau* PFTau::clone() const { return new PFTau(*this); }
00048 
00049 // Constituent getters and setters
00050 const PFJetRef& PFTau::jetRef() const {return jetRef_;}
00051 void PFTau::setjetRef(const PFJetRef& x) {jetRef_=x;}
00052 
00053 const PFTauTagInfoRef& PFTau::pfTauTagInfoRef() const {
00054 //  edm::LogWarning("DeprecatedPFTauMember")
00055 //      << "The PFTauTagInfoRef member is deprecated in the PFTau."
00056 //      << " For access to the underlying PFJet, please use the jetRef() method";
00057   return PFTauTagInfoRef_;
00058 }
00059 
00060 void PFTau::setpfTauTagInfoRef(const PFTauTagInfoRef x) {PFTauTagInfoRef_=x;}
00061 
00062 const PFCandidateRef& PFTau::leadPFChargedHadrCand() const {return leadPFChargedHadrCand_;}
00063 const PFCandidateRef& PFTau::leadPFNeutralCand() const {return leadPFNeutralCand_;}
00064 const PFCandidateRef& PFTau::leadPFCand() const {return leadPFCand_;}
00065 
00066 void PFTau::setleadPFChargedHadrCand(const PFCandidateRef& myLead) { leadPFChargedHadrCand_=myLead;}
00067 void PFTau::setleadPFNeutralCand(const PFCandidateRef& myLead) { leadPFNeutralCand_=myLead;}
00068 void PFTau::setleadPFCand(const PFCandidateRef& myLead) { leadPFCand_=myLead;}
00069 
00070 float PFTau::leadPFChargedHadrCandsignedSipt() const {return leadPFChargedHadrCandsignedSipt_;}
00071 void PFTau::setleadPFChargedHadrCandsignedSipt(const float& x){leadPFChargedHadrCandsignedSipt_=x;}
00072 
00073 const PFCandidateRefVector& PFTau::signalPFCands() const {return selectedSignalPFCands_;}
00074 void PFTau::setsignalPFCands(const PFCandidateRefVector& myParts)  { selectedSignalPFCands_ = myParts;}
00075 const PFCandidateRefVector& PFTau::signalPFChargedHadrCands() const {return selectedSignalPFChargedHadrCands_;}
00076 void PFTau::setsignalPFChargedHadrCands(const PFCandidateRefVector& myParts)  { selectedSignalPFChargedHadrCands_ = myParts;}
00077 const PFCandidateRefVector& PFTau::signalPFNeutrHadrCands() const {return selectedSignalPFNeutrHadrCands_;}
00078 void PFTau::setsignalPFNeutrHadrCands(const PFCandidateRefVector& myParts)  { selectedSignalPFNeutrHadrCands_ = myParts;}
00079 const PFCandidateRefVector& PFTau::signalPFGammaCands() const {return selectedSignalPFGammaCands_;}
00080 void PFTau::setsignalPFGammaCands(const PFCandidateRefVector& myParts)  { selectedSignalPFGammaCands_ = myParts;}
00081 
00082 const PFCandidateRefVector& PFTau::isolationPFCands() const {return selectedIsolationPFCands_;}
00083 void PFTau::setisolationPFCands(const PFCandidateRefVector& myParts)  { selectedIsolationPFCands_ = myParts;}
00084 const PFCandidateRefVector& PFTau::isolationPFChargedHadrCands() const {return selectedIsolationPFChargedHadrCands_;}
00085 void PFTau::setisolationPFChargedHadrCands(const PFCandidateRefVector& myParts)  { selectedIsolationPFChargedHadrCands_ = myParts;}
00086 const PFCandidateRefVector& PFTau::isolationPFNeutrHadrCands() const {return selectedIsolationPFNeutrHadrCands_;}
00087 void PFTau::setisolationPFNeutrHadrCands(const PFCandidateRefVector& myParts)  { selectedIsolationPFNeutrHadrCands_ = myParts;}
00088 const PFCandidateRefVector& PFTau::isolationPFGammaCands() const {return selectedIsolationPFGammaCands_;}
00089 void PFTau::setisolationPFGammaCands(const PFCandidateRefVector& myParts)  { selectedIsolationPFGammaCands_ = myParts;}
00090 
00091 // PiZero and decay mode information
00092 const std::vector<RecoTauPiZero>& PFTau::signalPiZeroCandidates() const {
00093    return signalPiZeroCandidates_;
00094 }
00095 void PFTau::setsignalPiZeroCandidates(const std::vector<RecoTauPiZero>& cands) {
00096    signalPiZeroCandidates_ = cands;
00097 }
00098 
00099 const std::vector<RecoTauPiZero>& PFTau::isolationPiZeroCandidates() const {
00100    return isolationPiZeroCandidates_;
00101 }
00102 void PFTau::setisolationPiZeroCandidates(const std::vector<RecoTauPiZero>& cands){
00103    signalPiZeroCandidates_ = cands;
00104 }
00105 
00106 PFTau::hadronicDecayMode PFTau::decayMode() const {
00107   unsigned int nCharged = signalPFChargedHadrCands().size();
00108   unsigned int nPiZeros = signalPiZeroCandidates().size();
00109   // If no tracks exist, this is definitely not a tau!
00110   if(!nCharged) return kNull;
00111   // Find the maximum number of PiZeros our parameterization can hold
00112   const unsigned int maxPiZeros = kOneProngNPiZero;
00113   // Determine our track index
00114   unsigned int trackIndex = (nCharged-1)*(maxPiZeros+1);
00115   // Check if we handle the given number of tracks
00116   if(trackIndex >= kRareDecayMode) return kRareDecayMode;
00117 
00118   nPiZeros = (nPiZeros <= maxPiZeros) ? nPiZeros : maxPiZeros;
00119   return static_cast<PFTau::hadronicDecayMode>(trackIndex + nPiZeros);
00120 }
00121 
00122 
00123 // Setting information about the isolation region
00124 float PFTau::isolationPFChargedHadrCandsPtSum() const {return isolationPFChargedHadrCandsPtSum_;}
00125 void PFTau::setisolationPFChargedHadrCandsPtSum(const float& x){isolationPFChargedHadrCandsPtSum_=x;}
00126 
00127 float PFTau::isolationPFGammaCandsEtSum() const {return isolationPFGammaCandsEtSum_;}
00128 void PFTau::setisolationPFGammaCandsEtSum(const float& x){isolationPFGammaCandsEtSum_=x;}
00129 
00130 float PFTau::maximumHCALPFClusterEt() const {return maximumHCALPFClusterEt_;}
00131 void PFTau::setmaximumHCALPFClusterEt(const float& x){maximumHCALPFClusterEt_=x;}
00132 
00133 // Electron variables
00134 float PFTau::emFraction() const {return emFraction_;}
00135 float PFTau::hcalTotOverPLead() const {return hcalTotOverPLead_;}
00136 float PFTau::hcalMaxOverPLead() const {return hcalMaxOverPLead_;}
00137 float PFTau::hcal3x3OverPLead() const {return hcal3x3OverPLead_;}
00138 float PFTau::ecalStripSumEOverPLead() const {return ecalStripSumEOverPLead_;}
00139 float PFTau::bremsRecoveryEOverPLead() const {return bremsRecoveryEOverPLead_;}
00140 reco::TrackRef PFTau::electronPreIDTrack() const {return electronPreIDTrack_;}
00141 float PFTau::electronPreIDOutput() const {return electronPreIDOutput_;}
00142 bool PFTau::electronPreIDDecision() const {return electronPreIDDecision_;}
00143 
00144 void PFTau::setemFraction(const float& x) {emFraction_ = x;}
00145 void PFTau::sethcalTotOverPLead(const float& x) {hcalTotOverPLead_ = x;}
00146 void PFTau::sethcalMaxOverPLead(const float& x) {hcalMaxOverPLead_ = x;}
00147 void PFTau::sethcal3x3OverPLead(const float& x) {hcal3x3OverPLead_ = x;}
00148 void PFTau::setecalStripSumEOverPLead(const float& x) {ecalStripSumEOverPLead_ = x;}
00149 void PFTau::setbremsRecoveryEOverPLead(const float& x) {bremsRecoveryEOverPLead_ = x;}
00150 void PFTau::setelectronPreIDTrack(const reco::TrackRef& x) {electronPreIDTrack_ = x;}
00151 void PFTau::setelectronPreIDOutput(const float& x) {electronPreIDOutput_ = x;}
00152 void PFTau::setelectronPreIDDecision(const bool& x) {electronPreIDDecision_ = x;}
00153 
00154 // Muon variables
00155 bool PFTau::hasMuonReference() const { // check if muon ref exists
00156     if( leadPFChargedHadrCand_.isNull() ) return false;
00157     else if( leadPFChargedHadrCand_.isNonnull() ){
00158         reco::MuonRef muonRef = leadPFChargedHadrCand_->muonRef();
00159         if( muonRef.isNull() )   return false;
00160         else if( muonRef.isNonnull() ) return  true;
00161     }
00162     return false;
00163 }
00164 
00165 float PFTau::caloComp() const {return caloComp_;}
00166 float PFTau::segComp() const {return segComp_;}
00167 bool  PFTau::muonDecision() const {return muonDecision_;}
00168 void PFTau::setCaloComp(const float& x) {caloComp_ = x;}
00169 void PFTau::setSegComp (const float& x) {segComp_  = x;}
00170 void PFTau::setMuonDecision(const bool& x) {muonDecision_ = x;}
00171 //
00172 
00173 
00174 CandidatePtr PFTau::sourceCandidatePtr( size_type i ) const {
00175     if( i!=0 ) return CandidatePtr();
00176     return  refToPtr( jetRef() );
00177 }
00178 
00179 
00180 bool PFTau::overlap(const Candidate& theCand) const {
00181     const RecoCandidate* theRecoCand=dynamic_cast<const RecoCandidate *>(&theCand);
00182     return (theRecoCand!=0 && (checkOverlap(track(),theRecoCand->track())));
00183 }
00184 
00185 void PFTau::dump(std::ostream& out) const {
00186 
00187     if(!out) return;
00188 
00189     if (pfTauTagInfoRef().isNonnull()) {
00190       out << "Its TauTagInfo constituents :"<<std::endl;
00191       out<<"# Tracks "<<pfTauTagInfoRef()->Tracks().size()<<std::endl;
00192       out<<"# PF charged hadr. cand's "<<pfTauTagInfoRef()->PFChargedHadrCands().size()<<std::endl;
00193       out<<"# PF neutral hadr. cand's "<<pfTauTagInfoRef()->PFNeutrHadrCands().size()<<std::endl;
00194       out<<"# PF gamma cand's "<<pfTauTagInfoRef()->PFGammaCands().size()<<std::endl;
00195     }
00196     if (jetRef().isNonnull()) {
00197       out << "Its constituents :"<< std::endl;
00198       out<<"# PF charged hadr. cand's "<< jetRef()->chargedHadronMultiplicity()<<std::endl;
00199       out<<"# PF neutral hadr. cand's "<< jetRef()->neutralHadronMultiplicity()<<std::endl;
00200       out<<"# PF gamma cand's "<< jetRef()->photonMultiplicity()<<std::endl;
00201       out<<"# Electron cand's "<< jetRef()->electronMultiplicity()<<std::endl;
00202     }
00203     out<<"in detail :"<<std::endl;
00204 
00205     out<<"Pt of the PFTau "<<pt()<<std::endl;
00206     PFCandidateRef theLeadPFCand = leadPFChargedHadrCand();
00207     if(!theLeadPFCand){
00208         out<<"No Lead PFCand "<<std::endl;
00209     }else{
00210         out<<"Lead PFCand Particle Id " << (*theLeadPFCand).particleId() << std::endl;
00211         out<<"Lead PFCand Pt "<<(*theLeadPFCand).pt()<<std::endl;
00212         out<<"Lead PFCand Charge "<<(*theLeadPFCand).charge()<<std::endl;
00213         out<<"Lead PFCand TrkRef "<<(*theLeadPFCand).trackRef().isNonnull()<<std::endl;
00214         out<<"Inner point position (x,y,z) of the PFTau ("<<vx()<<","<<vy()<<","<<vz()<<")"<<std::endl;
00215         out<<"Charge of the PFTau "<<charge()<<std::endl;
00216         out<<"Et of the highest Et HCAL PFCluster "<<maximumHCALPFClusterEt()<<std::endl;
00217         out<<"Number of SignalPFChargedHadrCands = "<<signalPFChargedHadrCands().size()<<std::endl;
00218         out<<"Number of SignalPFGammaCands = "<<signalPFGammaCands().size()<<std::endl;
00219         out<<"Number of IsolationPFChargedHadrCands = "<<isolationPFChargedHadrCands().size()<<std::endl;
00220         out<<"Number of IsolationPFGammaCands = "<<isolationPFGammaCands().size()<<std::endl;
00221         out<<"Sum of Pt of charged hadr. PFCandidates in isolation annulus around Lead PF = "<<isolationPFChargedHadrCandsPtSum()<<std::endl;
00222         out<<"Sum of Et of gamma PFCandidates in other isolation annulus around Lead PF = "<<isolationPFGammaCandsEtSum()<<std::endl;
00223 
00224     }
00225     // return out;
00226 }
00227 
00228 std::ostream& operator<<(std::ostream& out, const reco::PFTau& tau) {
00229 
00230    if(!out) return out;
00231 
00232    out << std::setprecision(3)
00233      <<"PFTau "
00234      << " charge: " << tau.charge() << " "
00235      << " pt:" <<tau.pt()<<" "
00236      << " eta:" <<tau.eta()<<" "
00237      << " phi:" <<tau.phi()<<" "
00238      << " mass:" << tau.mass() << " "
00239      << " dm: " << tau.decayMode() << " "
00240      <<tau.signalPFCands().size()<<","
00241      <<tau.signalPFChargedHadrCands().size()<<","
00242      <<tau.signalPFGammaCands().size()<<","
00243      <<tau.signalPiZeroCandidates().size()<<","
00244      <<tau.signalPFNeutrHadrCands().size()<<"  "
00245 
00246      <<tau.isolationPFCands().size()<<","
00247      <<tau.isolationPFChargedHadrCands().size()<<","
00248      <<tau.isolationPFGammaCands().size()<<","
00249      <<tau.isolationPiZeroCandidates().size()<<","
00250      <<tau.isolationPFNeutrHadrCands().size();
00251 
00252    return out;
00253 }
00254 
00255 }