CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoTauTag/RecoTau/src/HPSPFRecoTauAlgorithm.cc

Go to the documentation of this file.
00001 #include "RecoTauTag/RecoTau/interface/HPSPFRecoTauAlgorithm.h"
00002 #include "Math/GenVector/VectorUtil.h"
00003 using namespace reco;
00004 
00005 HPSPFRecoTauAlgorithm::HPSPFRecoTauAlgorithm():
00006     PFRecoTauAlgorithmBase()
00007 {
00008 }
00009 
00010 HPSPFRecoTauAlgorithm::HPSPFRecoTauAlgorithm(const edm::ParameterSet& config):
00011     PFRecoTauAlgorithmBase(config)
00012 {
00013   configure(config);
00014 }
00015 
00016 HPSPFRecoTauAlgorithm::~HPSPFRecoTauAlgorithm()
00017 {
00018   if(candidateMerger_ !=0 ) delete candidateMerger_;
00019 }
00020 
00021 PFTau
00022 HPSPFRecoTauAlgorithm::buildPFTau(const PFTauTagInfoRef& tagInfo,const Vertex& vertex)
00023 {
00024   PFTau pfTau;
00025 
00026   pfTaus_.clear();
00027   //make the strips globally.
00028   std::vector<PFCandidateRefVector> strips = candidateMerger_->mergeCandidates(tagInfo->PFCands());
00029 
00030 
00031   //Sort the hadrons globally and once!
00032   PFCandidateRefVector hadrons = tagInfo->PFChargedHadrCands();
00033   if(hadrons.size()>1)
00034     TauTagTools::sortRefVectorByPt(hadrons);
00035 
00036 
00037   //OK For this Tau Tag Info we should create all the possible taus
00038 
00039   //One Prongs
00040   if(doOneProngs_)
00041     buildOneProng(tagInfo,hadrons);
00042 
00043   //One Prong Strips
00044   if(doOneProngStrips_)
00045     buildOneProngStrip(tagInfo,strips,hadrons);
00046 
00047   //One Prong TwoStrips
00048   if(doOneProngTwoStrips_)
00049     buildOneProngTwoStrips(tagInfo,strips,hadrons);
00050 
00051   //Three Prong
00052   if(doThreeProngs_)
00053     buildThreeProngs(tagInfo,hadrons);
00054 
00055 
00056   //Lets see if we created any taus
00057   if(pfTaus_.size()>0) {
00058 
00059     pfTau = getBestTauCandidate(pfTaus_);
00060 
00061     //Set the IP for the leading track
00062     if(TransientTrackBuilder_!=0 &&pfTau.leadPFChargedHadrCand()->trackRef().isNonnull()) {
00063       const TransientTrack leadTrack=TransientTrackBuilder_->build(pfTau.leadPFChargedHadrCand()->trackRef());
00064       if(pfTau.pfTauTagInfoRef().isNonnull())
00065         if(pfTau.pfTauTagInfoRef()->pfjetRef().isNonnull()) {
00066           PFJetRef jet = pfTau.pfTauTagInfoRef()->pfjetRef();
00067           GlobalVector jetDir(jet->px(),jet->py(),jet->pz());
00068           if(IPTools::signedTransverseImpactParameter(leadTrack,jetDir,vertex).first)
00069             pfTau.setleadPFChargedHadrCandsignedSipt(IPTools::signedTransverseImpactParameter(leadTrack,jetDir,vertex).second.significance());
00070         }
00071     }
00072   }
00073   else {      //null PFTau
00074     //Simone asked that in case there is no tau returned make a tau
00075     //without refs and the LV of the jet
00076     pfTau.setpfTauTagInfoRef(tagInfo);
00077     pfTau.setP4(tagInfo->pfjetRef()->p4());
00078   }
00079 
00080   return pfTau;
00081 }
00082 
00083 
00084 
00085 //Create one prong tau
00086 void
00087 HPSPFRecoTauAlgorithm::buildOneProng(const reco::PFTauTagInfoRef& tagInfo,const reco::PFCandidateRefVector& hadrons)
00088 {
00089   PFTauCollection  taus;
00090 
00091   if(hadrons.size()>0)
00092     for(unsigned int h=0;h<hadrons.size();++h) {
00093       PFCandidateRef hadron = hadrons.at(h);
00094 
00095       //In the one prong case the lead Track pt should be above the tau Threshold!
00096       //since all the tau is just one track!
00097       if(hadron->pt()>tauThreshold_)
00098         if(hadron->pt()>leadPionThreshold_)
00099           //The track should be within the matching cone
00100           if(ROOT::Math::VectorUtil::DeltaR(hadron->p4(),tagInfo->pfjetRef()->p4())<matchingCone_) {
00101             //OK Lets create a Particle Flow Tau!
00102             PFTau tau = PFTau(hadron->charge(),hadron->p4(),hadron->vertex());
00103 
00104             //Associate the Tag Info to the tau
00105             tau.setpfTauTagInfoRef(tagInfo);
00106 
00107             //Put the Hadron in the signal Constituents
00108             PFCandidateRefVector signal;
00109             signal.push_back(hadron);
00110 
00111             //Set The signal candidates of the PF Tau
00112             tau.setsignalPFChargedHadrCands(signal);
00113             tau.setsignalPFCands(signal);
00114             tau.setleadPFChargedHadrCand(hadron);
00115             tau.setleadPFCand(hadron);
00116 
00117             //Fill isolation variables
00118             associateIsolationCandidates(tau,0.0);
00119 
00120             //Apply Muon rejection algorithms
00121             applyMuonRejection(tau);
00122             applyElectronRejection(tau,0.0);
00123 
00124             //Save this candidate
00125             taus.push_back(tau);
00126           }
00127     }
00128   if(taus.size()>0) {
00129     pfTaus_.push_back(getBestTauCandidate(taus));
00130   }
00131 
00132 }
00133 
00134 //Build one Prong + Strip
00135 
00136 void
00137 HPSPFRecoTauAlgorithm::buildOneProngStrip(const reco::PFTauTagInfoRef& tagInfo,const  std::vector<PFCandidateRefVector>& strips,const reco::PFCandidateRefVector & hadrons)
00138 
00139 {
00140   //Create output Collection
00141   PFTauCollection taus;
00142 
00143 
00144 
00145 
00146   //make taus like this only if there is at least one hadron+ 1 strip
00147   if(hadrons.size()>0&&strips.size()>0){
00148     //Combinatorics between strips and clusters
00149     for(std::vector<PFCandidateRefVector>::const_iterator candVector=strips.begin();candVector!=strips.end();++candVector)
00150       for(PFCandidateRefVector::const_iterator hadron=hadrons.begin();hadron!=hadrons.end();++hadron) {
00151 
00152         //First Cross cleaning ! If you asked to clusterize the candidates
00153         //with tracks too then you should not double count the track
00154         PFCandidateRefVector emConstituents = *candVector;
00155         removeCandidateFromRefVector(*hadron,emConstituents);
00156 
00157         //Create a LorentzVector for the strip
00158         math::XYZTLorentzVector strip = createMergedLorentzVector(emConstituents);
00159 
00160         //TEST: Apply Strip Constraint
00161         applyMassConstraint(strip,0.1349);
00162 
00163         //create the Particle Flow Tau: Hadron plus Strip
00164         PFTau tau((*hadron)->charge(),
00165                   (*hadron)->p4()+strip,
00166                   (*hadron)->vertex());
00167 
00168         //Check tau threshold,  mass, Matching Cone window
00169         if(tau.pt()>tauThreshold_&&strip.pt()>stripPtThreshold_)
00170           if(tau.mass()>oneProngStripMassWindow_.at(0)&&tau.mass()<oneProngStripMassWindow_.at(1))//Apply mass window
00171             if(ROOT::Math::VectorUtil::DeltaR(tau.p4(),tagInfo->pfjetRef()->p4())<matchingCone_) { //Apply matching cone
00172               //Set The Tag Infor ref
00173               tau.setpfTauTagInfoRef(tagInfo);
00174 
00175               //Create the signal vectors
00176               PFCandidateRefVector signal;
00177               PFCandidateRefVector signalH;
00178               PFCandidateRefVector signalG;
00179 
00180               //Store the hadron in the PFTau
00181               signalH.push_back(*hadron);
00182               signal.push_back(*hadron);
00183 
00184               //calculate the cone size : For the strip use it as one candidate !
00185               double tauCone=0.0;
00186               if(coneMetric_ =="angle")
00187                 tauCone=std::max(fabs(ROOT::Math::VectorUtil::Angle(tau.p4(),(*hadron)->p4())),
00188                                  fabs(ROOT::Math::VectorUtil::Angle(tau.p4(),strip)));
00189               else if(coneMetric_ == "DR")
00190                 tauCone=std::max(ROOT::Math::VectorUtil::DeltaR(tau.p4(),(*hadron)->p4()),
00191                                  ROOT::Math::VectorUtil::DeltaR(tau.p4(),strip));
00192 
00193               if(emConstituents.size()>0)
00194                 for(PFCandidateRefVector::const_iterator j=emConstituents.begin();j!=emConstituents.end();++j)  {
00195                   signal.push_back(*j);
00196                   signalG.push_back(*j);
00197                 }
00198 
00199               //Set the PFTau
00200               tau.setsignalPFChargedHadrCands(signalH);
00201               tau.setsignalPFGammaCands(signalG);
00202               tau.setsignalPFCands(signal);
00203               tau.setleadPFChargedHadrCand(*hadron);
00204               tau.setleadPFNeutralCand(emConstituents.at(0));
00205 
00206               //Set the lead Candidate->Can be the hadron or the leading PFGamma(When we clear the Dataformat we will put the strip)
00207               if((*hadron)->pt()>emConstituents.at(0)->pt())
00208                 tau.setleadPFCand(*hadron);
00209               else
00210                 tau.setleadPFCand(emConstituents.at(0));
00211 
00212               //Apply the signal cone size formula
00213               if(isNarrowTau(tau,tauCone)) {
00214                 //calculate the isolation Deposits
00215                 associateIsolationCandidates(tau,tauCone);
00216                 //Set Muon Rejection
00217                 applyMuonRejection(tau);
00218                 applyElectronRejection(tau,strip.energy());
00219 
00220                 taus.push_back(tau);
00221               }
00222             }
00223       }
00224   }
00225 
00226   if(taus.size()>0) {
00227     pfTaus_.push_back(getBestTauCandidate(taus));
00228   }
00229 
00230 }
00231 
00232 void
00233 HPSPFRecoTauAlgorithm::buildOneProngTwoStrips(const reco::PFTauTagInfoRef& tagInfo,const  std::vector<PFCandidateRefVector>& strips,const reco::PFCandidateRefVector& hadrons)
00234 {
00235 
00236 
00237   PFTauCollection taus;
00238 
00239   //make taus like this only if there is at least one hadron+ 2 strips
00240   if(hadrons.size()>0&&strips.size()>1){
00241     //Combinatorics between strips and clusters
00242     for(unsigned int Nstrip1=0;Nstrip1<strips.size()-1;++Nstrip1)
00243       for(unsigned int Nstrip2=Nstrip1+1;Nstrip2<strips.size();++Nstrip2)
00244         for(PFCandidateRefVector::const_iterator hadron=hadrons.begin();hadron!=hadrons.end();++hadron) {
00245 
00246 
00247 
00248           //Create the strips and the vectors .Again cross clean the track if associated
00249           PFCandidateRefVector emConstituents1 = strips.at(Nstrip1);
00250           PFCandidateRefVector emConstituents2 = strips.at(Nstrip2);
00251           removeCandidateFromRefVector(*hadron,emConstituents1);
00252           removeCandidateFromRefVector(*hadron,emConstituents2);
00253 
00254 
00255           //Create a LorentzVector for the strip
00256           math::XYZTLorentzVector strip1 = createMergedLorentzVector(emConstituents1);
00257           math::XYZTLorentzVector strip2 = createMergedLorentzVector(emConstituents2);
00258 
00259 
00260 
00261           //Apply Mass Constraints
00262           applyMassConstraint(strip1,0.0);
00263           applyMassConstraint(strip2,0.0);
00264 
00265 
00266           PFTau tau((*hadron)->charge(),
00267                     (*hadron)->p4()+strip1+strip2,
00268                     (*hadron)->vertex());
00269 
00270 
00271           if(tau.pt()>tauThreshold_&&strip1.pt()>stripPtThreshold_&&strip2.pt()>stripPtThreshold_)
00272             if((strip1+strip2).M() >oneProngTwoStripsPi0MassWindow_.at(0) &&(strip1+strip2).M() <oneProngTwoStripsPi0MassWindow_.at(1) )//pi0 conmstraint for two strips
00273               if(tau.mass()>oneProngTwoStripsMassWindow_.at(0)&&tau.mass()<oneProngTwoStripsMassWindow_.at(1))//Apply mass window
00274                 if(ROOT::Math::VectorUtil::DeltaR(tau.p4(),tagInfo->pfjetRef()->p4())<matchingCone_) { //Apply matching cone
00275                   //create the PFTau
00276                   tau.setpfTauTagInfoRef(tagInfo);
00277 
00278 
00279                   //Create the signal vectors
00280                   PFCandidateRefVector signal;
00281                   PFCandidateRefVector signalH;
00282                   PFCandidateRefVector signalG;
00283 
00284                   //Store the hadron in the PFTau
00285                   signalH.push_back(*hadron);
00286                   signal.push_back(*hadron);
00287 
00288                   //calculate the cone size from the reconstructed Objects
00289                   double tauCone=1000.0;
00290                   if(coneMetric_ =="angle") {
00291                     tauCone=std::max(std::max(fabs(ROOT::Math::VectorUtil::Angle(tau.p4(),(*hadron)->p4())),
00292                                               fabs(ROOT::Math::VectorUtil::Angle(tau.p4(),strip1))),
00293                                               fabs(ROOT::Math::VectorUtil::Angle(tau.p4(),strip2)));
00294                   }
00295                   else if(coneMetric_ =="DR") {
00296                     tauCone=std::max(std::max((ROOT::Math::VectorUtil::DeltaR(tau.p4(),(*hadron)->p4())),
00297                                               (ROOT::Math::VectorUtil::DeltaR(tau.p4(),strip1))),
00298                                               (ROOT::Math::VectorUtil::DeltaR(tau.p4(),strip2)));
00299 
00300                   }
00301 
00302                   for(PFCandidateRefVector::const_iterator j=emConstituents1.begin();j!=emConstituents1.end();++j)  {
00303                     signal.push_back(*j);
00304                     signalG.push_back(*j);
00305                   }
00306 
00307                   for(PFCandidateRefVector::const_iterator j=emConstituents2.begin();j!=emConstituents2.end();++j)  {
00308                     signal.push_back(*j);
00309                     signalG.push_back(*j);
00310                   }
00311 
00312                   //Set the PFTau
00313                   tau.setsignalPFChargedHadrCands(signalH);
00314                   tau.setsignalPFGammaCands(signalG);
00315                   tau.setsignalPFCands(signal);
00316                   tau.setleadPFChargedHadrCand(*hadron);
00317 
00318                   //Set the lead Candidate->Can be the hadron or the leading PFGamma(When we clear the Dataformat we will put the strip)
00319                   if((*hadron)->pt()>emConstituents1.at(0)->pt())
00320                     tau.setleadPFCand(*hadron);
00321                   else
00322                     tau.setleadPFCand(emConstituents1.at(0));
00323 
00324                   //Apply the cone size formula
00325                   if(isNarrowTau(tau,tauCone)) {
00326 
00327                     //calculate the isolation Deposits
00328                     associateIsolationCandidates(tau,tauCone);
00329 
00330                     applyMuonRejection(tau);
00331 
00332                     //For two strips take the nearest strip to the track
00333                     if(ROOT::Math::VectorUtil::DeltaR(strip1,(*hadron)->p4())<
00334                        ROOT::Math::VectorUtil::DeltaR(strip2,(*hadron)->p4()))
00335                       applyElectronRejection(tau,strip1.energy());
00336                     else
00337                       applyElectronRejection(tau,strip2.energy());
00338 
00339                     taus.push_back(tau);
00340                   }
00341                 }
00342         }
00343   }
00344 
00345   if(taus.size()>0) {
00346     pfTaus_.push_back(getBestTauCandidate(taus));
00347   }
00348 }
00349 
00350 
00351 
00352 void
00353 HPSPFRecoTauAlgorithm::buildThreeProngs(const reco::PFTauTagInfoRef& tagInfo,const reco::PFCandidateRefVector& hadrons)
00354 {
00355   PFTauCollection taus;
00356   //get Hadrons
00357 
00358 
00359   //Require at least three hadrons
00360   if(hadrons.size()>2)
00361     for(unsigned int a=0;a<hadrons.size()-2;++a)
00362       for(unsigned int b=a+1;b<hadrons.size()-1;++b)
00363         for(unsigned int c=b+1;c<hadrons.size();++c) {
00364 
00365           PFCandidateRef h1 = hadrons.at(a);
00366           PFCandidateRef h2 = hadrons.at(b);
00367           PFCandidateRef h3 = hadrons.at(c);
00368 
00369           //check charge Compatibility and lead track
00370           int charge=h1->charge()+h2->charge()+h3->charge();
00371           if(std::abs(charge)==1 && h1->pt()>leadPionThreshold_)
00372             //check the track refs
00373             if(h1->trackRef()!=h2->trackRef()&&h1->trackRef()!=h3->trackRef()&&h2->trackRef()!=h3->trackRef())
00374             {
00375 
00376               //create the tau
00377               PFTau tau = PFTau(charge,h1->p4()+h2->p4()+h3->p4(),h1->vertex());
00378               tau.setpfTauTagInfoRef(tagInfo);
00379 
00380               if(tau.pt()>tauThreshold_)//Threshold
00381                 if(ROOT::Math::VectorUtil::DeltaR(tau.p4(),tagInfo->pfjetRef()->p4())<matchingCone_) {//Matching Cone
00382 
00383                   PFCandidateRefVector signal;
00384                   signal.push_back(h1);
00385                   signal.push_back(h2);
00386                   signal.push_back(h3);
00387                   //calculate the tau cone by getting the maximum distance
00388                   double tauCone = 10000.0;
00389                   if(coneMetric_=="DR")
00390                   {
00391                     tauCone = std::max(ROOT::Math::VectorUtil::DeltaR(tau.p4(),h1->p4()),
00392                               std::max(ROOT::Math::VectorUtil::DeltaR(tau.p4(),h2->p4()),
00393                                        ROOT::Math::VectorUtil::DeltaR(tau.p4(),h3->p4())));
00394                   }
00395                   else if(coneMetric_=="angle")
00396                   {
00397                     tauCone =std::max(fabs(ROOT::Math::VectorUtil::Angle(tau.p4(),h1->p4())),
00398                              std::max(fabs(ROOT::Math::VectorUtil::Angle(tau.p4(),h2->p4())),
00399                                       fabs(ROOT::Math::VectorUtil::Angle(tau.p4(),h3->p4()))));
00400                   }
00401 
00402                   //Set The PFTau
00403                   tau.setsignalPFChargedHadrCands(signal);
00404                   tau.setsignalPFCands(signal);
00405                   tau.setleadPFChargedHadrCand(h1);
00406                   tau.setleadPFCand(h1);
00407 
00408                   if(isNarrowTau(tau,tauCone)) {
00409                     //calculate the isolation Deposits
00410                     associateIsolationCandidates(tau,tauCone);
00411                     applyMuonRejection(tau);
00412                     applyElectronRejection(tau,0.0);
00413                     taus.push_back(tau);
00414 
00415                   }
00416                 }
00417             }
00418         }
00419 
00420   if(taus.size()>0) {
00421     PFTau bestTau  = getBestTauCandidate(taus);
00422     if(refitThreeProng(bestTau))
00423       //Apply mass constraint
00424       if(bestTau.mass()>threeProngMassWindow_.at(0)&&bestTau.mass()<threeProngMassWindow_.at(1))//MassWindow
00425         pfTaus_.push_back(bestTau);
00426   }
00427 
00428 }
00429 
00430 
00431 bool
00432 HPSPFRecoTauAlgorithm::isNarrowTau(const reco::PFTau& tau,double cone)
00433 {
00434   double allowedConeSize=coneSizeFormula.Eval(tau.energy(),tau.et());
00435   if (allowedConeSize<minSignalCone_) allowedConeSize=minSignalCone_;
00436   if (allowedConeSize>maxSignalCone_) allowedConeSize=maxSignalCone_;
00437 
00438   if(cone<allowedConeSize)
00439     return true;
00440   else
00441     return false;
00442 }
00443 
00444 
00445 void
00446 HPSPFRecoTauAlgorithm::associateIsolationCandidates(reco::PFTau& tau,
00447                                                     double tauCone)
00448 {
00449 
00450 
00451   using namespace reco;
00452 
00453   //Information to get filled
00454   double sumPT=0;
00455   double sumET=0;
00456 
00457   if(tau.pfTauTagInfoRef().isNull()) return;
00458 
00459   PFCandidateRefVector hadrons;
00460   PFCandidateRefVector gammas;
00461   PFCandidateRefVector neutral;
00462 
00463   //Remove candidates outside the cones
00464   if(useIsolationAnnulus_)
00465   {
00466 
00467     for(unsigned int i=0;i<tau.pfTauTagInfoRef()->PFChargedHadrCands().size();++i) {
00468       double DR = ROOT::Math::VectorUtil::DeltaR(tau.p4(),tau.pfTauTagInfoRef()->PFChargedHadrCands().at(i)->p4());
00469 
00470       if(DR>tauCone && DR<chargeIsolationCone_)
00471         hadrons.push_back(tau.pfTauTagInfoRef()->PFChargedHadrCands().at(i));
00472     }
00473 
00474     for(unsigned int i=0;i<tau.pfTauTagInfoRef()->PFGammaCands().size();++i) {
00475       double DR = ROOT::Math::VectorUtil::DeltaR(tau.p4(),tau.pfTauTagInfoRef()->PFGammaCands().at(i)->p4());
00476 
00477       if(DR>tauCone && DR<gammaIsolationCone_)
00478         gammas.push_back(tau.pfTauTagInfoRef()->PFGammaCands().at(i));
00479     }
00480     for(unsigned int i=0;i<tau.pfTauTagInfoRef()->PFNeutrHadrCands().size();++i) {
00481       double DR = ROOT::Math::VectorUtil::DeltaR(tau.p4(),tau.pfTauTagInfoRef()->PFNeutrHadrCands().at(i)->p4());
00482       if(DR>tauCone && DR <neutrHadrIsolationCone_)
00483         neutral.push_back(tau.pfTauTagInfoRef()->PFNeutrHadrCands().at(i));
00484     }
00485   }
00486   else
00487   {
00488 
00489     for(unsigned int i=0;i<tau.pfTauTagInfoRef()->PFChargedHadrCands().size();++i) {
00490       double DR = ROOT::Math::VectorUtil::DeltaR(tau.p4(),tau.pfTauTagInfoRef()->PFChargedHadrCands().at(i)->p4());
00491 
00492       if(DR<chargeIsolationCone_)
00493         hadrons.push_back(tau.pfTauTagInfoRef()->PFChargedHadrCands().at(i));
00494     }
00495 
00496     for(unsigned int i=0;i<tau.pfTauTagInfoRef()->PFGammaCands().size();++i) {
00497       double DR = ROOT::Math::VectorUtil::DeltaR(tau.p4(),tau.pfTauTagInfoRef()->PFGammaCands().at(i)->p4());
00498 
00499       if(DR<gammaIsolationCone_)
00500         gammas.push_back(tau.pfTauTagInfoRef()->PFGammaCands().at(i));
00501     }
00502     for(unsigned int i=0;i<tau.pfTauTagInfoRef()->PFNeutrHadrCands().size();++i) {
00503       double DR = ROOT::Math::VectorUtil::DeltaR(tau.p4(),tau.pfTauTagInfoRef()->PFNeutrHadrCands().at(i)->p4());
00504       if(DR <neutrHadrIsolationCone_)
00505         neutral.push_back(tau.pfTauTagInfoRef()->PFNeutrHadrCands().at(i));
00506     }
00507 
00508   }
00509 
00510   //remove the signal Constituents from the collections
00511   for(PFCandidateRefVector::const_iterator i=tau.signalPFChargedHadrCands().begin();i!=tau.signalPFChargedHadrCands().end();++i)
00512   {
00513     removeCandidateFromRefVector(*i,hadrons);
00514   }
00515 
00516   for(PFCandidateRefVector::const_iterator i=tau.signalPFGammaCands().begin();i!=tau.signalPFGammaCands().end();++i)
00517   {
00518     removeCandidateFromRefVector(*i,gammas);
00519     removeCandidateFromRefVector(*i,hadrons);//special case where we included a hadron if the strip!
00520   }
00521 
00522 
00523   //calculate isolation deposits
00524   for(unsigned int i=0;i<hadrons.size();++i)
00525   {
00526     sumPT+=hadrons.at(i)->pt();
00527   }
00528 
00529   for(unsigned int i=0;i<gammas.size();++i)
00530   {
00531     sumET+=gammas.at(i)->pt();
00532   }
00533 
00534 
00535   tau.setisolationPFChargedHadrCandsPtSum(sumPT);
00536   tau.setisolationPFGammaCandsEtSum(sumET);
00537   tau.setisolationPFChargedHadrCands(hadrons);
00538   tau.setisolationPFNeutrHadrCands(neutral);
00539   tau.setisolationPFGammaCands(gammas);
00540 
00541   PFCandidateRefVector isoAll = hadrons;
00542   for(unsigned int i=0;i<gammas.size();++i)
00543     isoAll.push_back(gammas.at(i));
00544 
00545   for(unsigned int i=0;i<neutral.size();++i)
00546     isoAll.push_back(neutral.at(i));
00547 
00548   tau.setisolationPFCands(isoAll);
00549 }
00550 
00551 void
00552 HPSPFRecoTauAlgorithm::applyMuonRejection(reco::PFTau& tau)
00553 {
00554 
00555   // Require that no signal track has segment matches
00556 
00557   //Also:
00558   //The segment compatibility is the number of matched Muon Segments
00559   //the old available does not exist in the muons anymore so i will fill the data format with that
00560   bool decision=true;
00561   float caloComp=0.0;
00562   float segComp=0.0;
00563 
00564   if(tau.leadPFChargedHadrCand().isNonnull()) {
00565     MuonRef mu =tau.leadPFChargedHadrCand()->muonRef();
00566     if(mu.isNonnull()){
00567       segComp=(float)(mu->matches().size());
00568       if(mu->caloCompatibility()>caloComp)
00569         caloComp = mu->caloCompatibility();
00570 
00571       if(segComp<1.0)
00572         decision=false;
00573 
00574       tau.setCaloComp(caloComp);
00575       tau.setSegComp(segComp);
00576       tau.setMuonDecision(decision);
00577     }
00578 
00579   }
00580 }
00581 
00582 
00583 void
00584 HPSPFRecoTauAlgorithm::applyElectronRejection(reco::PFTau& tau,double stripEnergy )
00585 {
00586   //Here we apply the common electron rejection variables.
00587   //The only not common is the E/P that is applied in the decay mode
00588   //construction
00589 
00590 
00591   if(tau.leadPFChargedHadrCand().isNonnull()) {
00592     PFCandidateRef leadCharged = tau.leadPFChargedHadrCand();
00593     math::XYZVector caloDir(leadCharged->positionAtECALEntrance().x(),
00594                             leadCharged->positionAtECALEntrance().y(),
00595                             leadCharged->positionAtECALEntrance().z());
00596 
00597     tau.setmaximumHCALPFClusterEt(leadCharged->hcalEnergy()*sin(caloDir.theta()));
00598 
00599 
00600 
00601     if(leadCharged->trackRef().isNonnull()) {
00602       TrackRef track = leadCharged->trackRef();
00603       tau.setemFraction(leadCharged->ecalEnergy()/(leadCharged->ecalEnergy()+leadCharged->hcalEnergy()));
00604       //For H/P trust particle Flow ! :Just take HCAL energy of the candidate
00605       //end of story
00606       tau.sethcalTotOverPLead(leadCharged->hcalEnergy()/track->p());
00607       tau.sethcalMaxOverPLead(leadCharged->hcalEnergy()/track->p());
00608       tau.sethcal3x3OverPLead(leadCharged->hcalEnergy()/track->p());
00609       tau.setelectronPreIDTrack(track);
00610       tau.setelectronPreIDOutput(leadCharged->mva_e_pi());
00611       //Since PF uses brem recovery we will store the default ecal energy here
00612       tau.setbremsRecoveryEOverPLead(leadCharged->ecalEnergy()/track->p());
00613       tau.setecalStripSumEOverPLead((leadCharged->ecalEnergy()-stripEnergy)/track->p());
00614       bool electronDecision;
00615       if(std::abs(leadCharged->pdgId())==11)
00616         electronDecision=true;
00617       else
00618         electronDecision=false;
00619       tau.setelectronPreIDDecision(electronDecision);
00620     }
00621   }
00622 }
00623 
00624 
00625 
00626 
00627 void
00628 HPSPFRecoTauAlgorithm::configure(const edm::ParameterSet& p)
00629 {
00630   emMerger_                      = p.getParameter<std::string>("emMergingAlgorithm");
00631   overlapCriterion_              = p.getParameter<std::string>("candOverlapCriterion");
00632   doOneProngs_                   = p.getParameter<bool>("doOneProng");
00633   doOneProngStrips_              = p.getParameter<bool>("doOneProngStrip");
00634   doOneProngTwoStrips_           = p.getParameter<bool>("doOneProngTwoStrips");
00635   doThreeProngs_                 = p.getParameter<bool>("doThreeProng");
00636   tauThreshold_                  = p.getParameter<double>("tauPtThreshold");
00637   leadPionThreshold_             = p.getParameter<double>("leadPionThreshold");
00638   stripPtThreshold_               = p.getParameter<double>("stripPtThreshold");
00639   chargeIsolationCone_           = p.getParameter<double>("chargeHadrIsolationConeSize");
00640   gammaIsolationCone_            = p.getParameter<double>("gammaIsolationConeSize");
00641   neutrHadrIsolationCone_        = p.getParameter<double>("neutrHadrIsolationConeSize");
00642   useIsolationAnnulus_           = p.getParameter<bool>("useIsolationAnnulus");
00643   oneProngStripMassWindow_       = p.getParameter<std::vector<double> >("oneProngStripMassWindow");
00644   oneProngTwoStripsMassWindow_   = p.getParameter<std::vector<double> >("oneProngTwoStripsMassWindow");
00645   oneProngTwoStripsPi0MassWindow_= p.getParameter<std::vector<double> >("oneProngTwoStripsPi0MassWindow");
00646   threeProngMassWindow_          = p.getParameter<std::vector<double> >("threeProngMassWindow");
00647   matchingCone_                  = p.getParameter<double>("matchingCone");
00648   coneMetric_                    = p.getParameter<std::string>("coneMetric");
00649   coneSizeFormula_               = p.getParameter<std::string>("coneSizeFormula");
00650   minSignalCone_                 = p.getParameter<double>("minimumSignalCone");
00651   maxSignalCone_                 = p.getParameter<double>("maximumSignalCone");
00652 
00653 
00654 
00655   //Initialize The Merging Algorithm!
00656   if(emMerger_ =="StripBased")
00657     candidateMerger_ =  new PFCandidateStripMerger(p);
00658   //Add the Pi0 Merger from Evan here
00659 
00660 
00661   if(oneProngStripMassWindow_.size()!=2)
00662     throw cms::Exception("") << "OneProngStripMassWindow must be a vector of size 2 [min,max] " << std::endl;
00663   if(oneProngTwoStripsMassWindow_.size()!=2)
00664     throw cms::Exception("") << "OneProngTwoStripsMassWindow must be a vector of size 2 [min,max] " << std::endl;
00665   if(threeProngMassWindow_.size()!=2)
00666     throw cms::Exception("") << "ThreeProngMassWindow must be a vector of size 2 [min,max] " << std::endl;
00667   if(coneMetric_!= "angle" && coneMetric_ != "DR")
00668     throw cms::Exception("") << "Cone Metric should be angle or DR " << std::endl;
00669 
00670   coneSizeFormula = TauTagTools::computeConeSizeTFormula(coneSizeFormula_,"Signal cone size Formula");
00671 
00672 
00673 }
00674 
00675 
00676 
00677 math::XYZTLorentzVector
00678 HPSPFRecoTauAlgorithm::createMergedLorentzVector(const reco::PFCandidateRefVector& cands)
00679 {
00680   math::XYZTLorentzVector sum;
00681   for(unsigned int i=0;i<cands.size();++i) {
00682     sum+=cands.at(i)->p4();
00683   }
00684   return sum;
00685 }
00686 
00687 void
00688 HPSPFRecoTauAlgorithm::removeCandidateFromRefVector(const reco::PFCandidateRef& cand,reco::PFCandidateRefVector& vec)
00689 {
00690   PFCandidateRefVector newVec;
00691 
00692   PFCandidateRefVector::iterator it;
00693   it = std::find(vec.begin(),vec.end(),cand);
00694   if(it!=vec.end())
00695     vec.erase(it);
00696 }
00697 
00698 void
00699 HPSPFRecoTauAlgorithm::applyMassConstraint(math::XYZTLorentzVector& vec,double mass)
00700 {
00701   double factor = sqrt(vec.energy()*vec.energy()-mass*mass)/vec.P();
00702   vec.SetCoordinates(vec.px()*factor,vec.py()*factor,vec.pz()*factor,vec.energy());
00703 }
00704 
00705 
00706 
00707 
00708 bool
00709 HPSPFRecoTauAlgorithm::refitThreeProng(reco::PFTau& tau)
00710 {
00711   bool response=false;
00712   //Get Hadrons
00713   reco::PFCandidateRefVector hadrons = tau.signalPFChargedHadrCands();
00714   PFCandidateRef  h1  = hadrons.at(0);
00715   PFCandidateRef  h2  = hadrons.at(1);
00716   PFCandidateRef  h3  = hadrons.at(2);
00717 
00718   //Make transient tracks
00719   std::vector<TransientTrack> transientTracks;
00720   transientTracks.push_back(TransientTrackBuilder_->build(h1->trackRef()));
00721   transientTracks.push_back(TransientTrackBuilder_->build(h2->trackRef()));
00722   transientTracks.push_back(TransientTrackBuilder_->build(h3->trackRef()));
00723 
00724   //Apply the Vertex Fit
00725   KalmanVertexFitter fitter(true);
00726   TransientVertex myVertex = fitter.vertex(transientTracks);
00727 
00728   //Just require a valid vertex+ 3 refitted tracks
00729   if(myVertex.isValid()&&
00730      myVertex.hasRefittedTracks()&&
00731      myVertex.refittedTracks().size()==3) {
00732 
00733     response=true;
00734     math::XYZPoint vtx(myVertex.position().x(),myVertex.position().y(),myVertex.position().z());
00735 
00736     //Create a LV for each refitted track
00737     math::XYZTLorentzVector p1(myVertex.refittedTracks().at(0).track().px(),
00738                                myVertex.refittedTracks().at(0).track().py(),
00739                                myVertex.refittedTracks().at(0).track().pz(),
00740                                sqrt(myVertex.refittedTracks().at(0).track().momentum().mag2() +0.139*0.139));
00741 
00742     math::XYZTLorentzVector p2(myVertex.refittedTracks().at(1).track().px(),
00743                                myVertex.refittedTracks().at(1).track().py(),
00744                                myVertex.refittedTracks().at(1).track().pz(),
00745                                sqrt(myVertex.refittedTracks().at(1).track().momentum().mag2() +0.139*0.139));
00746 
00747     math::XYZTLorentzVector p3(myVertex.refittedTracks().at(2).track().px(),
00748                                myVertex.refittedTracks().at(2).track().py(),
00749                                myVertex.refittedTracks().at(2).track().pz(),
00750                                sqrt(myVertex.refittedTracks().at(2).track().momentum().mag2() +0.139*0.139));
00751 
00752     //Update the tau p4
00753     tau.setP4(p1+p2+p3);
00754     //Update the vertex
00755     tau.setVertex(vtx);
00756   }
00757   return response;
00758 
00759 }
00760 
00761 
00762 reco::PFTau
00763 HPSPFRecoTauAlgorithm::getBestTauCandidate(reco::PFTauCollection& taus)
00764 {
00765   reco::PFTauCollection::iterator it;
00766   if(overlapCriterion_ =="Isolation"){
00767     HPSTauIsolationSorter sorter;
00768     it = std::min_element(taus.begin(),taus.end(),sorter);
00769 
00770   }
00771   else if(overlapCriterion_ =="Pt"){
00772     HPSTauPtSorter sorter;
00773     it = std::min_element(taus.begin(),taus.end(),sorter);
00774   }
00775 
00776   return *it;
00777 }
00778