CMS 3D CMS Logo

MonitorElement.h
Go to the documentation of this file.
1 #ifndef DQMSERVICES_CORE_MONITOR_ELEMENT_H
2 # define DQMSERVICES_CORE_MONITOR_ELEMENT_H
3 
6 # include "TF1.h"
7 # include "TH1F.h"
8 # include "TH1S.h"
9 # include "TH1D.h"
10 # include "TH2F.h"
11 # include "TH2S.h"
12 # include "TH2D.h"
13 # include "TH3F.h"
14 # include "TProfile.h"
15 # include "TProfile2D.h"
16 # include "TObjString.h"
17 # include "TAxis.h"
18 # include <sys/time.h>
19 # include <string>
20 # include <set>
21 # include <map>
22 # include <sstream>
23 # include <iomanip>
24 # include <cassert>
25 # include <cstdint>
26 
27 # ifndef DQM_ROOT_METHODS
28 # define DQM_ROOT_METHODS 1
29 # endif
30 
31 class QCriterion;
32 
33 // tag for a special constructor, see below
35 
38 {
39  friend class DQMStore;
40  friend class DQMService;
41 public:
42  struct Scalar
43  {
44  int64_t num;
45  double real;
47  };
48 
49  enum Kind
50  {
51  DQM_KIND_INVALID = DQMNet::DQM_PROP_TYPE_INVALID,
52  DQM_KIND_INT = DQMNet::DQM_PROP_TYPE_INT,
53  DQM_KIND_REAL = DQMNet::DQM_PROP_TYPE_REAL,
54  DQM_KIND_STRING = DQMNet::DQM_PROP_TYPE_STRING,
55  DQM_KIND_TH1F = DQMNet::DQM_PROP_TYPE_TH1F,
56  DQM_KIND_TH1S = DQMNet::DQM_PROP_TYPE_TH1S,
57  DQM_KIND_TH1D = DQMNet::DQM_PROP_TYPE_TH1D,
58  DQM_KIND_TH2F = DQMNet::DQM_PROP_TYPE_TH2F,
59  DQM_KIND_TH2S = DQMNet::DQM_PROP_TYPE_TH2S,
60  DQM_KIND_TH2D = DQMNet::DQM_PROP_TYPE_TH2D,
61  DQM_KIND_TH3F = DQMNet::DQM_PROP_TYPE_TH3F,
62  DQM_KIND_TPROFILE = DQMNet::DQM_PROP_TYPE_TPROF,
63  DQM_KIND_TPROFILE2D = DQMNet::DQM_PROP_TYPE_TPROF2D
64  };
65 
66 private:
67  DQMNet::CoreObject data_; //< Core object information.
68  Scalar scalar_; //< Current scalar value.
69  TH1 *object_; //< Current ROOT object value.
70  TH1 *reference_; //< Current ROOT reference object.
71  TH1 *refvalue_; //< Soft reference if any.
72  std::vector<QReport> qreports_; //< QReports associated to this object.
73 
74  MonitorElement *initialise(Kind kind);
75  MonitorElement *initialise(Kind kind, TH1 *rootobj);
76  MonitorElement *initialise(Kind kind, const std::string &value);
77  void globalize() {
78  data_.moduleId = 0;
79  }
80  void setLumi(uint32_t ls) {
81  data_.lumi = ls;
82  }
83 
84 public:
87  MonitorElement(const std::string *path, const std::string &name, uint32_t run, uint32_t moduleId);
91  MonitorElement &operator=(const MonitorElement &) = delete;
92  MonitorElement &operator=(MonitorElement &&) = delete;
93  ~MonitorElement();
94 
96  bool operator<(const MonitorElement &x) const
97  {
98  return DQMNet::setOrder(data_, x.data_);
99  }
100 
102  static bool CheckBinLabels(const TAxis* a1, const TAxis * a2);
103 
105  Kind kind() const
106  { return Kind(data_.flags & DQMNet::DQM_PROP_TYPE_MASK); }
107 
109  uint32_t flags() const
110  { return data_.flags; }
111 
113  const std::string &getName() const
114  { return data_.objname; }
115 
117  const std::string &getPathname() const
118  { return *data_.dirname; }
119 
121  const std::string getFullname() const
122  {
124  path.reserve(data_.dirname->size() + data_.objname.size() + 2);
125  path += *data_.dirname;
126  if (! data_.dirname->empty())
127  path += '/';
128  path += data_.objname;
129  return path;
130  }
131 
133  bool wasUpdated() const
134  { return data_.flags & DQMNet::DQM_PROP_NEW; }
135 
137  void update()
138  { data_.flags |= DQMNet::DQM_PROP_NEW; }
139 
142  void setResetMe(bool /* flag */)
143  { data_.flags |= DQMNet::DQM_PROP_RESET; }
144 
146  bool getLumiFlag() const
147  { return data_.flags & DQMNet::DQM_PROP_LUMI; }
148 
150  void setLumiFlag()
151  { data_.flags |= DQMNet::DQM_PROP_LUMI; }
152 
157 
158  // A static assert to check that T actually fits in
159  // int64_t.
160  template <typename T>
162  {
163  int checkArray[sizeof(int64_t) - sizeof(T) + 1];
164  };
165 
166  void Fill(long long x) { fits_in_int64_t<long long>(); doFill(static_cast<int64_t>(x)); }
167  void Fill(unsigned long long x) { fits_in_int64_t<unsigned long long>(); doFill(static_cast<int64_t>(x)); }
168  void Fill(unsigned long x) { fits_in_int64_t<unsigned long>(); doFill(static_cast<int64_t>(x)); }
169  void Fill(long x) { fits_in_int64_t<long>(); doFill(static_cast<int64_t>(x)); }
170  void Fill(unsigned int x) { fits_in_int64_t<unsigned int>(); doFill(static_cast<int64_t>(x)); }
171  void Fill(int x) { fits_in_int64_t<int>(); doFill(static_cast<int64_t>(x)); }
172  void Fill(short x) { fits_in_int64_t<short>(); doFill(static_cast<int64_t>(x)); }
173  void Fill(unsigned short x) { fits_in_int64_t<unsigned short>(); doFill(static_cast<int64_t>(x)); }
174  void Fill(char x) { fits_in_int64_t<char>(); doFill(static_cast<int64_t>(x)); }
175  void Fill(unsigned char x) { fits_in_int64_t<unsigned char>(); doFill(static_cast<int64_t>(x)); }
176 
177  void Fill(float x) { Fill(static_cast<double>(x)); }
178  void Fill(double x);
179  void Fill(std::string &value);
180 
181  void Fill(double x, double yw);
182  void Fill(double x, double y, double zw);
183  void Fill(double x, double y, double z, double w);
184  void ShiftFillLast(double y, double ye = 0., int32_t xscale = 1);
185  void Reset();
186 
187  std::string valueString() const;
188  std::string tagString() const;
189  std::string tagLabelString() const;
190  std::string effLabelString() const;
191  std::string qualityTagString(const DQMNet::QValue &qv) const;
192  void packScalarData(std::string &into, const char *prefix) const;
193  void packQualityData(std::string &into) const;
194 
196  bool hasError() const
197  { return data_.flags & DQMNet::DQM_PROP_REPORT_ERROR; }
198 
200  bool hasWarning() const
201  { return data_.flags & DQMNet::DQM_PROP_REPORT_WARN; }
202 
204  bool hasOtherReport() const
205  { return data_.flags & DQMNet::DQM_PROP_REPORT_OTHER; }
206 
209  bool isEfficiency() const
210  { return data_.flags & DQMNet::DQM_PROP_EFFICIENCY_PLOT; }
211 
213  const QReport *getQReport(const std::string &qtname) const;
214 
216  std::vector<QReport *> getQReports() const;
217 
219  std::vector<QReport *> getQWarnings() const;
220 
222  std::vector<QReport *> getQErrors() const;
223 
226  std::vector<QReport *> getQOthers() const;
227 
229  void runQTests();
230 
231 private:
232  void doFill(int64_t x);
233  void incompatible(const char *func) const;
234  TH1 *accessRootObject(const char *func, int reqdim) const;
235 
236 public:
237 #if DQM_ROOT_METHODS
238  double getMean(int axis = 1) const;
239  double getMeanError(int axis = 1) const;
240  double getRMS(int axis = 1) const;
241  double getRMSError(int axis = 1) const;
242  int getNbinsX() const;
243  int getNbinsY() const;
244  int getNbinsZ() const;
245  double getBinContent(int binx) const;
246  double getBinContent(int binx, int biny) const;
247  double getBinContent(int binx, int biny, int binz) const;
248  double getBinError(int binx) const;
249  double getBinError(int binx, int biny) const;
250  double getBinError(int binx, int biny, int binz) const;
251  double getEntries() const;
252  double getBinEntries(int bin) const;
253 
254 private:
255  double getYmin() const;
256  double getYmax() const;
257 
258 public:
259  std::string getAxisTitle(int axis = 1) const;
260  std::string getTitle() const;
261  void setBinContent(int binx, double content);
262  void setBinContent(int binx, int biny, double content);
263  void setBinContent(int binx, int biny, int binz, double content);
264  void setBinError(int binx, double error);
265  void setBinError(int binx, int biny, double error);
266  void setBinError(int binx, int biny, int binz, double error);
267  void setBinEntries(int bin, double nentries);
268  void setEntries(double nentries);
269  void setBinLabel(int bin, const std::string &label, int axis = 1);
270  void setAxisRange(double xmin, double xmax, int axis = 1);
271  void setAxisTitle(const std::string &title, int axis = 1);
272  void setAxisTimeDisplay(int value, int axis = 1);
273  void setAxisTimeFormat(const char *format = "", int axis = 1);
274 
275 private:
276  void setAxisTimeOffset(double toffset, const char *option="local", int axis = 1);
277 
278 public:
279  void setTitle(const std::string &title);
280 #endif // DQM_ROOT_METHODS
281 
282 private:
284  bool isSoftResetEnabled() const
285  { return refvalue_ != nullptr; }
286 
288  bool isAccumulateEnabled() const
289  { return data_.flags & DQMNet::DQM_PROP_ACCUMULATE; }
290 
292  bool markedToDelete() const
293  { return data_.flags & DQMNet::DQM_PROP_MARKTODELETE; }
294 
299 
300 private:
302  void resetUpdate()
303  { data_.flags &= ~DQMNet::DQM_PROP_NEW; }
304 
306  bool resetMe() const
307  { return data_.flags & DQMNet::DQM_PROP_RESET; }
308 
311  void setAccumulate(bool /* flag */)
312  { data_.flags |= DQMNet::DQM_PROP_ACCUMULATE; }
313 
314  TAxis *getAxis(const char *func, int axis) const;
315 
316  // ------------ Operations for MEs that are normally never reset ---------
317 public:
318  void softReset();
319 private:
320  void disableSoftReset();
321  void addProfiles(TProfile *h1, TProfile *h2, TProfile *sum, float c1, float c2);
322  void addProfiles(TProfile2D *h1, TProfile2D *h2, TProfile2D *sum, float c1, float c2);
323  void copyFunctions(TH1 *from, TH1 *to);
324  void copyFrom(TH1 *from);
325 
326 
327  // --- Operations on MEs that are normally reset at end of monitoring cycle ---
328  void getQReport(bool create, const std::string &qtname, QReport *&qr, DQMNet::QValue *&qv);
329  void addQReport(const DQMNet::QValue &desc, QCriterion *qc);
330  void addQReport(QCriterion *qc);
331  void updateQReportStats();
332 
333 public:
334  TObject *getRootObject() const;
335  TH1 *getTH1() const;
336  TH1F *getTH1F() const;
337  TH1S *getTH1S() const;
338  TH1D *getTH1D() const;
339  TH2F *getTH2F() const;
340  TH2S *getTH2S() const;
341  TH2D *getTH2D() const;
342  TH3F *getTH3F() const;
343  TProfile *getTProfile() const;
344  TProfile2D *getTProfile2D() const;
345 
346  TObject *getRefRootObject() const;
347  TH1 *getRefTH1() const;
348  TH1F *getRefTH1F() const;
349  TH1S *getRefTH1S() const;
350  TH1D *getRefTH1D() const;
351  TH2F *getRefTH2F() const;
352  TH2S *getRefTH2S() const;
353  TH2D *getRefTH2D() const;
354  TH3F *getRefTH3F() const;
355  TProfile *getRefTProfile() const;
356  TProfile2D *getRefTProfile2D() const;
357 
358  int64_t getIntValue() const
359  {
360  assert(kind() == DQM_KIND_INT);
361  return scalar_.num;
362  }
363 
364  double getFloatValue() const
365  {
366  assert(kind() == DQM_KIND_REAL);
367  return scalar_.real;
368  }
369 
371  {
372  assert(kind() == DQM_KIND_STRING);
373  return scalar_.str;
374  }
375 
376  DQMNet::TagList getTags() const // DEPRECATED
377  {
378  DQMNet::TagList tags;
379  if (data_.flags & DQMNet::DQM_PROP_TAGGED)
380  tags.push_back(data_.tag);
381  return tags;
382  }
383 
384  const uint32_t getTag() const
385  { return data_.tag; }
386 
387  const uint32_t run() const {return data_.run;}
388  const uint32_t lumi() const {return data_.lumi;}
389  const uint32_t moduleId() const {return data_.moduleId;}
390 };
391 
392 #endif // DQMSERVICES_CORE_MONITOR_ELEMENT_H
static const uint32_t DQM_PROP_REPORT_WARN
Definition: DQMNet.h:48
uint32_t moduleId
Definition: DQMNet.h:105
int64_t getIntValue() const
static const uint32_t DQM_PROP_TYPE_TH1S
Definition: DQMNet.h:33
static const uint32_t DQM_PROP_TYPE_TPROF
Definition: DQMNet.h:41
const double w
Definition: UKUtility.cc:23
def create(alignables, pedeDump, additionalData, outputFile, config)
static const uint32_t DQM_PROP_TYPE_TH2D
Definition: DQMNet.h:37
std::vector< QReport > qreports_
void Fill(unsigned short x)
const std::string & getPathname() const
get pathname of parent folder
double getFloatValue() const
Definition: DQMNet.h:23
uint32_t flags() const
Get the object flags.
void setLumi(uint32_t ls)
bool isSoftResetEnabled() const
whether soft-reset is enabled; default is false
uint32_t flags
Definition: DQMNet.h:99
bool hasError() const
true if at least of one of the quality tests returned an error
static const uint32_t DQM_PROP_TAGGED
Definition: DQMNet.h:55
const std::string & getName() const
get name of ME
void Fill(float x)
static const uint32_t DQM_PROP_EFFICIENCY_PLOT
Definition: DQMNet.h:64
bool hasOtherReport() const
true if at least of one of the tests returned some other (non-ok) status
static const uint32_t DQM_PROP_TYPE_TH3F
Definition: DQMNet.h:38
static const uint32_t DQM_PROP_RESET
Definition: DQMNet.h:57
uint32_t tag
Definition: DQMNet.h:100
const std::string * dirname
Definition: DQMNet.h:106
static const uint32_t DQM_PROP_TYPE_TH1F
Definition: DQMNet.h:32
void Fill(short x)
static const uint32_t DQM_PROP_MARKTODELETE
Definition: DQMNet.h:65
bool operator<(const MonitorElement &x) const
Compare monitor elements, for ordering in sets.
uint32_t run
Definition: DQMNet.h:102
void Fill(unsigned long long x)
bool wasUpdated() const
true if ME was updated in last monitoring cycle
void Fill(long long x)
static const uint32_t DQM_PROP_ACCUMULATE
Definition: DQMNet.h:56
void update()
Mark the object updated.
void Fill(int x)
const uint32_t getTag() const
void setEfficiencyFlag()
void setLumiFlag()
this ME is meant to be stored for each luminosity section
static const uint32_t DQM_PROP_TYPE_INT
Definition: DQMNet.h:29
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
bool hasWarning() const
true if at least of one of the quality tests returned a warning
void setAccumulate(bool)
uint32_t lumi
Definition: DQMNet.h:103
bool resetMe() const
true if ME should be reset at end of monitoring cycle
void resetUpdate()
reset "was updated" flag
bool markedToDelete() const
true if ME is marked for deletion
static const uint32_t DQM_PROP_REPORT_ERROR
Definition: DQMNet.h:47
Definition: value.py:1
static const uint32_t DQM_PROP_REPORT_OTHER
Definition: DQMNet.h:49
static const uint32_t DQM_PROP_TYPE_TH1D
Definition: DQMNet.h:34
void Fill(unsigned int x)
void Fill(char x)
std::string objname
Definition: DQMNet.h:107
const std::string & getStringValue() const
DQMNet::CoreObject data_
bin
set the eta bin as selection string.
const std::string getFullname() const
get full name of ME including Pathname
doFill
Definition: cuy.py:575
const uint32_t lumi() const
void Fill(long x)
def ls(path, rec=False)
Definition: eostools.py:348
void Fill(unsigned long x)
void setResetMe(bool)
bool isEfficiency() const
static bool setOrder(const CoreObject &a, const CoreObject &b)
Definition: DQMNet.h:180
const uint32_t moduleId() const
DQMNet::TagList getTags() const
static const uint32_t DQM_PROP_TYPE_TH2S
Definition: DQMNet.h:36
void Fill(unsigned char x)
std::vector< uint32_t > TagList
Definition: DQMNet.h:85
auto zw(V v) -> Vec2< typename std::remove_reference< decltype(v[0])>::type >
Definition: ExtVec.h:75
bool getLumiFlag() const
true if ME is meant to be stored for each luminosity section
static const uint32_t DQM_PROP_TYPE_STRING
Definition: DQMNet.h:31
void Reset(std::vector< TH2F > &depth)
static const uint32_t DQM_PROP_NEW
Definition: DQMNet.h:59
static const uint32_t DQM_PROP_TYPE_MASK
Definition: DQMNet.h:26
long double T
static const uint32_t DQM_PROP_LUMI
Definition: DQMNet.h:61
static const uint32_t DQM_PROP_TYPE_REAL
Definition: DQMNet.h:30
Kind kind() const
Get the type of the monitor element.
static const uint32_t DQM_PROP_TYPE_INVALID
Definition: DQMNet.h:28
static const uint32_t DQM_PROP_TYPE_TPROF2D
Definition: DQMNet.h:42
static const uint32_t DQM_PROP_TYPE_TH2F
Definition: DQMNet.h:35
bool isAccumulateEnabled() const
whether ME contents should be accumulated over multiple monitoring periods; default: false ...
const uint32_t run() const