Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:10:11

0001 #define __STDC_FORMAT_MACROS 1
0002 #include "DQMServices/Core/interface/MonitorElement.h"
0003 #include "TClass.h"
0004 #include "TMath.h"
0005 #include "TList.h"
0006 #include "THashList.h"
0007 #include <iostream>
0008 #include <cassert>
0009 #include <cfloat>
0010 #include <cinttypes>
0011 
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 
0014 namespace dqm::impl {
0015 
0016   static TH1 *checkRootObject(const std::string &name, TObject *tobj, const char *func, int reqdim) {
0017     if (!tobj)
0018       throw cms::Exception("MonitorElementError") << "Method '" << func
0019                                                   << "' cannot be invoked on monitor"
0020                                                      " element '"
0021                                                   << name << "' because it is not a ROOT object.";
0022 
0023     auto *h = static_cast<TH1 *>(tobj);
0024     int ndim = h->GetDimension();
0025     if (reqdim < 0 || reqdim > ndim)
0026       throw cms::Exception("MonitorElementError") << "Method '" << func
0027                                                   << "' cannot be invoked on monitor"
0028                                                      " element '"
0029                                                   << name << "' because it requires " << reqdim
0030                                                   << " dimensions; this"
0031                                                      " object of type '"
0032                                                   << typeid(*h).name() << "' has " << ndim << " dimensions";
0033 
0034     return h;
0035   }
0036 
0037   MonitorElement::MonitorElement(MonitorElementData &&data) {
0038     this->mutable_ = std::make_shared<MutableMonitorElementData>();
0039     this->mutable_->data_ = std::move(data);
0040     syncCoreObject();
0041   }
0042   MonitorElement::MonitorElement(std::shared_ptr<MutableMonitorElementData> data) { switchData(std::move(data)); }
0043   MonitorElement::MonitorElement(MonitorElement *me) { switchData(me); }
0044 
0045   MonitorElementData MonitorElement::cloneMEData() {
0046     MonitorElementData out;
0047     auto access = this->access();
0048     out.key_ = access.key;
0049     out.value_.scalar_ = access.value.scalar_;
0050     if (access.value.object_) {
0051       out.value_.object_ = std::unique_ptr<TH1>(static_cast<TH1 *>(access.value.object_->Clone()));
0052     }
0053     return out;
0054   }
0055 
0056   std::shared_ptr<MutableMonitorElementData> MonitorElement::release() {
0057     auto data = this->mutable_;
0058     this->mutable_.reset();
0059     return data;
0060   }
0061 
0062   void MonitorElement::switchData(MonitorElement *other) {
0063     assert(other);
0064     this->mutable_ = other->mutable_;
0065     syncCoreObject();
0066   }
0067 
0068   void MonitorElement::switchData(std::shared_ptr<MutableMonitorElementData> data) {
0069     this->mutable_ = std::move(data);
0070     syncCoreObject();
0071   }
0072 
0073   void MonitorElement::switchObject(std::unique_ptr<TH1> &&newobject) {
0074     auto access = this->accessMut();
0075     // Assume kind etc. matches.
0076     // This should free the old object.
0077     access.value.object_ = std::move(newobject);
0078   }
0079 
0080   void MonitorElement::syncCoreObject() {
0081     auto access = this->accessMut();
0082     syncCoreObject(access);
0083   }
0084 
0085   void MonitorElement::syncCoreObject(AccessMut &access) {
0086     data_.flags &= ~DQMNet::DQM_PROP_TYPE_MASK;
0087     data_.flags |= (int)access.key.kind_;
0088 
0089     // mark as updated.
0090     data_.flags |= DQMNet::DQM_PROP_NEW;
0091 
0092     // lumi flag is approximately equivalent to Scope::LUMI.
0093     data_.flags &= ~DQMNet::DQM_PROP_LUMI;
0094     if (access.key.scope_ == MonitorElementData::Scope::LUMI) {
0095       data_.flags |= DQMNet::DQM_PROP_LUMI;
0096     }
0097 
0098     // these are unsupported and always off.
0099     data_.flags &= ~DQMNet::DQM_PROP_HAS_REFERENCE;
0100     data_.flags &= ~DQMNet::DQM_PROP_TAGGED;
0101     data_.flags &= ~DQMNet::DQM_PROP_RESET;
0102     data_.flags &= ~DQMNet::DQM_PROP_ACCUMULATE;
0103 
0104     // we use ROOT's internal efficiency flag as the truth
0105     data_.flags &= ~DQMNet::DQM_PROP_EFFICIENCY_PLOT;
0106     if (access.value.object_ && access.value.object_->TestBit(TH1::kIsAverage)) {
0107       data_.flags |= DQMNet::DQM_PROP_EFFICIENCY_PLOT;
0108     }
0109 
0110     data_.tag = 0;
0111 
0112     // don't touch version (a timestamp).
0113 
0114     // we could set proper values here, but nobody should use them.
0115     data_.run = 0;
0116     data_.lumi = 0;
0117 
0118     // these are relics from the threaded migration and should not be used anywhere.
0119     data_.streamId = 0;
0120     data_.moduleId = 0;
0121 
0122     // leaking a pointer here, but that should be fine.
0123     data_.dirname = access.key.path_.getDirname();
0124 
0125     data_.objname = access.key.path_.getObjectname();
0126 
0127     data_.flags &= ~DQMNet::DQM_PROP_REPORT_ALARM;
0128     data_.qreports.clear();
0129     for (QReport const &qr : access.value.qreports_) {
0130       data_.qreports.push_back(qr.getValue());
0131       switch (qr.getStatus()) {
0132         case dqm::qstatus::STATUS_OK:
0133           break;
0134         case dqm::qstatus::WARNING:
0135           data_.flags |= DQMNet::DQM_PROP_REPORT_WARN;
0136           break;
0137         case dqm::qstatus::ERROR:
0138           data_.flags |= DQMNet::DQM_PROP_REPORT_ERROR;
0139           break;
0140         default:
0141           data_.flags |= DQMNet::DQM_PROP_REPORT_OTHER;
0142           break;
0143       }
0144     }
0145   }
0146 
0147   MonitorElement::~MonitorElement() {}
0148 
0149   //utility function to check the consistency of the axis labels
0150   //taken from TH1::CheckBinLabels which is not public
0151   bool MonitorElement::CheckBinLabels(const TAxis *a1, const TAxis *a2) {
0152     // check that axis have same labels
0153     THashList *l1 = (const_cast<TAxis *>(a1))->GetLabels();
0154     THashList *l2 = (const_cast<TAxis *>(a2))->GetLabels();
0155 
0156     if (!l1 && !l2)
0157       return true;
0158     if (!l1 || !l2) {
0159       return false;
0160     }
0161     // check now labels sizes  are the same
0162     if (l1->GetSize() != l2->GetSize()) {
0163       return false;
0164     }
0165     for (int i = 1; i <= a1->GetNbins(); ++i) {
0166       TString label1 = a1->GetBinLabel(i);
0167       TString label2 = a2->GetBinLabel(i);
0168       if (label1 != label2) {
0169         return false;
0170       }
0171     }
0172     return true;
0173   }
0174 
0175   /// "Fill" ME methods for string
0176   void MonitorElement::Fill(std::string &value) {
0177     auto access = this->accessMut();
0178     update();
0179     if (kind() == Kind::STRING) {
0180       access.value.scalar_.str = value;
0181     } else {
0182       incompatible(__PRETTY_FUNCTION__);
0183     }
0184   }
0185 
0186   /// "Fill" ME methods for double
0187   void MonitorElement::Fill(double x) {
0188     auto access = this->accessMut();
0189     update();
0190     if (kind() == Kind::INT)
0191       access.value.scalar_.num = static_cast<int64_t>(x);
0192     else if (kind() == Kind::REAL)
0193       access.value.scalar_.real = x;
0194     else if (kind() == Kind::TH1F)
0195       accessRootObject(access, __PRETTY_FUNCTION__, 1)->Fill(x, 1);
0196     else if (kind() == Kind::TH1S)
0197       accessRootObject(access, __PRETTY_FUNCTION__, 1)->Fill(x, 1);
0198     else if (kind() == Kind::TH1I)
0199       accessRootObject(access, __PRETTY_FUNCTION__, 1)->Fill(x, 1);
0200     else if (kind() == Kind::TH1D)
0201       accessRootObject(access, __PRETTY_FUNCTION__, 1)->Fill(x, 1);
0202     else
0203       incompatible(__PRETTY_FUNCTION__);
0204   }
0205 
0206   /// "Fill" ME method for int64_t
0207   void MonitorElement::doFill(int64_t x) {
0208     auto access = this->accessMut();
0209     update();
0210     if (kind() == Kind::INT)
0211       access.value.scalar_.num = static_cast<int64_t>(x);
0212     else if (kind() == Kind::REAL)
0213       access.value.scalar_.real = static_cast<double>(x);
0214     else if (kind() == Kind::TH1F)
0215       accessRootObject(access, __PRETTY_FUNCTION__, 1)->Fill(static_cast<double>(x), 1);
0216     else if (kind() == Kind::TH1S)
0217       accessRootObject(access, __PRETTY_FUNCTION__, 1)->Fill(static_cast<double>(x), 1);
0218     else if (kind() == Kind::TH1I)
0219       accessRootObject(access, __PRETTY_FUNCTION__, 1)->Fill(static_cast<double>(x), 1);
0220     else if (kind() == Kind::TH1D)
0221       accessRootObject(access, __PRETTY_FUNCTION__, 1)->Fill(static_cast<double>(x), 1);
0222     else
0223       incompatible(__PRETTY_FUNCTION__);
0224   }
0225 
0226   /// can be used with 2D (x,y) or 1D (x, w) histograms
0227   void MonitorElement::Fill(double x, double yw) {
0228     auto access = this->accessMut();
0229     update();
0230     if (kind() == Kind::TH1F)
0231       accessRootObject(access, __PRETTY_FUNCTION__, 1)->Fill(x, yw);
0232     else if (kind() == Kind::TH1S)
0233       accessRootObject(access, __PRETTY_FUNCTION__, 1)->Fill(x, yw);
0234     else if (kind() == Kind::TH1D)
0235       accessRootObject(access, __PRETTY_FUNCTION__, 1)->Fill(x, yw);
0236     else if (kind() == Kind::TH1I)
0237       accessRootObject(access, __PRETTY_FUNCTION__, 1)->Fill(x, yw);
0238     else if (kind() == Kind::TH2F)
0239       static_cast<TH2F *>(accessRootObject(access, __PRETTY_FUNCTION__, 2))->Fill(x, yw, 1);
0240     else if (kind() == Kind::TH2S)
0241       static_cast<TH2S *>(accessRootObject(access, __PRETTY_FUNCTION__, 2))->Fill(x, yw, 1);
0242     else if (kind() == Kind::TH2D)
0243       static_cast<TH2D *>(accessRootObject(access, __PRETTY_FUNCTION__, 2))->Fill(x, yw, 1);
0244     else if (kind() == Kind::TH2I)
0245       static_cast<TH2I *>(accessRootObject(access, __PRETTY_FUNCTION__, 2))->Fill(x, yw, 1);
0246     else if (kind() == Kind::TPROFILE)
0247       static_cast<TProfile *>(accessRootObject(access, __PRETTY_FUNCTION__, 1))->Fill(x, yw, 1);
0248     else
0249       incompatible(__PRETTY_FUNCTION__);
0250   }
0251 
0252   /// shift bin to the left and fill last bin with new entry
0253   /// 1st argument is y value, 2nd argument is y error (default 0)
0254   /// can be used with 1D or profile histograms only
0255   void MonitorElement::ShiftFillLast(double y, double ye, int xscale) {
0256     // TODO: this should take the lock only once to be actually safe.
0257     // But since it is not const, we don't even claim it is thread-safe.
0258     update();
0259     if (kind() == Kind::TH1F || kind() == Kind::TH1S || kind() == Kind::TH1D || kind() == Kind::TH1I) {
0260       int nbins = getNbinsX();
0261       auto entries = (int)getEntries();
0262       // first fill bins from left to right
0263       int index = entries + 1;
0264       int xlow = 2;
0265       int xup = nbins;
0266       // if more entries than bins then start shifting
0267       if (entries >= nbins) {
0268         index = nbins;
0269         xlow = entries - nbins + 3;
0270         xup = entries + 1;
0271         // average first bin
0272         double y1 = getBinContent(1);
0273         double y2 = getBinContent(2);
0274         double y1err = getBinError(1);
0275         double y2err = getBinError(2);
0276         double N = entries - nbins + 1.;
0277         if (ye == 0. || y1err == 0. || y2err == 0.) {
0278           // for errors zero calculate unweighted mean and its error
0279           double sum = N * y1 + y2;
0280           y1 = sum / (N + 1.);
0281           // FIXME check if correct
0282           double s = (N + 1.) * (N * y1 * y1 + y2 * y2) - sum * sum;
0283           if (s >= 0.)
0284             y1err = sqrt(s) / (N + 1.);
0285           else
0286             y1err = 0.;
0287         } else {
0288           // for errors non-zero calculate weighted mean and its error
0289           double denom = (1. / y1err + 1. / y2err);
0290           double mean = (y1 / y1err + y2 / y2err) / denom;
0291           // FIXME check if correct
0292           y1err = sqrt(((y1 - mean) * (y1 - mean) / y1err + (y2 - mean) * (y2 - mean) / y2err) / denom / 2.);
0293           y1 = mean;  // set y1 to mean for filling below
0294         }
0295         setBinContent(1, y1);
0296         setBinError(1, y1err);
0297         // shift remaining bins to the left
0298         for (int i = 3; i <= nbins; i++) {
0299           setBinContent(i - 1, getBinContent(i));
0300           setBinError(i - 1, getBinError(i));
0301         }
0302       }
0303       // fill last bin with new values
0304       setBinContent(index, y);
0305       setBinError(index, ye);
0306       // set entries
0307       setEntries(entries + 1);
0308       // set axis labels and reset drawing option
0309       char buffer[10];
0310       sprintf(buffer, "%d", xlow * xscale);
0311       std::string a(buffer);
0312       setBinLabel(2, a);
0313       sprintf(buffer, "%d", xup * xscale);
0314       std::string b(buffer);
0315       setBinLabel(nbins, b);
0316       setBinLabel(1, "av.");
0317     } else
0318       incompatible(__PRETTY_FUNCTION__);
0319   }
0320   /// can be used with 3D (x, y, z) or 2D (x, y, w) histograms
0321   void MonitorElement::Fill(double x, double y, double zw) {
0322     auto access = this->accessMut();
0323     update();
0324     if (kind() == Kind::TH2F)
0325       static_cast<TH2F *>(accessRootObject(access, __PRETTY_FUNCTION__, 2))->Fill(x, y, zw);
0326     else if (kind() == Kind::TH2S)
0327       static_cast<TH2S *>(accessRootObject(access, __PRETTY_FUNCTION__, 2))->Fill(x, y, zw);
0328     else if (kind() == Kind::TH2D)
0329       static_cast<TH2D *>(accessRootObject(access, __PRETTY_FUNCTION__, 2))->Fill(x, y, zw);
0330     else if (kind() == Kind::TH2I)
0331       static_cast<TH2I *>(accessRootObject(access, __PRETTY_FUNCTION__, 2))->Fill(x, y, zw);
0332     else if (kind() == Kind::TH3F)
0333       static_cast<TH3F *>(accessRootObject(access, __PRETTY_FUNCTION__, 2))->Fill(x, y, zw, 1);
0334     else if (kind() == Kind::TPROFILE)
0335       static_cast<TProfile *>(accessRootObject(access, __PRETTY_FUNCTION__, 2))->Fill(x, y, zw);
0336     else if (kind() == Kind::TPROFILE2D)
0337       static_cast<TProfile2D *>(accessRootObject(access, __PRETTY_FUNCTION__, 2))->Fill(x, y, zw, 1);
0338     else
0339       incompatible(__PRETTY_FUNCTION__);
0340   }
0341 
0342   /// can be used with 3D (x, y, z, w) histograms
0343   void MonitorElement::Fill(double x, double y, double z, double w) {
0344     auto access = this->accessMut();
0345     update();
0346     if (kind() == Kind::TH3F)
0347       static_cast<TH3F *>(accessRootObject(access, __PRETTY_FUNCTION__, 2))->Fill(x, y, z, w);
0348     else if (kind() == Kind::TPROFILE2D)
0349       static_cast<TProfile2D *>(accessRootObject(access, __PRETTY_FUNCTION__, 2))->Fill(x, y, z, w);
0350     else
0351       incompatible(__PRETTY_FUNCTION__);
0352   }
0353 
0354   /// reset ME (ie. contents, errors, etc)
0355   void MonitorElement::Reset() {
0356     auto access = this->accessMut();
0357     update();
0358     if (kind() == Kind::INT)
0359       access.value.scalar_.num = 0;
0360     else if (kind() == Kind::REAL)
0361       access.value.scalar_.real = 0;
0362     else if (kind() == Kind::STRING)
0363       access.value.scalar_.str.clear();
0364     else
0365       return accessRootObject(access, __PRETTY_FUNCTION__, 1)->Reset();
0366   }
0367 
0368   /// convert scalar data into a string.
0369   void MonitorElement::packScalarData(std::string &into, const char *prefix) const {
0370     auto access = this->access();
0371     char buf[64];
0372     if (kind() == Kind::INT) {
0373       snprintf(buf, sizeof(buf), "%s%" PRId64, prefix, access.value.scalar_.num);
0374       into = buf;
0375     } else if (kind() == Kind::REAL) {
0376       snprintf(buf, sizeof(buf), "%s%.*g", prefix, DBL_DIG + 2, access.value.scalar_.real);
0377       into = buf;
0378     } else if (kind() == Kind::STRING) {
0379       into.reserve(strlen(prefix) + access.value.scalar_.str.size());
0380       into += prefix;
0381       into += access.value.scalar_.str;
0382     } else
0383       incompatible(__PRETTY_FUNCTION__);
0384   }
0385 
0386   /// serialise quality report information into a string.
0387   void MonitorElement::packQualityData(std::string &into) const { DQMNet::packQualityData(into, data_.qreports); }
0388 
0389   /// returns value of ME in string format (eg. "f = 3.14151926" for double numbers);
0390   /// relevant only for scalar or string MEs
0391   std::string MonitorElement::valueString() const {
0392     std::string result;
0393     if (kind() == Kind::INT)
0394       packScalarData(result, "i=");
0395     else if (kind() == Kind::REAL)
0396       packScalarData(result, "f=");
0397     else if (kind() == Kind::STRING)
0398       packScalarData(result, "s=");
0399     else
0400       incompatible(__PRETTY_FUNCTION__);
0401 
0402     return result;
0403   }
0404 
0405   /// return tagged value of ME in string format
0406   /// (eg. <name>f=3.14151926</name> for double numbers);
0407   /// relevant only for sending scalar or string MEs over TSocket
0408   std::string MonitorElement::tagString() const {
0409     std::string result;
0410     std::string val(valueString());
0411     result.reserve(6 + 2 * data_.objname.size() + val.size());
0412     result += '<';
0413     result += data_.objname;
0414     result += '>';
0415     result += val;
0416     result += '<';
0417     result += '/';
0418     result += data_.objname;
0419     result += '>';
0420     return result;
0421   }
0422 
0423   /// return label string for the monitor element tag (eg. <name>t=12345</name>)
0424   std::string MonitorElement::tagLabelString() const {
0425     char buf[32];
0426     std::string result;
0427     size_t len = sprintf(buf, "t=%" PRIu32, data_.tag);
0428 
0429     result.reserve(6 + 2 * data_.objname.size() + len);
0430     result += '<';
0431     result += data_.objname;
0432     result += '>';
0433     result += buf;
0434     result += '<';
0435     result += '/';
0436     result += data_.objname;
0437     result += '>';
0438     return result;
0439   }
0440 
0441   /// return label string for the monitor element tag (eg. <name>t=12345</name>)
0442   std::string MonitorElement::effLabelString() const {
0443     std::string result;
0444 
0445     result.reserve(6 + 2 * data_.objname.size() + 3);
0446     result += '<';
0447     result += data_.objname;
0448     result += '>';
0449     result += "e=1";
0450     result += '<';
0451     result += '/';
0452     result += data_.objname;
0453     result += '>';
0454     return result;
0455   }
0456 
0457   std::string MonitorElement::qualityTagString(const DQMNet::QValue &qv) const {
0458     char buf[64];
0459     std::string result;
0460     size_t titlelen = data_.objname.size() + qv.qtname.size() + 1;
0461     size_t buflen = sprintf(buf, "qr=st:%d:%.*g:", qv.code, DBL_DIG + 2, qv.qtresult);
0462 
0463     result.reserve(7 + 2 * titlelen + buflen + qv.algorithm.size() + qv.message.size());
0464     result += '<';
0465     result += data_.objname;
0466     result += '.';
0467     result += qv.qtname;
0468     result += '>';
0469     result += buf;
0470     result += qv.algorithm;
0471     result += ':';
0472     result += qv.message;
0473     result += '<';
0474     result += '/';
0475     result += data_.objname;
0476     result += '.';
0477     result += qv.qtname;
0478     result += '>';
0479     return result;
0480   }
0481 
0482   const MonitorElementData::QReport *MonitorElement::getQReport(const std::string &qtname) const {
0483     MonitorElementData::MonitorElementData::QReport *qr;
0484     DQMNet::QValue *qv;
0485     const_cast<MonitorElement *>(this)->getQReport(false, qtname, qr, qv);
0486     return qr;
0487   }
0488 
0489   template <typename FILTER>
0490   std::vector<MonitorElementData::QReport *> MonitorElement::filterQReports(FILTER filter) const {
0491     auto access = this->access();
0492     std::vector<MonitorElementData::QReport *> result;
0493     for (MonitorElementData::QReport const &qr : access.value.qreports_) {
0494       if (filter(qr)) {
0495         // const_cast here because this API always violated cons'ness. Should
0496         // make the result type const and fix all usages.
0497         result.push_back(const_cast<MonitorElementData::QReport *>(&qr));
0498       }
0499     }
0500     return result;
0501   }
0502 
0503   std::vector<MonitorElementData::QReport *> MonitorElement::getQReports() const {
0504     return filterQReports([](MonitorElementData::QReport const &qr) { return true; });
0505   }
0506 
0507   std::vector<MonitorElementData::QReport *> MonitorElement::getQWarnings() const {
0508     return filterQReports(
0509         [](MonitorElementData::QReport const &qr) { return qr.getStatus() == dqm::qstatus::WARNING; });
0510   }
0511 
0512   std::vector<MonitorElementData::QReport *> MonitorElement::getQErrors() const {
0513     return filterQReports([](MonitorElementData::QReport const &qr) { return qr.getStatus() == dqm::qstatus::ERROR; });
0514   }
0515 
0516   std::vector<MonitorElementData::QReport *> MonitorElement::getQOthers() const {
0517     return filterQReports([](MonitorElementData::QReport const &qr) {
0518       return qr.getStatus() != dqm::qstatus::STATUS_OK && qr.getStatus() != dqm::qstatus::WARNING &&
0519              qr.getStatus() != dqm::qstatus::ERROR;
0520     });
0521   }
0522 
0523   void MonitorElement::incompatible(const char *func) const {
0524     throw cms::Exception("MonitorElementError") << "Method '" << func
0525                                                 << "' cannot be invoked on monitor"
0526                                                    " element '"
0527                                                 << data_.objname << "'";
0528   }
0529 
0530   TH1 const *MonitorElement::accessRootObject(Access const &access, const char *func, int reqdim) const {
0531     if (kind() < Kind::TH1F)
0532       throw cms::Exception("MonitorElement") << "Method '" << func
0533                                              << "' cannot be invoked on monitor"
0534                                                 " element '"
0535                                              << data_.objname << "' because it is not a root object";
0536     return access.value.object_.get();
0537   }
0538   TH1 *MonitorElement::accessRootObject(AccessMut const &access, const char *func, int reqdim) const {
0539     if (kind() < Kind::TH1F)
0540       throw cms::Exception("MonitorElement") << "Method '" << func
0541                                              << "' cannot be invoked on monitor"
0542                                                 " element '"
0543                                              << data_.objname << "' because it is not a root object";
0544     return checkRootObject(data_.objname, access.value.object_.get(), func, reqdim);
0545   }
0546 
0547   /*** getter methods (wrapper around ROOT methods) ****/
0548   //
0549   /// get mean value of histogram along x, y or z axis (axis=1, 2, 3 respectively)
0550   double MonitorElement::getMean(int axis /* = 1 */) const {
0551     auto access = this->access();
0552     return accessRootObject(access, __PRETTY_FUNCTION__, axis - 1)->GetMean(axis);
0553   }
0554 
0555   /// get mean value uncertainty of histogram along x, y or z axis
0556   /// (axis=1, 2, 3 respectively)
0557   double MonitorElement::getMeanError(int axis /* = 1 */) const {
0558     auto access = this->access();
0559     return accessRootObject(access, __PRETTY_FUNCTION__, axis - 1)->GetMeanError(axis);
0560   }
0561 
0562   /// get RMS of histogram along x, y or z axis (axis=1, 2, 3 respectively)
0563   double MonitorElement::getRMS(int axis /* = 1 */) const {
0564     auto access = this->access();
0565     return accessRootObject(access, __PRETTY_FUNCTION__, axis - 1)->GetRMS(axis);
0566   }
0567 
0568   /// get RMS uncertainty of histogram along x, y or z axis(axis=1,2,3 respectively)
0569   double MonitorElement::getRMSError(int axis /* = 1 */) const {
0570     auto access = this->access();
0571     return accessRootObject(access, __PRETTY_FUNCTION__, axis - 1)->GetRMSError(axis);
0572   }
0573 
0574   /// get # of bins in X-axis
0575   int MonitorElement::getNbinsX() const {
0576     auto access = this->access();
0577     return accessRootObject(access, __PRETTY_FUNCTION__, 1)->GetNbinsX();
0578   }
0579 
0580   /// get # of bins in Y-axis
0581   int MonitorElement::getNbinsY() const {
0582     auto access = this->access();
0583     return accessRootObject(access, __PRETTY_FUNCTION__, 2)->GetNbinsY();
0584   }
0585 
0586   /// get # of bins in Z-axis
0587   int MonitorElement::getNbinsZ() const {
0588     auto access = this->access();
0589     return accessRootObject(access, __PRETTY_FUNCTION__, 3)->GetNbinsZ();
0590   }
0591 
0592   /// get content of bin (1-D)
0593   double MonitorElement::getBinContent(int binx) const {
0594     auto access = this->access();
0595     return accessRootObject(access, __PRETTY_FUNCTION__, 1)->GetBinContent(binx);
0596   }
0597 
0598   /// get content of bin (2-D)
0599   double MonitorElement::getBinContent(int binx, int biny) const {
0600     auto access = this->access();
0601     return accessRootObject(access, __PRETTY_FUNCTION__, 2)->GetBinContent(binx, biny);
0602   }
0603 
0604   /// get content of bin (3-D)
0605   double MonitorElement::getBinContent(int binx, int biny, int binz) const {
0606     auto access = this->access();
0607     return accessRootObject(access, __PRETTY_FUNCTION__, 3)->GetBinContent(binx, biny, binz);
0608   }
0609 
0610   /// get uncertainty on content of bin (1-D) - See TH1::GetBinError for details
0611   double MonitorElement::getBinError(int binx) const {
0612     auto access = this->access();
0613     return accessRootObject(access, __PRETTY_FUNCTION__, 1)->GetBinError(binx);
0614   }
0615 
0616   /// get uncertainty on content of bin (2-D) - See TH1::GetBinError for details
0617   double MonitorElement::getBinError(int binx, int biny) const {
0618     auto access = this->access();
0619     return accessRootObject(access, __PRETTY_FUNCTION__, 2)->GetBinError(binx, biny);
0620   }
0621 
0622   /// get uncertainty on content of bin (3-D) - See TH1::GetBinError for details
0623   double MonitorElement::getBinError(int binx, int biny, int binz) const {
0624     auto access = this->access();
0625     return accessRootObject(access, __PRETTY_FUNCTION__, 3)->GetBinError(binx, biny, binz);
0626   }
0627 
0628   /// get # of entries
0629   double MonitorElement::getEntries() const {
0630     auto access = this->access();
0631     return accessRootObject(access, __PRETTY_FUNCTION__, 1)->GetEntries();
0632   }
0633 
0634   /// get global bin number (for 2-D profiles)
0635   int MonitorElement::getBin(int binx, int biny) const {
0636     auto access = this->access();
0637     if (kind() == Kind::TPROFILE2D)
0638       return static_cast<TProfile2D const *>(accessRootObject(access, __PRETTY_FUNCTION__, 1))->GetBin(binx, biny);
0639     else {
0640       incompatible(__PRETTY_FUNCTION__);
0641       return 0;
0642     }
0643   }
0644 
0645   /// get # of bin entries (for profiles)
0646   double MonitorElement::getBinEntries(int bin) const {
0647     auto access = this->access();
0648     if (kind() == Kind::TPROFILE)
0649       return static_cast<TProfile const *>(accessRootObject(access, __PRETTY_FUNCTION__, 1))->GetBinEntries(bin);
0650     else if (kind() == Kind::TPROFILE2D)
0651       return static_cast<TProfile2D const *>(accessRootObject(access, __PRETTY_FUNCTION__, 1))->GetBinEntries(bin);
0652     else {
0653       incompatible(__PRETTY_FUNCTION__);
0654       return 0;
0655     }
0656   }
0657 
0658   /// get # of bin entries (for 2-D profiles)
0659   double MonitorElement::getBinEntries(int binx, int biny) const {
0660     auto access = this->access();
0661     if (kind() == Kind::TPROFILE2D) {
0662       int globBin =
0663           static_cast<TProfile2D const *>(accessRootObject(access, __PRETTY_FUNCTION__, 1))->GetBin(binx, biny);
0664       return static_cast<TProfile2D const *>(accessRootObject(access, __PRETTY_FUNCTION__, 1))->GetBinEntries(globBin);
0665     } else {
0666       incompatible(__PRETTY_FUNCTION__);
0667       return 0;
0668     }
0669   }
0670 
0671   /// get integral of bins
0672   double MonitorElement::integral() const {
0673     auto access = this->access();
0674     return accessRootObject(access, __PRETTY_FUNCTION__, 1)->Integral();
0675   }
0676 
0677   /// get x-, y- or z-axis title (axis=1, 2, 3 respectively)
0678   std::string MonitorElement::getAxisTitle(int axis /* = 1 */) const {
0679     auto access = this->access();
0680     return getAxis(access, __PRETTY_FUNCTION__, axis)->GetTitle();
0681   }
0682 
0683   /// get MonitorElement title
0684   std::string MonitorElement::getTitle() const {
0685     auto access = this->access();
0686     return accessRootObject(access, __PRETTY_FUNCTION__, 1)->GetTitle();
0687   }
0688 
0689   /*** setter methods (wrapper around ROOT methods) ****/
0690   //
0691   /// set content of bin (1-D)
0692   void MonitorElement::setBinContent(int binx, double content) {
0693     auto access = this->accessMut();
0694     accessRootObject(access, __PRETTY_FUNCTION__, 1)->SetBinContent(binx, content);
0695   }
0696 
0697   /// set content of bin (2-D)
0698   void MonitorElement::setBinContent(int binx, int biny, double content) {
0699     auto access = this->accessMut();
0700     accessRootObject(access, __PRETTY_FUNCTION__, 2)->SetBinContent(binx, biny, content);
0701   }
0702 
0703   /// set content of bin (3-D)
0704   void MonitorElement::setBinContent(int binx, int biny, int binz, double content) {
0705     auto access = this->accessMut();
0706     accessRootObject(access, __PRETTY_FUNCTION__, 3)->SetBinContent(binx, biny, binz, content);
0707   }
0708 
0709   /// set uncertainty on content of bin (1-D)
0710   void MonitorElement::setBinError(int binx, double error) {
0711     auto access = this->accessMut();
0712     accessRootObject(access, __PRETTY_FUNCTION__, 1)->SetBinError(binx, error);
0713   }
0714 
0715   /// set uncertainty on content of bin (2-D)
0716   void MonitorElement::setBinError(int binx, int biny, double error) {
0717     auto access = this->accessMut();
0718     accessRootObject(access, __PRETTY_FUNCTION__, 2)->SetBinError(binx, biny, error);
0719   }
0720 
0721   /// set uncertainty on content of bin (3-D)
0722   void MonitorElement::setBinError(int binx, int biny, int binz, double error) {
0723     auto access = this->accessMut();
0724     accessRootObject(access, __PRETTY_FUNCTION__, 3)->SetBinError(binx, biny, binz, error);
0725   }
0726 
0727   /// set # of bin entries (to be used for profiles)
0728   void MonitorElement::setBinEntries(int bin, double nentries) {
0729     auto access = this->accessMut();
0730     if (kind() == Kind::TPROFILE)
0731       static_cast<TProfile *>(accessRootObject(access, __PRETTY_FUNCTION__, 1))->SetBinEntries(bin, nentries);
0732     else if (kind() == Kind::TPROFILE2D)
0733       static_cast<TProfile2D *>(accessRootObject(access, __PRETTY_FUNCTION__, 1))->SetBinEntries(bin, nentries);
0734     else
0735       incompatible(__PRETTY_FUNCTION__);
0736   }
0737 
0738   /// set # of entries
0739   void MonitorElement::setEntries(double nentries) {
0740     auto access = this->accessMut();
0741     accessRootObject(access, __PRETTY_FUNCTION__, 1)->SetEntries(nentries);
0742   }
0743 
0744   /// Replace entries with results of dividing num by denom
0745   void MonitorElement::divide(
0746       const MonitorElement *num, const MonitorElement *denom, double c1, double c2, const char *options) {
0747     if (num->kind() < Kind::TH1F)
0748       num->incompatible(__PRETTY_FUNCTION__);
0749     if (denom->kind() < Kind::TH1F)
0750       denom->incompatible(__PRETTY_FUNCTION__);
0751 
0752     TH1 const *numH = static_cast<TH1 const *>(num->getRootObject());
0753     TH1 const *denomH = static_cast<TH1 const *>(denom->getRootObject());
0754     TH1 *thisH = getTH1();
0755 
0756     //Need to take locks in a consistent order to avoid deadlocks. Use pointer value order of underlying ROOT object..
0757     //This is known as the monitor pattern.
0758     std::array<const MonitorElement *, 3> order{{this, num, denom}};
0759     std::sort(order.begin(), order.end(), [](auto const *lhs, auto const *rhs) {
0760       return lhs->mutable_->data_.value_.object_.get() < rhs->mutable_->data_.value_.object_.get();
0761     });
0762 
0763     auto a0 = order[0]->access();
0764     auto a1 = order[1]->access();
0765     auto a2 = order[2]->access();
0766 
0767     //Have ROOT do check that the types are compatible
0768     thisH->Divide(numH, denomH, c1, c2, options);
0769   }
0770 
0771   /// set bin label for x, y or z axis (axis=1, 2, 3 respectively)
0772   void MonitorElement::setBinLabel(int bin, const std::string &label, int axis /* = 1 */) {
0773     bool fail = false;
0774     {
0775       auto access = this->accessMut();
0776       update();
0777       if (getAxis(access, __PRETTY_FUNCTION__, axis)->GetNbins() >= bin) {
0778         getAxis(access, __PRETTY_FUNCTION__, axis)->SetBinLabel(bin, label.c_str());
0779       } else {
0780         fail = true;
0781       }
0782     }
0783     // do this with the ME lock released to prevent a deadlock
0784     if (fail) {
0785       // this also takes the lock, make sure to release it before going to edm
0786       // (which might take more locks)
0787       auto name = getFullname();
0788       edm::LogWarning("MonitorElement") << "*** MonitorElement: WARNING:"
0789                                         << "setBinLabel: attempting to set label of non-existent bin number for ME: "
0790                                         << name << " \n";
0791     }
0792   }
0793 
0794   /// set x-, y- or z-axis range (axis=1, 2, 3 respectively)
0795   void MonitorElement::setAxisRange(double xmin, double xmax, int axis /* = 1 */) {
0796     auto access = this->accessMut();
0797     getAxis(access, __PRETTY_FUNCTION__, axis)->SetRangeUser(xmin, xmax);
0798   }
0799 
0800   /// set x-, y- or z-axis title (axis=1, 2, 3 respectively)
0801   void MonitorElement::setAxisTitle(const std::string &title, int axis /* = 1 */) {
0802     auto access = this->accessMut();
0803     getAxis(access, __PRETTY_FUNCTION__, axis)->SetTitle(title.c_str());
0804   }
0805 
0806   /// set x-, y-, or z-axis to display time values
0807   void MonitorElement::setAxisTimeDisplay(int value, int axis /* = 1 */) {
0808     auto access = this->accessMut();
0809     getAxis(access, __PRETTY_FUNCTION__, axis)->SetTimeDisplay(value);
0810   }
0811 
0812   /// set the format of the time values that are displayed on an axis
0813   void MonitorElement::setAxisTimeFormat(const char *format /* = "" */, int axis /* = 1 */) {
0814     auto access = this->accessMut();
0815     getAxis(access, __PRETTY_FUNCTION__, axis)->SetTimeFormat(format);
0816   }
0817 
0818   /// set (ie. change) histogram/profile title
0819   void MonitorElement::setTitle(const std::string &title) {
0820     auto access = this->accessMut();
0821     accessRootObject(access, __PRETTY_FUNCTION__, 1)->SetTitle(title.c_str());
0822   }
0823 
0824   TAxis *MonitorElement::getAxis(AccessMut const &access, const char *func, int axis) const {
0825     TH1 *h = accessRootObject(access, func, axis - 1);
0826     TAxis *a = nullptr;
0827     if (axis == 1)
0828       a = h->GetXaxis();
0829     else if (axis == 2)
0830       a = h->GetYaxis();
0831     else if (axis == 3)
0832       a = h->GetZaxis();
0833 
0834     if (!a)
0835       throw cms::Exception("MonitorElementError") << "No such axis " << axis
0836                                                   << " in monitor element"
0837                                                      " '"
0838                                                   << data_.objname << "' of type '" << typeid(*h).name() << "'";
0839 
0840     return a;
0841   }
0842 
0843   TAxis const *MonitorElement::getAxis(Access const &access, const char *func, int axis) const {
0844     TH1 const *h = accessRootObject(access, func, axis - 1);
0845     TAxis const *a = nullptr;
0846     if (axis == 1)
0847       a = h->GetXaxis();
0848     else if (axis == 2)
0849       a = h->GetYaxis();
0850     else if (axis == 3)
0851       a = h->GetZaxis();
0852 
0853     if (!a)
0854       throw cms::Exception("MonitorElementError") << "No such axis " << axis
0855                                                   << " in monitor element"
0856                                                      " '"
0857                                                   << data_.objname << "' of type '" << typeid(*h).name() << "'";
0858 
0859     return a;
0860   }
0861 
0862   void MonitorElement::setXTitle(std::string const &title) {
0863     auto access = this->accessMut();
0864     update();
0865     access.value.object_->SetXTitle(title.c_str());
0866   }
0867   void MonitorElement::setYTitle(std::string const &title) {
0868     auto access = this->accessMut();
0869     update();
0870     access.value.object_->SetYTitle(title.c_str());
0871   }
0872 
0873   void MonitorElement::enableSumw2() {
0874     auto access = this->accessMut();
0875     update();
0876     if (access.value.object_->GetSumw2() == nullptr) {
0877       access.value.object_->Sumw2();
0878     }
0879   }
0880 
0881   void MonitorElement::disableAlphanumeric() {
0882     auto access = this->accessMut();
0883     update();
0884     access.value.object_->GetXaxis()->SetNoAlphanumeric(false);
0885     access.value.object_->GetYaxis()->SetNoAlphanumeric(false);
0886   }
0887 
0888   void MonitorElement::setOption(const char *option) {
0889     auto access = this->accessMut();
0890     update();
0891     access.value.object_->SetOption(option);
0892   }
0893   double MonitorElement::getAxisMin(int axis) const {
0894     auto access = this->access();
0895     return getAxis(access, __PRETTY_FUNCTION__, axis)->GetXmin();
0896   }
0897 
0898   double MonitorElement::getAxisMax(int axis) const {
0899     auto access = this->access();
0900     return getAxis(access, __PRETTY_FUNCTION__, axis)->GetXmax();
0901   }
0902 
0903   void MonitorElement::setCanExtend(unsigned int value) {
0904     auto access = this->accessMut();
0905     access.value.object_->SetCanExtend(value);
0906   }
0907 
0908   void MonitorElement::setStatOverflows(bool value) {
0909     auto access = this->accessMut();
0910     if (value == kTRUE)
0911       access.value.object_->SetStatOverflows(TH1::kConsider);
0912     else
0913       access.value.object_->SetStatOverflows(TH1::kIgnore);
0914   }
0915 
0916   bool MonitorElement::getStatOverflows() {
0917     auto access = this->accessMut();
0918     auto value = access.value.object_->GetStatOverflows();
0919     if (value == TH1::kConsider)
0920       return true;
0921     else
0922       return false;
0923   }
0924 
0925   int64_t MonitorElement::getIntValue() const {
0926     assert(kind() == Kind::INT);
0927     auto access = this->access();
0928     return access.value.scalar_.num;
0929   }
0930   double MonitorElement::getFloatValue() const {
0931     assert(kind() == Kind::REAL);
0932     auto access = this->access();
0933     return access.value.scalar_.real;
0934   }
0935   const std::string &MonitorElement::getStringValue() const {
0936     assert(kind() == Kind::STRING);
0937     auto access = this->access();
0938     return access.value.scalar_.str;
0939   }
0940 
0941   void MonitorElement::getQReport(bool create,
0942                                   const std::string &qtname,
0943                                   MonitorElementData::QReport *&qr,
0944                                   DQMNet::QValue *&qv) {
0945     auto access = this->accessMut();
0946 
0947     syncCoreObject(access);
0948 
0949     assert(access.value.qreports_.size() == data_.qreports.size());
0950 
0951     qr = nullptr;
0952     qv = nullptr;
0953 
0954     size_t pos = 0, end = access.value.qreports_.size();
0955     while (pos < end && data_.qreports[pos].qtname != qtname)
0956       ++pos;
0957 
0958     if (pos == end && !create)
0959       return;
0960     else if (pos == end) {
0961       DQMNet::QValue q;
0962       q.code = dqm::qstatus::DID_NOT_RUN;
0963       q.qtresult = 0;
0964       q.qtname = qtname;
0965       q.message = "NO_MESSAGE_ASSIGNED";
0966       q.algorithm = "UNKNOWN_ALGORITHM";
0967       access.value.qreports_.push_back(MonitorElementData::QReport(q));
0968       syncCoreObject(access);
0969     }
0970 
0971     qr = &access.value.qreports_[pos];
0972     qv = &(qr->getValue());
0973   }
0974 
0975   // -------------------------------------------------------------------
0976   // TODO: all of these are UNSAFE and have to be NON-const.
0977   TObject const *MonitorElement::getRootObject() const {
0978     auto access = this->access();
0979     return access.value.object_.get();
0980   }
0981 
0982   TH1 *MonitorElement::getTH1() {
0983     auto access = this->accessMut();
0984     return accessRootObject(access, __PRETTY_FUNCTION__, 0);
0985   }
0986 
0987   TH1F *MonitorElement::getTH1F() {
0988     auto access = this->accessMut();
0989     assert(kind() == Kind::TH1F);
0990     return static_cast<TH1F *>(accessRootObject(access, __PRETTY_FUNCTION__, 1));
0991   }
0992 
0993   TH1S *MonitorElement::getTH1S() {
0994     auto access = this->accessMut();
0995     assert(kind() == Kind::TH1S);
0996     return static_cast<TH1S *>(accessRootObject(access, __PRETTY_FUNCTION__, 1));
0997   }
0998 
0999   TH1I *MonitorElement::getTH1I() {
1000     auto access = this->accessMut();
1001     assert(kind() == Kind::TH1I);
1002     return static_cast<TH1I *>(accessRootObject(access, __PRETTY_FUNCTION__, 1));
1003   }
1004 
1005   TH1D *MonitorElement::getTH1D() {
1006     auto access = this->accessMut();
1007     assert(kind() == Kind::TH1D);
1008     return static_cast<TH1D *>(accessRootObject(access, __PRETTY_FUNCTION__, 1));
1009   }
1010 
1011   TH2F *MonitorElement::getTH2F() {
1012     auto access = this->accessMut();
1013     assert(kind() == Kind::TH2F);
1014     return static_cast<TH2F *>(accessRootObject(access, __PRETTY_FUNCTION__, 2));
1015   }
1016 
1017   TH2S *MonitorElement::getTH2S() {
1018     auto access = this->accessMut();
1019     assert(kind() == Kind::TH2S);
1020     return static_cast<TH2S *>(accessRootObject(access, __PRETTY_FUNCTION__, 2));
1021   }
1022 
1023   TH2I *MonitorElement::getTH2I() {
1024     auto access = this->accessMut();
1025     assert(kind() == Kind::TH2I);
1026     return static_cast<TH2I *>(accessRootObject(access, __PRETTY_FUNCTION__, 2));
1027   }
1028 
1029   TH2D *MonitorElement::getTH2D() {
1030     auto access = this->accessMut();
1031     assert(kind() == Kind::TH2D);
1032     return static_cast<TH2D *>(accessRootObject(access, __PRETTY_FUNCTION__, 2));
1033   }
1034 
1035   TH3F *MonitorElement::getTH3F() {
1036     auto access = this->accessMut();
1037     assert(kind() == Kind::TH3F);
1038     return static_cast<TH3F *>(accessRootObject(access, __PRETTY_FUNCTION__, 3));
1039   }
1040 
1041   TProfile *MonitorElement::getTProfile() {
1042     auto access = this->accessMut();
1043     assert(kind() == Kind::TPROFILE);
1044     return static_cast<TProfile *>(accessRootObject(access, __PRETTY_FUNCTION__, 1));
1045   }
1046 
1047   TProfile2D *MonitorElement::getTProfile2D() {
1048     auto access = this->accessMut();
1049     assert(kind() == Kind::TPROFILE2D);
1050     return static_cast<TProfile2D *>(accessRootObject(access, __PRETTY_FUNCTION__, 2));
1051   }
1052 
1053 }  // namespace dqm::impl