Go to the documentation of this file.00001 #include "DataFormats/Math/interface/LorentzVector.h"
00002 #include "RecoMET/METAlgorithms/interface/CaloSpecificAlgo.h"
00003 #include "DataFormats/CaloTowers/interface/CaloTowerDetId.h"
00004 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
00005 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
00006 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
00007
00008 using namespace reco;
00009 using namespace std;
00010
00011
00012
00013
00014
00015
00016
00017
00018 reco::CaloMET CaloSpecificAlgo::addInfo(edm::Handle<edm::View<Candidate> > towers, CommonMETData met, bool noHF, double globalThreshold)
00019 {
00020
00021 SpecificCaloMETData specific;
00022
00023 specific.MaxEtInEmTowers = 0.0;
00024 specific.MaxEtInHadTowers = 0.0;
00025 specific.HadEtInHO = 0.0;
00026 specific.HadEtInHB = 0.0;
00027 specific.HadEtInHF = 0.0;
00028 specific.HadEtInHE = 0.0;
00029 specific.EmEtInEB = 0.0;
00030 specific.EmEtInEE = 0.0;
00031 specific.EmEtInHF = 0.0;
00032 specific.EtFractionHadronic = 0.0;
00033 specific.EtFractionEm = 0.0;
00034 specific.CaloSETInpHF = 0.0;
00035 specific.CaloSETInmHF = 0.0;
00036 specific.CaloMETInpHF = 0.0;
00037 specific.CaloMETInmHF = 0.0;
00038 specific.CaloMETPhiInpHF = 0.0;
00039 specific.CaloMETPhiInmHF = 0.0;
00040 specific.METSignificance = 0.0;
00041
00042 double totalEt = 0.0;
00043 double totalEm = 0.0;
00044 double totalHad = 0.0;
00045 double MaxTowerEm = 0.0;
00046 double MaxTowerHad = 0.0;
00047 double sumEtInpHF = 0.0;
00048 double sumEtInmHF = 0.0;
00049 double MExInpHF = 0.0;
00050 double MEyInpHF = 0.0;
00051 double MExInmHF = 0.0;
00052 double MEyInmHF = 0.0;
00053
00054 if( towers->size() == 0 )
00055 {
00056
00057 const LorentzVector p4( met.mex, met.mey, 0.0, met.met );
00058 const Point vtx( 0.0, 0.0, 0.0 );
00059 CaloMET specificmet( specific, met.sumet, p4, vtx );
00060 return specificmet;
00061 }
00062
00063 edm::View<Candidate>::const_iterator towerCand = towers->begin();
00064 for( ; towerCand != towers->end(); towerCand++ )
00065 {
00066 const Candidate* candidate = &(*towerCand);
00067 if (candidate) {
00068 const CaloTower* calotower = dynamic_cast<const CaloTower*> (candidate);
00069 if (calotower)
00070 {
00071 double caloTowerEt=calotower->et();
00072 double caloTowerHadEt=calotower->hadEt();
00073 double caloTowerEmEt=calotower->emEt();
00074 if(caloTowerEt < globalThreshold) continue;
00075 totalEt += caloTowerEt;
00076 totalEm += caloTowerEmEt;
00077
00078
00079
00080 bool hadIsDone = false;
00081 bool emIsDone = false;
00082 int cell = calotower->constituentsSize();
00083 while ( --cell >= 0 && (!hadIsDone || !emIsDone) )
00084 {
00085 DetId id = calotower->constituent( cell );
00086 if( !hadIsDone && id.det() == DetId::Hcal )
00087 {
00088 HcalSubdetector subdet = HcalDetId(id).subdet();
00089 if( subdet == HcalBarrel || subdet == HcalOuter )
00090 {
00091 if( caloTowerHadEt > MaxTowerHad ) MaxTowerHad = caloTowerHadEt;
00092 specific.HadEtInHB += caloTowerHadEt;
00093 specific.HadEtInHO += calotower->outerEt();
00094 }
00095 else if( subdet == HcalEndcap )
00096 {
00097 if( caloTowerHadEt > MaxTowerHad ) MaxTowerHad = caloTowerHadEt;
00098 specific.HadEtInHE += caloTowerHadEt;
00099 }
00100 else if( subdet == HcalForward )
00101 {
00102 if (!noHF)
00103 {
00104 if( caloTowerHadEt > MaxTowerHad ) MaxTowerHad = caloTowerHadEt;
00105 if( caloTowerEmEt > MaxTowerEm ) MaxTowerEm = caloTowerEmEt;
00106
00107 specific.HadEtInHF += caloTowerHadEt;
00108 specific.EmEtInHF += caloTowerEmEt;
00109 }
00110 else
00111 {
00112
00113
00114 totalEm -= caloTowerEmEt;
00115 totalEt -= caloTowerEt;
00116 }
00117
00118
00119 if (calotower->eta()>=0)
00120 {
00121 sumEtInpHF += caloTowerEt;
00122 MExInpHF -= (caloTowerEt * cos(calotower->phi()));
00123 MEyInpHF -= (caloTowerEt * sin(calotower->phi()));
00124 }
00125 else
00126 {
00127 sumEtInmHF += caloTowerEt;
00128 MExInmHF -= (caloTowerEt * cos(calotower->phi()));
00129 MEyInmHF -= (caloTowerEt * sin(calotower->phi()));
00130 }
00131 }
00132 hadIsDone = true;
00133 }
00134 else if( !emIsDone && id.det() == DetId::Ecal )
00135 {
00136 EcalSubdetector subdet = EcalSubdetector( id.subdetId() );
00137 if( caloTowerEmEt > MaxTowerEm ) MaxTowerEm = caloTowerEmEt;
00138 if( subdet == EcalBarrel )
00139 {
00140 specific.EmEtInEB += caloTowerEmEt;
00141 }
00142 else if( subdet == EcalEndcap )
00143 {
00144 specific.EmEtInEE += caloTowerEmEt;
00145 }
00146 emIsDone = true;
00147 }
00148 }
00149 }
00150 }
00151 }
00152
00153
00154
00155 totalHad += (totalEt - totalEm);
00156
00157 if(!noHF)
00158 {
00159 LorentzVector METpHF(MExInpHF, MEyInpHF, 0, sqrt(MExInpHF*MExInpHF + MEyInpHF*MEyInpHF));
00160 LorentzVector METmHF(MExInmHF, MEyInmHF, 0, sqrt(MExInmHF*MExInmHF + MEyInmHF*MEyInmHF));
00161 specific.CaloMETInpHF = METpHF.pt();
00162 specific.CaloMETInmHF = METmHF.pt();
00163 specific.CaloMETPhiInpHF = METpHF.Phi();
00164 specific.CaloMETPhiInmHF = METmHF.Phi();
00165 specific.CaloSETInpHF = sumEtInpHF;
00166 specific.CaloSETInmHF = sumEtInmHF;
00167 }
00168 else
00169 {
00170 met.mex -= (MExInmHF + MExInpHF);
00171 met.mey -= (MEyInmHF + MEyInpHF);
00172 met.sumet -= (sumEtInpHF + sumEtInmHF);
00173 met.met = sqrt(met.mex*met.mex + met.mey*met.mey);
00174 }
00175
00176 specific.MaxEtInEmTowers = MaxTowerEm;
00177 specific.MaxEtInHadTowers = MaxTowerHad;
00178 specific.EtFractionHadronic = totalHad / totalEt;
00179 specific.EtFractionEm = totalEm / totalEt;
00180
00181 const LorentzVector p4( met.mex, met.mey, 0.0, met.met );
00182 const Point vtx( 0.0, 0.0, 0.0 );
00183
00184
00185 CaloMET specificmet( specific, met.sumet, p4, vtx );
00186 return specificmet;
00187 }
00188