00001 #include "RecoTauTag/TauTagTools/interface/PFCandidateStripMerger.h" 00002 #include "Math/GenVector/VectorUtil.h" 00003 00004 using namespace reco; 00005 using namespace std; 00006 PFCandidateStripMerger::PFCandidateStripMerger(): 00007 PFCandidateMergerBase() 00008 {} 00009 00010 PFCandidateStripMerger::PFCandidateStripMerger(const edm::ParameterSet& config): 00011 PFCandidateMergerBase(config) 00012 { 00013 inputPdgIds_ = config.getParameter<std::vector<int> >("stripCandidatesPdgIds"); 00014 etaAssociationDistance_ = config.getParameter<double>("stripEtaAssociationDistance"); 00015 phiAssociationDistance_ = config.getParameter<double>("stripPhiAssociationDistance"); 00016 00017 } 00018 00019 00020 PFCandidateStripMerger::~PFCandidateStripMerger() 00021 {} 00022 00023 00024 bool 00025 PFCandidateStripMerger::candidateMatches(const reco::PFCandidateRef& cand) 00026 { 00027 bool matches = false; 00028 for(unsigned int i=0; i < inputPdgIds_.size(); ++i) { 00029 if(std::abs(cand->pdgId()) == inputPdgIds_.at(i)) { 00030 matches = true; 00031 continue; 00032 } 00033 } 00034 00035 return matches; 00036 } 00037 00038 00039 00040 00041 00042 vector<PFCandidateRefVector> 00043 PFCandidateStripMerger::mergeCandidates(const PFCandidateRefVector& candidates) 00044 { 00045 00046 //Copy the input getting the relevant candidates and sort by pt 00047 PFCandidateRefVector cands; 00048 for(unsigned int i=0;i<candidates.size();++i) 00049 if(candidateMatches(candidates.at(i))) 00050 cands.push_back(candidates.at(i)); 00051 00052 if(cands.size()>1) 00053 TauTagTools::sortRefVectorByPt(cands); 00054 00055 std::vector<PFCandidateRefVector> strips; 00056 00057 00058 //Repeat while there are still unclusterized gammas 00059 while(cands.size()>0) { 00060 00061 //save the non associated candidates to a different collection 00062 PFCandidateRefVector notAssociated; 00063 00064 //Create a cluster from the Seed Photon 00065 PFCandidateRefVector strip; 00066 math::XYZTLorentzVector stripVector; 00067 strip.push_back(cands.at(0)); 00068 stripVector=cands.at(0)->p4(); 00069 00070 //Loop and associate 00071 for(unsigned int i=1;i<cands.size();++i) { 00072 if(fabs(cands.at(i)->eta()-stripVector.eta())<etaAssociationDistance_ && 00073 fabs(ROOT::Math::VectorUtil::DeltaPhi(cands.at(i)->p4(),stripVector))<phiAssociationDistance_) { 00074 strip.push_back(cands.at(i)); 00075 stripVector+=cands.at(i)->p4(); 00076 } 00077 else { 00078 notAssociated.push_back(cands.at(i)); 00079 } 00080 } 00081 //Save the strip 00082 strips.push_back(strip); 00083 00084 //Swap the candidate vector with the non associated vector 00085 cands.swap(notAssociated); 00086 //Clear 00087 notAssociated.clear(); 00088 } 00089 00090 return strips; 00091 }