CMS 3D CMS Logo

HcalTemplateAnalysis.cc

Go to the documentation of this file.
00001 #include "DQM/HcalMonitorTasks/interface/HcalTemplateAnalysis.h"
00002 
00003 HcalTemplateAnalysis::HcalTemplateAnalysis() {
00004 }
00005 
00006 HcalTemplateAnalysis::~HcalTemplateAnalysis() {
00007 }
00008 
00009 void HcalTemplateAnalysis::reset(){}
00010 
00011 void HcalTemplateAnalysis::setup(const edm::ParameterSet& ps){
00012   
00013   outputFile_ = ps.getUntrackedParameter<string>("analysisFile", "");
00014   if ( outputFile_.size() != 0 ) 
00015     cout << "Hcal DQM analysis histograms will be saved to " << outputFile_.c_str() << endl;    
00016   else outputFile_ = "Hcal_DQM_Analysis.root";
00017   
00018   etaMax_ = ps.getUntrackedParameter<double>("MaxEta", 29.5);
00019   etaMin_ = ps.getUntrackedParameter<double>("MinEta", -29.5);
00020   etaBins_ = (int)(etaMax_ - etaMin_);
00021   
00022   phiMax_ = ps.getUntrackedParameter<double>("MaxPhi", 73);
00023   phiMin_ = ps.getUntrackedParameter<double>("MinPhi", 0);
00024   phiBins_ = (int)(phiMax_ - phiMin_);
00025 
00026   rechitEnergy_HB = new TH1F("HB RecHit Energies","HB RecHit Energies",250,0,250);
00027   rechitTime_HB = new TH1F("HB RecHit Times","HB RecHit Times",250,0,250);
00028   rechitEnergy_HF = new TH1F("HF RecHit Energies","HF RecHit Energies",250,0,250);
00029   rechitTime_HF = new TH1F("HF RecHit Times","HF RecHit Times",250,0,250);
00030   digiShape = new TH1F("HB Digi Shape","HB Digi Shape",20,0,19);
00031   digiOccupancy = new TH2F("HB Digi Occupancy","HB Digi Occupancy",etaBins_,etaMin_,etaMax_,phiBins_,phiMin_,phiMax_);
00032   return;
00033 }
00034 
00035 void HcalTemplateAnalysis::processEvent(const HBHEDigiCollection& hbhe,
00036                                         const HODigiCollection& ho,
00037                                         const HFDigiCollection& hf,
00038                                         const HBHERecHitCollection& hbHits, 
00039                                         const HORecHitCollection& hoHits,
00040                                         const HFRecHitCollection& hfHits,
00041                                         const LTCDigiCollection& ltc,
00042                                         const HcalDbService& cond){
00043 
00044 
00045   try{
00046     for (HBHEDigiCollection::const_iterator j=hbhe.begin(); j!=hbhe.end(); j++){
00047       const HBHEDataFrame digi = (const HBHEDataFrame)(*j);           
00048       for (int i=0; i<digi.size(); i++) {
00049         digiShape->Fill(i,digi.sample(i).adc());        
00050       }
00051       digiOccupancy->Fill(digi.id().ieta(),digi.id().iphi());
00052     }
00053   } catch (...) {
00054     printf("HcalTemplateAnalysis::processEvent  No HBHE Digis.\n");
00055   }
00056 
00057 
00058 
00059 
00060   HBHERecHitCollection::const_iterator _ib;
00061   //  HORecHitCollection::const_iterator _io;
00062   HFRecHitCollection::const_iterator _if;
00063   if(hbHits.size()>0){
00064     for (_ib=hbHits.begin(); _ib!=hbHits.end(); _ib++) { // loop over all hits
00065       if(_ib->energy()>0.0){
00066         if((HcalSubdetector)(_ib->id().subdet())==HcalBarrel){
00067           rechitEnergy_HB->Fill(_ib->energy());
00068           rechitTime_HB->Fill(_ib->time());
00069           
00070         }
00071       }
00072     }
00073   }
00074   if(hfHits.size()>0){
00075     for (_if=hfHits.begin(); _if!=hfHits.end(); _if++) { // loop over all hits
00076       if(_if->energy()>0.0){
00077         rechitEnergy_HF->Fill(_if->energy());
00078         rechitTime_HF->Fill(_if->time());
00079       }
00080     }
00081   }
00082 
00083   return;
00084 }
00085 
00086 void HcalTemplateAnalysis::done(){
00087 
00088   TFile *writefile = new TFile(outputFile_.c_str(),"RECREATE");
00089   writefile->cd();
00090   
00091   rechitEnergy_HB->Write();
00092   rechitTime_HB->Write();
00093   rechitEnergy_HF->Write();
00094   rechitTime_HF->Write();
00095   digiShape->Write();
00096   digiOccupancy->Write();  
00097   
00098   writefile->Write();
00099   writefile->Close();
00100 
00101 }

Generated on Tue Jun 9 17:33:02 2009 for CMSSW by  doxygen 1.5.4