00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012 #include "RecoMET/METAlgorithms/interface/GenSpecificAlgo.h"
00013 #include "TMath.h"
00014
00015 #include <set>
00016
00017
00018 reco::GenMET GenSpecificAlgo::addInfo(edm::Handle<edm::View<reco::Candidate> > particles, CommonMETData *met, double globalThreshold, bool onlyFiducial, bool usePt)
00019 {
00020 fillCommonMETData(met, particles, globalThreshold, onlyFiducial, usePt);
00021
00022 SpecificGenMETData specific = mkSpecificGenMETData(particles, usePt);
00023
00024 const LorentzVector p4( met->mex, met->mey, met->mez, met->met );
00025 const Point vtx( 0.0, 0.0, 0.0 );
00026
00027 reco::GenMET genMet(specific, met->sumet, p4, vtx );
00028 return genMet;
00029 }
00030
00031
00032 void GenSpecificAlgo::fillCommonMETData(CommonMETData *met, edm::Handle<edm::View<reco::Candidate> >& particles, double globalThreshold, bool onlyFiducial, bool usePt)
00033 {
00034 double sum_et = 0.0;
00035 double sum_ex = 0.0;
00036 double sum_ey = 0.0;
00037 double sum_ez = 0.0;
00038
00039 for( edm::View<reco::Candidate>::const_iterator iParticle = (particles.product())->begin(); iParticle != (particles.product())->end(); ++iParticle)
00040 {
00041 if( (onlyFiducial && TMath::Abs(iParticle->eta()) >= 5.0)) continue;
00042
00043 if( iParticle->et() <= globalThreshold ) continue;
00044
00045 if(usePt)
00046 {
00047 double phi = iParticle->phi();
00048 double et = iParticle->pt();
00049 sum_ez += iParticle->pz();
00050 sum_et += et;
00051 sum_ex += et*cos(phi);
00052 sum_ey += et*sin(phi);
00053 }
00054 else
00055 {
00056 double phi = iParticle->phi();
00057 double theta = iParticle->theta();
00058 double e = iParticle->energy();
00059 double et = e*sin(theta);
00060 sum_ez += e*cos(theta);
00061 sum_et += et;
00062 sum_ex += et*cos(phi);
00063 sum_ey += et*sin(phi);
00064 }
00065 }
00066
00067 met->mex = -sum_ex;
00068 met->mey = -sum_ey;
00069 met->mez = -sum_ez;
00070 met->met = sqrt( sum_ex*sum_ex + sum_ey*sum_ey );
00071 met->sumet = sum_et;
00072 met->phi = atan2( -sum_ey, -sum_ex );
00073 }
00074
00075
00076 SpecificGenMETData GenSpecificAlgo::mkSpecificGenMETData(edm::Handle<edm::View<reco::Candidate> >& particles, bool usePt)
00077 {
00078 const static int neutralEMpdgId[] = { 22 };
00079 const static std::set<int> neutralEMpdgIdSet(neutralEMpdgId, neutralEMpdgId + sizeof(neutralEMpdgId)/sizeof(int));
00080
00081 const static int chargedEMpdgId[] = { 11 };
00082 const static std::set<int> chargedEMpdgIdSet(chargedEMpdgId, chargedEMpdgId + sizeof(chargedEMpdgId)/sizeof(int));
00083
00084 const static int muonpdgId[] = { 13 };
00085 const static std::set<int> muonpdgIdSet(muonpdgId, muonpdgId + sizeof(muonpdgId)/sizeof(int));
00086
00087 const static int neutralHADpdgId[] = {
00088 130 ,
00089 310 ,
00090 3122 ,
00091 2112 ,
00092 3222
00093 };
00094 const static std::set<int> neutralHADpdgIdSet(neutralHADpdgId, neutralHADpdgId + sizeof(neutralHADpdgId)/sizeof(int));
00095
00096 const static int chargedHADpdgId[] = {
00097 211 ,
00098 321 ,
00099 2212 ,
00100 3312 ,
00101 3112 ,
00102 3322 ,
00103 3334
00104 };
00105 const static std::set<int> chargedHADpdgIdSet(chargedHADpdgId, chargedHADpdgId + sizeof(chargedHADpdgId)/sizeof(int));
00106
00107 const static int invisiblepdgId[] = {
00108 12 ,
00109 14 ,
00110 16 ,
00111 1000022 ,
00112 1000012 ,
00113 1000014 ,
00114 1000016 ,
00115 2000012 ,
00116 2000014 ,
00117 2000016 ,
00118 39 ,
00119 1000039 ,
00120 5100039 ,
00121 4000012 ,
00122 4000014 ,
00123 4000016 ,
00124 9900012 ,
00125 9900014 ,
00126 9900016 ,
00127 };
00128 const static std::set<int> invisiblepdgIdSet(invisiblepdgId, invisiblepdgId + sizeof(invisiblepdgId)/sizeof(int));
00129
00130 SpecificGenMETData specific = SpecificGenMETData();
00131 specific.NeutralEMEtFraction = 0.0;
00132 specific.NeutralHadEtFraction = 0.0;
00133 specific.ChargedEMEtFraction = 0.0;
00134 specific.ChargedHadEtFraction = 0.0;
00135 specific.MuonEtFraction = 0.0;
00136 specific.InvisibleEtFraction = 0.0;
00137 double Et_unclassified = 0.0;
00138
00139 for(edm::View<reco::Candidate>::const_iterator iParticle = (particles.product())->begin(); iParticle != (particles.product())->end(); ++iParticle)
00140 {
00141 int pdgId = TMath::Abs( iParticle->pdgId() ) ;
00142 double pt = (usePt) ? iParticle->pt() : iParticle->et();
00143 if(neutralEMpdgIdSet.count(pdgId)) specific.NeutralEMEtFraction += pt;
00144 else if(chargedEMpdgIdSet.count(pdgId)) specific.ChargedEMEtFraction += pt;
00145 else if(muonpdgIdSet.count(pdgId)) specific.MuonEtFraction += pt;
00146 else if(neutralHADpdgIdSet.count(pdgId)) specific.NeutralHadEtFraction += pt;
00147 else if(chargedHADpdgIdSet.count(pdgId)) specific.ChargedHadEtFraction += pt;
00148 else if(invisiblepdgIdSet.count(pdgId)) specific.InvisibleEtFraction += pt;
00149 else Et_unclassified += pt;
00150 }
00151
00152 double Et_Total = specific.NeutralEMEtFraction + specific.NeutralHadEtFraction + specific.ChargedEMEtFraction +
00153 specific.ChargedHadEtFraction + specific.MuonEtFraction + specific.InvisibleEtFraction + Et_unclassified;
00154
00155 if(Et_Total)
00156 {
00157 specific.NeutralEMEtFraction /= Et_Total;
00158 specific.NeutralHadEtFraction /= Et_Total;
00159 specific.ChargedEMEtFraction /= Et_Total;
00160 specific.ChargedHadEtFraction /= Et_Total;
00161 specific.MuonEtFraction /= Et_Total;
00162 specific.InvisibleEtFraction /= Et_Total;
00163 }
00164
00165 return specific;
00166 }
00167
00168