CMS 3D CMS Logo

RPCDqmClient.cc
Go to the documentation of this file.
1 // Package: RPCDqmClient
2 // Original Author: Anna Cimmino
3 
7 //include client headers
13 //Geometry
17 //Framework
21 
23 
24  edm::LogVerbatim ("rpcdqmclient") << "[RPCDqmClient]: Constructor";
25 
26  offlineDQM_ = parameters_.getUntrackedParameter<bool> ("OfflineDQM",true);
27  useRollInfo_= parameters_.getUntrackedParameter<bool>("UseRollInfo", false);
28  //check enabling
29  enableDQMClients_ = parameters_.getUntrackedParameter<bool> ("EnableRPCDqmClient",true);
30  minimumEvents_= parameters_.getUntrackedParameter<int>("MinimumRPCEvents", 10000);
31 
32  std::string subsystemFolder = parameters_.getUntrackedParameter<std::string>("RPCFolder", "RPC");
33  std::string recHitTypeFolder= parameters_.getUntrackedParameter<std::string>("RecHitTypeFolder", "AllHits");
34  std::string summaryFolder = parameters_.getUntrackedParameter<std::string>("SummaryFolder", "SummaryHistograms");
35 
36  prefixDir_ = subsystemFolder+ "/"+ recHitTypeFolder;
37  globalFolder_ = subsystemFolder + "/"+ recHitTypeFolder + "/"+ summaryFolder;
38 
39  //get prescale factor
40  prescaleGlobalFactor_ = parameters_.getUntrackedParameter<int>("DiagnosticGlobalPrescale", 5);
41 
42  //make default client list
43  clientList_.push_back("RPCMultiplicityTest");
44  clientList_.push_back("RPCDeadChannelTest");
45  clientList_.push_back("RPCClusterSizeTest");
46  clientList_= parameters_.getUntrackedParameter<std::vector<std::string> >("RPCDqmClientList",clientList_);
47 
48  //get all the possible RPC DQM clients
49  this->makeClientMap(parameters_);
50 
51  //clear counters
52  lumiCounter_ = 0;
53 }
54 
56 
58 
59  if (!enableDQMClients_) {return;} ;
60  edm::LogVerbatim ("rpcdqmclient") << "[RPCDqmClient]: Begin Job";
61 
62  //Do whatever the begin jobs of all client modules do
63  for(std::vector<RPCClient*>::iterator it = clientModules_.begin(); it!=clientModules_.end(); it++ ){
64  (*it)->beginJob( globalFolder_ );
65  }
66 
67 }
68 
69 
70 
72 
73  if (!enableDQMClients_ ) {return;}
74  edm::LogVerbatim ("rpcdqmclient") << "[RPCDqmClient]: End DQM LB";
75 
76  if( myDetIds_.empty() ) {
77  //Get RPCdetId...
78 
79 
80  this->getRPCdetId( c);
81 
82  //...book summary histograms
83  for (std::vector<RPCClient*>::iterator it = clientModules_.begin(); it!=clientModules_.end(); it++ ){
84  (*it)->myBooker( ibooker);
85  }
86  }
87 
88  if (!offlineDQM_){ //Do this only for the online
89 
90  if (lumiCounter_ == 0){ //only for the first lumi section do this...
91  // ...get chamber based histograms and pass them to the client modules
92  this->getMonitorElements(igetter);
93  }
94 
95  //Do not perform client oparations every lumi block
96  lumiCounter_++;
97  if (lumiCounter_%prescaleGlobalFactor_ != 0) {return;}
98 
99  //Check if there's enough statistics
100  float rpcevents = minimumEvents_;
101  if(RPCEvents_) {rpcevents = RPCEvents_->getBinContent(1);}
102  if( rpcevents < minimumEvents_) {return;}
103 
104  edm::LogVerbatim ("rpcdqmclient") <<"[RPCDqmClient]: Client operations";
105  for (std::vector<RPCClient*>::iterator it = clientModules_.begin(); it!=clientModules_.end(); it++ ){
106  (*it)->clientOperation();
107  }
108  }//end of online operations
109 
110 
111 }
112 
113 
114 
116 
117  if (!enableDQMClients_) {return;}
118 
119  edm::LogVerbatim ("rpcdqmclient") << "[RPCDqmClient]: End DQM Job";
120 
121  if(offlineDQM_){ // ...get chamber based histograms and pass them to the client modules
122  this->getMonitorElements( igetter);
123  }
124 
125  float rpcevents = minimumEvents_;
126  if(RPCEvents_) {
127  rpcevents = RPCEvents_ ->getBinContent(1);}
128  if(rpcevents < minimumEvents_) {return;}
129 
130  edm::LogVerbatim ("rpcdqmclient") <<"[RPCDqmClient]: Client operations";
131  for (std::vector<RPCClient*>::iterator it = clientModules_.begin(); it!=clientModules_.end(); it++ ){
132  (*it)->clientOperation();
133  }
134 
135 }
136 
137 
138 
140 
141  std::vector<MonitorElement *> myMeVect;
142  std::vector<RPCDetId> myDetIds;
143 
144  //dbe_->setCurrentFolder(prefixDir_);
146  MonitorElement * myMe = nullptr;
147  std::string rollName= "";
148 
149  //loop on all geometry and get all histos
150  for ( auto& detId : myDetIds_ ) {
151  //Get name
152  RPCGeomServ RPCname(detId);
153  rollName = useRollInfo_ ? RPCname.name() : RPCname.chambername();
154 
155  //loop on clients
156  for( unsigned int cl = 0, nCL = clientModules_.size(); cl < nCL; ++cl ){
157  myMe = igetter.get(prefixDir_ +"/"+ folderStr->folderStructure(detId)+"/"+clientHisto_[cl]+ "_"+rollName);
158  if (!myMe){continue;}
159 
160  // dbe_->tag(myMe, clientTag_[cl]);
161  myMeVect.push_back(myMe);
162  myDetIds.push_back(detId);
163 
164  }//end loop on clients
165  }//end loop on all geometry and get all histos
166 
167  RPCEvents_ = igetter.get(prefixDir_ +"/RPCEvents");
168  for ( unsigned int cl = 0; cl < clientModules_.size(); ++cl ) {
169  clientModules_[cl]->getMonitorElements(myMeVect, myDetIds, clientHisto_[cl]);
170  }
171 
172  delete folderStr;
173 }
174 
175 
176 
178 
179  myDetIds_.clear();
180 
182  eventSetup.get<MuonGeometryRecord>().get(rpcGeo);
183 
184  for ( auto& det : rpcGeo->dets() ) {
185  const RPCChamber* ch = dynamic_cast< const RPCChamber* >(det);
186  if ( !ch ) continue;
187 
188  //Loop on rolls in given chamber
189  for ( auto& r : ch->rolls() ) {
190  RPCDetId detId = r->id();
191  myDetIds_.push_back(detId);
192  }
193  }
194 
195 }
196 
197 
198 
200 
201  for(unsigned int i = 0; i<clientList_.size(); i++){
202 
203  if( clientList_[i] == "RPCMultiplicityTest" ) {
204  clientHisto_.push_back("Multiplicity");
205  // clientTag_.push_back(rpcdqm::MULTIPLICITY);
206  clientModules_.push_back( new RPCMultiplicityTest(parameters_));
207  } else if ( clientList_[i] == "RPCDeadChannelTest" ){
208  clientHisto_.push_back("Occupancy");
209  clientModules_.push_back( new RPCDeadChannelTest(parameters_));
210  // clientTag_.push_back(rpcdqm::OCCUPANCY);
211  } else if ( clientList_[i] == "RPCClusterSizeTest" ){
212  clientHisto_.push_back("ClusterSize");
213  clientModules_.push_back( new RPCClusterSizeTest(parameters_));
214  // clientTag_.push_back(rpcdqm::CLUSTERSIZE);
215  } else if ( clientList_[i] == "RPCOccupancyTest" ){
216  clientHisto_.push_back("Occupancy");
217  clientModules_.push_back( new RPCOccupancyTest(parameters_));
218  // clientTag_.push_back(rpcdqm::OCCUPANCY);
219  } else if ( clientList_[i] == "RPCNoisyStripTest" ){
220  clientHisto_.push_back("Occupancy");
221  clientModules_.push_back( new RPCNoisyStripTest(parameters_));
222  //clientTag_.push_back(rpcdqm::OCCUPANCY);
223  }
224  }
225 
226  return;
227 
228 }
void dqmEndLuminosityBlock(DQMStore::IBooker &, DQMStore::IGetter &, edm::LuminosityBlock const &, edm::EventSetup const &) override
Definition: RPCDqmClient.cc:71
virtual std::string chambername()
Definition: RPCGeomServ.cc:122
T getUntrackedParameter(std::string const &, T const &) const
MonitorElement * RPCEvents_
Definition: RPCDqmClient.h:43
void beginJob() override
Definition: RPCDqmClient.cc:57
std::string prefixDir_
Definition: RPCDqmClient.h:39
~RPCDqmClient() override
Destructor.
Definition: RPCDqmClient.cc:55
void getMonitorElements(DQMStore::IGetter &)
std::vector< std::string > clientList_
Definition: RPCDqmClient.h:41
virtual std::string name()
Definition: RPCGeomServ.cc:20
void makeClientMap(const edm::ParameterSet &parameters_)
std::vector< std::string > clientHisto_
Definition: RPCDqmClient.h:45
void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override
const std::vector< const RPCRoll * > & rolls() const
Return the Rolls.
Definition: RPCChamber.cc:68
bool useRollInfo_
Definition: RPCDqmClient.h:38
std::string globalFolder_
Definition: RPCDqmClient.h:40
MonitorElement * get(std::string const &path)
Definition: DQMStore.cc:303
void getRPCdetId(const edm::EventSetup &)
int minimumEvents_
Definition: RPCDqmClient.h:36
std::string folderStructure(RPCDetId detId)
const DetContainer & dets() const override
Returm a vector of all GeomDet (including all GeomDetUnits)
Definition: RPCGeometry.cc:33
std::vector< RPCDetId > myDetIds_
Definition: RPCDqmClient.h:44
double getBinContent(int binx) const
get content of bin (1-D)
bool enableDQMClients_
Definition: RPCDqmClient.h:38
T get() const
Definition: EventSetup.h:71
std::vector< RPCClient * > clientModules_
Definition: RPCDqmClient.h:46
RPCDqmClient(const edm::ParameterSet &ps)
Constructor.
Definition: RPCDqmClient.cc:22
int prescaleGlobalFactor_
Definition: RPCDqmClient.h:36
bool offlineDQM_
Definition: RPCDqmClient.h:35