CMS 3D CMS Logo

DQMSourceExample.cc
Go to the documentation of this file.
1 /*
2  * \file DQMSourceExample.cc
3  * \author C.Leonidopoulos
4  * Last Update:
5  *
6  * Description: Simple example showing how to create a DQM source creating and
7  * filling monitoring elements
8  */
9 
12 
13 #include "TRandom.h"
14 #include <cmath>
15 
16 using namespace std;
17 using namespace edm;
18 
19 //==================================================================//
20 //================= Constructor and Destructor =====================//
21 //==================================================================//
23  parameters_ = ps;
24  initialize();
25 }
26 
28 
29 //==================================================================//
30 //======================= Initialise ===============================//
31 //==================================================================//
34  counterEvt_ = 0;
35  counterLS_ = 0;
36 
39 
41  monitorName_ = parameters_.getUntrackedParameter<string>("monitorName", "YourSubsystemName");
42  cout << "DQMSourceExample: Monitor name = " << monitorName_ << endl;
43  if (!monitorName_.empty())
44  monitorName_ = monitorName_ + "/";
45 
47  prescaleLS_ = parameters_.getUntrackedParameter<int>("prescaleLS", -1);
48  cout << "DQMSourceExample: DQM lumi section prescale = " << prescaleLS_ << " lumi section(s)" << endl;
49  prescaleEvt_ = parameters_.getUntrackedParameter<int>("prescaleEvt", -1);
50  cout << "DQMSourceExample: DQM event prescale = " << prescaleEvt_ << " events(s)" << endl;
51 
52  // read in files (use DQMStore.collateHistograms = True for summing
53  // dbe_->load("ref.root");
54  // dbe_->load("ref.root");
55 }
56 
57 //==================================================================//
58 //========================= beginJob ===============================//
59 //==================================================================//
63 
65  dbe_->setCurrentFolder(monitorName_ + "DQMsource/Summary");
66  summ = dbe_->book1D("summary", "Run Summary", 100, 0, 100);
67 
68  //-------------------------------------
69  // testing of Quality Tests
70  //-------------------------------------
71 
73  dbe_->setCurrentFolder(monitorName_ + "DQMsource/QTests");
74 
76  NBINS = 40;
77  XMIN = 0.;
78  XMAX = 40.;
79 
83  xTrue = dbe_->book1D("XTrue", "X Range QTest", NBINS, XMIN, XMAX);
84  xFalse = dbe_->book1D("XFalse", "X Range QTest", NBINS, XMIN, XMAX);
85  yTrue = dbe_->book1D("YTrue", "Y Range QTest", NBINS, XMIN, XMAX);
86  yFalse = dbe_->book1D("YFalse", "Y Range QTest", NBINS, XMIN, XMAX);
87  wExpTrue = dbe_->book2D("WExpTrue", "Contents Within Expected QTest", NBINS, XMIN, XMAX, NBINS, XMIN, XMAX);
88  wExpFalse = dbe_->book2D("WExpFalse", "Contents Within Expected QTest", NBINS, XMIN, XMAX, NBINS, XMIN, XMAX);
89  meanTrue = dbe_->book1D("MeanTrue", "Mean Within Expected QTest", NBINS, XMIN, XMAX);
90  meanFalse = dbe_->book1D("MeanFalse", "Mean Within Expected QTest", NBINS, XMIN, XMAX);
91  deadTrue = dbe_->book1D("DeadTrue", "Dead Channel QTest", NBINS, XMIN, XMAX);
92  deadFalse = dbe_->book1D("DeadFalse", "Dead Channel QTest", NBINS, XMIN, XMAX);
93  noisyTrue = dbe_->book1D("NoisyTrue", "Noisy Channel QTest", NBINS, XMIN, XMAX);
94  noisyFalse = dbe_->book1D("NoisyFalse", "Noisy Channel QTest", NBINS, XMIN, XMAX);
95 
96  //-------------------------------------
97  // book several ME more
98  //-------------------------------------
99 
101  dbe_->setCurrentFolder(monitorName_ + "DQMsource/C1");
102  const int NBINS2 = 10;
103 
104  i1 = dbe_->bookInt("int1");
105  f1 = dbe_->bookFloat("float1");
106  s1 = dbe_->bookString("s1", "My string");
107  h1 = dbe_->book1D("h1f", "Example TH1F 1D histogram.", NBINS2, XMIN, XMAX);
108  h2 = dbe_->book1S("h1s", "Example TH1S histogram.", NBINS, XMIN, XMAX);
109  // h3 = dbe_->book1DD("h1d", "Example TH1D histogram.", NBINS, XMIN,
110  // XMAX); h4 = dbe_->book2DD("h2d", "Example TH2D histogram.", NBINS,
111  // XMIN, XMAX,NBINS, XMIN, XMAX);
112  p1 = dbe_->bookProfile("prof1", "My profile 1D", NBINS, XMIN, XMAX, NBINS, XMIN, XMAX, "");
113  p2 = dbe_->bookProfile2D("prof2", "My profile 2D", NBINS, XMIN, XMAX, NBINS, XMIN, XMAX, NBINS, XMIN, XMAX, "");
114  h1hist = dbe_->book1D("history 1D", "Example 1 1D history plot", 30, 0., 30.);
115 
116  // set labels for h1
117  char temp[1024];
118  for (int i = 1; i <= NBINS2; ++i) {
119  sprintf(temp, " bin no. %d", i);
120  h1->setBinLabel(i, temp);
121  }
122 
123  // assign tag to MEs h1
124  const unsigned int detector_id = 17;
125  dbe_->tag(h1, detector_id);
126 
127  // tag full directory
128  dbe_->tagContents(monitorName_ + "DQMsource/C1", detector_id);
129 
130  /*
131  // contents of h5 & h6 will be reset at end of monitoring cycle
132  h5->setResetMe(true);
133  h6->setResetMe(true);
134  dbe_->showDirStructure();
135  std::vector<std::string> tags;
136  for (size_t i = 0, e = tags.size(); i < e; ++i)
137  std::cout << "TAGS [" << i << "] = " << tags[i] << std::endl;
138  */
139 
140  dbe_->showDirStructure();
141 }
142 
143 //==================================================================//
144 //========================= beginRun ===============================//
145 //==================================================================//
146 void DQMSourceExample::beginRun(const edm::Run &r, const EventSetup &context) {}
147 
148 //==================================================================//
149 //==================== analyse (takes each event) ==================//
150 //==================================================================//
151 void DQMSourceExample::analyze(const Event &iEvent, const EventSetup &iSetup) {
152  counterEvt_++;
153  if (prescaleEvt_ < 1)
154  return;
155  if (prescaleEvt_ > 0 && counterEvt_ % prescaleEvt_ != 0)
156  return;
157  // cout << " processing conterEvt_: " << counterEvt_ <<endl;
158 
159  // fill integer and float
160  // number exceeding 32 bits
161  i1->Fill(400000000000000LL); // FIXME use double
162  f1->Fill(-3.14);
163 
164  //----------------------------------------
165  // Filling the histograms with random data
166  //----------------------------------------
167 
168  srand(0);
169  // fill summ histo
170  if (counterEvt_ % 1000 == 0) {
171  cout << " # of events = " << counterEvt_ << endl;
172  summ->Fill(counterEvt_ / 1000., counterEvt_);
173  }
174  // fill summ histo
175  if (counterEvt_ % 100 == 0) {
176  h1hist->ShiftFillLast(gRandom->Gaus(12, 1.), 1., 5);
177  }
178 
179  float z = gRandom->Uniform(XMAX);
180  xTrue->Fill(z, 1. / log(z + 1.));
181  xFalse->Fill(z + (XMAX / 2.), z);
182  yTrue->Fill(z, 1. / log(z + 1.));
183  yFalse->Fill(z, z);
184  meanTrue->Fill(gRandom->Gaus(10, 2), 1.);
185  meanFalse->Fill(gRandom->Gaus(12, 3), 1.);
186  wExpTrue->Fill(gRandom->Gaus(12, 1), gRandom->Gaus(12, 1), 1.);
187  wExpFalse->Fill(gRandom->Gaus(20, 2), gRandom->Gaus(20, 2), 1.);
188  deadTrue->Fill(gRandom->Gaus(20, 10), 2.);
189  deadFalse->Fill(gRandom->Gaus(20, 4), 1.);
190  h2->Fill(gRandom->Gaus(20, 4), 1.);
191  // h3->Fill( XMIN, 0xffff00000000LL);
192  // h4->Fill( XMIN, XMIN, 0xffff00000000LL);
193 
194  // h1hist->Print();
195  // h1hist->Print();
196 
197  for (int i = 0; i != 10; ++i) {
198  float w = gRandom->Uniform(XMAX);
199  noisyTrue->Fill(w, 1.);
200  noisyFalse->Fill(z, 1.);
201  float x = gRandom->Gaus(12, 1);
202  float y = gRandom->Gaus(20, 2);
203  p1->Fill(x, y);
204  p2->Fill(x, y, (x + y) / 2.);
205  h1->Fill(y, 1.);
206  }
207 
208  // usleep(100);
209 }
210 
211 //==================================================================//
212 //============================= endRun =============================//
213 //==================================================================//
214 void DQMSourceExample::endRun(const Run &r, const EventSetup &context) {}
215 
216 //==================================================================//
217 //============================= endJob =============================//
218 //==================================================================//
219 void DQMSourceExample::endJob() { std::cout << "DQMSourceExample::endJob()" << std::endl; }
static AlgebraicMatrix initialize()
const double w
Definition: UKUtility.cc:23
void endJob() override
~DQMSourceExample() override
DQMSourceExample(const edm::ParameterSet &)
void beginRun(const edm::Run &r, const edm::EventSetup &c) override
int iEvent
Definition: GenABIO.cc:224
double p2[4]
Definition: TauolaWrapper.h:90
void beginJob() override
DQMStore * dbe_
HLT enums.
double p1[4]
Definition: TauolaWrapper.h:89
void endRun(const edm::Run &r, const edm::EventSetup &c) override
void analyze(const edm::Event &e, const edm::EventSetup &c) override
const int NBINS
Definition: Run.h:45