CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_4/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 
00517 std::string
00518 MonitorElement::effLabelString(void) const
00519 {
00520   std::string result;
00521 
00522   result.reserve(6 + 2*data_.objname.size() + 3);
00523   result += '<'; result += data_.objname; result += '>';
00524   result += "e=1";
00525   result += '<'; result += '/'; result += data_.objname; result += '>';
00526   return result;
00527 }
00528 
00529 std::string
00530 MonitorElement::qualityTagString(const DQMNet::QValue &qv) const
00531 {
00532   char buf[64];
00533   std::string result;
00534   size_t titlelen = data_.objname.size() + qv.qtname.size() + 1;
00535   size_t buflen = sprintf(buf, "qr=st:%d:%.*g:", qv.code, DBL_DIG+2, qv.qtresult);
00536 
00537   result.reserve(7 + 2*titlelen + buflen + qv.algorithm.size() + qv.message.size());
00538   result += '<'; result += data_.objname; result += '.'; result += qv.qtname; result += '>';
00539   result += buf; result += qv.algorithm; result += ':'; result += qv.message;
00540   result += '<'; result += '/'; result += data_.objname; result += '.'; result += qv.qtname; result += '>';
00541   return result;
00542 }
00543 
00544 const QReport *
00545 MonitorElement::getQReport(const std::string &qtname) const
00546 {
00547   QReport *qr;
00548   DQMNet::QValue *qv;
00549   const_cast<MonitorElement *>(this)->getQReport(false, qtname, qr, qv);
00550   return qr;
00551 }
00552 
00553 std::vector<QReport *>
00554 MonitorElement::getQReports(void) const
00555 {
00556   std::vector<QReport *> result;
00557   result.reserve(qreports_.size());
00558   for (size_t i = 0, e = qreports_.size(); i != e; ++i)
00559   {
00560     const_cast<MonitorElement *>(this)->qreports_[i].qvalue_
00561       = const_cast<DQMNet::QValue *>(&data_.qreports[i]);
00562     result.push_back(const_cast<QReport *>(&qreports_[i]));
00563   }
00564   return result;
00565 }
00566 
00567 std::vector<QReport *>
00568 MonitorElement::getQWarnings(void) const
00569 {
00570   std::vector<QReport *> result;
00571   result.reserve(qreports_.size());
00572   for (size_t i = 0, e = qreports_.size(); i != e; ++i)
00573     if (data_.qreports[i].code == dqm::qstatus::WARNING)
00574     {
00575       const_cast<MonitorElement *>(this)->qreports_[i].qvalue_
00576         = const_cast<DQMNet::QValue *>(&data_.qreports[i]);
00577       result.push_back(const_cast<QReport *>(&qreports_[i]));
00578     }
00579   return result;
00580 }
00581 
00582 std::vector<QReport *>
00583 MonitorElement::getQErrors(void) const
00584 {
00585   std::vector<QReport *> result;
00586   result.reserve(qreports_.size());
00587   for (size_t i = 0, e = qreports_.size(); i != e; ++i)
00588     if (data_.qreports[i].code == dqm::qstatus::ERROR)
00589     {
00590       const_cast<MonitorElement *>(this)->qreports_[i].qvalue_
00591         = const_cast<DQMNet::QValue *>(&data_.qreports[i]);
00592       result.push_back(const_cast<QReport *>(&qreports_[i]));
00593     }
00594   return result;
00595 }
00596 
00597 std::vector<QReport *>
00598 MonitorElement::getQOthers(void) const
00599 {
00600   std::vector<QReport *> result;
00601   result.reserve(qreports_.size());
00602   for (size_t i = 0, e = qreports_.size(); i != e; ++i)
00603     if (data_.qreports[i].code != dqm::qstatus::STATUS_OK
00604         && data_.qreports[i].code != dqm::qstatus::WARNING
00605         && data_.qreports[i].code != dqm::qstatus::ERROR)
00606     {
00607       const_cast<MonitorElement *>(this)->qreports_[i].qvalue_
00608         = const_cast<DQMNet::QValue *>(&data_.qreports[i]);
00609       result.push_back(const_cast<QReport *>(&qreports_[i]));
00610     }
00611   return result;
00612 }
00613 
00615 void
00616 MonitorElement::runQTests(void)
00617 {
00618   assert(qreports_.size() == data_.qreports.size());
00619 
00620   // Rerun quality tests where the ME or the quality algorithm was modified.
00621   bool dirty = wasUpdated();
00622   for (size_t i = 0, e = data_.qreports.size(); i < e; ++i)
00623   {
00624     DQMNet::QValue &qv = data_.qreports[i];
00625     QReport &qr = qreports_[i];
00626     QCriterion *qc = qr.qcriterion_;
00627     qr.qvalue_ = &qv;
00628 
00629     // if (qc && (dirty || qc->wasModified()))  // removed for new QTest (abm-090503)
00630     if (qc && dirty)
00631     {
00632       assert(qc->getName() == qv.qtname);
00633       std::string oldMessage = qv.message;
00634       int oldStatus = qv.code;
00635 
00636       qc->runTest(this, qr, qv);
00637 
00638       if (oldStatus != qv.code || oldMessage != qv.message)
00639         update();
00640     }
00641   }
00642 
00643   // Update QReport statistics.
00644   updateQReportStats();
00645 }
00646 
00647 void
00648 MonitorElement::incompatible(const char *func) const
00649 {
00650   raiseDQMError("MonitorElement", "Method '%s' cannot be invoked on monitor"
00651                 " element '%s'", func, data_.objname.c_str());
00652 }
00653 
00654 TH1 *
00655 MonitorElement::accessRootObject(const char *func, int reqdim) const
00656 {
00657   if (kind() < DQM_KIND_TH1F)
00658     raiseDQMError("MonitorElement", "Method '%s' cannot be invoked on monitor"
00659                   " element '%s' because it is not a root object",
00660                   func, data_.objname.c_str());
00661 
00662   return checkRootObject(data_.objname, object_, func, reqdim);
00663 }
00664 
00665 /*** getter methods (wrapper around ROOT methods) ****/
00666 // 
00668 double
00669 MonitorElement::getMean(int axis /* = 1 */) const
00670 { return accessRootObject(__PRETTY_FUNCTION__, axis-1)
00671     ->GetMean(axis); }
00672 
00675 double
00676 MonitorElement::getMeanError(int axis /* = 1 */) const
00677 { return accessRootObject(__PRETTY_FUNCTION__, axis-1)
00678     ->GetMeanError(axis); }
00679 
00681 double
00682 MonitorElement::getRMS(int axis /* = 1 */) const
00683 { return accessRootObject(__PRETTY_FUNCTION__, axis-1)
00684     ->GetRMS(axis); }
00685 
00687 double
00688 MonitorElement::getRMSError(int axis /* = 1 */) const
00689 { return accessRootObject(__PRETTY_FUNCTION__, axis-1)
00690     ->GetRMSError(axis); }
00691 
00693 int
00694 MonitorElement::getNbinsX(void) const
00695 { return accessRootObject(__PRETTY_FUNCTION__, 1)
00696     ->GetNbinsX(); }
00697 
00699 int
00700 MonitorElement::getNbinsY(void) const
00701 { return accessRootObject(__PRETTY_FUNCTION__, 2)
00702     ->GetNbinsY(); }
00703 
00705 int
00706 MonitorElement::getNbinsZ(void) const
00707 { return accessRootObject(__PRETTY_FUNCTION__, 3)
00708     ->GetNbinsZ(); }
00709 
00711 double
00712 MonitorElement::getBinContent(int binx) const
00713 { return accessRootObject(__PRETTY_FUNCTION__, 1)
00714     ->GetBinContent(binx); }
00715 
00717 double
00718 MonitorElement::getBinContent(int binx, int biny) const
00719 { return accessRootObject(__PRETTY_FUNCTION__, 2)
00720     ->GetBinContent(binx, biny); }
00721 
00723 double
00724 MonitorElement::getBinContent(int binx, int biny, int binz) const
00725 { return accessRootObject(__PRETTY_FUNCTION__, 3)
00726     ->GetBinContent(binx, biny, binz); }
00727 
00729 double
00730 MonitorElement::getBinError(int binx) const
00731 { return accessRootObject(__PRETTY_FUNCTION__, 1)
00732     ->GetBinError(binx); }
00733 
00735 double
00736 MonitorElement::getBinError(int binx, int biny) const
00737 { return accessRootObject(__PRETTY_FUNCTION__, 2)
00738     ->GetBinError(binx, biny); }
00739 
00741 double
00742 MonitorElement::getBinError(int binx, int biny, int binz) const
00743 { return accessRootObject(__PRETTY_FUNCTION__, 3)
00744     ->GetBinError(binx, biny, binz); }
00745 
00747 double
00748 MonitorElement::getEntries(void) const
00749 { return accessRootObject(__PRETTY_FUNCTION__, 1)
00750     ->GetEntries(); }
00751 
00753 double
00754 MonitorElement::getBinEntries(int bin) const
00755 {
00756   if (kind() == DQM_KIND_TPROFILE)
00757     return static_cast<TProfile *>(accessRootObject(__PRETTY_FUNCTION__, 1))
00758       ->GetBinEntries(bin);
00759   else if (kind() == DQM_KIND_TPROFILE2D)
00760     return static_cast<TProfile2D *>(accessRootObject(__PRETTY_FUNCTION__, 1))
00761       ->GetBinEntries(bin);
00762   else
00763   {
00764     incompatible(__PRETTY_FUNCTION__);
00765     return 0;
00766   }
00767 }
00768 
00770 double
00771 MonitorElement::getYmin(void) const
00772 {
00773   if (kind() == DQM_KIND_TPROFILE)
00774     return static_cast<TProfile *>(accessRootObject(__PRETTY_FUNCTION__, 1))
00775       ->GetYmin();
00776   else
00777   {
00778     incompatible(__PRETTY_FUNCTION__);
00779     return 0;
00780   }
00781 }
00782 
00784 double
00785 MonitorElement::getYmax(void) const
00786 {
00787   if (kind() == DQM_KIND_TPROFILE)
00788     return static_cast<TProfile *>(accessRootObject(__PRETTY_FUNCTION__, 1))
00789       ->GetYmax();
00790   else
00791   {
00792     incompatible(__PRETTY_FUNCTION__);
00793     return 0;
00794   }
00795 }
00796 
00798 std::string
00799 MonitorElement::getAxisTitle(int axis /* = 1 */) const
00800 { return getAxis(__PRETTY_FUNCTION__, axis)
00801     ->GetTitle(); }
00802 
00804 std::string
00805 MonitorElement::getTitle(void) const
00806 { return accessRootObject(__PRETTY_FUNCTION__, 1)
00807     ->GetTitle(); }
00808 
00809 /*** setter methods (wrapper around ROOT methods) ****/
00810 // 
00812 void
00813 MonitorElement::setBinContent(int binx, double content)
00814 {
00815   update();
00816   accessRootObject(__PRETTY_FUNCTION__, 1)
00817     ->SetBinContent(binx, content);
00818 }
00819 
00821 void
00822 MonitorElement::setBinContent(int binx, int biny, double content)
00823 {
00824   update();
00825   accessRootObject(__PRETTY_FUNCTION__, 2)
00826     ->SetBinContent(binx, biny, content); }
00827 
00829 void
00830 MonitorElement::setBinContent(int binx, int biny, int binz, double content)
00831 {
00832   update();
00833   accessRootObject(__PRETTY_FUNCTION__, 3)
00834     ->SetBinContent(binx, biny, binz, content); }
00835 
00837 void
00838 MonitorElement::setBinError(int binx, double error)
00839 {
00840   update();
00841   accessRootObject(__PRETTY_FUNCTION__, 1)
00842     ->SetBinError(binx, error);
00843 }
00844 
00846 void
00847 MonitorElement::setBinError(int binx, int biny, double error)
00848 {
00849   update();
00850   accessRootObject(__PRETTY_FUNCTION__, 2)
00851     ->SetBinError(binx, biny, error);
00852 }
00853 
00855 void
00856 MonitorElement::setBinError(int binx, int biny, int binz, double error)
00857 {
00858   update();
00859   accessRootObject(__PRETTY_FUNCTION__, 3)
00860     ->SetBinError(binx, biny, binz, error);
00861 }
00862 
00864 void
00865 MonitorElement::setBinEntries(int bin, double nentries)
00866 {
00867   update();
00868   if (kind() == DQM_KIND_TPROFILE)
00869     static_cast<TProfile *>(accessRootObject(__PRETTY_FUNCTION__, 1))
00870       ->SetBinEntries(bin, nentries);
00871   else if (kind() == DQM_KIND_TPROFILE2D)
00872     static_cast<TProfile2D *>(accessRootObject(__PRETTY_FUNCTION__, 1))
00873       ->SetBinEntries(bin, nentries);
00874   else
00875     incompatible(__PRETTY_FUNCTION__);
00876 }
00877 
00879 void
00880 MonitorElement::setEntries(double nentries)
00881 {
00882   update();
00883   accessRootObject(__PRETTY_FUNCTION__, 1)
00884     ->SetEntries(nentries);
00885 }
00886 
00888 void
00889 MonitorElement::setBinLabel(int bin, const std::string &label, int axis /* = 1 */)
00890 {
00891   update();
00892   if ( getAxis(__PRETTY_FUNCTION__, axis)->GetNbins() >= bin ) 
00893   {
00894     getAxis(__PRETTY_FUNCTION__, axis)
00895       ->SetBinLabel(bin, label.c_str());
00896   }
00897   else
00898   {
00899     //  edm::LogWarning ("MonitorElement") 
00900     std::cout << "*** MonitorElement: WARNING:"
00901               <<"setBinLabel: attempting to set label of non-existent bin number \n";
00902   }
00903 }
00904 
00906 void
00907 MonitorElement::setAxisRange(double xmin, double xmax, int axis /* = 1 */)
00908 {
00909   update();
00910   getAxis(__PRETTY_FUNCTION__, axis)
00911     ->SetRangeUser(xmin, xmax);
00912 }
00913 
00915 void
00916 MonitorElement::setAxisTitle(const std::string &title, int axis /* = 1 */)
00917 {
00918   update();
00919   getAxis(__PRETTY_FUNCTION__, axis)
00920     ->SetTitle(title.c_str());
00921 }
00922 
00924 void
00925 MonitorElement::setAxisTimeDisplay(int value, int axis /* = 1 */)
00926 {
00927   update();
00928   getAxis(__PRETTY_FUNCTION__, axis)
00929     ->SetTimeDisplay(value);
00930 }
00931 
00933 void
00934 MonitorElement::setAxisTimeFormat(const char *format /* = "" */, int axis /* = 1 */)
00935 {
00936   update();
00937   getAxis(__PRETTY_FUNCTION__, axis)
00938     ->SetTimeFormat(format);
00939 }
00940 
00942 void
00943 MonitorElement::setAxisTimeOffset(double toffset, const char *option /* ="local" */, int axis /* = 1 */)
00944 {
00945   update();
00946   getAxis(__PRETTY_FUNCTION__, axis)
00947     ->SetTimeOffset(toffset, option);
00948 }
00949 
00951 void
00952 MonitorElement::setTitle(const std::string &title)
00953 {
00954   update();
00955   accessRootObject(__PRETTY_FUNCTION__, 1)
00956     ->SetTitle(title.c_str());
00957 }
00958 
00959 TAxis *
00960 MonitorElement::getAxis(const char *func, int axis) const
00961 {
00962   TH1 *h = accessRootObject(func, axis-1);
00963   TAxis *a = 0;
00964   if (axis == 1)
00965     a = h->GetXaxis();
00966   else if (axis == 2)
00967     a = h->GetYaxis();
00968   else if (axis == 3)
00969     a = h->GetZaxis();
00970 
00971   if (! a)
00972     raiseDQMError("MonitorElement", "No such axis %d in monitor element"
00973                   " '%s' of type '%s'", axis, data_.objname.c_str(),
00974                   typeid(*h).name());
00975 
00976   return a;
00977 }
00978 
00979 // ------------ Operations for MEs that are normally never reset ---------
00980 
00983 void
00984 MonitorElement::softReset(void)
00985 {
00986   update();
00987 
00988   // Create the reference object the first time this is called.
00989   // On subsequent calls accumulate the current value to the
00990   // reference, and then reset the current value.  This way the
00991   // future contents will have the reference "subtracted".
00992   if (kind() == DQM_KIND_TH1F)
00993   {
00994     TH1F *orig = static_cast<TH1F *>(object_);
00995     TH1F *r = static_cast<TH1F *>(refvalue_);
00996     if (! r)
00997     {
00998       refvalue_ = r = new TH1F((std::string(orig->GetName()) + "_ref").c_str(),
00999                                orig->GetTitle(),
01000                                orig->GetNbinsX(),
01001                                orig->GetXaxis()->GetXmin(),
01002                                orig->GetXaxis()->GetXmax());
01003       r->SetDirectory(0);
01004       r->Reset();
01005     }
01006 
01007     r->Add(orig);
01008     orig->Reset();
01009   }
01010   else if (kind() == DQM_KIND_TH1S)
01011   {
01012     TH1S *orig = static_cast<TH1S *>(object_);
01013     TH1S *r = static_cast<TH1S *>(refvalue_);
01014     if (! r)
01015     {
01016       refvalue_ = r = new TH1S((std::string(orig->GetName()) + "_ref").c_str(),
01017                                orig->GetTitle(),
01018                                orig->GetNbinsX(),
01019                                orig->GetXaxis()->GetXmin(),
01020                                orig->GetXaxis()->GetXmax());
01021       r->SetDirectory(0);
01022       r->Reset();
01023     }
01024 
01025     r->Add(orig);
01026     orig->Reset();
01027   }
01028   else if (kind() == DQM_KIND_TH1D)
01029   {
01030     TH1D *orig = static_cast<TH1D *>(object_);
01031     TH1D *r = static_cast<TH1D *>(refvalue_);
01032     if (! r)
01033     {
01034       refvalue_ = r = new TH1D((std::string(orig->GetName()) + "_ref").c_str(),
01035                                orig->GetTitle(),
01036                                orig->GetNbinsX(),
01037                                orig->GetXaxis()->GetXmin(),
01038                                orig->GetXaxis()->GetXmax());
01039       r->SetDirectory(0);
01040       r->Reset();
01041     }
01042 
01043     r->Add(orig);
01044     orig->Reset();
01045   }
01046   else if (kind() == DQM_KIND_TH2F)
01047   {
01048     TH2F *orig = static_cast<TH2F *>(object_);
01049     TH2F *r = static_cast<TH2F *>(refvalue_);
01050     if (! r)
01051     {
01052       refvalue_ = r = new TH2F((std::string(orig->GetName()) + "_ref").c_str(),
01053                                orig->GetTitle(),
01054                                orig->GetNbinsX(),
01055                                orig->GetXaxis()->GetXmin(),
01056                                orig->GetXaxis()->GetXmax(),
01057                                orig->GetNbinsY(),
01058                                orig->GetYaxis()->GetXmin(),
01059                                orig->GetYaxis()->GetXmax());
01060       r->SetDirectory(0);
01061       r->Reset();
01062     }
01063 
01064     r->Add(orig);
01065     orig->Reset();
01066   }
01067   else if (kind() == DQM_KIND_TH2S)
01068   {
01069     TH2S *orig = static_cast<TH2S *>(object_);
01070     TH2S *r = static_cast<TH2S *>(refvalue_);
01071     if (! r)
01072     {
01073       refvalue_ = r = new TH2S((std::string(orig->GetName()) + "_ref").c_str(),
01074                                orig->GetTitle(),
01075                                orig->GetNbinsX(),
01076                                orig->GetXaxis()->GetXmin(),
01077                                orig->GetXaxis()->GetXmax(),
01078                                orig->GetNbinsY(),
01079                                orig->GetYaxis()->GetXmin(),
01080                                orig->GetYaxis()->GetXmax());
01081       r->SetDirectory(0);
01082       r->Reset();
01083     }
01084 
01085     r->Add(orig);
01086     orig->Reset();
01087   }
01088   else if (kind() == DQM_KIND_TH2D)
01089   {
01090     TH2D *orig = static_cast<TH2D *>(object_);
01091     TH2D *r = static_cast<TH2D *>(refvalue_);
01092     if (! r)
01093     {
01094       refvalue_ = r = new TH2D((std::string(orig->GetName()) + "_ref").c_str(),
01095                                orig->GetTitle(),
01096                                orig->GetNbinsX(),
01097                                orig->GetXaxis()->GetXmin(),
01098                                orig->GetXaxis()->GetXmax(),
01099                                orig->GetNbinsY(),
01100                                orig->GetYaxis()->GetXmin(),
01101                                orig->GetYaxis()->GetXmax());
01102       r->SetDirectory(0);
01103       r->Reset();
01104     }
01105 
01106     r->Add(orig);
01107     orig->Reset();
01108   }
01109   else if (kind() == DQM_KIND_TH3F)
01110   {
01111     TH3F *orig = static_cast<TH3F *>(object_);
01112     TH3F *r = static_cast<TH3F *>(refvalue_);
01113     if (! r)
01114     {
01115       refvalue_ = r = new TH3F((std::string(orig->GetName()) + "_ref").c_str(),
01116                                orig->GetTitle(),
01117                                orig->GetNbinsX(),
01118                                orig->GetXaxis()->GetXmin(),
01119                                orig->GetXaxis()->GetXmax(),
01120                                orig->GetNbinsY(),
01121                                orig->GetYaxis()->GetXmin(),
01122                                orig->GetYaxis()->GetXmax(),
01123                                orig->GetNbinsZ(),
01124                                orig->GetZaxis()->GetXmin(),
01125                                orig->GetZaxis()->GetXmax());
01126       r->SetDirectory(0);
01127       r->Reset();
01128     }
01129 
01130     r->Add(orig);
01131     orig->Reset();
01132   }
01133   else if (kind() == DQM_KIND_TPROFILE)
01134   {
01135     TProfile *orig = static_cast<TProfile *>(object_);
01136     TProfile *r = static_cast<TProfile *>(refvalue_);
01137     if (! r)
01138     {
01139       refvalue_ = r = new TProfile((std::string(orig->GetName()) + "_ref").c_str(),
01140                                    orig->GetTitle(),
01141                                    orig->GetNbinsX(),
01142                                    orig->GetXaxis()->GetXmin(),
01143                                    orig->GetXaxis()->GetXmax(),
01144                                    orig->GetYaxis()->GetXmin(),
01145                                    orig->GetYaxis()->GetXmax(),
01146                                    orig->GetErrorOption());
01147       r->SetDirectory(0);
01148       r->Reset();
01149     }
01150 
01151     addProfiles(r, orig, r, 1, 1);
01152     orig->Reset();
01153   }
01154   else if (kind() == DQM_KIND_TPROFILE2D)
01155   {
01156     TProfile2D *orig = static_cast<TProfile2D *>(object_);
01157     TProfile2D *r = static_cast<TProfile2D *>(refvalue_);
01158     if (! r)
01159     {
01160       refvalue_ = r = new TProfile2D((std::string(orig->GetName()) + "_ref").c_str(),
01161                                      orig->GetTitle(),
01162                                      orig->GetNbinsX(),
01163                                      orig->GetXaxis()->GetXmin(),
01164                                      orig->GetXaxis()->GetXmax(),
01165                                      orig->GetNbinsY(),
01166                                      orig->GetYaxis()->GetXmin(),
01167                                      orig->GetYaxis()->GetXmax(),
01168                                      orig->GetZaxis()->GetXmin(),
01169                                      orig->GetZaxis()->GetXmax(),
01170                                      orig->GetErrorOption());
01171       r->SetDirectory(0);
01172       r->Reset();
01173     }
01174 
01175     addProfiles(r, orig, r, 1, 1);
01176     orig->Reset();
01177   }
01178   else
01179     incompatible(__PRETTY_FUNCTION__);
01180 }
01181 
01183 void
01184 MonitorElement::disableSoftReset(void)
01185 {
01186   if (refvalue_)
01187   {
01188     if (kind() == DQM_KIND_TH1F
01189         || kind() == DQM_KIND_TH1S
01190         || kind() == DQM_KIND_TH1D
01191         || kind() == DQM_KIND_TH2F
01192         || kind() == DQM_KIND_TH2S
01193         || kind() == DQM_KIND_TH2D
01194         || kind() == DQM_KIND_TH3F)
01195     {
01196       TH1 *orig = static_cast<TH1 *>(object_);
01197       orig->Add(refvalue_);
01198     }
01199     else if (kind() == DQM_KIND_TPROFILE)
01200     {
01201       TProfile *orig = static_cast<TProfile *>(object_);
01202       TProfile *r = static_cast<TProfile *>(refvalue_);
01203       addProfiles(orig, r, orig, 1, 1);
01204     }
01205     else if (kind() == DQM_KIND_TPROFILE2D)
01206     {
01207       TProfile2D *orig = static_cast<TProfile2D *>(object_);
01208       TProfile2D *r = static_cast<TProfile2D *>(refvalue_);
01209       addProfiles(orig, r, orig, 1, 1);
01210     }
01211     else
01212       incompatible(__PRETTY_FUNCTION__);
01213 
01214     delete refvalue_;
01215     refvalue_ = 0;
01216   }
01217 }
01218   
01219 // implementation: Giuseppe.Della-Ricca@ts.infn.it
01220 // Can be called with sum = h1 or sum = h2
01221 void
01222 MonitorElement::addProfiles(TProfile *h1, TProfile *h2, TProfile *sum, float c1, float c2)
01223 {
01224   assert(h1);
01225   assert(h2);
01226   assert(sum);
01227 
01228   static const Int_t NUM_STAT = 6;
01229   Double_t stats1[NUM_STAT];
01230   Double_t stats2[NUM_STAT];
01231   Double_t stats3[NUM_STAT];
01232 
01233   bool isRebinOn = sum->TestBit(TH1::kCanRebin);
01234   sum->ResetBit(TH1::kCanRebin);
01235 
01236   for (Int_t i = 0; i < NUM_STAT; ++i)
01237     stats1[i] = stats2[i] = stats3[i] = 0;
01238 
01239   h1->GetStats(stats1);
01240   h2->GetStats(stats2);
01241 
01242   for (Int_t i = 0; i < NUM_STAT; ++i)
01243     stats3[i] = c1*stats1[i] + c2*stats2[i];
01244 
01245   stats3[1] = c1*TMath::Abs(c1)*stats1[1]
01246               + c2*TMath::Abs(c2)*stats2[1];
01247 
01248   Double_t entries = c1*h1->GetEntries() + c2* h2->GetEntries();
01249   TArrayD* h1sumw2 = h1->GetSumw2();
01250   TArrayD* h2sumw2 = h2->GetSumw2();
01251   for (Int_t bin = 0, nbin = sum->GetNbinsX()+1; bin <= nbin; ++bin) 
01252   {
01253     Double_t entries = c1*h1->GetBinEntries(bin)
01254                        + c2*h2->GetBinEntries(bin);
01255     Double_t content = c1*h1->GetBinEntries(bin)*h1->GetBinContent(bin)
01256                        + c2*h2->GetBinEntries(bin)*h2->GetBinContent(bin);
01257     Double_t error = TMath::Sqrt(c1*TMath::Abs(c1)*h1sumw2->fArray[bin]
01258                                  + c2*TMath::Abs(c2)*h2sumw2->fArray[bin]);
01259     sum->SetBinContent(bin, content);
01260     sum->SetBinError(bin, error);
01261     sum->SetBinEntries(bin, entries);
01262   }
01263   
01264   sum->SetEntries(entries);
01265   sum->PutStats(stats3);
01266   if (isRebinOn) sum->SetBit(TH1::kCanRebin);
01267 }
01268 
01269 // implementation: Giuseppe.Della-Ricca@ts.infn.it
01270 // Can be called with sum = h1 or sum = h2
01271 void
01272 MonitorElement::addProfiles(TProfile2D *h1, TProfile2D *h2, TProfile2D *sum, float c1, float c2)
01273 {
01274   assert(h1);
01275   assert(h2);
01276   assert(sum);
01277 
01278   static const Int_t NUM_STAT = 9;
01279   Double_t stats1[NUM_STAT];
01280   Double_t stats2[NUM_STAT];
01281   Double_t stats3[NUM_STAT];
01282 
01283   bool isRebinOn = sum->TestBit(TH1::kCanRebin);
01284   sum->ResetBit(TH1::kCanRebin);
01285 
01286   for (Int_t i = 0; i < NUM_STAT; ++i)
01287     stats1[i] = stats2[i] = stats3[i] = 0;
01288 
01289   h1->GetStats(stats1);
01290   h2->GetStats(stats2);
01291 
01292   for (Int_t i = 0; i < NUM_STAT; i++)
01293     stats3[i] = c1*stats1[i] + c2*stats2[i];
01294 
01295   stats3[1] = c1*TMath::Abs(c1)*stats1[1]
01296               + c2*TMath::Abs(c2)*stats2[1];
01297 
01298   Double_t entries = c1*h1->GetEntries() + c2*h2->GetEntries();
01299   TArrayD *h1sumw2 = h1->GetSumw2();
01300   TArrayD *h2sumw2 = h2->GetSumw2();
01301   for (Int_t xbin = 0, nxbin = sum->GetNbinsX()+1; xbin <= nxbin; ++xbin)
01302     for (Int_t ybin = 0, nybin = sum->GetNbinsY()+1; ybin <= nybin; ++ybin)
01303     {
01304       Int_t bin = sum->GetBin(xbin, ybin);
01305       Double_t entries = c1*h1->GetBinEntries(bin)
01306                          + c2*h2->GetBinEntries(bin);
01307       Double_t content = c1*h1->GetBinEntries(bin)*h1->GetBinContent(bin)
01308                          + c2*h2->GetBinEntries(bin)*h2->GetBinContent(bin);
01309       Double_t error = TMath::Sqrt(c1*TMath::Abs(c1)*h1sumw2->fArray[bin]
01310                                    + c2*TMath::Abs(c2)*h2sumw2->fArray[bin]);
01311 
01312       sum->SetBinContent(bin, content);
01313       sum->SetBinError(bin, error);
01314       sum->SetBinEntries(bin, entries);
01315     }
01316   sum->SetEntries(entries);
01317   sum->PutStats(stats3);
01318   if (isRebinOn) sum->SetBit(TH1::kCanRebin);
01319 }
01320 
01321 void
01322 MonitorElement::copyFunctions(TH1 *from, TH1 *to)
01323 {
01324   // will copy functions only if local-copy and original-object are equal
01325   // (ie. no soft-resetting or accumulating is enabled)
01326   if (isSoftResetEnabled() || isAccumulateEnabled())
01327     return;
01328 
01329   update();
01330   TList *fromf = from->GetListOfFunctions();
01331   TList *tof   = to->GetListOfFunctions();
01332   for (int i = 0, nfuncs = fromf ? fromf->GetSize() : 0; i < nfuncs; ++i)
01333   {
01334     TObject *obj = fromf->At(i);
01335     // not interested in statistics
01336     if (!strcmp(obj->IsA()->GetName(), "TPaveStats"))
01337       continue;
01338 
01339     if (TF1 *fn = dynamic_cast<TF1 *>(obj))
01340       tof->Add(new TF1(*fn));
01341     //else if (dynamic_cast<TPaveStats *>(obj))
01342     //  ; // FIXME? tof->Add(new TPaveStats(*stats));
01343     else
01344       raiseDQMError("MonitorElement", "Cannot extract function '%s' of type"
01345                     " '%s' from monitor element '%s' for a copy",
01346                     obj->GetName(), obj->IsA()->GetName(), data_.objname.c_str());
01347   }
01348 }
01349 
01350 void
01351 MonitorElement::copyFrom(TH1 *from)
01352 {
01353   TH1 *orig = accessRootObject(__PRETTY_FUNCTION__, 1);
01354   if (orig->GetTitle() != from->GetTitle())
01355     orig->SetTitle(from->GetTitle());
01356 
01357   if (!isAccumulateEnabled())
01358     orig->Reset();
01359 
01360   if (isSoftResetEnabled())
01361   {
01362     if (kind() == DQM_KIND_TH1F
01363         || kind() == DQM_KIND_TH1S
01364         || kind() == DQM_KIND_TH1D
01365         || kind() == DQM_KIND_TH2F
01366         || kind() == DQM_KIND_TH2S
01367         || kind() == DQM_KIND_TH2D
01368         || kind() == DQM_KIND_TH3F)
01369       // subtract "reference"
01370       orig->Add(from, refvalue_, 1, -1);
01371     else if (kind() == DQM_KIND_TPROFILE)
01372       // subtract "reference"
01373       addProfiles(static_cast<TProfile *>(from),
01374                   static_cast<TProfile *>(refvalue_),
01375                   static_cast<TProfile *>(orig),
01376                   1, -1);
01377     else if (kind() == DQM_KIND_TPROFILE2D)
01378       // subtract "reference"
01379       addProfiles(static_cast<TProfile2D *>(from),
01380                   static_cast<TProfile2D *>(refvalue_),
01381                   static_cast<TProfile2D *>(orig),
01382                   1, -1);
01383     else
01384       incompatible(__PRETTY_FUNCTION__);
01385   }
01386   else
01387     orig->Add(from);
01388 
01389   copyFunctions(from, orig);
01390 }
01391 
01392 // --- Operations on MEs that are normally reset at end of monitoring cycle ---
01393 void
01394 MonitorElement::getQReport(bool create, const std::string &qtname, QReport *&qr, DQMNet::QValue *&qv)
01395 {
01396   assert(qreports_.size() == data_.qreports.size());
01397 
01398   qr = 0;
01399   qv = 0;
01400 
01401   size_t pos = 0, end = qreports_.size();
01402   while (pos < end && data_.qreports[pos].qtname != qtname)
01403     ++pos;
01404 
01405   if (pos == end && ! create)
01406     return;
01407   else if (pos == end)
01408   {
01409     data_.qreports.push_back(DQMNet::QValue());
01410     qreports_.push_back(QReport(0, 0));
01411 
01412     DQMNet::QValue &q = data_.qreports.back();
01413     q.code = dqm::qstatus::DID_NOT_RUN;
01414     q.qtresult = 0;
01415     q.qtname = qtname;
01416     q.message = "NO_MESSAGE_ASSIGNED";
01417     q.algorithm = "UNKNOWN_ALGORITHM";
01418   }
01419 
01420   qr = &qreports_[pos];
01421   qv = &data_.qreports[pos];
01422 }
01423   
01425 void
01426 MonitorElement::addQReport(const DQMNet::QValue &desc, QCriterion *qc)
01427 {
01428   QReport *qr;
01429   DQMNet::QValue *qv;
01430   getQReport(true, desc.qtname, qr, qv);
01431   qr->qcriterion_ = qc;
01432   *qv = desc;
01433   update();
01434 }
01435 
01436 void
01437 MonitorElement::addQReport(QCriterion *qc)
01438 {
01439   QReport *qr;
01440   DQMNet::QValue *qv;
01441   getQReport(true, qc->getName(), qr, qv);
01442   qv->code = dqm::qstatus::DID_NOT_RUN;
01443   qv->message = "NO_MESSAGE_ASSIGNED";
01444   qr->qcriterion_ = qc;
01445   update();
01446 }
01447 
01449 void
01450 MonitorElement::updateQReportStats(void)
01451 {
01452   data_.flags &= ~DQMNet::DQM_PROP_REPORT_ALARM;
01453   for (size_t i = 0, e = data_.qreports.size(); i < e; ++i)
01454     switch (data_.qreports[i].code)
01455     {
01456     case dqm::qstatus::STATUS_OK:
01457       break;
01458     case dqm::qstatus::WARNING:
01459       data_.flags |= DQMNet::DQM_PROP_REPORT_WARN;
01460       break;
01461     case dqm::qstatus::ERROR:
01462       data_.flags |= DQMNet::DQM_PROP_REPORT_ERROR;
01463       break;
01464     default:
01465       data_.flags |= DQMNet::DQM_PROP_REPORT_OTHER;
01466       break;
01467     }
01468 }
01469 
01470 // -------------------------------------------------------------------
01471 TObject *
01472 MonitorElement::getRootObject(void) const
01473 {
01474   const_cast<MonitorElement *>(this)->update();
01475   return object_;
01476 }
01477 
01478 TH1 *
01479 MonitorElement::getTH1(void) const
01480 {
01481   const_cast<MonitorElement *>(this)->update();
01482   return accessRootObject(__PRETTY_FUNCTION__, 0);
01483 }
01484 
01485 TH1F *
01486 MonitorElement::getTH1F(void) const
01487 {
01488   assert(kind() == DQM_KIND_TH1F);
01489   const_cast<MonitorElement *>(this)->update();
01490   return static_cast<TH1F *>(accessRootObject(__PRETTY_FUNCTION__, 1));
01491 }
01492 
01493 TH1S *
01494 MonitorElement::getTH1S(void) const
01495 {
01496   assert(kind() == DQM_KIND_TH1S);
01497   const_cast<MonitorElement *>(this)->update();
01498   return static_cast<TH1S *>(accessRootObject(__PRETTY_FUNCTION__, 1));
01499 }
01500 
01501 TH1D *
01502 MonitorElement::getTH1D(void) const
01503 {
01504   assert(kind() == DQM_KIND_TH1D);
01505   const_cast<MonitorElement *>(this)->update();
01506   return static_cast<TH1D *>(accessRootObject(__PRETTY_FUNCTION__, 1));
01507 }
01508 
01509 TH2F *
01510 MonitorElement::getTH2F(void) const
01511 {
01512   assert(kind() == DQM_KIND_TH2F);
01513   const_cast<MonitorElement *>(this)->update();
01514   return static_cast<TH2F *>(accessRootObject(__PRETTY_FUNCTION__, 2));
01515 }
01516 
01517 TH2S *
01518 MonitorElement::getTH2S(void) const
01519 {
01520   assert(kind() == DQM_KIND_TH2S);
01521   const_cast<MonitorElement *>(this)->update();
01522   return static_cast<TH2S *>(accessRootObject(__PRETTY_FUNCTION__, 2));
01523 }
01524 
01525 TH2D *
01526 MonitorElement::getTH2D(void) const
01527 {
01528   assert(kind() == DQM_KIND_TH2D);
01529   const_cast<MonitorElement *>(this)->update();
01530   return static_cast<TH2D *>(accessRootObject(__PRETTY_FUNCTION__, 2));
01531 }
01532 
01533 TH3F *
01534 MonitorElement::getTH3F(void) const
01535 {
01536   assert(kind() == DQM_KIND_TH3F);
01537   const_cast<MonitorElement *>(this)->update();
01538   return static_cast<TH3F *>(accessRootObject(__PRETTY_FUNCTION__, 3));
01539 }
01540 
01541 TProfile *
01542 MonitorElement::getTProfile(void) const
01543 {
01544   assert(kind() == DQM_KIND_TPROFILE);
01545   const_cast<MonitorElement *>(this)->update();
01546   return static_cast<TProfile *>(accessRootObject(__PRETTY_FUNCTION__, 1));
01547 }
01548 
01549 TProfile2D *
01550 MonitorElement::getTProfile2D(void) const
01551 {
01552   assert(kind() == DQM_KIND_TPROFILE2D);
01553   const_cast<MonitorElement *>(this)->update();
01554   return static_cast<TProfile2D *>(accessRootObject(__PRETTY_FUNCTION__, 2));
01555 }
01556 
01557 // -------------------------------------------------------------------
01558 TObject *
01559 MonitorElement::getRefRootObject(void) const
01560 {
01561   const_cast<MonitorElement *>(this)->update();
01562   return reference_;
01563 }
01564 
01565 TH1 *
01566 MonitorElement::getRefTH1(void) const
01567 {
01568   const_cast<MonitorElement *>(this)->update();
01569   return checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 0);
01570 }
01571 
01572 TH1F *
01573 MonitorElement::getRefTH1F(void) const
01574 {
01575   assert(kind() == DQM_KIND_TH1F);
01576   const_cast<MonitorElement *>(this)->update();
01577   return static_cast<TH1F *>
01578     (checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 1));
01579 }
01580 
01581 TH1S *
01582 MonitorElement::getRefTH1S(void) const
01583 {
01584   assert(kind() == DQM_KIND_TH1S);
01585   const_cast<MonitorElement *>(this)->update();
01586   return static_cast<TH1S *>
01587     (checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 1));
01588 }
01589 
01590 TH1D *
01591 MonitorElement::getRefTH1D(void) const
01592 {
01593   assert(kind() == DQM_KIND_TH1D);
01594   const_cast<MonitorElement *>(this)->update();
01595   return static_cast<TH1D *>
01596     (checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 1));
01597 }
01598 
01599 TH2F *
01600 MonitorElement::getRefTH2F(void) const
01601 {
01602   assert(kind() == DQM_KIND_TH2F);
01603   const_cast<MonitorElement *>(this)->update();
01604   return static_cast<TH2F *>
01605     (checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 2));
01606 }
01607 
01608 TH2S *
01609 MonitorElement::getRefTH2S(void) const
01610 {
01611   assert(kind() == DQM_KIND_TH2S);
01612   const_cast<MonitorElement *>(this)->update();
01613   return static_cast<TH2S *>
01614     (checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 2));
01615 }
01616 
01617 TH2D *
01618 MonitorElement::getRefTH2D(void) const
01619 {
01620   assert(kind() == DQM_KIND_TH2D);
01621   const_cast<MonitorElement *>(this)->update();
01622   return static_cast<TH2D *>
01623     (checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 2));
01624 }
01625 
01626 TH3F *
01627 MonitorElement::getRefTH3F(void) const
01628 {
01629   assert(kind() == DQM_KIND_TH3F);
01630   const_cast<MonitorElement *>(this)->update();
01631   return static_cast<TH3F *>
01632     (checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 3));
01633 }
01634 
01635 TProfile *
01636 MonitorElement::getRefTProfile(void) const
01637 {
01638   assert(kind() == DQM_KIND_TPROFILE);
01639   const_cast<MonitorElement *>(this)->update();
01640   return static_cast<TProfile *>
01641     (checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 1));
01642 }
01643 
01644 TProfile2D *
01645 MonitorElement::getRefTProfile2D(void) const
01646 {
01647   assert(kind() == DQM_KIND_TPROFILE2D);
01648   const_cast<MonitorElement *>(this)->update();
01649   return static_cast<TProfile2D *>
01650     (checkRootObject(data_.objname, reference_, __PRETTY_FUNCTION__, 2));
01651 }