CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RPCMonitorRaw.cc
Go to the documentation of this file.
2 
3 
8 
12 
17 
18 #include "TH1F.h"
19 #include "TH2F.h"
20 #include "TFile.h"
21 
22 #include <fstream>
23 #include <iostream>
24 #include <iomanip>
25 #include <bitset>
26 
27 typedef std::map< std::pair<int,int>, int>::const_iterator IT;
28 
30  : theConfig(cfg)
31 {
32  for (unsigned int i=0; i<10;i++) theWatchedErrorHistoPos[i]=0;
33  std::vector<int> algos = cfg.getUntrackedParameter<std::vector<int> >("watchedErrors");
34  for (std::vector<int>::const_iterator it=algos.begin();it!= algos.end(); ++it) {
35  unsigned int ialgo = *it;
36  if (ialgo < 10) theWatchedErrorHistoPos[ialgo]=1; // real position initialisain is in begin job. here mark just switched on.
37  }
38 }
39 
40 RPCMonitorRaw::~RPCMonitorRaw() { LogTrace("") << "RPCMonitorRaw destructor"; }
41 
43 {
44 
45 // Get DQM interface
46  DQMStore* theDMBE = edm::Service<DQMStore>().operator->();
47 
48  theDMBE->setCurrentFolder("RPC/LinkMonitor");
49 
50  me_t[0]=theDMBE->book1D("recordType_790",RPCRawDataCountsHistoMaker::emptyRecordTypeHisto(790));
51  me_t[1]=theDMBE->book1D("recordType_791",RPCRawDataCountsHistoMaker::emptyRecordTypeHisto(791));
52  me_t[2]=theDMBE->book1D("recordType_792",RPCRawDataCountsHistoMaker::emptyRecordTypeHisto(792));
53  for (int i=0;i<3;++i)me_t[i]->getTH1F()->SetStats(0);
54 
55  me_e[0]=theDMBE->book1D("readoutErrors_790",RPCRawDataCountsHistoMaker::emptyReadoutErrorHisto(790));
56  me_e[1]=theDMBE->book1D("readoutErrors_791",RPCRawDataCountsHistoMaker::emptyReadoutErrorHisto(791));
57  me_e[2]=theDMBE->book1D("readoutErrors_792",RPCRawDataCountsHistoMaker::emptyReadoutErrorHisto(792));
58  for (int i=0;i<3;++i)me_e[i]->getTH1F()->SetStats(0);
59 
60  me_mapGoodEvents=theDMBE->book2D("mapGoodRecords","mapGoodRecords",36,-0.5,35.5, 3, 789.5,792.5);
61  me_mapGoodEvents->getTH2F()->SetNdivisions(3,"y");
62  me_mapGoodEvents->getTH2F()->SetXTitle("rmb");
63  me_mapGoodEvents->getTH2F()->SetYTitle("fed");
64  me_mapGoodEvents->getTH2F()->SetStats(0);
65  me_mapBadEvents =theDMBE->book2D("mapErrorRecords", "mapErrorRecords", 36,-0.5,35.5, 3, 789.5,792.5);
66  me_mapBadEvents->getTH2F()->SetXTitle("fed");
67  me_mapBadEvents->getTH2F()->SetYTitle("rmb");
68  me_mapBadEvents->getTH2F()->SetNdivisions(3,"y");
69  me_mapBadEvents->getTH2F()->SetStats(0);
70 
71  for (unsigned int i=0; i<=9; ++i) {
73  for (unsigned int fed=790; fed <=792; ++fed) {
75  MonitorElement* watched = theDMBE->book2D(histo->GetName(),histo);
76  theWatchedErrorHistos[fed-790].push_back(watched);
78  }
79  }
80  }
81 
82 }
84 {
85 bool writeHistos = theConfig.getUntrackedParameter<bool>("writeHistograms", false);
86  if (writeHistos) {
87  std::string histoFile = theConfig.getUntrackedParameter<std::string>("histoFileName");
88  TFile f(histoFile.c_str(),"RECREATE");
89  for (int i=0; i<3; ++i) {
90  me_t[i]->getTH1F()->Write();
91  me_e[i]->getTH1F()->Write();
92  std::vector<MonitorElement* > & wh = theWatchedErrorHistos[i];
93  for (std::vector<MonitorElement* >::const_iterator it=wh.begin(); it != wh.end(); ++it)(*it)->getTH2F()->Write();
94  }
95  me_mapGoodEvents->getTH2F()->Write();
96  me_mapBadEvents->getTH2F()->Write();
97  edm::LogInfo(" END JOB, histos saved!");
98  f.Close();
99  }
100 }
101 
102 
104 {
105 
107  ev.getByType( rawCounts);
108  const RPCRawDataCounts & counts = *rawCounts.product();
109 
110  //
111  // record type
112  //
113  for (IT it=counts.theRecordTypes.begin(); it != counts.theRecordTypes.end(); ++it)
114  me_t[it->first.first-790]->Fill(it->first.second,it->second);
115 
116  //
117  // good events topology
118  //
119  for (IT it = counts.theGoodEvents.begin(); it != counts.theGoodEvents.end(); ++it)
120  me_mapGoodEvents->Fill(it->first.second, it->first.first, it->second);
121 
122  //
123  // bad events topology
124  //
125  for (IT it = counts.theBadEvents.begin(); it != counts.theBadEvents.end(); ++it)
126  me_mapBadEvents->Fill(it->first.second, it->first.first, it->second);
127 
128 
129  //
130  // readout errors
131  //
132  for (IT it=counts.theReadoutErrors.begin(); it != counts.theReadoutErrors.end(); ++it) {
133  rpcrawtodigi::ReadoutError error(it->first.second);
134  LinkBoardElectronicIndex ele = error.where();
136 
137  int fed = it->first.first;
138  me_e[fed-790]->Fill(type, it->second);
139 
140  //
141  // in addition fill location map for selected errors
142  //
143  int idx = theWatchedErrorHistoPos[type]-1;
144  if ( idx >= 0) {
145  std::vector<MonitorElement* > & wh = theWatchedErrorHistos[fed-790];
146  MonitorElement* me = wh[idx];
147  me->Fill(ele.dccInputChannelNum, ele.tbLinkInputNum, it->second);
148  }
149  }
150 
151 
152 // for (int i=0; i<3; ++i) {
153 // me_t[i]->update();
154 // me_e[i]->update();
155 // std::vector<MonitorElement* > & wh = theWatchedErrorHistos[i];
156 // for (std::vector<MonitorElement* >::iterator it=wh.begin(); it != wh.end(); ++it) (*it)->update();
157 // }
158 // me_mapGoodEvents->update();
159 // me_mapBadEvents->update();
160 
161 
162 
163 }
std::map< std::pair< int, int >, int > theReadoutErrors
LinkBoardElectronicIndex where() const
Definition: ReadoutError.cc:20
type
Definition: HCALResponse.h:22
TFile * histoFile
Definition: hcalCalib.cc:45
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:514
MonitorElement * me_mapBadEvents
Definition: RPCMonitorRaw.h:36
std::vector< MonitorElement * > theWatchedErrorHistos[3]
Definition: RPCMonitorRaw.h:39
bool getByType(Handle< PROD > &result) const
Definition: Event.h:397
tuple histo
Definition: trackerHits.py:12
static TH1F * emptyReadoutErrorHisto(int fedId)
ReadoutErrorType type() const
Definition: ReadoutError.cc:14
void Fill(long long x)
MonitorElement * me_mapGoodEvents
Definition: RPCMonitorRaw.h:35
MonitorElement * me_t[3]
Definition: RPCMonitorRaw.h:33
virtual void endJob()
static TH2F * emptyReadoutErrorMapHisto(int fedId, int type)
virtual ~RPCMonitorRaw()
edm::ParameterSet theConfig
Definition: RPCMonitorRaw.h:38
virtual void beginJob()
RPCMonitorRaw(const edm::ParameterSet &cfg)
double f[11][100]
TH1F * getTH1F(std::string name, std::string process, std::string rootfolder, DQMStore *dbe_, bool verb, bool clone)
std::map< std::pair< int, int >, int > theGoodEvents
std::vector< LinkConnSpec >::const_iterator IT
unsigned int theWatchedErrorHistoPos[10]
Definition: RPCMonitorRaw.h:42
#define LogTrace(id)
std::map< std::pair< int, int >, int > theRecordTypes
virtual void analyze(const edm::Event &, const edm::EventSetup &)
get data, convert to digis attach againe to Event
TH1F * getTH1F(void) const
T const * product() const
Definition: Handle.h:74
static TH1F * emptyRecordTypeHisto(int fedId)
std::map< std::pair< int, int >, int > theBadEvents
TH2F * getTH2F(void) const
MonitorElement * book2D(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Book 2D histogram.
Definition: DQMStore.cc:642
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:232
MonitorElement * me_e[3]
Definition: RPCMonitorRaw.h:34