CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/DQMServices/Core/src/MonitorElement.cc

Go to the documentation of this file.
00001 #define __STDC_FORMAT_MACROS 1
00002 #define DQM_ROOT_METHODS 1
00003 #include "DQMServices/Core/interface/MonitorElement.h"
00004 #include "DQMServices/Core/interface/QTest.h"
00005 #include "DQMServices/Core/src/DQMError.h"
00006 #include "TClass.h"
00007 #include "TMath.h"
00008 #include "TList.h"
00009 #include <iostream>
00010 #include <cassert>
00011 #include <cfloat>
00012 #include <inttypes.h>
00013 
00014 static TH1 *
00015 checkRootObject(const std::string &name, TObject *tobj, const char *func, int reqdim)
00016 {
00017   if (! tobj)
00018     raiseDQMError("MonitorElement", "Method '%s' cannot be invoked on monitor"
00019                   " element '%s' because it is not a ROOT object.",
00020                   func, name.c_str());
00021 
00022   TH1 *h = static_cast<TH1 *>(tobj);
00023   int ndim = h->GetDimension();
00024   if (reqdim < 0 || reqdim > ndim)
00025     raiseDQMError("MonitorElement", "Method '%s' cannot be invoked on monitor"
00026                   " element '%s' because it requires %d dimensions; this"
00027                   " object of type '%s' has %d dimensions",
00028                   func, name.c_str(), reqdim, typeid(*h).name(), ndim);
00029 
00030   return h;
00031 }
00032 
00033 MonitorElement *
00034 MonitorElement::initialise(Kind kind)
00035 {
00036   switch (kind)
00037   {
00038   case DQM_KIND_INT:
00039   case DQM_KIND_REAL:
00040   case DQM_KIND_STRING:
00041   case DQM_KIND_TH1F:
00042   case DQM_KIND_TH1S:
00043   case DQM_KIND_TH1D:
00044   case DQM_KIND_TH2F:
00045   case DQM_KIND_TH2S:
00046   case DQM_KIND_TH2D:
00047   case DQM_KIND_TH3F:
00048   case DQM_KIND_TPROFILE:
00049   case DQM_KIND_TPROFILE2D:
00050     data_.flags &= ~DQMNet::DQM_PROP_TYPE_MASK;
00051     data_.flags |= kind;
00052     break;
00053 
00054   default:
00055     raiseDQMError("MonitorElement", "cannot initialise monitor element"
00056                   " to invalid type %d", (int) kind);
00057   }
00058 
00059   return this;
00060 }
00061 
00062 MonitorElement *
00063 MonitorElement::initialise(Kind kind, TH1 *rootobj)
00064 {
00065   initialise(kind);
00066   switch (kind)
00067   {
00068   case DQM_KIND_TH1F:
00069     assert(dynamic_cast<TH1F *>(rootobj));
00070     assert(! reference_ || dynamic_cast<TH1F *>(reference_));
00071     object_ = rootobj;
00072     break;
00073 
00074   case DQM_KIND_TH1S:
00075     assert(dynamic_cast<TH1S *>(rootobj));
00076     assert(! reference_ || dynamic_cast<TH1S *>(reference_));
00077     object_ = rootobj;
00078     break;
00079 
00080   case DQM_KIND_TH1D:
00081     assert(dynamic_cast<TH1D *>(rootobj));
00082     assert(! reference_ || dynamic_cast<TH1D *>(reference_));
00083     object_ = rootobj;
00084     break;
00085 
00086   case DQM_KIND_TH2F:
00087     assert(dynamic_cast<TH2F *>(rootobj));
00088     assert(! reference_ || dynamic_cast<TH2F *>(reference_));
00089     object_ = rootobj;
00090     break;
00091 
00092   case DQM_KIND_TH2S:
00093     assert(dynamic_cast<TH2S *>(rootobj));
00094     assert(! reference_ || dynamic_cast<TH2S *>(reference_));
00095     object_ = rootobj;
00096     break;
00097 
00098   case DQM_KIND_TH2D:
00099     assert(dynamic_cast<TH2D *>(rootobj));
00100     assert(! reference_ || dynamic_cast<TH1D *>(reference_));
00101     object_ = rootobj;
00102     break;
00103 
00104   case DQM_KIND_TH3F:
00105     assert(dynamic_cast<TH3F *>(rootobj));
00106     assert(! reference_ || dynamic_cast<TH3F *>(reference_));
00107     object_ = rootobj;
00108     break;
00109 
00110   case DQM_KIND_TPROFILE:
00111     assert(dynamic_cast<TProfile *>(rootobj));
00112     assert(! reference_ || dynamic_cast<TProfile *>(reference_));
00113     object_ = rootobj;
00114     break;
00115 
00116   case DQM_KIND_TPROFILE2D:
00117     assert(dynamic_cast<TProfile2D *>(rootobj));
00118     assert(! reference_ || dynamic_cast<TProfile2D *>(reference_));
00119     object_ = rootobj;
00120     break;
00121 
00122   default:
00123     raiseDQMError("MonitorElement", "cannot initialise monitor element"
00124                   " as a root object with type %d", (int) kind);
00125   }
00126 
00127   if (reference_)
00128     data_.flags |= DQMNet::DQM_PROP_HAS_REFERENCE;
00129 
00130   return this;
00131 }
00132 
00133 MonitorElement *
00134 MonitorElement::initialise(Kind kind, const std::string &value)
00135 {
00136   initialise(kind);
00137   if (kind == DQM_KIND_STRING)
00138     scalar_.str = value;
00139   else
00140     raiseDQMError("MonitorElement", "cannot initialise monitor element"
00141                   " as a string with type %d", (int) kind);
00142 
00143   return this;
00144 }
00145 
00146 MonitorElement::MonitorElement(void)
00147   : object_(0),
00148     reference_(0),
00149     refvalue_(0)
00150 {
00151   data_.version = 0;
00152   data_.dirname = 0;
00153   data_.tag = 0;
00154   data_.flags = DQM_KIND_INVALID | DQMNet::DQM_PROP_NEW;
00155   scalar_.num = 0;
00156   scalar_.real = 0;
00157 }
00158 
00159 MonitorElement::MonitorElement(const std::string *path, const std::string &name)
00160   : object_(0),
00161     reference_(0),
00162     refvalue_(0)
00163 {
00164   data_.version = 0;
00165   data_.dirname = path;
00166   data_.objname = name;
00167   data_.tag = 0;
00168   data_.flags = DQM_KIND_INVALID | DQMNet::DQM_PROP_NEW;
00169   scalar_.num = 0;
00170   scalar_.real = 0;
00171 }
00172 
00173 MonitorElement::MonitorElement(const MonitorElement &x)
00174   : data_(x.data_),
00175     scalar_(x.scalar_),
00176     object_(x.object_),
00177     reference_(x.reference_),
00178     refvalue_(x.refvalue_),
00179     qreports_(x.qreports_)
00180 {
00181   if (object_)
00182     object_ = static_cast<TH1 *>(object_->Clone());
00183 
00184   if (refvalue_)
00185     refvalue_ = static_cast<TH1 *>(refvalue_->Clone());
00186 }
00187 
00188 MonitorElement &
00189 MonitorElement::operator=(const MonitorElement &x)
00190 {
00191   if (this != &x)
00192   {
00193     delete object_;
00194     delete refvalue_;
00195 
00196     data_ = x.data_;
00197     scalar_ = x.scalar_;
00198     object_ = x.object_;
00199     reference_ = x.reference_;
00200     refvalue_ = x.refvalue_;
00201     qreports_ = x.qreports_;
00202 
00203     if (object_)
00204       object_ = static_cast<TH1 *>(object_->Clone());
00205 
00206     if (refvalue_)
00207       refvalue_ = static_cast<TH1 *>(refvalue_->Clone());
00208   }
00209 
00210   return *this;
00211 }
00212 
00213 MonitorElement::~MonitorElement(void)
00214 {
00215   delete object_;
00216   delete refvalue_;
00217 }
00218 
00220 void
00221 MonitorElement::Fill(std::string &value)
00222 {
00223   update();
00224   if (kind() == DQM_KIND_STRING)
00225     scalar_.str = value;
00226   else
00227     incompatible(__PRETTY_FUNCTION__);
00228 }
00229 
00231 void
00232 MonitorElement::Fill(double x)
00233 {
00234   update();
00235   if (kind() == DQM_KIND_INT)
00236     scalar_.num = static_cast<int64_t>(x);
00237   else if (kind() == DQM_KIND_REAL)
00238     scalar_.real = x;
00239   else if (kind() == DQM_KIND_TH1F)
00240     accessRootObject(__PRETTY_FUNCTION__, 1)
00241       ->Fill(x, 1);
00242   else if (kind() == DQM_KIND_TH1S)
00243     accessRootObject(__PRETTY_FUNCTION__, 1)
00244       ->Fill(x, 1);
00245   else if (kind() == DQM_KIND_TH1D)
00246     accessRootObject(__PRETTY_FUNCTION__, 1)
00247       ->Fill(x, 1);
00248   else
00249     incompatible(__PRETTY_FUNCTION__);
00250 }
00251 
00253 void
00254 MonitorElement::doFill(int64_t x)
00255 {
00256   update();
00257   if (kind() == DQM_KIND_INT)
00258     scalar_.num = static_cast<int64_t>(x);
00259   else if (kind() == DQM_KIND_REAL)
00260     scalar_.real = static_cast<double>(x);
00261   else if (kind() == DQM_KIND_TH1F)
00262     accessRootObject(__PRETTY_FUNCTION__, 1)
00263       ->Fill(static_cast<double>(x), 1);
00264   else if (kind() == DQM_KIND_TH1S)
00265     accessRootObject(__PRETTY_FUNCTION__, 1)
00266       ->Fill(static_cast<double>(x), 1);
00267   else if (kind() == DQM_KIND_TH1D)
00268     accessRootObject(__PRETTY_FUNCTION__, 1)
00269       ->Fill(static_cast<double>(x), 1);
00270   else
00271     incompatible(__PRETTY_FUNCTION__);
00272 }
00273 
00275 void
00276 MonitorElement::Fill(double x, double yw)
00277 {
00278   update();
00279   if (kind() == DQM_KIND_TH1F)
00280     accessRootObject(__PRETTY_FUNCTION__, 1)
00281       ->Fill(x, yw);
00282   else if (kind() == DQM_KIND_TH1S)
00283     accessRootObject(__PRETTY_FUNCTION__, 1)
00284       ->Fill(x, yw);
00285   else if (kind() == DQM_KIND_TH1D)
00286     accessRootObject(__PRETTY_FUNCTION__, 1)
00287       ->Fill(x, yw);
00288   else if (kind() == DQM_KIND_TH2F)
00289     static_cast<TH2F *>(accessRootObject(__PRETTY_FUNCTION__, 2))
00290       ->Fill(x, yw, 1);
00291   else if (kind() == DQM_KIND_TH2S)
00292     static_cast<TH2S *>(accessRootObject(__PRETTY_FUNCTION__, 2))
00293       ->Fill(x, yw, 1);
00294   else if (kind() == DQM_KIND_TH2D)
00295     static_cast<TH2D *>(accessRootObject(__PRETTY_FUNCTION__, 2))
00296       ->Fill(x, yw, 1);
00297   else if (kind() == DQM_KIND_TPROFILE)
00298     static_cast<TProfile *>(accessRootObject(__PRETTY_FUNCTION__, 1))
00299       ->Fill(x, yw, 1);
00300   else
00301     incompatible(__PRETTY_FUNCTION__);
00302 }
00303 
00307 void
00308 MonitorElement::ShiftFillLast(double y, double ye, int xscale)
00309 {
00310   update();
00311   if (kind() == DQM_KIND_TH1F 
00312       || kind() == DQM_KIND_TH1S 
00313       || kind() == DQM_KIND_TH1D) 
00314   {
00315     int nbins = getNbinsX();
00316     int entries = (int)getEntries();
00317     // first fill bins from left to right
00318     int index = entries + 1 ;
00319     int xlow = 2 ; int xup = nbins ;
00320     // if more entries than bins then start shifting
00321     if ( entries >= nbins ) 
00322     {
00323       index = nbins;
00324       xlow = entries - nbins + 3 ; xup = entries+1 ;
00325       // average first bin
00326       double y1 = getBinContent(1);
00327       double y2 = getBinContent(2);
00328       double y1err = getBinError(1);
00329       double y2err = getBinError(2);
00330       double N = entries - nbins + 1.;
00331       if ( ye == 0. || y1err == 0. || y2err == 0.) 
00332       {
00333         // for errors zero calculate unweighted mean and its error
00334         double sum = N*y1 + y2;
00335         y1 = sum/(N+1.) ;
00336         // FIXME check if correct
00337         double s=(N+1.)*(N*y1*y1 + y2*y2) - sum*sum;
00338         if (s>=0.) 
00339           y1err = sqrt(s)/(N+1.);  
00340         else
00341           y1err = 0.;
00342       }
00343       else 
00344       {
00345         // for errors non-zero calculate weighted mean and its error
00346         double denom = (1./y1err + 1./y2err);
00347         double mean = (y1/y1err + y2/y2err)/denom;
00348         // FIXME check if correct
00349         y1err = sqrt(((y1-mean)*(y1-mean)/y1err +
00350                       (y2-mean)*(y2-mean)/y2err)/denom/2.);
00351         y1 = mean; // set y1 to mean for filling below
00352       }
00353       setBinContent(1,y1);
00354       setBinError(1,y1err);
00355       // shift remaining bins to the left
00356       for ( int i = 3; i <= nbins ; i++) 
00357       {
00358         setBinContent(i-1,getBinContent(i));
00359         setBinError(i-1,getBinError(i));
00360       }
00361     }
00362     // fill last bin with new values
00363     setBinContent(index,y);
00364     setBinError(index,ye); 
00365     // set entries
00366     setEntries(entries+1);
00367     // set axis labels and reset drawing option
00368     char buffer [10];
00369     sprintf (buffer, "%d", xlow*xscale); 
00370     std::string a(buffer); setBinLabel(2,a);
00371     sprintf (buffer, "%d", xup*xscale); 
00372     std::string b(buffer); setBinLabel(nbins,b);
00373     setBinLabel(1,"av.");
00374   }
00375   else
00376     incompatible(__PRETTY_FUNCTION__);
00377 }
00379 void
00380 MonitorElement::Fill(double x, double y, double zw)
00381 {
00382   update();
00383   if (kind() == DQM_KIND_TH2F)
00384     static_cast<TH2F *>(accessRootObject(__PRETTY_FUNCTION__, 2))
00385       ->Fill(x, y, zw);
00386   else if (kind() == DQM_KIND_TH2S)
00387     static_cast<TH2S *>(accessRootObject(__PRETTY_FUNCTION__, 2))
00388       ->Fill(x, y, zw);
00389   else if (kind() == DQM_KIND_TH2D)
00390     static_cast<TH2D *>(accessRootObject(__PRETTY_FUNCTION__, 2))
00391       ->Fill(x, y, zw);
00392   else if (kind() == DQM_KIND_TH3F)
00393     static_cast<TH3F *>(accessRootObject(__PRETTY_FUNCTION__, 2))
00394       ->Fill(x, y, zw, 1);
00395   else if (kind() == DQM_KIND_TPROFILE)
00396     static_cast<TProfile *>(accessRootObject(__PRETTY_FUNCTION__, 2))
00397       ->Fill(x, y, zw);
00398   else if (kind() == DQM_KIND_TPROFILE2D)
00399     static_cast<TProfile2D *>(accessRootObject(__PRETTY_FUNCTION__, 2))
00400       ->Fill(x, y, zw, 1);
00401   else
00402     incompatible(__PRETTY_FUNCTION__);
00403 }
00404 
00406 void
00407 MonitorElement::Fill(double x, double y, double z, double w)
00408 {
00409   update();
00410   if (kind() == DQM_KIND_TH3F)
00411     static_cast<TH3F *>(accessRootObject(__PRETTY_FUNCTION__, 2))
00412       ->Fill(x, y, z, w);
00413   else if (kind() == DQM_KIND_TPROFILE2D)
00414     static_cast<TProfile2D *>(accessRootObject(__PRETTY_FUNCTION__, 2))
00415       ->Fill(x, y, z, w);
00416   else
00417     incompatible(__PRETTY_FUNCTION__);
00418 }
00419 
00421 void
00422 MonitorElement::Reset(void)
00423 {
00424   update();
00425   if (kind() == DQM_KIND_INT)
00426     scalar_.num = 0;
00427   else if (kind() == DQM_KIND_REAL)
00428     scalar_.real = 0;
00429   else if (kind() == DQM_KIND_STRING)
00430     scalar_.str.clear();
00431   else
00432     return accessRootObject(__PRETTY_FUNCTION__, 1)
00433       ->Reset();
00434 }
00435 
00437 void
00438 MonitorElement::packScalarData(std::string &into, const char *prefix) const
00439 {
00440   char buf[64];
00441   if (kind() == DQM_KIND_INT)
00442   {
00443     snprintf(buf, sizeof(buf), "%s%" PRId64, prefix, scalar_.num);
00444     into = buf;
00445   }
00446   else if (kind() == DQM_KIND_REAL)
00447   {
00448     snprintf(buf, sizeof(buf), "%s%.*g", prefix, DBL_DIG+2, scalar_.real);
00449     into = buf;
00450   }
00451   else if (kind() == DQM_KIND_STRING)
00452   {
00453     into.reserve(strlen(prefix) + scalar_.str.size());
00454     into += prefix;
00455     into += scalar_.str;
00456   }
00457   else
00458     incompatible(__PRETTY_FUNCTION__);
00459 }
00460 
00462 void
00463 MonitorElement::packQualityData(std::string &into) const
00464 {
00465   DQMNet::packQualityData(into, data_.qreports);
00466 }
00467 
00470 std::string
00471 MonitorElement::valueString(void) const
00472 {
00473   std::string result;
00474   if (kind() == DQM_KIND_INT)
00475     packScalarData(result, "i=");
00476   else if (kind() == DQM_KIND_REAL)
00477     packScalarData(result, "f=");
00478   else if (kind() == DQM_KIND_STRING)
00479     packScalarData(result, "s=");
00480   else
00481     incompatible(__PRETTY_FUNCTION__);
00482 
00483   return result;
00484 }
00485 
00489 std::string
00490 MonitorElement::tagString(void) const
00491 {
00492   std::string result;
00493   std::string val(valueString());
00494   result.reserve(6 + 2*data_.objname.size() + val.size());
00495   result += '<'; result += data_.objname; result += '>';
00496   result += val;
00497   result += '<'; result += '/'; result += data_.objname; result += '>';
00498   return result;
00499 }
00500 
00502 std::string
00503 MonitorElement::tagLabelString(void) const
00504 {
00505   char buf[32];
00506   std::string result;
00507   size_t len = sprintf(buf, "t=%" PRIu32, data_.tag);
00508 
00509   result.reserve(6 + 2*data_.objname.size() + len);
00510   result += '<'; result += data_.objname; result += '>';
00511   result += buf;
00512   result += '<'; result += '/'; result += data_.objname; result += '>';
00513   return result;
00514 }
00515 
00516 std::string
00517 MonitorElement::qualityTagString(const DQMNet::QValue &qv) const
00518 {
00519   char buf[64];
00520   std::string result;
00521   size_t titlelen = data_.objname.size() + qv.qtname.size() + 1;
00522   size_t buflen = sprintf(buf, "qr=st:%d:%.*g:", qv.code, DBL_DIG+2, qv.qtresult);
00523 
00524   result.reserve(7 + 2*titlelen + buflen + qv.algorithm.size() + qv.message.size());
00525   result += '<'; result += data_.objname; result += '.'; result += qv.qtname; result += '>';
00526   result += buf; result += qv.algorithm; result += ':'; result += qv.message;
00527   result += '<'; result += '/'; result += data_.objname; result += '.'; result += qv.qtname; result += '>';
00528   return result;
00529 }
00530 
00531 const QReport *
00532 MonitorElement::getQReport(const std::string &qtname) const
00533 {
00534   QReport *qr;
00535   DQMNet::QValue *qv;
00536   const_cast<MonitorElement *>(this)->getQReport(false, qtname, qr, qv);
00537   return qr;
00538 }
00539 
00540 std::vector<QReport *>
00541 MonitorElement::getQReports(void) const
00542 {
00543   std::vector<QReport *> result;
00544   result.reserve(qreports_.size());
00545   for (size_t i = 0, e = qreports_.size(); i != e; ++i)
00546   {
00547     const_cast<MonitorElement *>(this)->qreports_[i].qvalue_
00548       = const_cast<DQMNet::QValue *>(&data_.qreports[i]);
00549     result.push_back(const_cast<QReport *>(&qreports_[i]));
00550   }
00551   return result;
00552 }
00553 
00554 std::vector<QReport *>
00555 MonitorElement::getQWarnings(void) const
00556 {
00557   std::vector<QReport *> result;
00558   result.reserve(qreports_.size());
00559   for (size_t i = 0, e = qreports_.size(); i != e; ++i)
00560     if (data_.qreports[i].code == dqm::qstatus::WARNING)
00561     {
00562       const_cast<MonitorElement *>(this)->qreports_[i].qvalue_
00563         = const_cast<DQMNet::QValue *>(&data_.qreports[i]);
00564       result.push_back(const_cast<QReport *>(&qreports_[i]));
00565     }
00566   return result;
00567 }
00568 
00569 std::vector<QReport *>
00570 MonitorElement::getQErrors(void) const
00571 {
00572   std::vector<QReport *> result;
00573   result.reserve(qreports_.size());
00574   for (size_t i = 0, e = qreports_.size(); i != e; ++i)
00575     if (data_.qreports[i].code == dqm::qstatus::ERROR)
00576     {
00577       const_cast<MonitorElement *>(this)->qreports_[i].qvalue_
00578         = const_cast<DQMNet::QValue *>(&data_.qreports[i]);
00579       result.push_back(const_cast<QReport *>(&qreports_[i]));
00580     }
00581   return result;
00582 }
00583 
00584 std::vector<QReport *>
00585 MonitorElement::getQOthers(void) const
00586 {
00587   std::vector<QReport *> result;
00588   result.reserve(qreports_.size());
00589   for (size_t i = 0, e = qreports_.size(); i != e; ++i)
00590     if (data_.qreports[i].code != dqm::qstatus::STATUS_OK
00591         && data_.qreports[i].code != dqm::qstatus::WARNING
00592         && data_.qreports[i].code != dqm::qstatus::ERROR)
00593     {
00594       const_cast<MonitorElement *>(this)->qreports_[i].qvalue_
00595         = const_cast<DQMNet::QValue *>(&data_.qreports[i]);
00596       result.push_back(const_cast<QReport *>(&qreports_[i]));
00597     }
00598   return result;
00599 }
00600 
00602 void
00603 MonitorElement::runQTests(void)
00604 {
00605   assert(qreports_.size() == data_.qreports.size());
00606 
00607   // Rerun quality tests where the ME or the quality algorithm was modified.
00608   bool dirty = wasUpdated();
00609   for (size_t i = 0, e = data_.qreports.size(); i < e; ++i)
00610   {
00611     DQMNet::QValue &qv = data_.qreports[i];
00612     QReport &qr = qreports_[i];
00613     QCriterion *qc = qr.qcriterion_;
00614     qr.qvalue_ = &qv;
00615 
00616     // if (qc && (dirty || qc->wasModified()))  // removed for new QTest (abm-090503)
00617     if (qc && dirty)
00618     {
00619       assert(qc->getName() == qv.qtname);
00620       std::string oldMessage = qv.message;
00621       int oldStatus = qv.code;
00622 
00623       qc->runTest(this, qr, qv);
00624 
00625       if (oldStatus != qv.code || oldMessage != qv.message)
00626         update();
00627     }
00628   }
00629 
00630   // Update QReport statistics.
00631   updateQReportStats();
00632 }
00633 
00634 void
00635 MonitorElement::incompatible(const char *func) const
00636 {
00637   raiseDQMError("MonitorElement", "Method '%s' cannot be invoked on monitor"
00638                 " element '%s'", func, data_.objname.c_str());
00639 }
00640 
00641 TH1 *
00642 MonitorElement::accessRootObject(const char *func, int reqdim) const
00643 {
00644   if (kind() < DQM_KIND_TH1F)
00645     raiseDQMError("MonitorElement", "Method '%s' cannot be invoked on monitor"
00646                   " element '%s' because it is not a root object",
00647                   func, data_.objname.c_str());
00648 
00649   return checkRootObject(data_.objname, object_, func, reqdim);
00650 }
00651 
00652 /*** getter methods (wrapper around ROOT methods) ****/
00653 // 
00655 double
00656 MonitorElement::getMean(int axis /* = 1 */) const
00657 { return accessRootObject(__PRETTY_FUNCTION__, axis-1)
00658     ->GetMean(axis); }
00659 
00662 double
00663 MonitorElement::getMeanError(int axis /* = 1 */) const
00664 { return accessRootObject(__PRETTY_FUNCTION__, axis-1)
00665     ->GetMeanError(axis); }
00666 
00668 double
00669 MonitorElement::getRMS(int axis /* = 1 */) const
00670 { return accessRootObject(__PRETTY_FUNCTION__, axis-1)
00671     ->GetRMS(axis); }
00672 
00674 double
00675 MonitorElement::getRMSError(int axis /* = 1 */) const
00676 { return accessRootObject(__PRETTY_FUNCTION__, axis-1)
00677     ->GetRMSError(axis); }
00678 
00680 int
00681 MonitorElement::getNbinsX(void) const
00682 { return accessRootObject(__PRETTY_FUNCTION__, 1)
00683     ->GetNbinsX(); }
00684 
00686 int
00687 MonitorElement::getNbinsY(void) const
00688 { return accessRootObject(__PRETTY_FUNCTION__, 2)
00689     ->GetNbinsY(); }
00690 
00692 int
00693 MonitorElement::getNbinsZ(void) const
00694 { return accessRootObject(__PRETTY_FUNCTION__, 3)
00695     ->GetNbinsZ(); }
00696 
00698 double
00699 MonitorElement::getBinContent(int binx) const
00700 { return accessRootObject(__PRETTY_FUNCTION__, 1)
00701     ->GetBinContent(binx); }
00702 
00704 double
00705 MonitorElement::getBinContent(int binx, int biny) const
00706 { return accessRootObject(__PRETTY_FUNCTION__, 2)
00707     ->GetBinContent(binx, biny); }
00708 
00710 double
00711 MonitorElement::getBinContent(int binx, int biny, int binz) const
00712 { return accessRootObject(__PRETTY_FUNCTION__, 3)
00713     ->GetBinContent(binx, biny, binz); }
00714 
00716 double
00717 MonitorElement::getBinError(int binx) const
00718 { return accessRootObject(__PRETTY_FUNCTION__, 1)
00719     ->GetBinError(binx); }
00720 
00722 double
00723 MonitorElement::getBinError(int binx, int biny) const
00724 { return accessRootObject(__PRETTY_FUNCTION__, 2)
00725     ->GetBinError(binx, biny); }
00726 
00728 double
00729 MonitorElement::getBinError(int binx, int biny, int binz) const
00730 { return accessRootObject(__PRETTY_FUNCTION__, 3)
00731     ->GetBinError(binx, biny, binz); }
00732 
00734 double
00735 MonitorElement::getEntries(void) const
00736 { return accessRootObject(__PRETTY_FUNCTION__, 1)
00737     ->GetEntries(); }
00738 
00740 double
00741 MonitorElement::getBinEntries(int bin) const
00742 {
00743   if (kind() == DQM_KIND_TPROFILE)
00744     return static_cast<TProfile *>(accessRootObject(__PRETTY_FUNCTION__, 1))
00745       ->GetBinEntries(bin);
00746   else if (kind() == DQM_KIND_TPROFILE2D)
00747     return static_cast<TProfile2D *>(accessRootObject(__PRETTY_FUNCTION__, 1))
00748       ->GetBinEntries(bin);
00749   else
00750   {
00751     incompatible(__PRETTY_FUNCTION__);
00752     return 0;
00753   }
00754 }
00755 
00757 double
00758 MonitorElement::getYmin(void) const
00759 {
00760   if (kind() == DQM_KIND_TPROFILE)
00761     return static_cast<TProfile *>(accessRootObject(__PRETTY_FUNCTION__, 1))
00762       ->GetYmin();
00763   else
00764   {
00765     incompatible(__PRETTY_FUNCTION__);
00766     return 0;
00767   }
00768 }
00769 
00771 double
00772 MonitorElement::getYmax(void) const
00773 {
00774   if (kind() == DQM_KIND_TPROFILE)
00775     return static_cast<TProfile *>(accessRootObject(__PRETTY_FUNCTION__, 1))
00776       ->GetYmax();
00777   else
00778   {
00779     incompatible(__PRETTY_FUNCTION__);
00780     return 0;
00781   }
00782 }
00783 
00785 std::string
00786 MonitorElement::getAxisTitle(int axis /* = 1 */) const
00787 { return getAxis(__PRETTY_FUNCTION__, axis)
00788     ->GetTitle(); }
00789 
00791 std::string
00792 MonitorElement::getTitle(void) const
00793 { return accessRootObject(__PRETTY_FUNCTION__, 1)
00794     ->GetTitle(); }
00795 
00796 /*** setter methods (wrapper around ROOT methods) ****/
00797 // 
00799 void
00800 MonitorElement::setBinContent(int binx, double content)
00801 {
00802   update();
00803   accessRootObject(__PRETTY_FUNCTION__, 1)
00804     ->SetBinContent(binx, content);
00805 }
00806 
00808 void
00809 MonitorElement::setBinContent(int binx, int biny, double content)
00810 {
00811   update();
00812   accessRootObject(__PRETTY_FUNCTION__, 2)
00813     ->SetBinContent(binx, biny, content); }
00814 
00816 void
00817 MonitorElement::setBinContent(int binx, int biny, int binz, double content)
00818 {
00819   update();
00820   accessRootObject(__PRETTY_FUNCTION__, 3)
00821     ->SetBinContent(binx, biny, binz, content); }
00822 
00824 void
00825 MonitorElement::setBinError(int binx, double error)
00826 {
00827   update();
00828   accessRootObject(__PRETTY_FUNCTION__, 1)
00829     ->SetBinError(binx, error);
00830 }
00831 
00833 void
00834 MonitorElement::setBinError(int binx, int biny, double error)
00835 {
00836   update();
00837   accessRootObject(__PRETTY_FUNCTION__, 2)
00838     ->SetBinError(binx, biny, error);
00839 }
00840 
00842 void
00843 MonitorElement::setBinError(int binx, int biny, int binz, double error)
00844 {
00845   update();
00846   accessRootObject(__PRETTY_FUNCTION__, 3)
00847     ->SetBinError(binx, biny, binz, error);
00848 }
00849 
00851 void
00852 MonitorElement::setBinEntries(int bin, double nentries)
00853 {
00854   update();
00855   if (kind() == DQM_KIND_TPROFILE)
00856     static_cast<TProfile *>(accessRootObject(__PRETTY_FUNCTION__, 1))
00857       ->SetBinEntries(bin, nentries);
00858   else if (kind() == DQM_KIND_TPROFILE2D)
00859     static_cast<TProfile2D *>(accessRootObject(__PRETTY_FUNCTION__, 1))
00860       ->SetBinEntries(bin, nentries);
00861   else
00862     incompatible(__PRETTY_FUNCTION__);
00863 }
00864 
00866 void
00867 MonitorElement::setEntries(double nentries)
00868 {
00869   update();
00870   accessRootObject(__PRETTY_FUNCTION__, 1)
00871     ->SetEntries(nentries);
00872 }
00873 
00875 void
00876 MonitorElement::setBinLabel(int bin, const std::string &label, int axis /* = 1 */)
00877 {
00878   update();
00879   if ( getAxis(__PRETTY_FUNCTION__, axis)->GetNbins() >= bin ) 
00880   {
00881     getAxis(__PRETTY_FUNCTION__, axis)
00882       ->SetBinLabel(bin, label.c_str());
00883   }
00884   else
00885   {
00886     //  edm::LogWarning ("MonitorElement") 
00887     std::cout << "*** MonitorElement: WARNING:"
00888               <<"setBinLabel: attempting to set label of non-existent bin number \n";
00889   }
00890 }
00891 
00893 void
00894 MonitorElement::setAxisRange(double xmin, double xmax, int axis /* = 1 */)
00895 {
00896   update();
00897   getAxis(__PRETTY_FUNCTION__, axis)
00898     ->SetRangeUser(xmin, xmax);
00899 }
00900 
00902 void
00903 MonitorElement::setAxisTitle(const std::string &title, int axis /* = 1 */)
00904 {
00905   update();
00906   getAxis(__PRETTY_FUNCTION__, axis)
00907     ->SetTitle(title.c_str());
00908 }
00909 
00911 void
00912 MonitorElement::setAxisTimeDisplay(int value, int axis /* = 1 */)
00913 {
00914   update();
00915   getAxis(__PRETTY_FUNCTION__, axis)
00916     ->SetTimeDisplay(value);
00917 }
00918 
00920 void
00921 MonitorElement::setAxisTimeFormat(const char *format /* = "" */, int axis /* = 1 */)
00922 {
00923   update();
00924   getAxis(__PRETTY_FUNCTION__, axis)
00925     ->SetTimeFormat(format);
00926 }
00927 
00929 void
00930 MonitorElement::setAxisTimeOffset(double toffset, const char *option /* ="local" */, int axis /* = 1 */)
00931 {
00932   update();
00933   getAxis(__PRETTY_FUNCTION__, axis)
00934     ->SetTimeOffset(toffset, option);
00935 }
00936 
00938 void
00939 MonitorElement::setTitle(const std::string &title)
00940 {
00941   update();
00942   accessRootObject(__PRETTY_FUNCTION__, 1)
00943     ->SetTitle(title.c_str());
00944 }
00945 
00946 TAxis *
00947 MonitorElement::getAxis(const char *func, int axis) const
00948 {
00949   TH1 *h = accessRootObject(func, axis-1);
00950   TAxis *a = 0;
00951   if (axis == 1)
00952     a = h->GetXaxis();
00953   else if (axis == 2)
00954     a = h->GetYaxis();
00955   else if (axis == 3)
00956     a = h->GetZaxis();
00957 
00958   if (! a)
00959     raiseDQMError("MonitorElement", "No such axis %d in monitor element"
00960                   " '%s' of type '%s'", axis, data_.objname.c_str(),
00961                   typeid(*h).name());
00962 
00963   return a;
00964 }
00965 
00966 // ------------ Operations for MEs that are normally never reset ---------
00967 
00970 void
00971 MonitorElement::softReset(void)
00972 {
00973   update();
00974 
00975   // Create the reference object the first time this is called.
00976   // On subsequent calls accumulate the current value to the
00977   // reference, and then reset the current value.  This way the
00978   // future contents will have the reference "subtracted".
00979   if (kind() == DQM_KIND_TH1F)
00980   {
00981     TH1F *orig = static_cast<TH1F *>(object_);
00982     TH1F *r = static_cast<TH1F *>(refvalue_);
00983     if (! r)
00984     {
00985       refvalue_ = r = new TH1F((std::string(orig->GetName()) + "_ref").c_str(),
00986                                orig->GetTitle(),
00987                                orig->GetNbinsX(),
00988                                orig->GetXaxis()->GetXmin(),
00989                                orig->GetXaxis()->GetXmax());
00990       r->SetDirectory(0);
00991       r->Reset();
00992     }
00993 
00994     r->Add(orig);
00995     orig->Reset();
00996   }
00997   else if (kind() == DQM_KIND_TH1S)
00998   {
00999     TH1S *orig = static_cast<TH1S *>(object_);
01000     TH1S *r = static_cast<TH1S *>(refvalue_);
01001     if (! r)
01002     {
01003       refvalue_ = r = new TH1S((std::string(orig->GetName()) + "_ref").c_str(),
01004                                orig->GetTitle(),
01005                                orig->GetNbinsX(),
01006                                orig->GetXaxis()->GetXmin(),
01007                                orig->GetXaxis()->GetXmax());
01008       r->SetDirectory(0);
01009       r->Reset();
01010     }
01011 
01012     r->Add(orig);
01013     orig->Reset();
01014   }
01015   else if (kind() == DQM_KIND_TH1D)
01016   {
01017     TH1D *orig = static_cast<TH1D *>(object_);
01018     TH1D *r = static_cast<TH1D *>(refvalue_);
01019     if (! r)
01020     {
01021       refvalue_ = r = new TH1D((std::string(orig->GetName()) + "_ref").c_str(),
01022                                orig->GetTitle(),
01023                                orig->GetNbinsX(),
01024                                orig->GetXaxis()->GetXmin(),
01025                                orig->GetXaxis()->GetXmax());
01026       r->SetDirectory(0);
01027       r->Reset();
01028     }
01029 
01030     r->Add(orig);
01031     orig->Reset();
01032   }
01033   else if (kind() == DQM_KIND_TH2F)
01034   {
01035     TH2F *orig = static_cast<TH2F *>(object_);
01036     TH2F *r = static_cast<TH2F *>(refvalue_);
01037     if (! r)
01038     {
01039       refvalue_ = r = new TH2F((std::string(orig->GetName()) + "_ref").c_str(),
01040                                orig->GetTitle(),
01041                                orig->GetNbinsX(),
01042                                orig->GetXaxis()->GetXmin(),
01043                                orig->GetXaxis()->GetXmax(),
01044                                orig->GetNbinsY(),
01045                                orig->GetYaxis()->GetXmin(),
01046                                orig->GetYaxis()->GetXmax());
01047       r->SetDirectory(0);
01048       r->Reset();
01049     }
01050 
01051     r->Add(orig);
01052     orig->Reset();
01053   }
01054   else if (kind() == DQM_KIND_TH2S)
01055   {
01056     TH2S *orig = static_cast<TH2S *>(object_);
01057     TH2S *r = static_cast<TH2S *>(refvalue_);
01058     if (! r)
01059     {
01060       refvalue_ = r = new TH2S((std::string(orig->GetName()) + "_ref").c_str(),
01061                                orig->GetTitle(),
01062                                orig->GetNbinsX(),
01063                                orig->GetXaxis()->GetXmin(),
01064                                orig->GetXaxis()->GetXmax(),
01065                                orig->GetNbinsY(),
01066                                orig->GetYaxis()->GetXmin(),
01067                                orig->GetYaxis()->GetXmax());
01068       r->SetDirectory(0);
01069       r->Reset();
01070     }
01071 
01072     r->Add(orig);
01073     orig->Reset();
01074   }
01075   else if (kind() == DQM_KIND_TH2D)
01076   {
01077     TH2D *orig = static_cast<TH2D *>(object_);
01078     TH2D *r = static_cast<TH2D *>(refvalue_);
01079     if (! r)
01080     {
01081       refvalue_ = r = new TH2D((std::string(orig->GetName()) + "_ref").c_str(),
01082                                orig->GetTitle(),
01083                                orig->GetNbinsX(),
01084                                orig->GetXaxis()->GetXmin(),
01085                                orig->GetXaxis()->GetXmax(),
01086                                orig->GetNbinsY(),
01087                                orig->GetYaxis()->GetXmin(),
01088                                orig->GetYaxis()->GetXmax());
01089       r->SetDirectory(0);
01090       r->Reset();
01091     }
01092 
01093     r->Add(orig);
01094     orig->Reset();
01095   }
01096   else if (kind() == DQM_KIND_TH3F)
01097   {
01098     TH3F *orig = static_cast<TH3F *>(object_);
01099     TH3F *r = static_cast<TH3F *>(refvalue_);
01100     if (! r)
01101     {
01102       refvalue_ = r = new TH3F((std::string(orig->GetName()) + "_ref").c_str(),
01103                                orig->GetTitle(),
01104                                orig->GetNbinsX(),
01105                                orig->GetXaxis()->GetXmin(),
01106                                orig->GetXaxis()->GetXmax(),
01107                                orig->GetNbinsY(),
01108                                orig->GetYaxis()->GetXmin(),
01109                                orig->GetYaxis()->GetXmax(),
01110                                orig->GetNbinsZ(),
01111                                orig->GetZaxis()->GetXmin(),
01112                                orig->GetZaxis()->GetXmax());
01113       r->SetDirectory(0);
01114       r->Reset();
01115     }
01116 
01117     r->Add(orig);
01118     orig->Reset();
01119   }
01120   else if (kind() == DQM_KIND_TPROFILE)
01121   {
01122     TProfile *orig = static_cast<TProfile *>(object_);
01123     TProfile *r = static_cast<TProfile *>(refvalue_);
01124     if (! r)
01125     {
01126       refvalue_ = r = new TProfile((std::string(orig->GetName()) + "_ref").c_str(),
01127                                    orig->GetTitle(),
01128                                    orig->GetNbinsX(),
01129                                    orig->GetXaxis()->GetXmin(),
01130                                    orig->GetXaxis()->GetXmax(),
01131                                    orig->GetYaxis()->GetXmin(),
01132                                    orig->GetYaxis()->GetXmax(),
01133                                    orig->GetErrorOption());
01134       r->SetDirectory(0);
01135       r->Reset();
01136     }
01137 
01138     addProfiles(r, orig, r, 1, 1);
01139     orig->Reset();
01140   }
01141   else if (kind() == DQM_KIND_TPROFILE2D)
01142   {
01143     TProfile2D *orig = static_cast<TProfile2D *>(object_);
01144     TProfile2D *r = static_cast<TProfile2D *>(refvalue_);
01145     if (! r)
01146     {
01147       refvalue_ = r = new TProfile2D((std::string(orig->GetName()) + "_ref").c_str(),
01148                                      orig->GetTitle(),
01149                                      orig->GetNbinsX(),
01150                                      orig->GetXaxis()->GetXmin(),
01151                                      orig->GetXaxis()->GetXmax(),
01152                                      orig->GetNbinsY(),
01153                                      orig->GetYaxis()->GetXmin(),
01154                                      orig->GetYaxis()->GetXmax(),
01155                                      orig->GetZaxis()->GetXmin(),
01156                                      orig->GetZaxis()->GetXmax(),
01157                                      orig->GetErrorOption());
01158       r->SetDirectory(0);
01159       r->Reset();
01160     }
01161 
01162     addProfiles(r, orig, r, 1, 1);
01163     orig->Reset();
01164   }
01165   else
01166     incompatible(__PRETTY_FUNCTION__);
01167 }
01168 
01170 void
01171 MonitorElement::disableSoftReset(void)
01172 {
01173   if (refvalue_)
01174   {
01175     if (kind() == DQM_KIND_TH1F
01176         || kind() == DQM_KIND_TH1S
01177         || kind() == DQM_KIND_TH1D
01178         || kind() == DQM_KIND_TH2F
01179         || kind() == DQM_KIND_TH2S
01180         || kind() == DQM_KIND_TH2D
01181         || kind() == DQM_KIND_TH3F)
01182     {
01183       TH1 *orig = static_cast<TH1 *>(object_);
01184       orig->Add(refvalue_);
01185     }
01186     else if (kind() == DQM_KIND_TPROFILE)
01187     {
01188       TProfile *orig = static_cast<TProfile *>(object_);
01189       TProfile *r = static_cast<TProfile *>(refvalue_);
01190       addProfiles(orig, r, orig, 1, 1);
01191     }
01192     else if (kind() == DQM_KIND_TPROFILE2D)
01193     {
01194       TProfile2D *orig = static_cast<TProfile2D *>(object_);
01195       TProfile2D *r = static_cast<TProfile2D *>(refvalue_);
01196       addProfiles(orig, r, orig, 1, 1);
01197     }
01198     else
01199       incompatible(__PRETTY_FUNCTION__);
01200 
01201     delete refvalue_;
01202     refvalue_ = 0;
01203   }
01204 }
01205   
01206 // implementation: Giuseppe.Della-Ricca@ts.infn.it
01207 // Can be called with sum = h1 or sum = h2
01208 void
01209 MonitorElement::addProfiles(TProfile *h1, TProfile *h2, TProfile *sum, float c1, float c2)
01210 {
01211   assert(h1);
01212   assert(h2);
01213   assert(sum);
01214 
01215   static const Int_t NUM_STAT = 6;
01216   Double_t stats1[NUM_STAT];
01217   Double_t stats2[NUM_STAT];
01218   Double_t stats3[NUM_STAT];
01219 
01220   bool isRebinOn = sum->TestBit(TH1::kCanRebin);
01221   sum->ResetBit(TH1::kCanRebin);
01222 
01223   for (Int_t i = 0; i < NUM_STAT; ++i)
01224     stats1[i] = stats2[i] = stats3[i] = 0;
01225 
01226   h1->GetStats(stats1);
01227   h2->GetStats(stats2);
01228 
01229   for (Int_t i = 0; i < NUM_STAT; ++i)
01230     stats3[i] = c1*stats1[i] + c2*stats2[i];
01231 
01232   stats3[1] = c1*TMath::Abs(c1)*stats1[1]
01233               + c2*TMath::Abs(c2)*stats2[1];
01234 
01235   Double_t entries = c1*h1->GetEntries() + c2* h2->GetEntries();
01236   TArrayD* h1sumw2 = h1->GetSumw2();
01237   TArrayD* h2sumw2 = h2->GetSumw2();
01238   for (Int_t bin = 0, nbin = sum->GetNbinsX()+1; bin <= nbin; ++bin) 
01239   {
01240     Double_t entries = c1*h1->GetBinEntries(bin)
01241                        + c2*h2->GetBinEntries(bin);
01242     Double_t content = c1*h1->GetBinEntries(bin)*h1->GetBinContent(bin)
01243                        + c2*h2->GetBinEntries(bin)*h2->GetBinContent(bin);
01244     Double_t error = TMath::Sqrt(c1*TMath::Abs(c1)*h1sumw2->fArray[bin]
01245                                  + c2*TMath::Abs(c2)*h2sumw2->fArray[bin]);
01246     sum->SetBinContent(bin, content);
01247     sum->SetBinError(bin, error);
01248     sum->SetBinEntries(bin, entries);
01249   }
01250   
01251   sum->SetEntries(entries);
01252   sum->PutStats(stats3);
01253   if (isRebinOn) sum->SetBit(TH1::kCanRebin);
01254 }
01255 
01256 // implementation: Giuseppe.Della-Ricca@ts.infn.it
01257 // Can be called with sum = h1 or sum = h2
01258 void
01259 MonitorElement::addProfiles(TProfile2D *h1, TProfile2D *h2, TProfile2D *sum, float c1, float c2)
01260 {
01261   assert(h1);
01262   assert(h2);
01263   assert(sum);
01264 
01265   static const Int_t NUM_STAT = 9;
01266   Double_t stats1[NUM_STAT];
01267   Double_t stats2[NUM_STAT];
01268   Double_t stats3[NUM_STAT];
01269 
01270   bool isRebinOn = sum->TestBit(TH1::kCanRebin);
01271   sum->ResetBit(TH1::kCanRebin);
01272 
01273   for (Int_t i = 0; i < NUM_STAT; ++i)
01274     stats1[i] = stats2[i] = stats3[i] = 0;
01275 
01276   h1->GetStats(stats1);
01277   h2->GetStats(stats2);
01278 
01279   for (Int_t i = 0; i < NUM_STAT; i++)
01280     stats3[i] = c1*stats1[i] + c2*stats2[i];
01281 
01282   stats3[1] = c1*TMath::Abs(c1)*stats1[1]
01283               + c2*TMath::Abs(c2)*stats2[1];
01284 
01285   Double_t entries = c1*h1->GetEntries() + c2*h2->GetEntries();
01286   TArrayD *h1sumw2 = h1->GetSumw2();
01287   TArrayD *h2sumw2 = h2->GetSumw2();
01288   for (Int_t xbin = 0, nxbin = sum->GetNbinsX()+1; xbin <= nxbin; ++xbin)
01289     for (Int_t ybin = 0, nybin = sum->GetNbinsY()+1; ybin <= nybin; ++ybin)
01290     {
01291       Int_t bin = sum->GetBin(xbin, ybin);
01292       Double_t entries = c1*h1->GetBinEntries(bin)
01293                          + c2*h2->GetBinEntries(bin);
01294       Double_t content = c1*h1->GetBinEntries(bin)*h1->GetBinContent(bin)
01295                          + c2*h2->GetBinEntries(bin)*h2->GetBinContent(bin);
01296       Double_t error = TMath::Sqrt(c1*TMath::Abs(c1)*h1sumw2->fArray[bin]
01297                                    + c2*TMath::Abs(c2)*h2sumw2->fArray[bin]);
01298 
01299       sum->SetBinContent(bin, content);
01300       sum->SetBinError(bin, error);
01301       sum->SetBinEntries(bin, entries);
01302     }
01303   sum->SetEntries(entries);
01304   sum->PutStats(stats3);
01305   if (isRebinOn) sum->SetBit(TH1::kCanRebin);
01306 }
01307 
01308 void
01309 MonitorElement::copyFunctions(TH1 *from, TH1 *to)
01310 {
01311   // will copy functions only if local-copy and original-object are equal
01312   // (ie. no soft-resetting or accumulating is enabled)
01313   if (isSoftResetEnabled() || isAccumulateEnabled())
01314     return;
01315 
01316   update();
01317   TList *fromf = from->GetListOfFunctions();
01318   TList *tof   = to->GetListOfFunctions();
01319   for (int i = 0, nfuncs = fromf ? fromf->GetSize() : 0; i < nfuncs; ++i)
01320   {
01321     TObject *obj = fromf->At(i);
01322     // not interested in statistics
01323     if (!strcmp(obj->IsA()->GetName(), "TPaveStats"))
01324       continue;
01325 
01326     if (TF1 *fn = dynamic_cast<TF1 *>(obj))
01327       tof->Add(new TF1(*fn));
01328     //else if (dynamic_cast<TPaveStats *>(obj))
01329     //  ; // FIXME? tof->Add(new TPaveStats(*stats));
01330     else
01331       raiseDQMError("MonitorElement", "Cannot extract function '%s' of type"
01332                     " '%s' from monitor element '%s' for a copy",
01333                     obj->GetName(), obj->IsA()->GetName(), data_.objname.c_str());
01334   }
01335 }
01336 
01337 void
01338 MonitorElement::copyFrom(TH1 *from)
01339 {
01340   TH1 *orig = accessRootObject(__PRETTY_FUNCTION__, 1);
01341   if (orig->GetTitle() != from->GetTitle())
01342     orig->SetTitle(from->GetTitle());
01343 
01344   if (!isAccumulateEnabled())
01345     orig->Reset();
01346 
01347   if (isSoftResetEnabled())
01348   {
01349     if (kind() == DQM_KIND_TH1F
01350         || kind() == DQM_KIND_TH1S
01351         || kind() == DQM_KIND_TH1D
01352         || kind() == DQM_KIND_TH2F
01353         || kind() == DQM_KIND_TH2S
01354         || kind() == DQM_KIND_TH2D
01355         || kind() == DQM_KIND_TH3F)
01356       // subtract "reference"
01357       orig->Add(from, refvalue_, 1, -1);
01358     else if (kind() == DQM_KIND_TPROFILE)
01359       // subtract "reference"
01360       addProfiles(static_cast<TProfile *>(from),
01361                   static_cast<TProfile *>(refvalue_),
01362                   static_cast<TProfile *>(orig),
01363                   1, -1);
01364     else if (kind() == DQM_KIND_TPROFILE2D)
01365       // subtract "reference"
01366       addProfiles(static_cast<TProfile2D *>(from),
01367                   static_cast<TProfile2D *>(refvalue_),
01368                   static_cast<TProfile2D *>(orig),
01369                   1, -1);
01370     else
01371       incompatible(__PRETTY_FUNCTION__);
01372   }
01373   else
01374     orig->Add(from);
01375 
01376   copyFunctions(from, orig);
01377 }
01378 
01379 // --- Operations on MEs that are normally reset at end of monitoring cycle ---
01380 void
01381 MonitorElement::getQReport(bool create, const std::string &qtname, QReport *&qr, DQMNet::QValue *&qv)
01382 {
01383   assert(qreports_.size() == data_.qreports.size());
01384 
01385   qr = 0;
01386   qv = 0;
01387 
01388   size_t pos = 0, end = qreports_.size();
01389   while (pos < end && data_.qreports[pos].qtname != qtname)
01390     ++pos;
01391 
01392   if (pos == end && ! create)
01393     return;
01394   else if (pos == end)
01395   {
01396     data_.qreports.push_back(DQMNet::QValue());
01397     qreports_.push_back(QReport(0, 0));
01398 
01399     DQMNet::QValue &q = data_.qreports.back();
01400     q.code = dqm::qstatus::DID_NOT_RUN;
01401     q.qtresult = 0;
01402     q.qtname = qtname;
01403     q.message = "NO_MESSAGE_ASSIGNED";
01404     q.algorithm = "UNKNOWN_ALGORITHM";
01405   }
01406 
01407   qr = &qreports_[pos];
01408   qv = &data_.qreports[pos];
01409 }
01410   
01412 void
01413 MonitorElement::addQReport(const DQMNet::QValue &desc, QCriterion *qc)
01414 {
01415   QReport *qr;
01416   DQMNet::QValue *qv;
01417   getQReport(true, desc.qtname, qr, qv);
01418   qr->qcriterion_ = qc;
01419   *qv = desc;
01420   update();
01421 }
01422 
01423 void
01424 MonitorElement::addQReport(QCriterion *qc)
01425 {
01426   QReport *qr;
01427   DQMNet::QValue *qv;
01428   getQReport(true, qc->getName(), qr, qv);
01429   qv->code = dqm::qstatus::DID_NOT_RUN;
01430   qv->message = "NO_MESSAGE_ASSIGNED";
01431   qr->qcriterion_ = qc;
01432   update();
01433 }
01434 
01436 void
01437 MonitorElement::updateQReportStats(void)
01438 {
01439   data_.flags &= ~DQMNet::DQM_PROP_REPORT_ALARM;
01440   for (size_t i = 0, e = data_.qreports.size(); i < e; ++i)
01441     switch (data_.qreports[i].code)
01442     {
01443     case dqm::qstatus::STATUS_OK:
01444       break;
01445     case dqm::qstatus::WARNING:
01446       data_.flags |= DQMNet::DQM_PROP_REPORT_WARN;
01447       break;
01448     case dqm::qstatus::ERROR:
01449       data_.flags |= DQMNet::DQM_PROP_REPORT_ERROR;
01450       break;
01451     default:
01452       data_.flags |= DQMNet::DQM_PROP_REPORT_OTHER;
01453       break;
01454     }
01455 }
01456 
01457 // -------------------------------------------------------------------
01458 TObject *
01459 MonitorElement::getRootObject(void) const
01460 {
01461   const_cast<MonitorElement *>(this)->update();
01462   return object_;
01463 }
01464 
01465 TH1 *
01466 MonitorElement::getTH1(void) const
01467 {
01468   const_cast<MonitorElement *>(this)->update();
01469   return accessRootObject(__PRETTY_FUNCTION__, 0);
01470 }
01471 
01472 TH1F *
01473 MonitorElement::getTH1F(void) const
01474 {
01475   assert(kind() == DQM_KIND_TH1F);
01476   const_cast<MonitorElement *>(this)->update();
01477   return static_cast<TH1F *>(accessRootObject(__PRETTY_FUNCTION__, 1));
01478 }
01479 
01480 TH1S *
01481 MonitorElement::getTH1S(void) const
01482 {
01483   assert(kind() == DQM_KIND_TH1S);
01484   const_cast<MonitorElement *>(this)->update();
01485   return static_cast<TH1S *>(accessRootObject(__PRETTY_FUNCTION__, 1));
01486 }
01487 
01488 TH1D *
01489 MonitorElement::getTH1D(void) const
01490 {
01491   assert(kind() == DQM_KIND_TH1D);
01492   const_cast<MonitorElement *>(this)->update();
01493   return static_cast<TH1D *>(accessRootObject(__PRETTY_FUNCTION__, 1));
01494 }
01495 
01496 TH2F *
01497 MonitorElement::getTH2F(void) const
01498 {
01499   assert(kind() == DQM_KIND_TH2F);
01500   const_cast<MonitorElement *>(this)->update();
01501   return static_cast<TH2F *>(accessRootObject(__PRETTY_FUNCTION__, 2));
01502 }
01503 
01504 TH2S *
01505 MonitorElement::getTH2S(void) const
01506 {
01507   assert(kind() == DQM_KIND_TH2S);
01508   const_cast<MonitorElement *>(this)->update();
01509   return static_cast<TH2S *>(accessRootObject(__PRETTY_FUNCTION__, 2));
01510 }
01511 
01512 TH2D *
01513 MonitorElement::getTH2D(void) const
01514 {
01515   assert(kind() == DQM_KIND_TH2D);
01516   const_cast<MonitorElement *>(this)->update();
01517   return static_cast<TH2D *>(accessRootObject(__PRETTY_FUNCTION__, 2));
01518 }
01519 
01520 TH3F *
01521 MonitorElement::getTH3F(void) const
01522 {
01523   assert(kind() == DQM_KIND_TH3F);
01524   const_cast<MonitorElement *>(this)->update();
01525   return static_cast<TH3F *>(accessRootObject(__PRETTY_FUNCTION__, 3));
01526 }
01527 
01528 TProfile *
01529 MonitorElement::getTProfile(void) const
01530 {
01531   assert(kind() == DQM_KIND_TPROFILE);
01532   const_cast<MonitorElement *>(this)->update();
01533   return static_cast<TProfile *>(accessRootObject(__PRETTY_FUNCTION__, 1));
01534 }
01535 
01536 TProfile2D *
01537 MonitorElement::getTProfile2D(void) const
01538 {
01539   assert(kind() == DQM_KIND_TPROFILE2D);
01540   const_cast<MonitorElement *>(this)->update();
01541   return static_cast<TProfile2D *>(accessRootObject(__PRETTY_FUNCTION__, 2));
01542 }
01543 
01544 // -------------------------------------------------------------------
01545 TObject *
01546 MonitorElement::getRefRootObject(void) const
01547 {
01548   const_cast<MonitorElement *>(this)->update();
01549   return reference_;
01550 }
01551 
01552 TH1 *
01553 MonitorElement::getRefTH1(void) const
01554 {
01555   const_cast<MonitorElement *>(this)->update();
01556   return checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 0);
01557 }
01558 
01559 TH1F *
01560 MonitorElement::getRefTH1F(void) const
01561 {
01562   assert(kind() == DQM_KIND_TH1F);
01563   const_cast<MonitorElement *>(this)->update();
01564   return static_cast<TH1F *>
01565     (checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 1));
01566 }
01567 
01568 TH1S *
01569 MonitorElement::getRefTH1S(void) const
01570 {
01571   assert(kind() == DQM_KIND_TH1S);
01572   const_cast<MonitorElement *>(this)->update();
01573   return static_cast<TH1S *>
01574     (checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 1));
01575 }
01576 
01577 TH1D *
01578 MonitorElement::getRefTH1D(void) const
01579 {
01580   assert(kind() == DQM_KIND_TH1D);
01581   const_cast<MonitorElement *>(this)->update();
01582   return static_cast<TH1D *>
01583     (checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 1));
01584 }
01585 
01586 TH2F *
01587 MonitorElement::getRefTH2F(void) const
01588 {
01589   assert(kind() == DQM_KIND_TH2F);
01590   const_cast<MonitorElement *>(this)->update();
01591   return static_cast<TH2F *>
01592     (checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 2));
01593 }
01594 
01595 TH2S *
01596 MonitorElement::getRefTH2S(void) const
01597 {
01598   assert(kind() == DQM_KIND_TH2S);
01599   const_cast<MonitorElement *>(this)->update();
01600   return static_cast<TH2S *>
01601     (checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 2));
01602 }
01603 
01604 TH2D *
01605 MonitorElement::getRefTH2D(void) const
01606 {
01607   assert(kind() == DQM_KIND_TH2D);
01608   const_cast<MonitorElement *>(this)->update();
01609   return static_cast<TH2D *>
01610     (checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 2));
01611 }
01612 
01613 TH3F *
01614 MonitorElement::getRefTH3F(void) const
01615 {
01616   assert(kind() == DQM_KIND_TH3F);
01617   const_cast<MonitorElement *>(this)->update();
01618   return static_cast<TH3F *>
01619     (checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 3));
01620 }
01621 
01622 TProfile *
01623 MonitorElement::getRefTProfile(void) const
01624 {
01625   assert(kind() == DQM_KIND_TPROFILE);
01626   const_cast<MonitorElement *>(this)->update();
01627   return static_cast<TProfile *>
01628     (checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 1));
01629 }
01630 
01631 TProfile2D *
01632 MonitorElement::getRefTProfile2D(void) const
01633 {
01634   assert(kind() == DQM_KIND_TPROFILE2D);
01635   const_cast<MonitorElement *>(this)->update();
01636   return static_cast<TProfile2D *>
01637     (checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 2));
01638 }