Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013 #include "RecoTauTag/RecoTau/interface/RecoTauBuilderPlugins.h"
00014 #include "DataFormats/TrackReco/interface/Track.h"
00015 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00016 #include "DataFormats/ParticleFlowReco/interface/PFBlockElement.h"
00017 #include "DataFormats/ParticleFlowReco/interface/PFBlock.h"
00018 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
00019 #include <Math/VectorUtil.h>
00020 #include <algorithm>
00021
00022 namespace reco { namespace tau {
00023
00024 class RecoTauElectronRejectionPlugin : public RecoTauModifierPlugin {
00025 public:
00026 explicit RecoTauElectronRejectionPlugin(const edm::ParameterSet& pset);
00027 virtual ~RecoTauElectronRejectionPlugin() {}
00028 void operator()(PFTau&) const;
00029 private:
00030 double ElecPreIDLeadTkMatch_maxDR_;
00031 double EcalStripSumE_minClusEnergy_;
00032 double EcalStripSumE_deltaEta_;
00033 double EcalStripSumE_deltaPhiOverQ_minValue_;
00034 double EcalStripSumE_deltaPhiOverQ_maxValue_;
00035 double maximumForElectrionPreIDOutput_;
00036 std::string DataType_;
00037 };
00038
00039 RecoTauElectronRejectionPlugin::RecoTauElectronRejectionPlugin(
00040 const edm::ParameterSet& pset):RecoTauModifierPlugin(pset) {
00041
00042 ElecPreIDLeadTkMatch_maxDR_ =
00043 pset.getParameter<double>("ElecPreIDLeadTkMatch_maxDR");
00044 EcalStripSumE_minClusEnergy_ =
00045 pset.getParameter<double>("EcalStripSumE_minClusEnergy");
00046 EcalStripSumE_deltaEta_ =
00047 pset.getParameter<double>("EcalStripSumE_deltaEta");
00048 EcalStripSumE_deltaPhiOverQ_minValue_ =
00049 pset.getParameter<double>("EcalStripSumE_deltaPhiOverQ_minValue");
00050 EcalStripSumE_deltaPhiOverQ_maxValue_ =
00051 pset.getParameter<double>("EcalStripSumE_deltaPhiOverQ_maxValue");
00052 maximumForElectrionPreIDOutput_ =
00053 pset.getParameter<double>("maximumForElectrionPreIDOutput");
00054 DataType_ = pset.getParameter<std::string>("DataType");
00055 }
00056
00057 namespace {
00058 bool checkPos(std::vector<math::XYZPoint> CalPos,math::XYZPoint CandPos) {
00059 bool flag = false;
00060 for (unsigned int i=0;i<CalPos.size();i++) {
00061 if (CalPos[i] == CandPos) {
00062 flag = true;
00063 break;
00064 }
00065 }
00066 return flag;
00067 }
00068 }
00069
00070 void RecoTauElectronRejectionPlugin::operator()(PFTau& tau) const {
00071
00072 double myECALenergy = 0.;
00073 double myHCALenergy = 0.;
00074 double myHCALenergy3x3 = 0.;
00075 double myMaximumHCALPFClusterE = 0.;
00076 double myMaximumHCALPFClusterEt = 0.;
00077 double myStripClusterE = 0.;
00078 double myEmfrac = -1.;
00079 double myElectronPreIDOutput = -1111.;
00080 bool myElecPreid = false;
00081 reco::TrackRef myElecTrk;
00082
00083 typedef std::pair<reco::PFBlockRef, unsigned> ElementInBlock;
00084 typedef std::vector< ElementInBlock > ElementsInBlocks;
00085
00086 PFCandidateRef myleadPFChargedCand = tau.leadPFChargedHadrCand();
00087
00088 PFCandidateRefVector myPFCands;
00089 myPFCands.reserve(tau.isolationPFCands().size()+tau.signalPFCands().size());
00090
00091 std::copy(tau.isolationPFCands().begin(), tau.isolationPFCands().end(),
00092 std::back_inserter(myPFCands));
00093 std::copy(tau.signalPFCands().begin(), tau.signalPFCands().end(),
00094 std::back_inserter(myPFCands));
00095
00096
00097 if(myleadPFChargedCand.isNonnull()){
00098 myElectronPreIDOutput = myleadPFChargedCand->mva_e_pi();
00099
00100 math::XYZPointF myElecTrkEcalPos = myleadPFChargedCand->positionAtECALEntrance();
00101 myElecTrk = myleadPFChargedCand->trackRef();
00102
00103 if(myElecTrk.isNonnull()) {
00104
00105 if(DataType_ == "AOD"){
00106
00107 for(int i=0;i<(int)myPFCands.size();i++){
00108 myHCALenergy += myPFCands[i]->hcalEnergy();
00109 myECALenergy += myPFCands[i]->ecalEnergy();
00110
00111 math::XYZPointF candPos;
00112 if (myPFCands[i]->particleId()==1 || myPFCands[i]->particleId()==2)
00113 candPos = myPFCands[i]->positionAtECALEntrance();
00114 else
00115 candPos = math::XYZPointF(myPFCands[i]->px(),myPFCands[i]->py(),myPFCands[i]->pz());
00116
00117 double deltaR = ROOT::Math::VectorUtil::DeltaR(myElecTrkEcalPos,candPos);
00118 double deltaPhi = ROOT::Math::VectorUtil::DeltaPhi(myElecTrkEcalPos,candPos);
00119 double deltaEta = abs(myElecTrkEcalPos.eta()-candPos.eta());
00120 double deltaPhiOverQ = deltaPhi/(double)myElecTrk->charge();
00121
00122 if (myPFCands[i]->ecalEnergy() >= EcalStripSumE_minClusEnergy_ && deltaEta < EcalStripSumE_deltaEta_ &&
00123 deltaPhiOverQ > EcalStripSumE_deltaPhiOverQ_minValue_ && deltaPhiOverQ < EcalStripSumE_deltaPhiOverQ_maxValue_) {
00124 myStripClusterE += myPFCands[i]->ecalEnergy();
00125 }
00126 if (deltaR<0.184) {
00127 myHCALenergy3x3 += myPFCands[i]->hcalEnergy();
00128 }
00129 if (myPFCands[i]->hcalEnergy()>myMaximumHCALPFClusterE) {
00130 myMaximumHCALPFClusterE = myPFCands[i]->hcalEnergy();
00131 }
00132 if ((myPFCands[i]->hcalEnergy()*fabs(sin(candPos.Theta())))>myMaximumHCALPFClusterEt) {
00133 myMaximumHCALPFClusterEt = (myPFCands[i]->hcalEnergy()*fabs(sin(candPos.Theta())));
00134 }
00135 }
00136
00137 } else if(DataType_ == "RECO"){
00138
00139 std::vector<math::XYZPoint> hcalPosV; hcalPosV.clear();
00140 std::vector<math::XYZPoint> ecalPosV; ecalPosV.clear();
00141 for(int i=0;i<(int)myPFCands.size();i++){
00142 const ElementsInBlocks& elts = myPFCands[i]->elementsInBlocks();
00143 for(ElementsInBlocks::const_iterator it=elts.begin(); it!=elts.end(); ++it) {
00144 const reco::PFBlock& block = *(it->first);
00145 unsigned indexOfElementInBlock = it->second;
00146 const edm::OwnVector< reco::PFBlockElement >& elements = block.elements();
00147 assert(indexOfElementInBlock<elements.size());
00148
00149 const reco::PFBlockElement& element = elements[indexOfElementInBlock];
00150
00151 if(element.type()==reco::PFBlockElement::HCAL) {
00152 math::XYZPoint clusPos = element.clusterRef()->position();
00153 double en = (double)element.clusterRef()->energy();
00154 double et = (double)element.clusterRef()->energy()*fabs(sin(clusPos.Theta()));
00155 if (en>myMaximumHCALPFClusterE) {
00156 myMaximumHCALPFClusterE = en;
00157 }
00158 if (et>myMaximumHCALPFClusterEt) {
00159 myMaximumHCALPFClusterEt = et;
00160 }
00161 if (!checkPos(hcalPosV,clusPos)) {
00162 hcalPosV.push_back(clusPos);
00163 myHCALenergy += en;
00164 double deltaR = ROOT::Math::VectorUtil::DeltaR(myElecTrkEcalPos,clusPos);
00165 if (deltaR<0.184) {
00166 myHCALenergy3x3 += en;
00167 }
00168 }
00169 } else if(element.type()==reco::PFBlockElement::ECAL) {
00170 double en = (double)element.clusterRef()->energy();
00171 math::XYZPoint clusPos = element.clusterRef()->position();
00172 if (!checkPos(ecalPosV,clusPos)) {
00173 ecalPosV.push_back(clusPos);
00174 myECALenergy += en;
00175 double deltaPhi = ROOT::Math::VectorUtil::DeltaPhi(myElecTrkEcalPos,clusPos);
00176 double deltaEta = abs(myElecTrkEcalPos.eta()-clusPos.eta());
00177 double deltaPhiOverQ = deltaPhi/(double)myElecTrk->charge();
00178 if (en >= EcalStripSumE_minClusEnergy_ && deltaEta<EcalStripSumE_deltaEta_ && deltaPhiOverQ > EcalStripSumE_deltaPhiOverQ_minValue_ && deltaPhiOverQ < EcalStripSumE_deltaPhiOverQ_maxValue_) {
00179 myStripClusterE += en;
00180 }
00181 }
00182 }
00183 }
00184 }
00185 }
00186 }
00187 }
00188
00189 if ((myHCALenergy+myECALenergy)>0.)
00190 myEmfrac = myECALenergy/(myHCALenergy+myECALenergy);
00191 tau.setemFraction((float)myEmfrac);
00192
00193
00194 if (myElecTrk.isNonnull())
00195 {
00196 float myElectronMomentum = (float)myElecTrk->p();
00197 if (myElectronMomentum > 0.)
00198 {
00199 myHCALenergy /= myElectronMomentum;
00200 myMaximumHCALPFClusterE /= myElectronMomentum;
00201 myHCALenergy3x3 /= myElectronMomentum;
00202 myStripClusterE /= myElectronMomentum;
00203 }
00204 }
00205 tau.sethcalTotOverPLead((float)myHCALenergy);
00206 tau.sethcalMaxOverPLead((float)myMaximumHCALPFClusterE);
00207 tau.sethcal3x3OverPLead((float)myHCALenergy3x3);
00208 tau.setecalStripSumEOverPLead((float)myStripClusterE);
00209 tau.setmaximumHCALPFClusterEt(myMaximumHCALPFClusterEt);
00210 tau.setelectronPreIDOutput(myElectronPreIDOutput);
00211 if (myElecTrk.isNonnull())
00212 tau.setelectronPreIDTrack(myElecTrk);
00213 if (myElectronPreIDOutput > maximumForElectrionPreIDOutput_)
00214 myElecPreid = true;
00215 tau.setelectronPreIDDecision(myElecPreid);
00216
00217
00218
00219
00220
00221 }
00222 }}
00223 #include "FWCore/Framework/interface/MakerMacros.h"
00224 DEFINE_EDM_PLUGIN(RecoTauModifierPluginFactory,
00225 reco::tau::RecoTauElectronRejectionPlugin,
00226 "RecoTauElectronRejectionPlugin");