CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
CastorDigiMonitor Class Reference

#include <CastorDigiMonitor.h>

Inheritance diagram for CastorDigiMonitor:
CastorBaseMonitor

Public Member Functions

void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &)
 
 CastorDigiMonitor (const edm::ParameterSet &ps)
 
int ModSecToIndex (int module, int sector)
 
void processEvent (const CastorDigiCollection &cast, const CastorDbService &cond)
 
void setup (const edm::ParameterSet &ps)
 
 ~CastorDigiMonitor ()
 
- Public Member Functions inherited from CastorBaseMonitor
 CastorBaseMonitor ()
 
virtual ~CastorBaseMonitor ()
 

Private Attributes

int fVerbosity
 
MonitorElementh2digierr
 
MonitorElementh2QmeanMap
 
MonitorElementh2QmeantsvsCh
 
MonitorElementh2QrmsTSvsCh
 
MonitorElementh2QtsvsCh
 
MonitorElementh2reportMap
 
MonitorElementh2status
 
MonitorElementh2TSratio
 
MonitorElementhdigisize
 
MonitorElementhModule
 
MonitorElementhQIErms [10]
 
MonitorElementhReport
 
MonitorElementhSector
 
MonitorElementhTSratio
 
int ievt_
 
float QIEerrThreshold = 0.0001
 
double QmeanTS [224][10]
 
float Qrms_DEAD
 
double QrmsTS [224][10]
 
float RatioThresh1 = 0.
 
std::string subsystemname_
 
int TS_MAX = 10
 
const int TSped = 0
 

Additional Inherited Members

- Protected Attributes inherited from CastorBaseMonitor
std::string baseFolder_
 
edm::CPUTimer cpu_timer
 
int fVerbosity
 
std::string rootFolder_
 
bool showTiming
 

Detailed Description

Definition at line 16 of file CastorDigiMonitor.h.

Constructor & Destructor Documentation

CastorDigiMonitor::CastorDigiMonitor ( const edm::ParameterSet ps)

Definition at line 18 of file CastorDigiMonitor.cc.

References fVerbosity, edm::ParameterSet::getUntrackedParameter(), Qrms_DEAD, RatioThresh1, AlCaHLTBitMon_QueryRunRegistry::string, subsystemname_, and TS_MAX.

19 {
20  fVerbosity = ps.getUntrackedParameter<int>("debug",0);
21 subsystemname_=ps.getUntrackedParameter<std::string>("subSystemFolder","Castor");
22  RatioThresh1 = ps.getUntrackedParameter<double>("ratioThreshold",0.9);
23  Qrms_DEAD = ps.getUntrackedParameter<double>("QrmsDead",0.01); //fC
25  TS_MAX = ps.getUntrackedParameter<double>("qieTSmax",6);
26 }
T getUntrackedParameter(std::string const &, T const &) const
std::string subsystemname_
CastorDigiMonitor::~CastorDigiMonitor ( )

Definition at line 29 of file CastorDigiMonitor.cc.

29 { }

Member Function Documentation

void CastorDigiMonitor::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  iRun,
edm::EventSetup const &  iSetup 
)

Definition at line 42 of file CastorDigiMonitor.cc.

References DQMStore::IBooker::book1D(), DQMStore::IBooker::book2D(), DQMStore::IBooker::bookFloat(), DQMStore::IBooker::bookProfile2D(), gather_cfg::cout, fVerbosity, MonitorElement::getTH1F(), MonitorElement::getTH2F(), MonitorElement::getTProfile2D(), h2digierr, h2QmeanMap, h2QmeantsvsCh, h2QrmsTSvsCh, h2QtsvsCh, h2reportMap, h2status, h2TSratio, hdigisize, hModule, hQIErms, hReport, hSector, hTSratio, ievt_, QIEerrThreshold, QmeanTS, QrmsTS, alignCSCRings::s, indexGen::s2, DQMStore::IBooker::setCurrentFolder(), AlCaHLTBitMon_QueryRunRegistry::string, and subsystemname_.

Referenced by CastorMonitorModule::bookHistograms().

44 {
45  char s[60];
46  if(fVerbosity>0) std::cout << "CastorDigiMonitor::beginRun (start)" << std::endl;
47  char sTileIndex[50];
48  sprintf(sTileIndex,"Tile(=module*16+sector)");
49 
50  ievt_=0;
51 
52  ibooker.setCurrentFolder(subsystemname_ + "/CastorDigiMonitor");
53 
54  std::string s2 = "CASTOR QIE_capID+er+dv";
55  h2digierr=ibooker.bookProfile2D(s2,s2,14,0.,14., 16,0.,16.,100,0,1.e10,"");
56  h2digierr->getTProfile2D()->GetXaxis()->SetTitle("ModuleZ");
57  h2digierr->getTProfile2D()->GetYaxis()->SetTitle("SectorPhi");
58  h2digierr->getTProfile2D()->SetMaximum(1.);
59  h2digierr->getTProfile2D()->SetMinimum(QIEerrThreshold);
60  h2digierr->getTProfile2D()->SetOption("colz");
61 
62  sprintf(s,"CASTOR DeadChannelsMap");
63  h2status = ibooker.book2D(s,s,14,0.,14., 16,0.,16.);
64  h2status->getTH2F()->GetXaxis()->SetTitle("ModuleZ");
65  h2status->getTH2F()->GetYaxis()->SetTitle("SectorPhi");
66  h2status->getTH2F()->SetOption("colz");
67 
68  sprintf(s,"CASTOR AverageToMaxRatioMap");
69  h2TSratio = ibooker.book2D(s,s,14,0.,14., 16,0.,16.);
70  h2TSratio->getTH2F()->GetXaxis()->SetTitle("ModuleZ");
71  h2TSratio->getTH2F()->GetYaxis()->SetTitle("SectorPhi");
72  h2TSratio->getTH2F()->SetOption("colz");
73 
74  sprintf(s,"CASTOR AverageToMaxRatio");
75  hTSratio = ibooker.book1D(s,s,100,0.,1.);
76 
77  sprintf(s,"DigiSize");
78  hdigisize = ibooker.book1D(s,s,20,0.,20.);
79  sprintf(s,"Module(fC)_allTS");
80  hModule = ibooker.book1D(s,s,14,0.,14.);
81  hModule->getTH1F()->GetXaxis()->SetTitle("ModuleZ");
82  hModule->getTH1F()->GetYaxis()->SetTitle("QIE(fC)");
83  sprintf(s,"Sector(fC)_allTS");
84  hSector = ibooker.book1D(s,s,16,0.,16.);
85  hSector->getTH1F()->GetXaxis()->SetTitle("SectorPhi");
86  hSector->getTH1F()->GetYaxis()->SetTitle("QIE(fC)");
87 
88  sprintf(s,"QfC=f(x=Tile y=TS) (cumulative)");
89  h2QtsvsCh = ibooker.book2D(s,s,224,0.,224., 10,0.,10.);
90  h2QtsvsCh->getTH2F()->GetXaxis()->SetTitle(sTileIndex);
91  h2QtsvsCh->getTH2F()->GetYaxis()->SetTitle("TS");
92  h2QtsvsCh->getTH2F()->SetOption("colz");
93 
94  sprintf(s,"QmeanfC=f(Tile TS)");
95  h2QmeantsvsCh = ibooker.book2D(s,s,224,0.,224., 10,0.,10.);
96  h2QmeantsvsCh->getTH2F()->GetXaxis()->SetTitle(sTileIndex);
97  h2QmeantsvsCh->getTH2F()->GetYaxis()->SetTitle("Time Slice");
98  h2QmeantsvsCh->getTH2F()->SetOption("colz");
99 
100  sprintf(s,"QrmsfC=f(Tile TS)");
101  h2QrmsTSvsCh = ibooker.book2D(s,s,224,0.,224., 10,0.,10.);
102  h2QrmsTSvsCh->getTH2F()->GetXaxis()->SetTitle(sTileIndex);
103  h2QrmsTSvsCh->getTH2F()->GetYaxis()->SetTitle("TS");
104  h2QrmsTSvsCh->getTH2F()->SetOption("colz");
105 
106  sprintf(s,"CASTORreportSummaryMap");
107  h2reportMap = ibooker.book2D(s,s,14, 0,14, 16, 0,16);
108  h2reportMap->getTH2F()->GetXaxis()->SetTitle("moduleZ");
109  h2reportMap->getTH2F()->GetYaxis()->SetTitle("sectorPhi");
110  h2reportMap->getTH2F()->SetOption("colz");
111 
112  hReport = ibooker.bookFloat("CASTOR reportSummary");
113 
114  sprintf(s,"QmeanfC_map(allTS)");
115  h2QmeanMap = ibooker.book2D(s,s,14,0.,14., 16,0.,16.);
116  h2QmeanMap->getTH2F()->GetXaxis()->SetTitle("ModuleZ");
117  h2QmeanMap->getTH2F()->GetYaxis()->SetTitle("SectorPhi");
118  h2QmeanMap->getTH2F()->SetOption("textcolz");
119 
120  for(int ts=0; ts<=1; ts++) {
121  sprintf(s,"QIErms_TS=%d",ts);
122  hQIErms[ts] = ibooker.book1D(s,s,1000,0.,100.);
123  hQIErms[ts]->getTH1F()->GetXaxis()->SetTitle("QIErms(fC)");
124  }
125 
126  for(int ind=0; ind<224; ind++) for(int ts=0; ts<10; ts++)
127  QrmsTS[ind][ts] = QmeanTS[ind][ts]= 0.;
128 
129  if(fVerbosity>0) std::cout<<"CastorDigiMonitor::beginRun(end)"<<std::endl;
130  return;
131 }
MonitorElement * h2TSratio
double QrmsTS[224][10]
TProfile2D * getTProfile2D(void) const
MonitorElement * h2reportMap
MonitorElement * hTSratio
tuple s2
Definition: indexGen.py:106
MonitorElement * h2status
MonitorElement * h2digierr
MonitorElement * hModule
MonitorElement * bookProfile2D(Args &&...args)
Definition: DQMStore.h:163
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
MonitorElement * h2QmeanMap
MonitorElement * h2QtsvsCh
std::string subsystemname_
MonitorElement * h2QrmsTSvsCh
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:276
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
TH1F * getTH1F(void) const
MonitorElement * hdigisize
MonitorElement * hSector
double QmeanTS[224][10]
tuple cout
Definition: gather_cfg.py:145
MonitorElement * bookFloat(Args &&...args)
Definition: DQMStore.h:109
MonitorElement * hQIErms[10]
TH2F * getTH2F(void) const
MonitorElement * hReport
MonitorElement * h2QmeantsvsCh
int CastorDigiMonitor::ModSecToIndex ( int  module,
int  sector 
)

Definition at line 253 of file CastorDigiMonitor.cc.

Referenced by processEvent().

253  {
254  int ind = sector + module*16;
255  if(ind>223) ind=223;
256  return(ind);
257 }
Definition: vlib.h:208
void CastorDigiMonitor::processEvent ( const CastorDigiCollection cast,
const CastorDbService cond 
)

Definition at line 135 of file CastorDigiMonitor.cc.

References a, HcalQIESample::adc(), edm::SortedCollection< T, SORT >::begin(), HcalQIESample::capid(), gather_cfg::cout, HcalQIESample::dv(), edm::SortedCollection< T, SORT >::end(), HcalQIESample::er(), MonitorElement::Fill(), fVerbosity, MonitorElement::getTH1F(), MonitorElement::getTH2F(), MonitorElement::getTProfile2D(), h2digierr, h2QmeanMap, h2QmeantsvsCh, h2QrmsTSvsCh, h2QtsvsCh, h2reportMap, h2status, h2TSratio, hdigisize, hModule, hQIErms, hReport, hSector, hTSratio, i, CastorDataFrame::id(), ievt_, isum, j, LedMonAdc2fc, visualization-live-secondInstance_cfg::m, mod(), ModSecToIndex(), HcalCastorDetId::module(), lumiQueryAPI::q, QIEerrThreshold, QmeanTS, Qrms_DEAD, QrmsTS, alignCSCRings::r, RatioThresh1, alignCSCRings::s, CastorDataFrame::sample(), HcalCastorDetId::sector(), CastorDataFrame::size(), edm::SortedCollection< T, SORT >::size(), mathSSE::sqrt(), TS_MAX, and TSped.

Referenced by CastorMonitorModule::analyze().

136  {
137  if(fVerbosity>0) std::cout << "CastorDigiMonitor::processEvent (begin)"<< std::endl;
138 
139  if(castorDigis.size() <= 0) {
140  if(fVerbosity>0) std::cout<<"CastorPSMonitor::processEvent NO Castor Digis"<<std::endl;
141  return;
142  }
143 
145  j!=castorDigis.end(); j++)
146  {
147  const CastorDataFrame digi = (const CastorDataFrame)(*j);
148 
149  int capid1 = digi.sample(0).capid();
150  hdigisize->Fill(digi.size());
151  double sum = 0.;
152  for (int i=0; i<digi.size(); i++) {
153  int module = digi.id().module()-1;
154  int sector = digi.id().sector()-1;
155  int capid = digi.sample(i).capid();
156  int dv = digi.sample(i).dv();
157  int er = digi.sample(i).er();
158  int rawd = digi.sample(i).adc();
159  rawd = rawd&0x7F;
160  int err = (capid != capid1) | er<<1 | (!dv)<<2; // =0
161  h2digierr->Fill(module,sector,err);
162 // if(err !=0) continue;
163  int ind = ModSecToIndex(module,sector);
164  h2QtsvsCh->Fill(ind,i,LedMonAdc2fc[rawd]);
165  float q = LedMonAdc2fc[rawd];
166  sum += q; // sum += LedMonAdc2fc[rawd];
167  QrmsTS[ind][i] += (q*q);
168  QmeanTS[ind][i] += q;
169  if(err != 0 && fVerbosity>0)
170  std::cout<<"event/idigi=" <<ievt_<<"/"<<i<<" cap=cap1_dv_er_err: "<<
171  capid <<"="<< capid1 <<" "<< dv <<" "<< er<<" "<< err << std::endl;
172  if(capid1 < 3) capid1 = capid+1;
173  else capid1 = 0;
174  }
175 // hBunchOcc->Fill(iBunch,sum);
176  } //end for(CastorDigiCollection::const_iterator ...
177 
178  ievt_++;
179 
180  const float repChanBAD = 0.9;
181  const float repChanWarning = 0.95;
182  if(ievt_ %100 != 0) return;
183  float ModuleSum[14], SectorSum[16];
184  for(int m=0; m<14; m++) ModuleSum[m]=0.;
185  for(int s=0; s<16; s++) SectorSum[s]=0.;
186  for(int mod=0; mod<14; mod++) for(int sec=0; sec<16; sec++) {
187  for(int ts=0; ts<=1; ts++) {
188  int ind = ModSecToIndex(mod,sec);
189  double Qmean = QmeanTS[ind][ts]/ievt_;
190  double Qrms = sqrt(QrmsTS[ind][ts]/ievt_ - Qmean*Qmean);
191  hQIErms[ts]->Fill(Qrms);
192  }
193 
194  double sum=0.;
195  for(int ts=1; ts<=TS_MAX; ts++) {
196  int ind = ModSecToIndex(mod,sec) + 1;
197  double a=h2QtsvsCh->getTH2F()->GetBinContent(ind,ts);
198  h2QmeantsvsCh->getTH2F()->SetBinContent(ind,ts,a/double(ievt_));
199  sum += a;
200  double Qmean = QmeanTS[ind-1][ts-1]/ievt_;
201  double Qrms = QrmsTS[ind-1][ts-1]/ievt_ - Qmean*Qmean;
202  h2QrmsTSvsCh->getTH2F()->SetBinContent(ind,ts,sqrt(Qrms));
203  }
204  sum /= double(ievt_);
205  ModuleSum[mod] += sum;
206  SectorSum[sec] += sum;
207  float isum = float(int(sum*10.+0.5))/10.;
208  h2QmeanMap->getTH2F()->SetBinContent(mod+1,sec+1,isum);
209  } // end for(int mod=0; mod<14; mod++) for(int sec=0; sec<16; sec++)
210 
211  for(int mod=0; mod<14; mod++)
212  hModule->getTH1F()->SetBinContent(mod+1,ModuleSum[mod]);
213  for(int sec=0; sec<16; sec++)
214  hSector->getTH1F()->SetBinContent(sec+1,SectorSum[sec]);
215 
216  int nGoodCh = 0;
217  for(int mod=0; mod<14; mod++) for(int sec=0; sec<16;sec++) {
218  int ind = ModSecToIndex(mod,sec);
219  double Qmean = QmeanTS[ind][TSped]/ievt_;
220  double Qrms = QrmsTS[ind][TSped]/ievt_ - Qmean*Qmean;
221  float ChanStatus = 0.;
222  if(Qrms < Qrms_DEAD) ChanStatus = 1.;
223  h2status->getTH2F()->SetBinContent(mod+1,sec+1,ChanStatus);
224 
225  int tsm = 0;
226  float am=0.;
227  for(int ts=0; ts<TS_MAX; ts++) {
228  float a = h2QmeantsvsCh->getTH2F()->GetBinContent(ind+1,ts+1);
229  if(am < a) {am = a; tsm = ts;}
230  }
231 
232  double sum = 0.;
233  for(int ts=0; ts<TS_MAX; ts++) if(ts != tsm)
234  sum += h2QmeantsvsCh->getTH2F()->GetBinContent(ind+1,ts+1);
235  float r = 0.;
236  if(am > 0.) r = sum/(TS_MAX-1)/am;
237  h2TSratio->getTH2F()->SetBinContent(mod+1,sec+1,r);
238  hTSratio->Fill(r);
239  float statusTS = 1.0;
240  if(r > RatioThresh1) statusTS = repChanWarning;
241  else if(r > 0.99) statusTS = repChanBAD;
242  float gChanStatus = statusTS;
243  if(ChanStatus > 0.) gChanStatus = repChanBAD; // RMS
244  if(h2digierr->getTProfile2D()->GetBinContent(mod+1,sec+1)>QIEerrThreshold)
245  gChanStatus = repChanBAD;
246  h2reportMap->getTH2F()->SetBinContent(mod+1,sec+1,gChanStatus);
247  if(gChanStatus > repChanBAD) ++nGoodCh;
248  }
249  hReport->Fill(float(nGoodCh)/224.);
250  return;
251  }
int i
Definition: DBlmapReader.cc:9
MonitorElement * h2TSratio
double QrmsTS[224][10]
TProfile2D * getTProfile2D(void) const
int sector() const
get the sector (1-16)
MonitorElement * h2reportMap
MonitorElement * hTSratio
int adc() const
get the ADC sample
Definition: HcalQIESample.h:22
const HcalQIESample & sample(int i) const
access a sample
std::vector< CastorDataFrame >::const_iterator const_iterator
bool dv() const
is the Data Valid bit set?
Definition: HcalQIESample.h:28
int module() const
get the module (1-2 for EM, 1-12 for HAD)
double isum
void Fill(long long x)
MonitorElement * h2status
MonitorElement * h2digierr
MonitorElement * hModule
T sqrt(T t)
Definition: SSEVec.h:18
int j
Definition: DBlmapReader.cc:9
MonitorElement * h2QmeanMap
MonitorElement * h2QtsvsCh
MonitorElement * h2QrmsTSvsCh
static const float LedMonAdc2fc[128]
int capid() const
get the Capacitor id
Definition: HcalQIESample.h:26
int ModSecToIndex(int module, int sector)
TH1F * getTH1F(void) const
MonitorElement * hdigisize
MonitorElement * hSector
double a
Definition: hdecay.h:121
double QmeanTS[224][10]
tuple cout
Definition: gather_cfg.py:145
MonitorElement * hQIErms[10]
TH2F * getTH2F(void) const
const HcalCastorDetId & id() const
T mod(const T &a, const T &b)
Definition: ecalDccMap.h:4
int size() const
total number of samples in the digi
Definition: vlib.h:208
bool er() const
is the error bit set?
Definition: HcalQIESample.h:30
MonitorElement * hReport
MonitorElement * h2QmeantsvsCh
void CastorDigiMonitor::setup ( const edm::ParameterSet ps)
virtual

Reimplemented from CastorBaseMonitor.

Definition at line 32 of file CastorDigiMonitor.cc.

References gather_cfg::cout, fVerbosity, edm::ParameterSet::getUntrackedParameter(), CastorBaseMonitor::setup(), AlCaHLTBitMon_QueryRunRegistry::string, and subsystemname_.

Referenced by CastorMonitorModule::CastorMonitorModule().

33 {
35  ps.getUntrackedParameter<std::string>("subSystemFolder","Castor");
37  if(fVerbosity>0) std::cout << "CastorDigiMonitor::setup (end)"<<std::endl;
38  return;
39 }
T getUntrackedParameter(std::string const &, T const &) const
std::string subsystemname_
virtual void setup(const edm::ParameterSet &ps)
tuple cout
Definition: gather_cfg.py:145

Member Data Documentation

int CastorDigiMonitor::fVerbosity
private

Definition at line 28 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), CastorDigiMonitor(), processEvent(), and setup().

MonitorElement* CastorDigiMonitor::h2digierr
private

Definition at line 37 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), and processEvent().

MonitorElement* CastorDigiMonitor::h2QmeanMap
private

Definition at line 42 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), and processEvent().

MonitorElement* CastorDigiMonitor::h2QmeantsvsCh
private

Definition at line 41 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), and processEvent().

MonitorElement* CastorDigiMonitor::h2QrmsTSvsCh
private

Definition at line 32 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), and processEvent().

MonitorElement* CastorDigiMonitor::h2QtsvsCh
private

Definition at line 40 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), and processEvent().

MonitorElement* CastorDigiMonitor::h2reportMap
private

Definition at line 38 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), and processEvent().

MonitorElement* CastorDigiMonitor::h2status
private

Definition at line 36 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), and processEvent().

MonitorElement* CastorDigiMonitor::h2TSratio
private

Definition at line 35 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), and processEvent().

MonitorElement* CastorDigiMonitor::hdigisize
private

Definition at line 45 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), and processEvent().

MonitorElement* CastorDigiMonitor::hModule
private

Definition at line 43 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), and processEvent().

MonitorElement* CastorDigiMonitor::hQIErms[10]
private

Definition at line 33 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), and processEvent().

MonitorElement* CastorDigiMonitor::hReport
private

Definition at line 39 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), and processEvent().

MonitorElement* CastorDigiMonitor::hSector
private

Definition at line 44 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), and processEvent().

MonitorElement* CastorDigiMonitor::hTSratio
private

Definition at line 34 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), and processEvent().

int CastorDigiMonitor::ievt_
private

Definition at line 29 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), and processEvent().

float CastorDigiMonitor::QIEerrThreshold = 0.0001
private

Definition at line 50 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), and processEvent().

double CastorDigiMonitor::QmeanTS[224][10]
private

Definition at line 51 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), and processEvent().

float CastorDigiMonitor::Qrms_DEAD
private

Definition at line 30 of file CastorDigiMonitor.h.

Referenced by CastorDigiMonitor(), and processEvent().

double CastorDigiMonitor::QrmsTS[224][10]
private

Definition at line 51 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), and processEvent().

float CastorDigiMonitor::RatioThresh1 = 0.
private

Definition at line 49 of file CastorDigiMonitor.h.

Referenced by CastorDigiMonitor(), and processEvent().

std::string CastorDigiMonitor::subsystemname_
private

Definition at line 27 of file CastorDigiMonitor.h.

Referenced by bookHistograms(), CastorDigiMonitor(), and setup().

int CastorDigiMonitor::TS_MAX = 10
private

Definition at line 48 of file CastorDigiMonitor.h.

Referenced by CastorDigiMonitor(), and processEvent().

const int CastorDigiMonitor::TSped = 0
private

Definition at line 52 of file CastorDigiMonitor.h.

Referenced by processEvent().