00001 #include <iostream>
00002 #include <sstream>
00003 #include <istream>
00004 #include <fstream>
00005 #include <iomanip>
00006 #include <string>
00007 #include <cmath>
00008 #include <functional>
00009 #include <stdlib.h>
00010 #include <string.h>
00011
00012 #include "HLTrigger/HLTanalyzers/interface/HLTJets.h"
00013
00014 HLTJets::HLTJets() {
00015 evtCounter=0;
00016
00017
00018 _Monte=false;
00019 _Debug=false;
00020 _CalJetMin=0.;
00021 _GenJetMin=0.;
00022 }
00023
00024
00025 void HLTJets::setup(const edm::ParameterSet& pSet, TTree* HltTree) {
00026
00027 edm::ParameterSet myJetParams = pSet.getParameter<edm::ParameterSet>("RunParameters") ;
00028 vector<std::string> parameterNames = myJetParams.getParameterNames() ;
00029
00030 for ( vector<std::string>::iterator iParam = parameterNames.begin();
00031 iParam != parameterNames.end(); iParam++ ){
00032 if ( (*iParam) == "Monte" ) _Monte = myJetParams.getParameter<bool>( *iParam );
00033 else if ( (*iParam) == "Debug" ) _Debug = myJetParams.getParameter<bool>( *iParam );
00034 else if ( (*iParam) == "CalJetMin" ) _CalJetMin = myJetParams.getParameter<double>( *iParam );
00035 else if ( (*iParam) == "GenJetMin" ) _GenJetMin = myJetParams.getParameter<double>( *iParam );
00036 }
00037
00038 const int kMaxJetCal = 10000;
00039 jcalpt = new float[kMaxJetCal];
00040 jcalphi = new float[kMaxJetCal];
00041 jcaleta = new float[kMaxJetCal];
00042 jcale = new float[kMaxJetCal];
00043 const int kMaxJetgen = 10000;
00044 jgenpt = new float[kMaxJetgen];
00045 jgenphi = new float[kMaxJetgen];
00046 jgeneta = new float[kMaxJetgen];
00047 jgene = new float[kMaxJetgen];
00048 const int kMaxTower = 10000;
00049 towet = new float[kMaxTower];
00050 toweta = new float[kMaxTower];
00051 towphi = new float[kMaxTower];
00052 towen = new float[kMaxTower];
00053 towem = new float[kMaxTower];
00054 towhd = new float[kMaxTower];
00055 towoe = new float[kMaxTower];
00056 const int kMaxTau = 500;
00057 l2tauemiso = new float[kMaxTau];
00058 l25tauPt = new float[kMaxTau];
00059 l25tautckiso = new int[kMaxTau];
00060 l3tauPt = new float[kMaxTau];
00061 l3tautckiso = new int[kMaxTau];
00062 tauEta = new float[kMaxTau];
00063 tauPt = new float[kMaxTau];
00064 tauPhi = new float[kMaxTau];
00065
00066
00067
00068 HltTree->Branch("NrecoJetCal",&njetcal,"NrecoJetCal/I");
00069 HltTree->Branch("NrecoJetGen",&njetgen,"NrecoJetGen/I");
00070 HltTree->Branch("NrecoTowCal",&ntowcal,"NrecoTowCal/I");
00071 HltTree->Branch("recoJetCalPt",jcalpt,"recoJetCalPt[NrecoJetCal]/F");
00072 HltTree->Branch("recoJetCalPhi",jcalphi,"recoJetCalPhi[NrecoJetCal]/F");
00073 HltTree->Branch("recoJetCalEta",jcaleta,"recoJetCalEta[NrecoJetCal]/F");
00074 HltTree->Branch("recoJetCalE",jcale,"recoJetCalE[NrecoJetCal]/F");
00075 HltTree->Branch("recoJetGenPt",jgenpt,"recoJetGenPt[NrecoJetGen]/F");
00076 HltTree->Branch("recoJetGenPhi",jgenphi,"recoJetGenPhi[NrecoJetGen]/F");
00077 HltTree->Branch("recoJetGenEta",jgeneta,"recoJetGenEta[NrecoJetGen]/F");
00078 HltTree->Branch("recoJetGenE",jgene,"recoJetGenE[NrecoJetGen]/F");
00079 HltTree->Branch("recoTowEt",towet,"recoTowEt[NrecoTowCal]/F");
00080 HltTree->Branch("recoTowEta",toweta,"recoTowEta[NrecoTowCal]/F");
00081 HltTree->Branch("recoTowPhi",towphi,"recoTowPhi[NrecoTowCal]/F");
00082 HltTree->Branch("recoTowE",towen,"recoTowE[NrecoTowCal]/F");
00083 HltTree->Branch("recoTowEm",towem,"recoTowEm[NrecoTowCal]/F");
00084 HltTree->Branch("recoTowHad",towhd,"recoTowHad[NrecoTowCal]/F");
00085 HltTree->Branch("recoTowOE",towoe,"recoTowOE[NrecoTowCal]/F");
00086 HltTree->Branch("recoMetCal",&mcalmet,"recoMetCal/F");
00087 HltTree->Branch("recoMetCalPhi",&mcalphi,"recoMetCalPhi/F");
00088 HltTree->Branch("recoMetCalSum",&mcalsum,"recoMetCalSum/F");
00089 HltTree->Branch("recoMetGen",&mgenmet,"recoMetGen/F");
00090 HltTree->Branch("recoMetGenPhi",&mgenphi,"recoMetGenPhi/F");
00091 HltTree->Branch("recoMetGenSum",&mgensum,"recoMetGenSum/F");
00092 HltTree->Branch("recoHTCal",&htcalet,"recoHTCal/F");
00093 HltTree->Branch("recoHTCalPhi",&htcalphi,"recoHTCalPhi/F");
00094 HltTree->Branch("recoHTCalSum",&htcalsum,"recoHTCalSum/F");
00095
00096
00097
00098
00099 HltTree->Branch("NohTau",&nohtau,"NohTau/I");
00100 HltTree->Branch("ohTauEta",tauEta,"ohTauEta[NohTau]/F");
00101 HltTree->Branch("ohTauPhi",tauPhi,"ohTauPhi[NohTau]/F");
00102 HltTree->Branch("ohTauPt",tauPt,"ohTauPt[NohTau]/F");
00103 HltTree->Branch("ohTauEiso",l2tauemiso,"ohTauEiso[NohTau]/F");
00104 HltTree->Branch("ohTauL25Tpt",l25tauPt,"ohTauL25Tpt[NohTau]/F");
00105 HltTree->Branch("ohTauL25Tiso",l25tautckiso,"ohTauL25Tiso[NohTau]/I");
00106 HltTree->Branch("ohTauL3Tpt",l3tauPt,"ohTauL3Tpt[NohTau]/F");
00107 HltTree->Branch("ohTauL3Tiso",l3tautckiso,"ohTauL3Tiso[NohTau]/I");
00108
00109 }
00110
00111
00112 void HLTJets::analyze(const edm::Handle<CaloJetCollection> & calojets,
00113 const edm::Handle<GenJetCollection> & genjets,
00114 const edm::Handle<CaloMETCollection> & recmets,
00115 const edm::Handle<GenMETCollection> & genmets,
00116 const edm::Handle<METCollection> & ht,
00117 const edm::Handle<reco::HLTTauCollection> & taujets,
00118 const edm::Handle<CaloTowerCollection> & caloTowers,
00119 TTree * HltTree) {
00120
00121 if (_Debug) std::cout << " Beginning HLTJets " << std::endl;
00122
00123
00124 njetcal=0; njetgen=0;ntowcal=0;
00125 mcalmet=0.; mcalphi=0.;
00126 mgenmet=0.; mgenphi=0.;
00127 htcalet=0.,htcalphi=0.,htcalsum=0.;
00128
00129 if (calojets.isValid()) {
00130 CaloJetCollection mycalojets;
00131 mycalojets=*calojets;
00132 std::sort(mycalojets.begin(),mycalojets.end(),PtGreater());
00133 typedef CaloJetCollection::const_iterator cjiter;
00134 int jcal=0;
00135 for ( cjiter i=mycalojets.begin(); i!=mycalojets.end(); i++) {
00136
00137 if (i->pt()>_CalJetMin){
00138 jcalpt[jcal] = i->pt();
00139 jcalphi[jcal] = i->phi();
00140 jcaleta[jcal] = i->eta();
00141 jcale[jcal] = i->energy();
00142 jcal++;
00143 }
00144
00145 }
00146 njetcal = jcal;
00147 }
00148 else {njetcal = 0;}
00149
00150 if (caloTowers.isValid()) {
00151 ntowcal = caloTowers->size();
00152 int jtow = 0;
00153 for ( CaloTowerCollection::const_iterator tower=caloTowers->begin(); tower!=caloTowers->end(); tower++) {
00154 towet[jtow] = tower->et();
00155 toweta[jtow] = tower->eta();
00156 towphi[jtow] = tower->phi();
00157 towen[jtow] = tower->energy();
00158 towem[jtow] = tower->emEnergy();
00159 towhd[jtow] = tower->hadEnergy();
00160 towoe[jtow] = tower->outerEnergy();
00161 jtow++;
00162 }
00163 }
00164 else {ntowcal = 0;}
00165
00166 if (recmets.isValid()) {
00167 typedef CaloMETCollection::const_iterator cmiter;
00168 for ( cmiter i=recmets->begin(); i!=recmets->end(); i++) {
00169 mcalmet = i->pt();
00170 mcalphi = i->phi();
00171 mcalsum = i->sumEt();
00172 }
00173 }
00174
00175 if (ht.isValid()) {
00176 typedef METCollection::const_iterator iter;
00177 for ( iter i=ht->begin(); i!=ht->end(); i++) {
00178 htcalet = i->pt();
00179 htcalphi = i->phi();
00180 htcalsum = i->sumEt();
00181 }
00182 }
00183
00184 if (_Monte){
00185
00186 if (genjets.isValid()) {
00187 GenJetCollection mygenjets;
00188 mygenjets=*genjets;
00189 std::sort(mygenjets.begin(),mygenjets.end(),PtGreater());
00190 typedef GenJetCollection::const_iterator gjiter;
00191 int jgen=0;
00192 for ( gjiter i=mygenjets.begin(); i!=mygenjets.end(); i++) {
00193
00194 if (i->pt()>_GenJetMin){
00195 jgenpt[jgen] = i->pt();
00196 jgenphi[jgen] = i->phi();
00197 jgeneta[jgen] = i->eta();
00198 jgene[jgen] = i->energy();
00199 jgen++;
00200 }
00201
00202 }
00203 njetgen = jgen;
00204 }
00205 else {njetgen = 0;}
00206
00207 if (genmets.isValid()) {
00208 typedef GenMETCollection::const_iterator gmiter;
00209 for ( gmiter i=genmets->begin(); i!=genmets->end(); i++) {
00210 mgenmet = i->pt();
00211 mgenphi = i->phi();
00212 mgensum = i->sumEt();
00213 }
00214 }
00215
00216 }
00217
00218
00220
00221 if (taujets.isValid()) {
00222 nohtau = taujets->size();
00223 HLTTauCollection mytaujets;
00224 mytaujets=*taujets;
00225 std::sort(mytaujets.begin(),mytaujets.end(),GetPtGreater());
00226 typedef HLTTauCollection::const_iterator tauit;
00227 int itau=0;
00228 for(tauit i=mytaujets.begin(); i!=mytaujets.end(); i++){
00229
00230 tauEta[itau] = i->getEta();
00231 tauPhi[itau] = i->getPhi();
00232 tauPt[itau] = i->getPt();
00233
00234 l2tauemiso[itau] = i->getEMIsolationValue();
00235
00236 l25tauPt[itau] = i->getL25LeadTrackPtValue();
00237
00238 l25tautckiso[itau] = i->getL25TrackIsolationResponse();
00239
00240 l3tauPt[itau] = i->getL3LeadTrackPtValue();
00241
00242 l3tautckiso[itau] = i->getL3TrackIsolationResponse();
00243
00244 itau++;
00245 }
00246 }
00247 else {nohtau = 0;}
00248
00249 }