#include <RecoTBCalo/HcalPlotter/src/HcalQLPlotHistoMgr.h>
Public Types | |
enum | EventType { UNKNOWN = 0, PEDESTAL = 1, LED = 2, LASER = 3, BEAM = 4 } |
enum | HistType { ENERGY = 0, TIME = 1, PULSE = 2, ADC = 3 } |
Public Member Functions | |
TH1 * | GetAHistogram (const HcalCalibDetId &id, const HcalElectronicsId &eid, HistType ht, EventType et) |
TH1 * | GetAHistogram (const HcalDetId &id, const HcalElectronicsId &eid, HistType ht, EventType et) |
HcalQLPlotHistoMgr (TDirectory *parent, const edm::ParameterSet &histoParams) | |
Static Public Member Functions | |
static std::string | nameForEvent (EventType et) |
static std::string | nameForFlavor (HistType ht) |
Private Member Functions | |
TH1 * | GetAHistogramImpl (const char *name, HistType ht, EventType et) |
Private Attributes | |
TDirectory * | beamHistDir |
TDirectory * | ctrHistDir |
edm::ParameterSet | histoParams_ |
TDirectory * | laserHistDir |
TDirectory * | ledHistDir |
TDirectory * | otherHistDir |
TDirectory * | pedHistDir |
Definition at line 12 of file HcalQLPlotHistoMgr.h.
HcalQLPlotHistoMgr::HcalQLPlotHistoMgr | ( | TDirectory * | parent, | |
const edm::ParameterSet & | histoParams | |||
) |
Definition at line 15 of file HcalQLPlotHistoMgr.cc.
References beamHistDir, histoParams_, laserHistDir, ledHistDir, otherHistDir, and pedHistDir.
00016 { 00017 pedHistDir=parent->mkdir("PEDESTAL"); 00018 ledHistDir=parent->mkdir("LED"); 00019 laserHistDir=parent->mkdir("LASER"); 00020 beamHistDir=parent->mkdir("BEAM"); 00021 otherHistDir=parent->mkdir("OTHER"); 00022 histoParams_ = histoParams; 00023 }
TH1 * HcalQLPlotHistoMgr::GetAHistogram | ( | const HcalCalibDetId & | id, | |
const HcalElectronicsId & | eid, | |||
HistType | ht, | |||
EventType | et | |||
) |
Definition at line 70 of file HcalQLPlotHistoMgr.cc.
References HcalElectronicsId::dccid(), HcalElectronicsId::fiberChanId(), HcalElectronicsId::fiberIndex(), GetAHistogramImpl(), HcalBarrel, HcalEndcap, HcalForward, HcalOuter, hcalSubdet(), HcalElectronicsId::htrSlot(), HcalElectronicsId::htrTopBottom(), name, nameForFlavor(), HcalElectronicsId::readoutVMECrateId(), and HcalElectronicsId::spigot().
00073 { 00074 std::string flavor=nameForFlavor(ht); 00075 00076 char name[120]; 00077 00078 std::string subdetStr; 00079 switch (id.hcalSubdet()) { 00080 case (HcalBarrel) : subdetStr="HB"; break; 00081 case (HcalEndcap) : subdetStr="HE"; break; 00082 case (HcalOuter) : subdetStr="HO"; break; 00083 case (HcalForward) : subdetStr="HF"; break; 00084 default: subdetStr="Other"; break; 00085 } 00086 00087 std::string chanstring = id.cboxChannelString(); 00088 if (!chanstring.size()) { 00089 chanstring = "Unknown"; 00090 edm::LogInfo("HcalQLPlotHistoMgr::GetAHistogram") << "Unknown calibration channel " << id.cboxChannel(); 00091 } 00092 00093 sprintf(name,"%s_CALIB_%s_%d_%d_chan=%s_eid=%d_%d_%d_%d_HTR_%d:%d%c", 00094 flavor.c_str(),subdetStr.c_str(),id.ieta(),id.iphi(),chanstring.c_str(), 00095 eid.dccid(),eid.spigot(), eid.fiberIndex(), eid.fiberChanId(), 00096 eid.readoutVMECrateId(), eid.htrSlot(),(eid.htrTopBottom()==1)?('t'):('b') ); 00097 00098 return GetAHistogramImpl(name,ht,et); 00099 }
TH1 * HcalQLPlotHistoMgr::GetAHistogram | ( | const HcalDetId & | id, | |
const HcalElectronicsId & | eid, | |||
HistType | ht, | |||
EventType | et | |||
) |
Definition at line 45 of file HcalQLPlotHistoMgr.cc.
References HcalElectronicsId::dccid(), HcalElectronicsId::fiberChanId(), HcalElectronicsId::fiberIndex(), GetAHistogramImpl(), HcalBarrel, HcalEndcap, HcalForward, HcalOuter, HcalElectronicsId::htrSlot(), HcalElectronicsId::htrTopBottom(), name, nameForFlavor(), HcalElectronicsId::readoutVMECrateId(), and HcalElectronicsId::spigot().
Referenced by HcalQLPlotAnalAlgos::processDigi(), and HcalQLPlotAnalAlgos::processRH().
00048 { 00049 std::string flavor=nameForFlavor(ht); 00050 00051 char name[120]; 00052 00053 std::string subdetStr; 00054 switch (id.subdet()) { 00055 case (HcalBarrel) : subdetStr="HB"; break; 00056 case (HcalEndcap) : subdetStr="HE"; break; 00057 case (HcalOuter) : subdetStr="HO"; break; 00058 case (HcalForward) : subdetStr="HF"; break; 00059 default: subdetStr="Other"; break; 00060 } 00061 00062 sprintf(name,"%s_%s_%d_%d_%d_eid=%d_%d_%d_%d_HTR_%d:%d%c", 00063 flavor.c_str(),subdetStr.c_str(),id.ieta(),id.iphi(),id.depth(), 00064 eid.dccid(),eid.spigot(), eid.fiberIndex(), eid.fiberChanId(), 00065 eid.readoutVMECrateId(), eid.htrSlot(),(eid.htrTopBottom()==1)?('t'):('b') ); 00066 00067 return GetAHistogramImpl(name,ht,et); 00068 }
TH1 * HcalQLPlotHistoMgr::GetAHistogramImpl | ( | const char * | name, | |
HistType | ht, | |||
EventType | et | |||
) | [private] |
Definition at line 101 of file HcalQLPlotHistoMgr.cc.
References BEAM, BEAM_BINS, beamHistDir, e, ENERGY, exception, edm::ParameterSet::getParameter(), histoParams_, LASER, LASER_BINS, laserHistDir, LED, LED_BINS, ledHistDir, OTHER_BINS, otherHistDir, PED_BINS, PEDESTAL, pedHistDir, PULSE, PULSE_BINS, getDQMSummary::td, TIME, TIME_BINS, and UNKNOWN.
Referenced by GetAHistogram().
00103 { 00104 TDirectory* td; 00105 00106 switch (et) { 00107 case(PEDESTAL): td=pedHistDir; break; 00108 case(LED): td=ledHistDir; break; 00109 case(LASER): td=laserHistDir; break; 00110 case(BEAM): td=beamHistDir; break; 00111 case(UNKNOWN): td=otherHistDir; break; 00112 default: td=0; break; 00113 } 00114 00115 if (td==0) { 00116 printf("Unknown %d !\n", et); 00117 return 0; 00118 } 00119 00120 TH1* retval=0; 00121 00122 retval=(TH1*)td->Get(name); 00123 int bins=0; double lo=0, hi=0; 00124 00125 // If the histogram doesn't exist and we are authorized, 00126 // create it! 00127 // 00128 if (retval==0) { 00129 td->cd(); 00130 switch (ht) { 00131 case(ENERGY): { 00132 switch (et) { 00133 case(PEDESTAL): 00134 bins=PED_BINS; 00135 try { 00136 lo=histoParams_.getParameter<double>("pedGeVlo"); 00137 hi=histoParams_.getParameter<double>("pedGeVhi"); 00138 } catch (std::exception& e) { // can't find it! 00139 edm::LogError("HcalQLPlotHistoMgr::GetAHistogram") << "Parameter(s) pedGeVlo/hi not found."; 00140 throw e; 00141 } 00142 break; 00143 case(LED): 00144 bins=LED_BINS; 00145 try { 00146 lo=histoParams_.getParameter<double>("ledGeVlo"); 00147 hi=histoParams_.getParameter<double>("ledGeVhi"); 00148 } catch (std::exception& e) { // can't find it! 00149 edm::LogError("HcalQLPlotHistoMgr::GetAHistogram") << "Parameter(s) ledGeVlo/hi not found."; 00150 throw e; 00151 } 00152 break; 00153 case(LASER): 00154 bins=LASER_BINS; 00155 try { 00156 lo=histoParams_.getParameter<double>("laserGeVlo"); 00157 hi=histoParams_.getParameter<double>("laserGeVhi"); 00158 } catch (std::exception& e) { // can't find it! 00159 edm::LogError("HcalQLPlotHistoMgr::GetAHistogram") << "Parameter(s) laserGeVlo/hi not found."; 00160 throw e; 00161 } 00162 break; 00163 case(BEAM): 00164 bins=BEAM_BINS; 00165 try { 00166 lo=histoParams_.getParameter<double>("beamGeVlo"); 00167 hi=histoParams_.getParameter<double>("beamGeVhi"); 00168 } catch (std::exception& e) { // can't find it! 00169 edm::LogError("HcalQLPlotHistoMgr::GetAHistogram") << "Parameter(s) beamGeVlo/hi not found."; 00170 throw e; 00171 } 00172 break; 00173 case(UNKNOWN): 00174 bins=OTHER_BINS; 00175 try { 00176 lo=histoParams_.getParameter<double>("otherGeVlo"); 00177 hi=histoParams_.getParameter<double>("otherGeVhi"); 00178 } catch (std::exception& e) { // can't find it! 00179 edm::LogError("HcalQLPlotHistoMgr::GetAHistogram") << "Parameter(s) otherGeVlo/hi not found."; 00180 throw e; 00181 } 00182 break; 00183 default: break; 00184 }; 00185 } 00186 break; 00187 case(TIME): 00188 bins=TIME_BINS; 00189 try { 00190 lo=histoParams_.getParameter<double>("timeNSlo"); 00191 hi=histoParams_.getParameter<double>("timeNShi"); 00192 } catch (std::exception& e) { // can't find it! 00193 edm::LogError("HcalQLPlotHistoMgr::GetAHistogram") << "Parameter(s) timeNSlo/hi not found."; 00194 throw e; 00195 } 00196 break; 00197 case(PULSE): 00198 bins=PULSE_BINS; 00199 lo=-0.5; 00200 hi=PULSE_BINS-0.5; 00201 break; 00202 case(ADC): 00203 bins=PED_BINS; 00204 try { 00205 lo=histoParams_.getParameter<double>("pedADClo"); 00206 hi=histoParams_.getParameter<double>("pedADChi"); 00207 } catch (std::exception& e) { // can't find it! 00208 edm::LogError("HcalQLPlotHistoMgr::GetAHistogram") << "Parameter(s) pedADClo/hi not found."; 00209 throw e; 00210 } 00211 break; 00212 } 00213 00214 if (bins>0){ 00215 if (ht==PULSE){ 00216 retval=new TProfile(name,name,bins,lo,hi); 00217 retval->GetXaxis()->SetTitle("TimeSlice(25ns)"); 00218 retval->GetYaxis()->SetTitle("fC"); 00219 } 00220 else if (ht==TIME){ 00221 retval=new TH1F(name,name,bins,lo,hi); 00222 retval->GetXaxis()->SetTitle("Timing(ns)"); 00223 } 00224 else if (ht==ENERGY){ 00225 retval=new TH1F(name,name,bins,lo,hi); 00226 retval->GetXaxis()->SetTitle("Energy(GeV)"); 00227 } 00228 else if (ht==ADC){ 00229 retval=new TH1F(name,name,bins,lo,hi); 00230 retval->GetXaxis()->SetTitle("ADC Counts"); 00231 } 00232 } 00233 } 00234 00235 return retval; 00236 }
std::string HcalQLPlotHistoMgr::nameForEvent | ( | EventType | et | ) | [static] |
Definition at line 35 of file HcalQLPlotHistoMgr.cc.
References BEAM, LASER, LED, and PEDESTAL.
00035 { 00036 switch(et) { 00037 case(PEDESTAL): return "Pedestal"; break; 00038 case(LED): return "LED"; break; 00039 case(LASER): return "Laser"; break; 00040 case(BEAM): return "Beam"; break; 00041 default: return "Other"; break; 00042 } 00043 }
std::string HcalQLPlotHistoMgr::nameForFlavor | ( | HistType | ht | ) | [static] |
Definition at line 25 of file HcalQLPlotHistoMgr.cc.
References ENERGY, PULSE, and TIME.
Referenced by GetAHistogram().
00025 { 00026 switch (ht) { 00027 case(ENERGY): return "Energy"; break; 00028 case(TIME): return "Time"; break; 00029 case(PULSE): return "Pulse"; break; 00030 case(ADC): return "ADC"; break; 00031 default: return ""; break; 00032 } 00033 }
TDirectory* HcalQLPlotHistoMgr::beamHistDir [private] |
Definition at line 34 of file HcalQLPlotHistoMgr.h.
Referenced by GetAHistogramImpl(), and HcalQLPlotHistoMgr().
TDirectory* HcalQLPlotHistoMgr::ctrHistDir [private] |
Definition at line 35 of file HcalQLPlotHistoMgr.h.
Definition at line 37 of file HcalQLPlotHistoMgr.h.
Referenced by GetAHistogramImpl(), and HcalQLPlotHistoMgr().
TDirectory* HcalQLPlotHistoMgr::laserHistDir [private] |
Definition at line 33 of file HcalQLPlotHistoMgr.h.
Referenced by GetAHistogramImpl(), and HcalQLPlotHistoMgr().
TDirectory* HcalQLPlotHistoMgr::ledHistDir [private] |
Definition at line 32 of file HcalQLPlotHistoMgr.h.
Referenced by GetAHistogramImpl(), and HcalQLPlotHistoMgr().
TDirectory* HcalQLPlotHistoMgr::otherHistDir [private] |
Definition at line 36 of file HcalQLPlotHistoMgr.h.
Referenced by GetAHistogramImpl(), and HcalQLPlotHistoMgr().
TDirectory* HcalQLPlotHistoMgr::pedHistDir [private] |
Definition at line 31 of file HcalQLPlotHistoMgr.h.
Referenced by GetAHistogramImpl(), and HcalQLPlotHistoMgr().