CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 <stdint.h>
26 
27 # ifndef DQM_ROOT_METHODS
28 # define DQM_ROOT_METHODS 1
29 # endif
30 
31 class QCriterion;
32 
35 {
36  friend class DQMStore;
37  friend class DQMService;
38 public:
39  struct Scalar
40  {
41  int64_t num;
42  double real;
43  std::string str;
44  };
45 
46  enum Kind
47  {
61  };
62 
63  typedef std::vector<QReport>::const_iterator QReportIterator;
64 
65 private:
66  DQMNet::CoreObject data_; //< Core object information.
67  Scalar scalar_; //< Current scalar value.
68  TH1 *object_; //< Current ROOT object value.
69  TH1 *reference_; //< Current ROOT reference object.
70  TH1 *refvalue_; //< Soft reference if any.
71  std::vector<QReport> qreports_; //< QReports associated to this object.
72 
74  MonitorElement *initialise(Kind kind, TH1 *rootobj);
75  MonitorElement *initialise(Kind kind, const std::string &value);
76 
77 public:
78  MonitorElement(void);
79  MonitorElement(const std::string *path, const std::string &name);
82  ~MonitorElement(void);
83 
85  bool operator<(const MonitorElement &x) const
86  {
87  return DQMNet::setOrder(data_, x.data_);
88  }
89 
91  Kind kind(void) const
93 
95  uint32_t flags(void) const
96  { return data_.flags; }
97 
99  const std::string &getName(void) const
100  { return data_.objname; }
101 
103  const std::string &getPathname(void) const
104  { return *data_.dirname; }
105 
107  const std::string getFullname(void) const
108  {
109  std::string path;
110  path.reserve(data_.dirname->size() + data_.objname.size() + 2);
111  path += *data_.dirname;
112  if (! data_.dirname->empty())
113  path += '/';
114  path += data_.objname;
115  return path;
116  }
117 
119  bool wasUpdated(void) const
120  { return data_.flags & DQMNet::DQM_PROP_NEW; }
121 
123  void update(void)
125 
128  void setResetMe(bool flag)
130 
132  bool getLumiFlag(void) const
133  { return data_.flags & DQMNet::DQM_PROP_LUMI; }
134 
136  void setLumiFlag(void)
138 
139  // A static assert to check that T actually fits in
140  // int64_t.
141  template <typename T>
143  {
144  int checkArray[sizeof(int64_t) - sizeof(T) + 1];
145  };
146 
147  void Fill(long long x) { fits_in_int64_t<long long>(); doFill(static_cast<int64_t>(x)); }
148  void Fill(unsigned long long x) { fits_in_int64_t<unsigned long long>(); doFill(static_cast<int64_t>(x)); }
149  void Fill(unsigned long x) { fits_in_int64_t<unsigned long>(); doFill(static_cast<int64_t>(x)); }
150  void Fill(long x) { fits_in_int64_t<long>(); doFill(static_cast<int64_t>(x)); }
151  void Fill(unsigned int x) { fits_in_int64_t<unsigned int>(); doFill(static_cast<int64_t>(x)); }
152  void Fill(int x) { fits_in_int64_t<int>(); doFill(static_cast<int64_t>(x)); }
153  void Fill(short x) { fits_in_int64_t<short>(); doFill(static_cast<int64_t>(x)); }
154  void Fill(unsigned short x) { fits_in_int64_t<unsigned short>(); doFill(static_cast<int64_t>(x)); }
155  void Fill(char x) { fits_in_int64_t<char>(); doFill(static_cast<int64_t>(x)); }
156  void Fill(unsigned char x) { fits_in_int64_t<unsigned char>(); doFill(static_cast<int64_t>(x)); }
157 
158  void Fill(float x) { Fill(static_cast<double>(x)); }
159  void Fill(double x);
160  void Fill(std::string &value);
161 
162  void Fill(double x, double yw);
163  void Fill(double x, double y, double zw);
164  void Fill(double x, double y, double z, double w);
165  void ShiftFillLast(double y, double ye = 0., int32_t xscale = 1);
166  void Reset(void);
167 
168  std::string valueString(void) const;
169  std::string tagString(void) const;
170  std::string tagLabelString(void) const;
171  std::string qualityTagString(const DQMNet::QValue &qv) const;
172  void packScalarData(std::string &into, const char *prefix) const;
173  void packQualityData(std::string &into) const;
174 
176  bool hasError(void) const
178 
180  bool hasWarning(void) const
182 
184  bool hasOtherReport(void) const
186 
188  const QReport *getQReport(const std::string &qtname) const;
189 
191  std::vector<QReport *> getQReports(void) const;
192 
194  std::vector<QReport *> getQWarnings(void) const;
195 
197  std::vector<QReport *> getQErrors(void) const;
198 
201  std::vector<QReport *> getQOthers(void) const;
202 
204  void runQTests(void);
205 
206 private:
207  void doFill(int64_t x);
208  void incompatible(const char *func) const;
209  TH1 *accessRootObject(const char *func, int reqdim) const;
210 
211 public:
212 #if DQM_ROOT_METHODS
213  double getMean(int axis = 1) const;
214  double getMeanError(int axis = 1) const;
215  double getRMS(int axis = 1) const;
216  double getRMSError(int axis = 1) const;
217  int getNbinsX(void) const;
218  int getNbinsY(void) const;
219  int getNbinsZ(void) const;
220  double getBinContent(int binx) const;
221  double getBinContent(int binx, int biny) const;
222  double getBinContent(int binx, int biny, int binz) const;
223  double getBinError(int binx) const;
224  double getBinError(int binx, int biny) const;
225  double getBinError(int binx, int biny, int binz) const;
226  double getEntries(void) const;
227  double getBinEntries(int bin) const;
228 
229 private:
230  double getYmin(void) const;
231  double getYmax(void) const;
232 
233 public:
234  std::string getAxisTitle(int axis = 1) const;
235  std::string getTitle(void) const;
236  void setBinContent(int binx, double content);
237  void setBinContent(int binx, int biny, double content);
238  void setBinContent(int binx, int biny, int binz, double content);
239  void setBinError(int binx, double error);
240  void setBinError(int binx, int biny, double error);
241  void setBinError(int binx, int biny, int binz, double error);
242  void setBinEntries(int bin, double nentries);
243  void setEntries(double nentries);
244  void setBinLabel(int bin, const std::string &label, int axis = 1);
245  void setAxisRange(double xmin, double xmax, int axis = 1);
246  void setAxisTitle(const std::string &title, int axis = 1);
247  void setAxisTimeDisplay(int value, int axis = 1);
248  void setAxisTimeFormat(const char *format = "", int axis = 1);
249 
250 private:
251  void setAxisTimeOffset(double toffset, const char *option="local", int axis = 1);
252 
253 public:
254  void setTitle(const std::string &title);
255 #endif // DQM_ROOT_METHODS
256 
257 private:
259  bool isSoftResetEnabled(void) const
260  { return refvalue_ != 0; }
261 
263  bool isAccumulateEnabled(void) const
265 
266 private:
268  void resetUpdate(void)
269  { data_.flags &= ~DQMNet::DQM_PROP_NEW; }
270 
272  bool resetMe(void) const
273  { return data_.flags & DQMNet::DQM_PROP_RESET; }
274 
277  void setAccumulate(bool flag)
279 
280  TAxis *getAxis(const char *func, int axis) const;
281 
282  // ------------ Operations for MEs that are normally never reset ---------
283  void softReset(void);
284  void disableSoftReset(void);
285  void addProfiles(TProfile *h1, TProfile *h2, TProfile *sum, float c1, float c2);
286  void addProfiles(TProfile2D *h1, TProfile2D *h2, TProfile2D *sum, float c1, float c2);
287  void copyFunctions(TH1 *from, TH1 *to);
288  void copyFrom(TH1 *from);
289 
290 
291  // --- Operations on MEs that are normally reset at end of monitoring cycle ---
292  void getQReport(bool create, const std::string &qtname, QReport *&qr, DQMNet::QValue *&qv);
293  void addQReport(const DQMNet::QValue &desc, QCriterion *qc);
294  void addQReport(QCriterion *qc);
295  void updateQReportStats(void);
296 
297 public:
298  TObject *getRootObject(void) const;
299  TH1 *getTH1(void) const;
300  TH1F *getTH1F(void) const;
301  TH1S *getTH1S(void) const;
302  TH1D *getTH1D(void) const;
303  TH2F *getTH2F(void) const;
304  TH2S *getTH2S(void) const;
305  TH2D *getTH2D(void) const;
306  TH3F *getTH3F(void) const;
307  TProfile *getTProfile(void) const;
308  TProfile2D *getTProfile2D(void) const;
309 
310  TObject *getRefRootObject(void) const;
311  TH1 *getRefTH1(void) const;
312  TH1F *getRefTH1F(void) const;
313  TH1S *getRefTH1S(void) const;
314  TH1D *getRefTH1D(void) const;
315  TH2F *getRefTH2F(void) const;
316  TH2S *getRefTH2S(void) const;
317  TH2D *getRefTH2D(void) const;
318  TH3F *getRefTH3F(void) const;
319  TProfile *getRefTProfile(void) const;
320  TProfile2D *getRefTProfile2D(void) const;
321 
322  int64_t getIntValue(void) const
323  {
324  assert(kind() == DQM_KIND_INT);
325  return scalar_.num;
326  }
327 
328  double getFloatValue(void) const
329  {
330  assert(kind() == DQM_KIND_REAL);
331  return scalar_.real;
332  }
333 
334  const std::string &getStringValue(void) const
335  {
336  assert(kind() == DQM_KIND_STRING);
337  return scalar_.str;
338  }
339 
340  DQMNet::TagList getTags(void) const // DEPRECATED
341  {
344  tags.push_back(data_.tag);
345  return tags;
346  }
347 
348  const uint32_t getTag(void) const
349  { return data_.tag; }
350 };
351 
352 #endif // DQMSERVICES_CORE_MONITOR_ELEMENT_H
TH1F * getRefTH1F(void) const
TH2S * getTH2S(void) const
TH1S * getTH1S(void) const
static const uint32_t DQM_PROP_REPORT_WARN
Definition: DQMNet.h:47
const std::string & getName(void) const
get name of ME
void incompatible(const char *func) const
const QReport * getQReport(const std::string &qtname) const
get QReport corresponding to &lt;qtname&gt; (null pointer if QReport does not exist)
void resetUpdate(void)
reset &quot;was updated&quot; flag
void setBinContent(int binx, double content)
set content of bin (1-D)
static const uint32_t DQM_PROP_TYPE_TH1S
Definition: DQMNet.h:32
long int flag
Definition: mlp_lapack.h:47
TH1 * accessRootObject(const char *func, int reqdim) const
void copyFrom(TH1 *from)
MonitorElement * initialise(Kind kind)
static const uint32_t DQM_PROP_TYPE_TPROF
Definition: DQMNet.h:40
const std::string & label
Definition: MVAComputer.cc:186
MonitorElement & operator=(const MonitorElement &)
TProfile2D * getTProfile2D(void) const
TH1 * getRefTH1(void) const
void setAxisRange(double xmin, double xmax, int axis=1)
set x-, y- or z-axis range (axis=1, 2, 3 respectively)
static const uint32_t DQM_PROP_TYPE_TH2D
Definition: DQMNet.h:36
void setAxisTimeFormat(const char *format="", int axis=1)
set the format of the time values that are displayed on an axis
void updateQReportStats(void)
Refresh QReport stats, usually after MEs were read in from a file.
std::vector< QReport > qreports_
void Fill(unsigned short x)
std::string tagLabelString(void) const
return label string for the monitor element tag (eg. &lt;name&gt;t=12345&lt;/name&gt;)
Definition: DQMNet.h:22
std::string qualityTagString(const DQMNet::QValue &qv) const
std::string getAxisTitle(int axis=1) const
get x-, y- or z-axis title (axis=1, 2, 3 respectively)
TProfile2D * getRefTProfile2D(void) const
uint32_t flags
Definition: DQMNet.h:96
void disableSoftReset(void)
reverts action of softReset
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
bool hasError(void) const
true if at least of one of the quality tests returned an error
TH3F * getTH3F(void) const
static const uint32_t DQM_PROP_TAGGED
Definition: DQMNet.h:54
TH1D * getTH1D(void) const
void runQTests(void)
run all quality tests
void Fill(float x)
void softReset(void)
static const uint32_t DQM_PROP_TYPE_TH3F
Definition: DQMNet.h:37
static const uint32_t DQM_PROP_RESET
Definition: DQMNet.h:56
TH2D * getTH2D(void) const
uint32_t tag
Definition: DQMNet.h:97
const std::string * dirname
Definition: DQMNet.h:99
static const uint32_t DQM_PROP_TYPE_TH1F
Definition: DQMNet.h:31
void update(void)
Mark the object updated.
void Fill(short x)
double getEntries(void) const
get # of entries
double getMean(int axis=1) const
get mean value of histogram along x, y or z axis (axis=1, 2, 3 respectively)
bool operator<(const MonitorElement &x) const
Compare monitor elements, for ordering in sets.
double double double z
void Fill(unsigned long long x)
double getMeanError(int axis=1) const
void packScalarData(std::string &into, const char *prefix) const
convert scalar data into a string.
const uint32_t getTag(void) const
int getNbinsY(void) const
get # of bins in Y-axis
void Fill(long long x)
std::vector< QReport * > getQErrors(void) const
get errors from last set of quality tests
int path() const
Definition: HLTadd.h:3
bool wasUpdated(void) const
true if ME was updated in last monitoring cycle
static const uint32_t DQM_PROP_ACCUMULATE
Definition: DQMNet.h:55
TH3F * getRefTH3F(void) const
bool getLumiFlag(void) const
true if ME is meant to be stored for each luminosity section
void Fill(int x)
void ShiftFillLast(double y, double ye=0., int32_t xscale=1)
void packQualityData(std::string &into) const
serialise quality report information into a string.
const std::string & getPathname(void) const
get pathname of parent folder
TH2F * getRefTH2F(void) const
bool hasWarning(void) const
true if at least of one of the quality tests returned a warning
void setAxisTimeDisplay(int value, int axis=1)
set x-, y-, or z-axis to display time values
static const uint32_t DQM_PROP_TYPE_INT
Definition: DQMNet.h:28
void doFill(int64_t x)
&quot;Fill&quot; ME method for int64_t
DQMNet::TagList getTags(void) const
int getNbinsZ(void) const
get # of bins in Z-axis
void addProfiles(TProfile *h1, TProfile *h2, TProfile *sum, float c1, float c2)
std::string getTitle(void) const
get MonitorElement title
double getFloatValue(void) const
int checkArray[sizeof(int64_t)-sizeof(T)+1]
TH1 * getTH1(void) const
static const uint32_t DQM_PROP_REPORT_ERROR
Definition: DQMNet.h:46
void setBinError(int binx, double error)
set uncertainty on content of bin (1-D)
std::vector< QReport >::const_iterator QReportIterator
static const uint32_t DQM_PROP_REPORT_OTHER
Definition: DQMNet.h:48
Kind kind(void) const
Get the type of the monitor element.
static const uint32_t DQM_PROP_TYPE_TH1D
Definition: DQMNet.h:33
void setEntries(double nentries)
set # of entries
uint32_t flags(void) const
Get the object flags.
const std::string getFullname(void) const
get full name of ME including Pathname
const std::string & getStringValue(void) const
bool isAccumulateEnabled(void) const
whether ME contents should be accumulated over multiple monitoring periods; default: false ...
void Fill(unsigned int x)
void Fill(char x)
std::string objname
Definition: DQMNet.h:100
std::string valueString(void) const
TH2D * getRefTH2D(void) const
DQMNet::CoreObject data_
TAxis * getAxis(const char *func, int axis) const
void setTitle(const std::string &title)
set (ie. change) histogram/profile title
tuple tags
Definition: o2o.py:248
static std::string from(" from ")
void Fill(long x)
double getRMSError(int axis=1) const
get RMS uncertainty of histogram along x, y or z axis(axis=1,2,3 respectively)
std::vector< uint32_t > TagList
Definition: DQMNet.h:82
TObject * getRootObject(void) const
std::vector< QReport * > getQReports(void) const
get map of QReports
void setAxisTimeOffset(double toffset, const char *option="local", int axis=1)
set the time offset, if option = &quot;gmt&quot; then the offset is treated as a GMT time
void Fill(unsigned long x)
std::vector< QReport * > getQOthers(void) const
double getBinError(int binx) const
get uncertainty on content of bin (1-D) - See TH1::GetBinError for details
int64_t getIntValue(void) const
TH1F * getTH1F(void) const
std::string tagString(void) const
string const
Definition: compareJSON.py:14
bool hasOtherReport(void) const
true if at least of one of the tests returned some other (non-ok) status
TProfile * getRefTProfile(void) const
void copyFunctions(TH1 *from, TH1 *to)
bool resetMe(void) const
true if ME should be reset at end of monitoring cycle
double getBinContent(int binx) const
get content of bin (1-D)
double getRMS(int axis=1) const
get RMS of histogram along x, y or z axis (axis=1, 2, 3 respectively)
double getYmax(void) const
get max Y value (for profiles)
static bool setOrder(const CoreObject &a, const CoreObject &b)
Definition: DQMNet.h:173
TProfile * getTProfile(void) const
TH2S * getRefTH2S(void) const
bool isSoftResetEnabled(void) const
whether soft-reset is enabled; default is false
double getBinEntries(int bin) const
get # of bin entries (for profiles)
std::vector< QReport * > getQWarnings(void) const
get warnings from last set of quality tests
static const uint32_t DQM_PROP_TYPE_TH2S
Definition: DQMNet.h:35
int getNbinsX(void) const
get # of bins in X-axis
void Fill(unsigned char x)
void setResetMe(bool flag)
void setAccumulate(bool flag)
static const uint32_t DQM_PROP_TYPE_STRING
Definition: DQMNet.h:30
static const uint32_t DQM_PROP_NEW
Definition: DQMNet.h:58
Definition: DDAxes.h:10
void setLumiFlag(void)
this ME is meant to be stored for each luminosity section
static const uint32_t DQM_PROP_TYPE_MASK
Definition: DQMNet.h:25
TH2F * getTH2F(void) const
long double T
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
void Reset(void)
reset ME (ie. contents, errors, etc)
TObject * getRefRootObject(void) const
static const uint32_t DQM_PROP_LUMI
Definition: DQMNet.h:60
void setBinEntries(int bin, double nentries)
set # of bin entries (to be used for profiles)
TH1S * getRefTH1S(void) const
static const uint32_t DQM_PROP_TYPE_REAL
Definition: DQMNet.h:29
SurfaceDeformation * create(int type, const std::vector< double > &params)
void addQReport(const DQMNet::QValue &desc, QCriterion *qc)
Add quality report, from DQMStore.
static const uint32_t DQM_PROP_TYPE_INVALID
Definition: DQMNet.h:27
static const uint32_t DQM_PROP_TYPE_TPROF2D
Definition: DQMNet.h:41
double getYmin(void) const
get min Y value (for profiles)
static const uint32_t DQM_PROP_TYPE_TH2F
Definition: DQMNet.h:34
TH1D * getRefTH1D(void) const