#include <RecoTBCalo/HcalPlotter/src/HcalQLPlotAnalAlgos.h>
Public Member Functions | |
void | end (void) |
HcalQLPlotAnalAlgos (const char *outputFilename, edm::ParameterSet histoParams) | |
void | processDigi (const HcalCalibDigiCollection &calibdigic, double calibFC2GeV) |
void | processDigi (const HFDigiCollection &hfdigic) |
void | processDigi (const HODigiCollection &hodigic) |
void | processDigi (const HBHEDigiCollection &hbhedigic) |
void | processRH (const HFRecHitCollection &hfrhc, const HFDigiCollection &hfdgc) |
void | processRH (const HORecHitCollection &horhc, const HODigiCollection &hodgc) |
void | processRH (const HBHERecHitCollection &hbherhc, const HBHEDigiCollection &hbhedgc) |
void | SetEventType (const HcalTBTriggerData &trigd) |
Private Member Functions | |
HcalCalibRecHit | recoCalib (const HcalCalibDataFrame &cdigi, double calibFC2GeV) |
Private Attributes | |
HcalQLPlotHistoMgr * | histos_ |
TFile * | mf_ |
HcalQLPlotHistoMgr::EventType | triggerID_ |
Definition at line 17 of file HcalQLPlotAnalAlgos.h.
HcalQLPlotAnalAlgos::HcalQLPlotAnalAlgos | ( | const char * | outputFilename, | |
edm::ParameterSet | histoParams | |||
) |
Definition at line 44 of file HcalQLPlotAnalAlgos.cc.
References histos_, mf_, triggerID_, and HcalQLPlotHistoMgr::UNKNOWN.
00046 { 00047 triggerID_=HcalQLPlotHistoMgr::UNKNOWN; 00048 00049 mf_ = new TFile(outputFilename,"RECREATE"); 00050 histos_ = new HcalQLPlotHistoMgr(mf_,histoParams); 00051 }
Definition at line 58 of file HcalQLPlotAnalAlgos.cc.
References mf_.
Referenced by HcalQLPlotAnal::endJob().
00059 { 00060 mf_->Write(); 00061 }
void HcalQLPlotAnalAlgos::processDigi | ( | const HcalCalibDigiCollection & | calibdigic, | |
double | calibFC2GeV | |||
) |
Definition at line 279 of file HcalQLPlotAnalAlgos.cc.
References HcalCalibRecHit::amplitude(), edm::SortedCollection< T, SORT >::begin(), edm::SortedCollection< T, SORT >::end(), HcalQLPlotHistoMgr::ENERGY, HcalQLPlotHistoMgr::GetAHistogram(), histos_, id, HcalQLPlotHistoMgr::PULSE, recoCalib(), HcalQLPlotHistoMgr::TIME, HcalCalibRecHit::time(), and triggerID_.
00281 { 00282 HcalCalibDigiCollection::const_iterator it; 00283 00284 for (it = calibdigic.begin(); 00285 it != calibdigic.end(); 00286 it++) { 00287 HcalCalibDetId id (it->id()); 00288 HcalElectronicsId eid (it->elecId()); 00289 00290 TH1* phist=histos_->GetAHistogram(id,eid,HcalQLPlotHistoMgr::PULSE,triggerID_); 00291 if (phist){ 00292 for (int bin=0; bin<it->size(); bin++) 00293 phist->Fill(bin*1.0,(*it)[bin].nominal_fC()); 00294 } 00295 00296 // HACK-reco the calib digi into a rechit: 00297 // 00298 HcalCalibRecHit rh = recoCalib(*it, calibFC2GeV); 00299 00300 TH1* ehist=histos_->GetAHistogram(id,eid,HcalQLPlotHistoMgr::ENERGY,triggerID_); 00301 if (ehist){ 00302 ehist->Fill(rh.amplitude()); 00303 } 00304 00305 TH1* thist=histos_->GetAHistogram(id,eid,HcalQLPlotHistoMgr::TIME,triggerID_); 00306 if (thist){ 00307 thist->Fill(rh.time()); 00308 } 00309 } 00310 }
void HcalQLPlotAnalAlgos::processDigi | ( | const HFDigiCollection & | hfdigic | ) |
void HcalQLPlotAnalAlgos::processDigi | ( | const HODigiCollection & | hodigic | ) |
void HcalQLPlotAnalAlgos::processDigi | ( | const HBHEDigiCollection & | hbhedigic | ) |
Definition at line 174 of file HcalQLPlotAnalAlgos.cc.
References HcalQLPlotHistoMgr::ADC, edm::SortedCollection< T, SORT >::begin(), edm::SortedCollection< T, SORT >::end(), HcalQLPlotHistoMgr::GetAHistogram(), histos_, id, HcalQLPlotHistoMgr::PEDESTAL, HcalQLPlotHistoMgr::PULSE, and triggerID_.
Referenced by HcalQLPlotAnal::analyze().
00175 { 00176 HBHEDigiCollection::const_iterator it; 00177 00178 for (it = hbhedigic.begin(); 00179 it != hbhedigic.end(); 00180 it++) { 00181 HcalDetId id (it->id()); 00182 HcalElectronicsId eid (it->elecId()); 00183 00184 TH1* phist=histos_->GetAHistogram(id,eid,HcalQLPlotHistoMgr::PULSE,triggerID_); 00185 if (phist){ 00186 for (int bin=0; bin<it->size(); bin++) 00187 phist->Fill(bin*1.0,(*it)[bin].nominal_fC()); 00188 } 00189 00190 if (triggerID_ == HcalQLPlotHistoMgr::PEDESTAL) { 00191 phist=histos_->GetAHistogram(id,eid,HcalQLPlotHistoMgr::ADC,triggerID_); 00192 if (phist){ 00193 for (int bin=0; bin<it->size(); bin++) 00194 phist->Fill((*it)[bin].adc()); 00195 } 00196 } 00197 } 00198 }
void HcalQLPlotAnalAlgos::processRH | ( | const HFRecHitCollection & | hfrhc, | |
const HFDigiCollection & | hfdgc | |||
) |
void HcalQLPlotAnalAlgos::processRH | ( | const HORecHitCollection & | horhc, | |
const HODigiCollection & | hodgc | |||
) |
void HcalQLPlotAnalAlgos::processRH | ( | const HBHERecHitCollection & | hbherhc, | |
const HBHEDigiCollection & | hbhedgc | |||
) |
Definition at line 81 of file HcalQLPlotAnalAlgos.cc.
References edm::SortedCollection< T, SORT >::begin(), edm::SortedCollection< T, SORT >::end(), HcalQLPlotHistoMgr::ENERGY, edm::SortedCollection< T, SORT >::find(), HcalQLPlotHistoMgr::GetAHistogram(), histos_, id, it, HcalQLPlotHistoMgr::TIME, and triggerID_.
Referenced by HcalQLPlotAnal::analyze().
00083 { 00084 HBHERecHitCollection::const_iterator it; 00085 00086 for (it = hbherhc.begin(); 00087 it != hbherhc.end(); 00088 it++) { 00089 HcalDetId id (it->id()); 00090 HcalElectronicsId eid; 00091 HBHEDigiCollection::const_iterator dit = hbhedgc.find(id); 00092 if (dit != hbhedgc.end()) 00093 eid = dit->elecId(); 00094 else { 00095 edm::LogWarning("HcalQLPlotAnalAlgos::processRH") << 00096 "No digi found for id" << id; 00097 continue; 00098 } 00099 00100 TH1* ehist=histos_->GetAHistogram(id,eid,HcalQLPlotHistoMgr::ENERGY,triggerID_); 00101 if (ehist){ 00102 ehist->Fill(it->energy()); 00103 } 00104 00105 TH1* thist=histos_->GetAHistogram(id,eid,HcalQLPlotHistoMgr::TIME,triggerID_); 00106 if (thist){ 00107 thist->Fill(it->time()); 00108 } 00109 } 00110 }
HcalCalibRecHit HcalQLPlotAnalAlgos::recoCalib | ( | const HcalCalibDataFrame & | cdigi, | |
double | calibFC2GeV | |||
) | [private] |
Definition at line 252 of file HcalQLPlotAnalAlgos.cc.
References i, HcalCalibDataFrame::id(), HcalCalibDataFrame::presamples(), and HcalCalibDataFrame::size().
Referenced by processDigi().
00254 { 00255 double nominal_ped = (cdigi[0].nominal_fC() + cdigi[1].nominal_fC())/2.0; 00256 00257 double totamp = 0.0; 00258 double maxA = -1e99; 00259 int maxI = -1; 00260 for (int i=0; i<cdigi.size(); i++) { 00261 double ampl = (cdigi[i].nominal_fC()-nominal_ped)*calibFC2GeV; 00262 totamp += ampl; 00263 00264 if (ampl > maxA) { 00265 maxA = ampl; 00266 maxI = i; 00267 } 00268 } 00269 00270 maxA=fabs(maxA); 00271 float t0 = (maxI > 0) ? (fabs((cdigi[maxI-1].nominal_fC()-nominal_ped))*calibFC2GeV):0.0; 00272 float t2 = fabs((cdigi[maxI+1].nominal_fC()-nominal_ped)*calibFC2GeV); 00273 float wpksamp = (maxA + 2.0*t2) / (t0 + maxA + t2); 00274 float time = (maxI - cdigi.presamples() + wpksamp)*25.0; 00275 00276 return HcalCalibRecHit(cdigi.id(),totamp,time); 00277 }
void HcalQLPlotAnalAlgos::SetEventType | ( | const HcalTBTriggerData & | trigd | ) |
Definition at line 63 of file HcalQLPlotAnalAlgos.cc.
References HcalQLPlotHistoMgr::BEAM, e, exception, HcalQLPlotHistoMgr::LASER, HcalQLPlotHistoMgr::LED, HcalQLPlotHistoMgr::PEDESTAL, triggerID_, HcalQLPlotHistoMgr::UNKNOWN, HcalTBTriggerData::wasBeamTrigger(), HcalTBTriggerData::wasInSpillPedestalTrigger(), HcalTBTriggerData::wasLaserTrigger(), HcalTBTriggerData::wasLEDTrigger(), HcalTBTriggerData::wasOutSpillPedestalTrigger(), and HcalTBTriggerData::wasSpillIgnorantPedestalTrigger().
Referenced by HcalQLPlotAnal::analyze().
00064 { 00065 if( trigd.wasInSpillPedestalTrigger() || 00066 trigd.wasOutSpillPedestalTrigger() || 00067 trigd.wasSpillIgnorantPedestalTrigger() ) 00068 triggerID_=HcalQLPlotHistoMgr::PEDESTAL; 00069 if( trigd.wasLEDTrigger() ) triggerID_=HcalQLPlotHistoMgr::LED; 00070 if( trigd.wasLaserTrigger() ) triggerID_=HcalQLPlotHistoMgr::LASER; 00071 if( trigd.wasBeamTrigger() ) triggerID_=HcalQLPlotHistoMgr::BEAM; 00072 00073 if( triggerID_ == HcalQLPlotHistoMgr::UNKNOWN) { 00074 edm::LogError("HcalQLPlotAnalAlgos::begin") << 00075 "Trigger Type unrecognized, aborting"; 00076 std::exception e; 00077 throw e; 00078 } 00079 }
HcalQLPlotHistoMgr* HcalQLPlotAnalAlgos::histos_ [private] |
Definition at line 41 of file HcalQLPlotAnalAlgos.h.
Referenced by HcalQLPlotAnalAlgos(), processDigi(), and processRH().
TFile* HcalQLPlotAnalAlgos::mf_ [private] |
Definition at line 42 of file HcalQLPlotAnalAlgos.h.
Referenced by end(), and HcalQLPlotAnalAlgos().
Definition at line 40 of file HcalQLPlotAnalAlgos.h.
Referenced by HcalQLPlotAnalAlgos(), processDigi(), processRH(), and SetEventType().