00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015 #include "RecoJets/JetAnalyzers/interface/JetToDigiDump.h"
00016 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00017 #include "DataFormats/JetReco/interface/CaloJet.h"
00018 #include "DataFormats/CaloTowers/interface/CaloTowerCollection.h"
00019 #include "DataFormats/DetId/interface/DetId.h"
00020 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
00021 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00022 #include "DataFormats/HcalRecHit/interface/HBHERecHit.h"
00023 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00024 #include "DataFormats/EcalDetId/interface/EEDetId.h"
00025 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00026 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
00027 #include "DataFormats/EcalDigi/interface/EBDataFrame.h"
00028 #include "FWCore/Framework/interface/Event.h"
00029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00030 #include <TROOT.h>
00031 #include <TSystem.h>
00032 #include <TFile.h>
00033 #include <TCanvas.h>
00034 #include <cmath>
00035 using namespace edm;
00036 using namespace reco;
00037 using namespace std;
00038
00039 JetToDigiDump::JetToDigiDump( const ParameterSet & cfg ) :
00040 DumpLevel( cfg.getParameter<string>( "DumpLevel" ) ),
00041 CaloJetAlg( cfg.getParameter<string>( "CaloJetAlg" ) )
00042 {
00043 }
00044
00045 void JetToDigiDump::beginJob( const EventSetup & ) {
00046 if(DumpLevel=="Jets")
00047 {
00048 cout << "Dump of Jets" << endl;
00049 Dump=1;
00050 }
00051 else if(DumpLevel=="Towers")
00052 {
00053 cout << "Dump of Jets and constituent CaloTowers" << endl;
00054 Dump=2;
00055 }
00056 else if(DumpLevel=="RecHits")
00057 {
00058 cout << "Dump of Jets, constituent CaloTowers, and constituent RecHits" << endl;
00059 Dump=3;
00060 }
00061 else if(DumpLevel=="Digis")
00062 {
00063 cout << "Dump of Jets, constituent CaloTowers, constituent RecHits and associated Digis" << endl;
00064 Dump=4;
00065 }
00066 cout << "Jet Algorithm being dumped is " << CaloJetAlg << endl;
00067
00068 evtCount = 0;
00069 }
00070
00071 void JetToDigiDump::analyze( const Event& evt, const EventSetup& es ) {
00072
00073 int jetInd;
00074 Handle<CaloJetCollection> caloJets;
00075 Handle<CaloTowerCollection> caloTowers;
00076 Handle<HBHERecHitCollection> HBHERecHits;
00077 Handle<HORecHitCollection> HORecHits;
00078 Handle<HFRecHitCollection> HFRecHits;
00079 Handle<EBRecHitCollection> EBRecHits;
00080 Handle<EERecHitCollection> EERecHits;
00081 Handle<HBHEDigiCollection> HBHEDigis;
00082 Handle<HODigiCollection> HODigis;
00083 Handle<HFDigiCollection> HFDigis;
00084 Handle<edm::SortedCollection<EBDataFrame> > EBDigis;
00085
00086
00087 evt.getByLabel( CaloJetAlg, caloJets );
00088 evt.getByLabel( "towerMaker", caloTowers );
00089 evt.getByLabel( "hbhereco", HBHERecHits );
00090 evt.getByLabel( "horeco", HORecHits );
00091 evt.getByLabel( "hfreco", HFRecHits );
00092 evt.getByLabel( "ecalRecHit", "EcalRecHitsEB", EBRecHits );
00093 evt.getByLabel( "ecalRecHit", "EcalRecHitsEE", EERecHits );
00094 evt.getByLabel( "hcaldigi", HBHEDigis );
00095 evt.getByLabel( "hcaldigi", HODigis );
00096 evt.getByLabel( "hcaldigi", HFDigis );
00097 evt.getByLabel( "ecalSelectiveReadout", "ebDigis", EBDigis );
00098
00099 cout << endl << "Evt: "<<evtCount <<", Num Jets=" <<caloJets->end() - caloJets->begin() << endl;
00100 if(Dump>=1)cout <<" *********************************************************" <<endl;
00101 jetInd = 0;
00102 if(Dump>=1)for( CaloJetCollection::const_iterator jet = caloJets->begin(); jet != caloJets->end(); ++ jet ) {
00103 cout <<" Jet: "<<jetInd<<", eta="<<jet->eta()<<", phi="<<jet->phi()<<", pt="<<jet->pt()<<\
00104 ",E="<<jet->energy()<<", EB E="<<jet->emEnergyInEB()<<" ,HB E="<<jet->hadEnergyInHB()<<\
00105 ", HO E="<<jet->hadEnergyInHO()<<" ,EE E="<< jet->emEnergyInEE()\
00106 <<", HE E="<<jet->hadEnergyInHE()<<", HF E="<<jet->hadEnergyInHF()+jet->emEnergyInHF()<<", Num Towers="<<jet->nConstituents()<<endl;
00107 if(Dump>=2)cout <<" ====================================================="<<endl;
00108 float sumTowerE = 0.0;
00109
00110 if(Dump>=2)for (int i = 0; i <jet->nConstituents(); i++) {
00111 const CaloTower& tower = *(jet->getCaloConstituent (i));
00112 int ietaTower = tower.id().ieta();
00113 int iphiTower = tower.id().iphi();
00114 sumTowerE += tower.energy();
00115 size_t numRecHits = tower.constituentsSize();
00116 cout << " Tower " << i <<": ieta=" << ietaTower << ", eta=" << tower.eta() <<", iphi=" << iphiTower << ", phi=" << tower.phi() << \
00117 ", energy=" << tower.energy() << ", EM=" << tower.emEnergy()<< ", HAD=" << tower.hadEnergy()\
00118 << ", HO=" << tower.outerEnergy() <<", Num Rec Hits =" << numRecHits << endl;
00119 if(Dump>=3)cout << " ------------------------------------------------"<<endl;
00120 float sumRecHitE = 0.0;
00121 if(Dump>=3)for(size_t j = 0; j <numRecHits ; j++) {
00122 DetId RecHitDetID=tower.constituent(j);
00123 DetId::Detector DetNum=RecHitDetID.det();
00124 if( DetNum == DetId::Hcal ){
00125
00126 HcalDetId HcalID = RecHitDetID;
00127 HcalSubdetector HcalNum = HcalID.subdet();
00128 if( HcalNum == HcalBarrel ){
00129 HBHERecHitCollection::const_iterator theRecHit=HBHERecHits->find(HcalID);
00130 sumRecHitE += theRecHit->energy();
00131 HBHEDigiCollection::const_iterator theDigis=HBHEDigis->find(HcalID);
00132 cout << " RecHit: " << j << ": HB, ieta=" << HcalID.ieta() << ", iphi=" << HcalID.iphi()<<
00133 ", depth=" << HcalID.depth() << ", energy=" << theRecHit->energy() << ", time=" <<\
00134 theRecHit->time() <<", All Digis=" << theDigis->size() << ", presamples =" <<\
00135 theDigis->presamples() <<endl;
00136
00137
00138
00139
00140
00141 float SumDigiCharge = 0.0;
00142 float EstimatedPedestal=0.0;
00143 int SamplesToAdd = 4;
00144 if(Dump>=4)cout << " ......................................"<<endl;
00145 if(Dump>=4)for(int k=0; k<theDigis->size(); k++){
00146 const HcalQIESample QIE = theDigis->sample(k);
00147 if(k>=theDigis->presamples()&&k<theDigis->presamples()+SamplesToAdd)SumDigiCharge+=QIE.nominal_fC();
00148 if(k<theDigis->presamples()-1)EstimatedPedestal+=QIE.nominal_fC()*SamplesToAdd/(theDigis->presamples()-1);
00149 cout << " Digi: " << k << ", cap ID = " << QIE.capid() << ": ADC Counts = " << QIE.adc() << ", nominal fC = " << QIE.nominal_fC() <<endl;
00150 }
00151 if(Dump>=4)cout << " 4 Digi fC ="<<SumDigiCharge<<", est. ped. fC="<<EstimatedPedestal<<", est. GeV/fc="<<theRecHit->energy()/(SumDigiCharge-EstimatedPedestal) << endl;
00152 if(Dump>=4)cout << " ......................................"<<endl;
00153 }
00154 else if( HcalNum == HcalEndcap ){
00155 HBHERecHitCollection::const_iterator theRecHit=HBHERecHits->find(HcalID);
00156 if( (abs(HcalID.ieta())==28||abs(HcalID.ieta())==29)&&HcalID.depth()==3){
00157 sumRecHitE += theRecHit->energy()/2;
00158 }
00159 else{
00160 sumRecHitE += theRecHit->energy();
00161 }
00162 HBHEDigiCollection::const_iterator theDigis=HBHEDigis->find(HcalID);
00163 cout << " RecHit: " << j << ": HE, ieta=" << HcalID.ieta() << ", iphi=" << HcalID.iphi()<<
00164 ", depth=" << HcalID.depth() << ", energy=" << theRecHit->energy() << ", time=" <<\
00165 theRecHit->time() <<", All Digis=" << theDigis->size() << ", presamples =" <<\
00166 theDigis->presamples() <<endl;
00167 float SumDigiCharge = 0.0;
00168 float EstimatedPedestal=0.0;
00169 int SamplesToAdd = 4;
00170 if(Dump>=4)cout << " ......................................"<<endl;
00171 if(Dump>=4)for(int k=0; k<theDigis->size(); k++){
00172 const HcalQIESample QIE = theDigis->sample(k);
00173 if(k>=theDigis->presamples()&&k<theDigis->presamples()+SamplesToAdd)SumDigiCharge+=QIE.nominal_fC();
00174 if(k<theDigis->presamples()-1)EstimatedPedestal+=QIE.nominal_fC()*SamplesToAdd/(theDigis->presamples()-1);
00175 cout << " Digi: " << k << ", cap ID = " << QIE.capid() << ": ADC Counts = " << QIE.adc() << ", nominal fC = " << QIE.nominal_fC() <<endl;
00176 }
00177 if(Dump>=4)cout << " 4 Digi fC ="<<SumDigiCharge<<", est. ped. fC="<<EstimatedPedestal<<", est. GeV/fc="<<theRecHit->energy()/(SumDigiCharge-EstimatedPedestal) << endl;
00178 if(Dump>=4)cout << " ......................................"<<endl;
00179 }
00180 else if( HcalNum == HcalOuter ){
00181 HORecHitCollection::const_iterator theRecHit=HORecHits->find(HcalID);
00182 sumRecHitE += theRecHit->energy();
00183 HODigiCollection::const_iterator theDigis=HODigis->find(HcalID);
00184 cout << " RecHit: " << j << ": HO, ieta=" << HcalID.ieta() << ", iphi=" << HcalID.iphi()<<
00185 ", depth=" << HcalID.depth() << ", energy=" << theRecHit->energy() << ", time=" <<\
00186 theRecHit->time() <<", All Digis=" << theDigis->size() << ", presamples =" <<\
00187 theDigis->presamples() <<endl;
00188 float SumDigiCharge = 0.0;
00189 float EstimatedPedestal=0.0;
00190 int SamplesToAdd = 4;
00191 if(Dump>=4)cout << " ......................................"<<endl;
00192 if(Dump>=4)for(int k=0; k<theDigis->size(); k++){
00193 const HcalQIESample QIE = theDigis->sample(k);
00194 if(k>=theDigis->presamples()&&k<theDigis->presamples()+SamplesToAdd)SumDigiCharge+=QIE.nominal_fC();
00195 if(k<theDigis->presamples()-1)EstimatedPedestal+=QIE.nominal_fC()*SamplesToAdd/(theDigis->presamples()-1);
00196 cout << " Digi: " << k << ", cap ID = " << QIE.capid() << ": ADC Counts = " << QIE.adc() << ", nominal fC = " << QIE.nominal_fC() <<endl;
00197 }
00198 if(Dump>=4)cout << " 4 Digi fC ="<<SumDigiCharge<<", est. ped. fC="<<EstimatedPedestal<<", est. GeV/fc="<<theRecHit->energy()/(SumDigiCharge-EstimatedPedestal) << endl;
00199 if(Dump>=4)cout << " ......................................"<<endl;
00200 }
00201 else if( HcalNum == HcalForward ){
00202 HFRecHitCollection::const_iterator theRecHit=HFRecHits->find(HcalID);
00203 sumRecHitE += theRecHit->energy();
00204 HFDigiCollection::const_iterator theDigis=HFDigis->find(HcalID);
00205 cout << " RecHit: " << j << ": HF, ieta=" << HcalID.ieta() << ", iphi=" << HcalID.iphi()<<
00206 ", depth=" << HcalID.depth() << ", energy=" << theRecHit->energy() << ", time=" <<\
00207 theRecHit->time() <<", All Digis=" << theDigis->size() << ", presamples =" <<\
00208 theDigis->presamples() <<endl;
00209 float SumDigiCharge = 0.0;
00210 float EstimatedPedestal=0.0;
00211 int SamplesToAdd = 1;
00212 if(Dump>=4)cout << " ......................................"<<endl;
00213 if(Dump>=4)for(int k=0; k<theDigis->size(); k++){
00214 const HcalQIESample QIE = theDigis->sample(k);
00215 if(k>=theDigis->presamples()&&k<theDigis->presamples()+SamplesToAdd)SumDigiCharge+=QIE.nominal_fC();
00216 if(k<theDigis->presamples()-1)EstimatedPedestal+=QIE.nominal_fC()*SamplesToAdd/(theDigis->presamples()-1);
00217 cout << " Digi: " << k << ", cap ID = " << QIE.capid() << ": ADC Counts = " << QIE.adc() << ", nominal fC = " << QIE.nominal_fC() <<endl;
00218 }
00219 if(Dump>=4)cout << " 1 Digi fC ="<<SumDigiCharge<<", est. ped. fC="<<EstimatedPedestal<<", est. GeV/fc="<<theRecHit->energy()/(SumDigiCharge-EstimatedPedestal) << endl;
00220 if(Dump>=4)cout << " ......................................"<<endl;
00221 }
00222 }
00223 if( DetNum == DetId::Ecal ){
00224 int EcalNum = RecHitDetID.subdetId();
00225 if( EcalNum == 1 ){
00226 EBDetId EcalID = RecHitDetID;
00227 EBRecHitCollection::const_iterator theRecHit=EBRecHits->find(EcalID);
00228 edm::SortedCollection<EBDataFrame>::const_iterator theDigis=EBDigis->find(EcalID);
00229 sumRecHitE += theRecHit->energy();
00230 cout << " RecHit " << j << ": EB, ieta=" << EcalID.ieta() << ", iphi=" << EcalID.iphi() << ", SM=" << EcalID.ism() << ", energy=" << theRecHit->energy() <<", All Digis=" << theDigis->size()<< endl;
00231 if(Dump>=4)cout << " ......................................"<<endl;
00232 if(Dump>=4)for(int k=0; k<theDigis->size(); k++){
00233 const EcalMGPASample MGPA = theDigis->sample(k);
00234 cout << " Digi: " << k << ": ADC Sample = " << MGPA.adc() << ", Gain ID = "<< MGPA.gainId() <<endl;
00235 }
00236 if(Dump>=4)cout << " ......................................"<<endl;
00237 }
00238 else if( EcalNum == 2 ){
00239 EEDetId EcalID = RecHitDetID;
00240 EERecHitCollection::const_iterator theRecHit=EERecHits->find(EcalID);
00241 sumRecHitE += theRecHit->energy();
00242 cout << " RecHit " << j << ": EE, ix=" << EcalID.ix() << ", iy=" << EcalID.iy() << ", energy=" << theRecHit->energy() << endl;
00243 }
00244 }
00245 }
00246 if(Dump>=3){
00247 if( abs(ietaTower)==28||abs(ietaTower)==29){
00248 cout << " Splitted Sum of RecHit Energies=" << sumRecHitE <<", CaloTower energy=" << tower.energy() << endl;
00249 }
00250 else{
00251 cout << " Sum of RecHit Energies=" << sumRecHitE <<", CaloTower energy=" << tower.energy() << endl;
00252 }
00253 }
00254 if(Dump>=3)cout << " ------------------------------------------------"<<endl;
00255 }
00256 if(Dump>=2)cout << " Sum of tower energies=" << sumTowerE << ", CaloJet energy=" << jet->energy() << endl;
00257 jetInd++;
00258 if(Dump>=2)cout <<" ====================================================="<<endl;
00259 }
00260 evtCount++;
00261 if(Dump>=1)cout <<" *********************************************************" <<endl;
00262
00263 }
00264
00265 void JetToDigiDump::endJob() {
00266
00267
00268 }
00269 #include "FWCore/Framework/interface/MakerMacros.h"
00270 DEFINE_FWK_MODULE(JetToDigiDump);