CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CastorDigiMonitor.cc
Go to the documentation of this file.
5 
6 //****************************************************//
7 //********** CastorDigiMonitor: ******************//
8 //********** Author: Dmytro Volyanskyy *************//
9 //********** Date : 29.08.2008 (first version) ******//
12 //****************************************************//
13 //---- critical revision 26.06.2014 (Vladimir Popov)
14 //==================================================================//
15 //======================= Constructor ==============================//
17 {
18 subsystemname_=ps.getUntrackedParameter<std::string>("subSystemFolder","Castor");
19  fVerbosity = ps.getUntrackedParameter<int>("debug",0);
20 }
21 
22 //======================= Destructor ===============================//
24 
25 //=========================== setup ===============//
27 {
29  ps.getUntrackedParameter<std::string>("subSystemFolder","Castor");
31  if(fVerbosity>0) std::cout << "CastorDigiMonitor::setup (end)"<<std::endl;
32  return;
33 }
34 
35 //================= bookHistograms ===================//
37  const edm::Run& iRun, const edm::EventSetup& iSetup)
38 {
39  char s[60];
40  if(fVerbosity>0) std::cout << "CastorDigiMonitor::beginRun (start)" << std::endl;
41  ievt_=0;
42 // pset_.getUntrackedParameter<std::string>("subSystemFolder","Castor");
43 // baseFolder_ = rootFolder_+"CastorDigiMonitor";
44 
45  ibooker.setCurrentFolder(subsystemname_ + "/CastorDigiMonitor");
46  std::string s2 = "QIE_capID+er+dv";
47  h2digierr = ibooker.book2D(s2,s2,14,0.,14., 16,0.,16.);
48  h2digierr->getTH2F()->GetXaxis()->SetTitle("ModuleZ");
49  h2digierr->getTH2F()->GetYaxis()->SetTitle("SectorPhi");
50  h2digierr->getTH2F()->SetOption("colz");
51 
52  sprintf(s,"BunchOccupancy(fC)_all_TS");
53  hBunchOcc = ibooker.bookProfile(s,s,4000,0,4000, 100,0,1.e10,"");
54  hBunchOcc->getTProfile()->GetXaxis()->SetTitle("BX");
55  hBunchOcc->getTProfile()->GetYaxis()->SetTitle("QIE(fC)");
56 
57  sprintf(s,"DigiSize");
58  hdigisize = ibooker.book1D(s,s,20,0.,20.);
59  sprintf(s,"Module(fC)_allTS");
60  hModule = ibooker.book1D(s,s,14,0.,14.);
61  hModule->getTH1F()->GetXaxis()->SetTitle("ModuleZ");
62  hModule->getTH1F()->GetYaxis()->SetTitle("QIE(fC)");
63  sprintf(s,"Sector(fC)_allTS");
64  hSector = ibooker.book1D(s,s,16,0.,16.);
65  hSector->getTH1F()->GetXaxis()->SetTitle("SectorPhi");
66  hSector->getTH1F()->GetYaxis()->SetTitle("QIE(fC)");
67  sprintf(s,"QfC=f(Tile TS) (cumulative)");
68  h2QtsvsCh = ibooker.book2D(s,s,224,0.,224., 10,0.,10.);
69  h2QtsvsCh->getTH2F()->GetXaxis()->SetTitle("Tile(=sector*14+module)");
70  h2QtsvsCh->getTH2F()->GetYaxis()->SetTitle("TS");
71  h2QtsvsCh->getTH2F()->SetOption("colz");
72  sprintf(s,"QmeanfC=f(Tile TS)");
73  h2QmeantsvsCh = ibooker.book2D(s,s,224,0.,224., 10,0.,10.);
74  h2QmeantsvsCh->getTH2F()->GetXaxis()->SetTitle("Tile(=sector*14+module)");
75  h2QmeantsvsCh->getTH2F()->GetYaxis()->SetTitle("TS");
76  h2QmeantsvsCh->getTH2F()->SetOption("colz");
77  sprintf(s,"QmeanfC_map(allTS)");
78  h2QmeanMap = ibooker.book2D(s,s,14,0.,14., 16,0.,16.);
79  h2QmeanMap->getTH2F()->GetXaxis()->SetTitle("ModuleZ");
80  h2QmeanMap->getTH2F()->GetYaxis()->SetTitle("SectorPhi");
81  h2QmeanMap->getTH2F()->SetOption("textcolz");
82 
83  if(fVerbosity>0) std::cout<<"CastorDigiMonitor::beginRun(end)"<<std::endl;
84  return;
85 }
86 
87 
88 //=============== processEvent =========
90  const CastorDbService& cond, int iBunch) {
91  if(fVerbosity>0) std::cout << "CastorDigiMonitor::processEvent (begin)"<< std::endl;
92 
93 // CaloSamples tool;
94 
95  if(castorDigis.size() <= 0) {
96  if(fVerbosity>0) std::cout<<"CastorPSMonitor::processEvent NO Castor Digis"<<std::endl;
97  return;
98  }
99 
100  for(CastorDigiCollection::const_iterator j=castorDigis.begin();
101  j!=castorDigis.end(); j++)
102  {
103  const CastorDataFrame digi = (const CastorDataFrame)(*j);
104 
105  int capid1 = digi.sample(0).capid();
106  hdigisize->Fill(digi.size());
107  double sum = 0.;
108 // for (int i=1; i<digi.size(); i++) {
109  for (int i=0; i<digi.size(); i++) {
110  int module = digi.id().module()-1;
111  int sector = digi.id().sector()-1;
112  if(capid1 < 3) capid1++;
113  else capid1 = 0;
114  int capid = digi.sample(i).capid();
115  int dv = digi.sample(i).dv();
116  int er = digi.sample(i).er();
117  int rawd = digi.sample(i).adc();
118  rawd = rawd&0x7F;
119  int err = (capid != capid1) | er<<1 | (!dv)<<2; // =0
120  if(err !=0) h2digierr->Fill(module,sector);
121  int ind = sector*14 + module;
122  h2QtsvsCh->Fill(ind,i,LedMonAdc2fc[rawd]);
123  sum += LedMonAdc2fc[rawd];
124 // if(err != 0 && fVerbosity>0)
125 // std::cout<<"event/idigi=" <<ievt_<<"/" <<i<< " cap_cap1_dv_er: " <<
126 // capid <<"="<< capid1 <<" "<< dv <<" "<< er<<" "<< err << std::endl;
127  capid1 = capid;
128  }
129  hBunchOcc->Fill(iBunch,sum);
130  } //end for(CastorDigiCollection::const_iterator ...
131 
132  ievt_++;
133  if(ievt_ %100 == 0) {
134  float ModuleSum[14], SectorSum[16];
135  for(int m=0; m<14; m++) ModuleSum[m]=0.;
136  for(int s=0; s<16; s++) SectorSum[s]=0.;
137  for(int mod=0; mod<14; mod++) for(int sec=0; sec<16; sec++) {
138  double sum=0.;
139  for(int ts=1; ts<=10; ts++) {
140  int ind = sec*14 + mod +1;
141  double a=h2QtsvsCh->getTH2F()->GetBinContent(ind,ts);
142  h2QmeantsvsCh->getTH2F()->SetBinContent(ind,ts,a/double(ievt_));
143  sum += a;
144  }
145  sum /= double(ievt_);
146  ModuleSum[mod] += sum;
147  SectorSum[sec] += sum;
148  float isum = float(int(sum*10.+0.5))/10.;
149  h2QmeanMap->getTH2F()->SetBinContent(mod+1,sec+1,isum);
150  } // end for(int mod=0; mod<14; mod++) for(int sec=0; sec<16; sec++)
151 
152  for(int mod=0; mod<14; mod++)
153  hModule->getTH1F()->SetBinContent(mod+1,ModuleSum[mod]);
154  for(int sec=0; sec<16; sec++)
155  hSector->getTH1F()->SetBinContent(sec+1,SectorSum[sec]);
156  } //end if(ievt_ %100 == 0) {
157 
158  if(fVerbosity>0) std::cout << "CastorDigiMonitor::processEvent (end)"<< std::endl;
159  return;
160  }
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
int sector() const
get the sector (1-16)
MonitorElement * bookProfile(Args &&...args)
Definition: DQMStore.h:157
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
tuple s2
Definition: indexGen.py:106
void Fill(long long x)
CastorDigiMonitor(const edm::ParameterSet &ps)
MonitorElement * h2digierr
MonitorElement * hModule
void processEvent(const CastorDigiCollection &cast, const CastorDbService &cond, int bunch)
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &)
int j
Definition: DBlmapReader.cc:9
MonitorElement * h2QmeanMap
void setup(const edm::ParameterSet &ps)
MonitorElement * h2QtsvsCh
std::string subsystemname_
const_iterator end() const
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
static const float LedMonAdc2fc[128]
int capid() const
get the Capacitor id
Definition: HcalQIESample.h:26
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
TH1F * getTH1F(void) const
MonitorElement * hdigisize
virtual void setup(const edm::ParameterSet &ps)
MonitorElement * hSector
size_type size() const
double a
Definition: hdecay.h:121
TProfile * getTProfile(void) const
tuple cout
Definition: gather_cfg.py:121
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
MonitorElement * hBunchOcc
Definition: vlib.h:208
bool er() const
is the error bit set?
Definition: HcalQIESample.h:30
const_iterator begin() const
Definition: Run.h:41
MonitorElement * h2QmeantsvsCh