CMS 3D CMS Logo

Classes | Public Member Functions | Public Attributes | Private Member Functions | Private Attributes

HDQMInspector Class Reference

#include <HDQMInspector.h>

List of all members.

Classes

struct  DetIdItemList

Public Member Functions

void closeFile ()
void createTrend (const std::string ListItems, const std::string CanvasName="", const int logy=0, const std::string Conditions="", std::string const &Labels="", const unsigned int firstRun=1, const unsigned int lastRun=0xFFFFFFFE, int const UseYRange=0, double const &YMin=999999, double const &YMax=-999999)
void createTrendLastRuns (const std::string ListItems, const std::string CanvasName="", const int logy=0, const std::string Conditions="", std::string const &Labels="", const unsigned int nRuns=10, int const UseYRange=0, double const &YMin=999999, double const &YMax=-999999)
double findGraphMax (TGraphErrors *)
double findGraphMin (TGraphErrors *)
std::vector< std::string > getListItems ()
std::vector< unsigned int > getRuns ()
std::vector< float > getSummary ()
std::vector< unsigned int > getvDetId ()
 HDQMInspector ()
 HDQMInspector (const HDQMInspectorConfigBase *InConfig)
std::string readListFromFile (const std::string &listFileName)
void setBlackList (std::string const &ListItems)
void setDB (const std::string &DBName, const std::string &DBTag, const std::string &DBauth="")
void setDebug (int i)
void setDoStat (int i)
void setSeparator (std::string const in)
void setSkip0s (bool const in)
void setSkip99s (bool const in)
void setWhiteList (std::string const &ListItems)
virtual ~HDQMInspector ()

Public Attributes

TFile * fOutFile

Private Member Functions

void accessDB ()
bool ApplyConditions (std::string &, std::vector< DetIdItemList > &)
void InitializeIOVList ()
bool isListed (unsigned int run, std::vector< unsigned int > &vList)
void plot (size_t &nPads, std::string CanvasName, int logy=0, std::string const &Labels="", int const UseYRange=0, double const XMin=999999, double const YMin=-999999)
void setItems (std::string)
bool setRange (unsigned int &firstRun, unsigned int &lastRun)
void style ()
void unpackConditions (std::string &, std::vector< DetIdItemList > &)
size_t unpackItems (std::string &)

Private Attributes

std::vector< unsigned int > blackList
std::string DBauth_
std::string DBName_
std::string DBTag_
const HDQMInspectorConfigBasefHDQMInspectorConfig
std::string fSep
bool fSkip0s
bool fSkip99s
int iDebug
int iDoStat
std::vector< unsigned int > iovList
CondCachedIter< HDQMSummary > * Iterator
std::vector< unsigned int > vdetId_
std::vector< DetIdItemListvDetIdItemList_
std::vector< std::string > vlistItems_
std::vector< unsigned int > vRun_
std::vector< float > vSummary_
std::vector< unsigned int > whiteList

Detailed Description

Definition at line 28 of file HDQMInspector.h.


Constructor & Destructor Documentation

HDQMInspector::HDQMInspector ( ) [inline]

Definition at line 31 of file HDQMInspector.h.

                 :
    DBName_(""),
    DBTag_(""),
    DBauth_(""),
    Iterator(0),
    iDebug(0),
    iDoStat(0),
    fSkip99s(false),
    fSkip0s(false),
    fHDQMInspectorConfig(0x0),
    fSep("@")
    {
    };
HDQMInspector::HDQMInspector ( const HDQMInspectorConfigBase InConfig) [inline]

Definition at line 45 of file HDQMInspector.h.

                                                        :
    DBName_(""),
    DBTag_(""),
    DBauth_(""),
    Iterator(0),
    iDebug(0),
    iDoStat(0),
    fSkip99s(false),
    fSkip0s(false),
    fHDQMInspectorConfig(InConfig),
    fSep("@")
    {
    };
virtual HDQMInspector::~HDQMInspector ( ) [inline, virtual]

Definition at line 59 of file HDQMInspector.h.

References Iterator.

                           {
    delete Iterator;
  };

Member Function Documentation

void HDQMInspector::accessDB ( ) [private]

Definition at line 82 of file HDQMInspector.cc.

References gather_cfg::cout, DBauth_, DBName_, DBTag_, InitializeIOVList(), and Iterator.

Referenced by setDB().

{
  //double start, end;
  // start = clock();   
 
  if(Iterator!=0)
    delete Iterator;
  
  Iterator = new CondCachedIter<HDQMSummary>(); 

  std::cout << "creating connection" << std::endl;
  Iterator->create(DBName_,DBTag_,DBauth_);
  std::cout << "connection created" << std::endl;

  
  InitializeIOVList();
  //  end = clock();
  //  if(iDebug)
  //std::cout <<"Time Creation link with Database = " <<  ((double) (end - start)) << " (a.u.)" <<std::endl; 
}
bool HDQMInspector::ApplyConditions ( std::string &  Conditions,
std::vector< DetIdItemList > &  vdetIdItemList 
) [private]

Definition at line 736 of file HDQMInspector.cc.

References gather_cfg::cout, cond::rpcobgas::detid, fSep, iDebug, NULL, and makeHLTPrescaleTable::values.

Referenced by createTrend().

{
  double resultdbl=1;
  char cConditions[1024];
  char singleCondition[1024];
  char condCVal[1024];   

  sprintf(cConditions,"%s",Conditions.c_str());
  if (iDebug) {
    std::cout << "Conditions " << cConditions << std::endl;
  }
  for(size_t ic=0;ic<vdetIdItemList.size();++ic)
    for(size_t jc=0;jc<vdetIdItemList[ic].items.size();++jc){
      //scientific precision doesn't work in HDQMExpressionEvaluator...
      //sprintf(condCVal,"%g",vdetIdItemList[ic].values[jc]);
      sprintf(condCVal,"%f",vdetIdItemList[ic].values[jc]);
      sprintf(singleCondition,"%d%s%s",vdetIdItemList[ic].detid,fSep.c_str(),vdetIdItemList[ic].items[jc].c_str());
      //printf("dhidas %d  %s  %s\n",vdetIdItemList[ic].detid,fSep.c_str(),vdetIdItemList[ic].items[jc].c_str());
      //printf("dhidas %s %s\n", cConditions, singleCondition);
      char* fpos = strstr(cConditions,singleCondition);
      //printf("dhidas %s %s %i\n", fpos, condCVal, strlen(condCVal));
      strncpy(fpos,condCVal,strlen(condCVal));
      memset(fpos+strlen(condCVal),' ',strlen(singleCondition)-strlen(condCVal));
      //std::cout << "fpos " << fpos << " len condCVal " << strlen(condCVal) << " strlen(singleCondition) " << strlen(singleCondition) << " len cConditions " << strlen(cConditions)<<std::endl;
      //std::cout << "Conditions Replace: Condition " << singleCondition << " string changed in " << cConditions << std::endl;
    }

  std::string stringToEvaluate;
  char * pch;
  pch = strtok (cConditions," ");
  while (pch != NULL){
    stringToEvaluate.append(pch);
    pch = strtok (NULL, " ");
  } 
  //for(size_t i=0;i<strlen(cConditions);++i)
  // if(cConditions[i] != " ")
  //  stringToEvaluate.push_back(cConditions[i]);

  if(iDebug) {
    std::cout << "Conditions After SubStitution " << stringToEvaluate << std::endl;
  }
  TFormula Formula("condition", stringToEvaluate.c_str());
  resultdbl = Formula.Eval(0);
  if(iDebug) {
    std::cout << "Result " << resultdbl << std::endl;
  }
  if(!resultdbl) {
    return false;
  }
  return true;
}
void HDQMInspector::closeFile ( ) [inline]

Definition at line 86 of file HDQMInspector.h.

References fOutFile.

  { 
    if( fOutFile ) {
      fOutFile->Close();
    }
  }
void HDQMInspector::createTrend ( const std::string  ListItems,
const std::string  CanvasName = "",
const int  logy = 0,
const std::string  Conditions = "",
std::string const &  Labels = "",
const unsigned int  firstRun = 1,
const unsigned int  lastRun = 0xFFFFFFFE,
int const  UseYRange = 0,
double const &  YMin = 999999,
double const &  YMax = -999999 
)

Definition at line 320 of file HDQMInspector.cc.

References ApplyConditions(), blackList, gather_cfg::cout, cond::rpcobgas::detid, HDQMSummary::getRunNr(), HDQMSummary::getSummaryObj(), i, iDebug, isListed(), HDQMInspector::DetIdItemList::items, plot(), RecoTauValidation_cfi::reference, setRange(), unpackConditions(), unpackItems(), HDQMInspector::DetIdItemList::values, vDetIdItemList_, vRun_, vSummary_, and whiteList.

Referenced by createTrendLastRuns().

{
  std::cout << "\n****************\nCreateTrend\n****************\n" << std::endl;
  std::cout << "ListItems : " << ListItems << std::endl;
  std::cout << "Conditions : " << Conditions << std::endl;

  vRun_.clear();
  vSummary_.clear();
  vDetIdItemList_.clear();

  std::vector<DetIdItemList> vDetIdItemListCut;
  
  size_t nPads=unpackItems(ListItems);   

  unpackConditions(Conditions,vDetIdItemListCut);
 
  //   double start = clock(); 

  std::cout << "firstRun " << firstRun << " lastRun " << lastRun << std::endl;
  if(!setRange(firstRun,lastRun)){
    Iterator->rewind();
    return;
  }
  const HDQMSummary* reference;
  while((reference = Iterator->next())) { 


    // Check the run and black and white lists
    if(Iterator->getStartTime()<firstRun || Iterator->getStartTime()>lastRun || isListed(reference->getRunNr(), blackList)) {
      continue;
    }
    if (whiteList.size() > 0 && !isListed(reference->getRunNr(), whiteList)) {
      continue;
    }

    if(vDetIdItemListCut.size()){
      for(size_t ij=0;ij!=vDetIdItemListCut.size();++ij){
        vDetIdItemListCut[ij].values=reference->getSummaryObj(vDetIdItemListCut[ij].detid, vDetIdItemListCut[ij].items);
      }

      if(!ApplyConditions(Conditions,vDetIdItemListCut))
        continue;
    }

    vRun_.push_back(reference->getRunNr());

    for(size_t ij=0;ij!=vDetIdItemList_.size();++ij){
      vDetIdItemList_[ij].values=reference->getSummaryObj(vDetIdItemList_[ij].detid, vDetIdItemList_[ij].items);
     
      vSummary_.insert(vSummary_.end(),vDetIdItemList_[ij].values.begin(),vDetIdItemList_[ij].values.end());   
      if(iDebug){
        std::cout << ListItems  << " run " << vRun_.back() << " values \n" ;
        DetIdItemList detiditemlist=vDetIdItemList_[ij];
        for(size_t i=0;i<detiditemlist.items.size();++i) {
          std::cout << "\t" << detiditemlist.items[i] << " " << detiditemlist.values[i] <<" " << i << " \n";
        }
        std::cout << "\n" << std::endl;
      }
    }
  }

  if(vRun_.size()) {
    plot(nPads, CanvasName, logy, Labels, UseYRange, YMin, YMax);
  }
   
     
  std::cout << "\n****** Ignore this error *****\n" << std::endl;
  Iterator->rewind();
  std::cout << "\n******************************\n" << std::endl;
}
void HDQMInspector::createTrendLastRuns ( const std::string  ListItems,
const std::string  CanvasName = "",
const int  logy = 0,
const std::string  Conditions = "",
std::string const &  Labels = "",
const unsigned int  nRuns = 10,
int const  UseYRange = 0,
double const &  YMin = 999999,
double const &  YMax = -999999 
)

Definition at line 297 of file HDQMInspector.cc.

References createTrend(), first, iovList, and prof2calltree::last.

{
  unsigned int first,last;
  unsigned int iovListSize = iovList.size();

  if (iovListSize>0) 
  { 
    last = iovList.back();

    if (iovListSize>=nRuns) {
      first = iovList.at(iovListSize-nRuns);
    } else {
      first = *iovList.begin();
    }
  }
  else return;

  createTrend(ListItems,CanvasName,logy,Conditions,Labels,first,last, UseYRange, YMin, YMax);

  return;
}
double HDQMInspector::findGraphMax ( TGraphErrors *  g)

Definition at line 832 of file HDQMInspector.cc.

References n, and detailsBasic3DVector::y.

Referenced by plot().

{
  // GetMaximum() doesn't work ....
  int n = g->GetN();
  double* y = g->GetY();
  int locmax = TMath::LocMax(n,y);
  return  y[locmax];
}
double HDQMInspector::findGraphMin ( TGraphErrors *  g)

Definition at line 842 of file HDQMInspector.cc.

References n, and detailsBasic3DVector::y.

Referenced by plot().

{
  // GetMinimum() doesn't work ....
  int n = g->GetN();
  double* y = g->GetY();
  int locmin = TMath::LocMin(n,y);
  return  y[locmin];
}
std::vector<std::string> HDQMInspector::getListItems ( ) [inline]

Definition at line 103 of file HDQMInspector.h.

References vlistItems_.

{ return vlistItems_;}
std::vector<unsigned int> HDQMInspector::getRuns ( ) [inline]

Definition at line 100 of file HDQMInspector.h.

References vRun_.

{ return vRun_;}
std::vector<float> HDQMInspector::getSummary ( ) [inline]

Definition at line 101 of file HDQMInspector.h.

References vSummary_.

{ return vSummary_;}
std::vector<unsigned int> HDQMInspector::getvDetId ( ) [inline]

Definition at line 104 of file HDQMInspector.h.

References vdetId_.

{ return vdetId_;}
void HDQMInspector::InitializeIOVList ( ) [private]

Definition at line 248 of file HDQMInspector.cc.

References gather_cfg::cout, iDebug, iovList, and RecoTauValidation_cfi::reference.

Referenced by accessDB().

{
  const HDQMSummary* reference;
  while((reference = Iterator->next())) {
    iovList.push_back(Iterator->getStartTime());
    if (iDebug) {
      std::cout << "iovList " << iovList.back() << std::endl;
    }
  } 
  Iterator->rewind();
}
bool HDQMInspector::isListed ( unsigned int  run,
std::vector< unsigned int > &  vList 
) [private]

Definition at line 230 of file HDQMInspector.cc.

References gather_cfg::cout, iDebug, and DTTTrigCorrFirst::run.

Referenced by createTrend().

{
  // This routine expectes a sorted list and returns true if the run is in the list,
  // false otherwise

  // Binary search is much faster, but you MUST give it a sorted list.
  if (std::binary_search(vList.begin(), vList.end(), run)) {
    if(iDebug) {
      std::cout << "\n Run "<< run << " is listed !!\n" << std::endl;
    }
    return true;
  }

  return false;

}
void HDQMInspector::plot ( size_t &  nPads,
std::string  CanvasName,
int  logy = 0,
std::string const &  Labels = "",
int const  UseYRange = 0,
double const  XMin = 999999,
double const  YMin = -999999 
) [private]

Definition at line 391 of file HDQMInspector.cc.

References funct::C, HDQMInspectorConfigBase::computeIntegral(), gather_cfg::cout, fHDQMInspectorConfig, spr::find(), findGraphMax(), findGraphMin(), fSep, fSkip0s, fSkip99s, i, iDebug, getHLTprescales::index, j, max(), min, mergeVDriftHistosByStation::name, mathSSE::sqrt(), style(), indexGen::title, HDQMInspectorConfigBase::translateDetId(), vdetId_, vDetIdItemList_, vlistItems_, vRun_, vSummary_, and X.

Referenced by createTrend().

{
  std::cout << "\n********\nplot\n*****\n"<< std::endl;

  style();

  double *X, *Y, *EX, *EY, *YCumul;
  X=new double[vRun_.size()];
  Y=new double[vRun_.size()];
  EX=new double[vRun_.size()]; 
  EY=new double[vRun_.size()];  
  YCumul=new double[vRun_.size()];
  
  size_t index;
  TCanvas *C;
  TGraphErrors *graph;

  if(CanvasName==""){
    char name[128];
    sprintf(name,"%d",(int) clock());
    CanvasName=std::string(name);
  }
  
  std::string rootCName = CanvasName;
  rootCName.replace(rootCName.find("."),rootCName.size()-rootCName.find("."),"");
 
  C=new TCanvas(rootCName.c_str(),"");
  int ndiv=(int) sqrt(nPads);
  C->Divide(ndiv,nPads/ndiv+ (nPads%ndiv?1:0));
 
  int padCount=0;

  vlistItems_.clear();
  vdetId_.clear();

  for(size_t ic=0;ic<vDetIdItemList_.size();++ic){
    vlistItems_.insert(vlistItems_.end(),vDetIdItemList_[ic].items.begin(),vDetIdItemList_[ic].items.end());
    vdetId_.insert(vdetId_.end(),vDetIdItemList_[ic].items.size(),vDetIdItemList_[ic].detid);
  }

  // Vector of graphs in this request and DetNames which correspond to them
  std::vector<TGraphErrors*> VectorOfGraphs;
  std::vector<std::string> VectorOfDetNames;

  for(size_t i=0;i<vlistItems_.size();++i){
    std::cout <<  "TkRegion " << vdetId_[i] << " " << vlistItems_[i] << std::endl;

    if(vlistItems_.at(i).find("Summary")!= std::string::npos) vlistItems_.at(i).replace(vlistItems_.at(i).find("Summary_"),8,"");
    if(vlistItems_.at(i).find(fSep)!= std::string::npos) vlistItems_.at(i).replace(vlistItems_.at(i).find(fSep),fSep.size(),"_");
    
 
    std::stringstream ss;
    if (fHDQMInspectorConfig != 0x0) {
      ss << fHDQMInspectorConfig->translateDetId( vdetId_[i] ) << vlistItems_[i];
      VectorOfDetNames.push_back( fHDQMInspectorConfig->translateDetId( vdetId_[i] ));
    } else {
      ss << "Id " << vdetId_[i] << " " << vlistItems_[i];
      VectorOfDetNames.push_back( "???" );
    }

    
    bool const itemForIntegration = fHDQMInspectorConfig ? fHDQMInspectorConfig->computeIntegral(vlistItems_[i]) : false;
   
    int addShift=0;
    for(size_t j=0;j<vRun_.size();++j){
      index=j*vlistItems_.size()+i;
      X[j]=vRun_[j];
      EX[j]=0;
      Y[j]=vSummary_[index];
      //if (Y[j]==-10 || Y[j]==-9999 || Y[j] ==-99) {EY[j] = 0; Y[j] = 0;}
       
      // -9999 : existing HDQMSummary object in DB but part of the information not uploaded
      // -99   : HDQMSummary object not existing for this detId, informations are missing for all quantities 
      // -10 bad fit ?

      //std::cout << "dhidas " << vlistItems_[i] << "  " << vRun_[j] << "  " << vSummary_[index] << std::endl;
     
      if(vlistItems_[i].find("mean")!=std::string::npos){
        //if the quantity requested is mean, the error is evaluated as the error on the mean=rms/sqrt(entries)
        EY[j]=vSummary_[index+2]>0?vSummary_[index+1]/sqrt(vSummary_[index+2]):0;
        addShift=2;
      }else if (vlistItems_[i].find("entries")!=std::string::npos) {
        addShift=0;
      }else if (vlistItems_[i].find("landauPeak")!=std::string::npos){
        EY[j]=vSummary_[index+1];
        addShift=1;
      }
      else if (vlistItems_[i].find("gaussMean")!=std::string::npos){
        EY[j]=vSummary_[index+1];
        addShift=1;
      }
      else if (vlistItems_[i].find("Chi2NDF")!=std::string::npos || vlistItems_[i].find("rms")!=std::string::npos){
        EY[j]= 0.;
      }
      else {
        //EY[j]=vSummary_[index+1];
        EY[j]=0;// dhidas hack fix for now.  needs to be fixed
        addShift=1;
      }

      // integrate
      if (j == 0 ) YCumul[j] = Y[j]; 
      else         YCumul[j] = Y[j] + YCumul[j-1];

      // dhidas HACK for now
      EY[j] = 0;

      if(iDebug) {
        std::cout << index-j*vlistItems_.size() <<  " " << j  << " " << X[j]  << " " << Y[j] << " " << EY[j] << std::endl;
      }
    }

    C->cd(++padCount);
    gPad->SetLogy(logy);

    // Loop over all values and correct them for user defined range
    if (UseYRange != 0) {
      for (size_t iRun = 0; iRun != vRun_.size(); ++iRun) {
        if (UseYRange % 2 == 1 && Y[iRun] < YMin) {
          Y[iRun] = YMin;
          EY[iRun] = 0;
        }
        if (UseYRange >= 2  && Y[iRun] > YMax) {
          Y[iRun] = YMax;
          EY[iRun] = 0;
        }
      }
    }
    
    graph = new TGraphErrors((int) vRun_.size(),X,Y,EX,EY);
    if( fSkip99s || fSkip0s ) {
      int iptTGraph = 0;
      for (size_t ipt = 0; ipt != vRun_.size(); ++ipt) {
        // skip 99s or 0s when requested
        // std::cout << "point = " << Y[ipt] << std::endl;
        // if( Y[ipt] == 0 ) {
        //   std::cout << "fSkip0s = " << fSkip0s << std::endl;
        // }
        // if( (Y[ipt] == -10 || Y[ipt] == -9999 || Y[ipt] == -999 || Y[ipt] == -99) ) {
        //   std::cout << "fSkip99s = " << fSkip99s << std::endl;
        // }
        if( ((Y[ipt] == -10 || Y[ipt] == -9999 || Y[ipt] == -999 || Y[ipt] == -99) && fSkip99s) || (Y[ipt] == 0 && fSkip0s) ) {
          // std::cout << "removing point Y["<<ipt<<"] = " << Y[ipt] << ", when graph->GetN() = " << graph->GetN() << " and iptTGraph = " << iptTGraph << std::endl;
          // Int_t point = graph->RemovePoint(iptTGraph);
          // std::cout << "point removed = " << point << std::endl;
          graph->RemovePoint(iptTGraph);
        }
        else {
          // The TGraph is shrinked everytime a point is removed. We use another counter that
          // is increased only when not removing elements from the TGraph.
          ++iptTGraph;
        }
      }
    }
        
    graph->SetTitle(ss.str().c_str());
    if (UseYRange % 2 == 1) {
      graph->SetMinimum(YMin);
    }
    if (UseYRange >= 2) {
      graph->SetMaximum(YMax);
    }

    graph->Draw("Ap");
    graph->SetName(ss.str().c_str());
    graph->GetXaxis()->SetTitle("Run number");
    graph->Write();

    // put the graph in the vector eh.
    VectorOfGraphs.push_back(graph);

    // dhidas
    // Want to get some values into file... testing
    //for (int iDean = 0; iDean != graph.GetN(); ++iDean) {
    //  static std::ofstream OutFile("DeanOut.txt");
    //  fprintf("%9i %9i %12.3f\n", iDean, graph.GetX()[iDean], graph.GetY()[iDean]);
    //}

    if (itemForIntegration)
    {  
      std::stringstream ss2; std::stringstream ss3; std::stringstream ss4;
      std::string title =  vlistItems_.at(i);
     
      ss2 << title << "_Integral";
      ss3 << title << "_Integrated.gif";
      ss4 << title << "_Integrated.root";

      TCanvas* C2 = new TCanvas(ss2.str().c_str(),"");
      TGraphErrors* graph2 = new TGraphErrors((int) vRun_.size(),X,YCumul,EX,EX);
      graph2->SetTitle(ss2.str().c_str());
      graph2->SetMarkerColor(1);
      graph2->Draw("Ap");
      graph2->SetName(ss2.str().c_str());
      graph2->GetXaxis()->SetTitle("Run number");
      graph2->Write();
      C2->Write();
      C2->SaveAs(ss3.str().c_str());
      C2->SaveAs(ss4.str().c_str());
      // dhidas commented out below because it doesn't seem useful.  uncomment if you like, it was just annoying me.
      //C2->SaveAs(ss3.str().replace(ss3.str().find("."),ss3.str().size()-ss3.str().find("."),".C").c_str());
      }
    i+=addShift;
  }
  C->Write();
  C->SaveAs(CanvasName.c_str());
  // dhidas commented out below because it doesn't seem useful.  uncomment if you like, it was just annoying me.
  //C->SaveAs(CanvasName.replace(CanvasName.find("."),CanvasName.size()-CanvasName.find("."),".C").c_str());//savewith .C
  // dhidas commented out below because it doesn't seem useful.  uncomment if you like, it was just annoying me.
  //C->SaveAs(CanvasName.replace(CanvasName.find("."),CanvasName.size()-CanvasName.find("."),".C").c_str());//savewith .C


  // Okay, we wrote the first canvas, not let's try to overlay the graphs on another one..
  if (VectorOfGraphs.size() > 1) {

    // Create the legend for this overlay graph
    TLegend OverlayLegend(0.80,0.35,0.99,0.65);

    // Use for storing the global min/max.
    float max = -9999;
    float min =  9999;

    // Canvas we'll paint the overlays on
    TCanvas DeanCan("DeanCan", "DeanCan");
    TVirtualPad* VPad = DeanCan.cd();
    VPad->SetRightMargin(0.21);
    VPad->SetTopMargin(0.13);

    // Replace default legend names with labels if they exist
    TString const LNames = Labels;
    TObjArray* MyArrayPtr = LNames.Tokenize(",");
    if (MyArrayPtr) {
      MyArrayPtr->SetOwner(kTRUE);
      for( int i = 0; i <= MyArrayPtr->GetLast(); ++i ) {
        if( i < int(VectorOfDetNames.size()) ) {
          VectorOfDetNames[i] = ((TObjString*) MyArrayPtr->At(i) )->GetString().Data();
        }
      }
      MyArrayPtr->Delete();
    }
        

    // Let's loop over all graphs in this request
    for (size_t i = 0; i != VectorOfGraphs.size(); ++i) {

      // Strip off the det name in the i-th hist title
      TString MyTitle = VectorOfGraphs[i]->GetTitle();
      std::cout << "dhidas " << MyTitle << " : " << VectorOfDetNames[i] << std::endl;
      MyTitle.ReplaceAll(VectorOfDetNames[i]+"_", "");
      MyTitle.ReplaceAll("_"+VectorOfDetNames[i], "");
      MyTitle.ReplaceAll(VectorOfDetNames[i], "");
      std::cout << "dhidas " << MyTitle << std::endl;
      VectorOfGraphs[i]->SetTitle( MyTitle );

      // Add this to the legend, sure, good
      OverlayLegend.AddEntry(VectorOfGraphs[i], VectorOfDetNames[i].c_str(), "p");

      // You have to get the min and max by hand because root is completely retarded
      if (min > findGraphMin(VectorOfGraphs[i]) ) {
        min = findGraphMin(VectorOfGraphs[i]);
      }
      if (max < findGraphMax(VectorOfGraphs[i])) {
        max = findGraphMax(VectorOfGraphs[i]);
      }

      // let's use these colors and shapes for now
      VectorOfGraphs[i]->SetMarkerStyle(20+i);
      VectorOfGraphs[i]->SetMarkerColor(2+i);
    }
    // May as well set the min and max for first graph we'll draw
    VectorOfGraphs[0]->SetMinimum((min)-((max)-(min))/5.);
    VectorOfGraphs[0]->SetMaximum((max)+((max)-(min))/5.);
    if (UseYRange % 2 == 1) {
      VectorOfGraphs[0]->SetMinimum(YMin);
    }
    if (UseYRange >= 2) {
      VectorOfGraphs[0]->SetMaximum(YMax);
    }

    // Draw the first one with axis (A) and the rest just points (p), draw the legend, and save that canvas
    
    VectorOfGraphs[0]->Draw("Ap");
    for (size_t i = 1; i != VectorOfGraphs.size(); ++i) {
      VectorOfGraphs[i]->Draw("p");
    }
    OverlayLegend.Draw("same");
    //OverlayLegend.SetTextSize(1.5);
    DeanCan.SaveAs(CanvasName.replace(CanvasName.find("."),CanvasName.size()-CanvasName.find("."),"_Overlay.gif").c_str());
  }

  // While I'm here I may as well try deleting the graphs since people don't like to clean up after themselves
  for (size_t i = 0; i != VectorOfGraphs.size(); ++i) {
    delete VectorOfGraphs[i];
  }

  // Why do people never put a friggin return statement?
  return;

}
std::string HDQMInspector::readListFromFile ( const std::string &  listFileName)

Definition at line 204 of file HDQMInspector.cc.

References gather_cfg::cout, geometryCSVtoXML::line, and pos.

{
  std::ifstream listFile;
  listFile.open(listFileName.c_str());
  std::string listString;
  if( !listFile ) {
    std::cout << "Warning: list file" << listFileName << " not found" << std::endl;
    return listString;
  }
  while( !listFile.eof() ) {
    std::string line;
    listFile >> line;
    if( line != "" ) {
      listString += line;
      listString += ",";
    }
  }
  // Remove the last ","
  std::string::size_type pos = listString.find_last_of(",");
  if( pos != std::string::npos ) {
    listString.erase(pos);
  }
  std::cout << "whiteList = " << listString << std::endl;
  return listString;
}
void HDQMInspector::setBlackList ( std::string const &  ListItems)

Definition at line 104 of file HDQMInspector.cc.

References blackList, dtNoiseDBValidation_cfg::cerr, cmsRelvalreport::exit, i, and python::multivaluedict::sort().

{

  // Run over entire input string
  for (std::string::const_iterator Pos = ListItems.begin(); Pos != ListItems.end(); ) {

    // The rest of the string
    std::string Remainder(Pos, ListItems.end());

    // This entry will be from the beginning of the remainder to either a ","
    // or the end of the string
    std::string Entry = Remainder.substr(0, Remainder.find(","));

    // If we find a "-" we know it's a blacklist range
    if ( Entry.find("-") ) {

      // Get the first and last runs from this range
      int const FirstRun = atoi( Entry.substr(0, Entry.find("-")).c_str() );
      int const LastRun  = atoi( Entry.substr(Entry.find("-")+1).c_str() );

      // If you entered it stupidly we're going to stop here.
      if (FirstRun > LastRun) {
        std::cerr << "ERROR: FirstRun > LastRun in blackList" << std::endl;
        exit(1);
      }

      // For now the simplest thing to do is fill in gaps including each end
      for (int i = FirstRun; i <= LastRun; ++i) {
        blackList.push_back(i);
      }

    } else {
      // If we didn't see a "-" just add it to the list
      blackList.push_back( atoi(Entry.c_str()) );
    }

    // This is to make sure we are in the correct position as we go on.
    Pos += Entry.size();
    if (Pos != ListItems.end()) {
      Pos += 1;
    }

  }

  // sort the list for faster searching later
  std::sort(blackList.begin(), blackList.end());

  return;
}
void HDQMInspector::setDB ( const std::string &  DBName,
const std::string &  DBTag,
const std::string &  DBauth = "" 
)

Definition at line 58 of file HDQMInspector.cc.

References accessDB(), dtNoiseDBValidation_cfg::cerr, gather_cfg::cout, DBauth_, DBName_, DBTag_, cmsRelvalreport::exit, and fOutFile.

{
  if( DBName_==DBName && DBTag_==DBTag && DBauth_ == DBauth)
    return;

  DBName_= DBName;
  DBTag_= DBTag;
  DBauth_ = DBauth;

  std::cout << "Name of DB = "<< DBName << std::endl;
  std::cout << "DBTag = "<< DBTag << std::endl;
  std::cout << "DBauth = "<< DBauth << std::endl;
  std::cout <<std::endl;

  accessDB();
  
  fOutFile = new TFile( "historicDQM.root","RECREATE" );
  if (!fOutFile->IsOpen()) {
    std::cerr << "ERROR: cannot open output file" << std::endl;
    exit(1);
  }
  fOutFile->cd();
}
void HDQMInspector::setDebug ( int  i) [inline]

Definition at line 73 of file HDQMInspector.h.

References i, and iDebug.

{iDebug=i;}
void HDQMInspector::setDoStat ( int  i) [inline]

Definition at line 74 of file HDQMInspector.h.

References i, and iDoStat.

void HDQMInspector::setItems ( std::string  itemD) [private]

Definition at line 788 of file HDQMInspector.cc.

References gather_cfg::cout, HDQMInspector::DetIdItemList::detid, fSep, iDebug, HDQMInspector::DetIdItemList::items, and vDetIdItemList_.

Referenced by unpackItems().

{
  DetIdItemList detiditemlist;
  detiditemlist.detid=atol(itemD.substr(0,itemD.find(fSep)).c_str());

  std::string item=itemD.substr(itemD.find(fSep)+fSep.size());
  detiditemlist.items.push_back(item);

  if(iDebug)
    std::cout << "Found new item " << detiditemlist.items.back() << " for detid " << detiditemlist.detid << std::endl;

  if(item.find("mean")!=std::string::npos){
    detiditemlist.items.push_back(item.replace(item.find("mean"),4,"rms")); 
    if(iDebug)
      std::cout << "Found new item " << detiditemlist.items.back() << std::endl;
    detiditemlist.items.push_back(item.replace(item.find("rms"),3,"entries")); 
    if(iDebug)
      std::cout << "Found new item " << detiditemlist.items.back() << std::endl;
  }
  else if(item.find("landauPeak")!=std::string::npos){
    detiditemlist.items.push_back(item.replace(item.find("landauPeak"),10,"landauPeakErr")); 
    if(iDebug)
      std::cout << "Found new item " << detiditemlist.items.back() << std::endl;
  }
  else if(item.find("gaussMean")!=std::string::npos){
    detiditemlist.items.push_back(item.replace(item.find("gaussMean"),9,"gaussSigma")); 
    if(iDebug)
      std::cout << "Found new item " << detiditemlist.items.back() << std::endl;
  }

  if(vDetIdItemList_.size()) {
    if(vDetIdItemList_.back().detid==detiditemlist.detid) {
      vDetIdItemList_.back().items.insert(vDetIdItemList_.back().items.end(),detiditemlist.items.begin(),detiditemlist.items.end());
    } else {
      vDetIdItemList_.push_back(detiditemlist);
    }
  } else {
    vDetIdItemList_.push_back(detiditemlist);
  }

  return;
}
bool HDQMInspector::setRange ( unsigned int &  firstRun,
unsigned int &  lastRun 
) [private]

Definition at line 260 of file HDQMInspector.cc.

References gather_cfg::cout, first, dtTTrigAnalyzer_cfg::firstRun, i, iDebug, iovList, prof2calltree::last, and MergeJob_cfg::lastRun.

Referenced by createTrend().

{
  unsigned int first,last;

  for(size_t i=0;i<iovList.size();++i) {
    if (iDebug) {
      std::cout << iovList.at(i)<< std::endl;
    }
  }

  std::vector<unsigned int>::iterator iter;

  iter=std::lower_bound(iovList.begin(),iovList.end(),firstRun);
  if (iter!=iovList.end())
    first=*iter;
  else{
    std::cout << "firstRun (" << firstRun << ") > last iov ("<<iovList.back()<< ")"<<std::endl; 
    return false;
  }

  iter=std::lower_bound(iovList.begin(),iovList.end(),lastRun);
  if (iter!=iovList.end()){
    if (*iter>lastRun) last = *(iter-1);
    else last=*iter;
  }
  else{
    last=iovList.back();
  }
  
  firstRun=first;
  lastRun=last; 
  std::cout << "setting Range firstRun (" << first << ") - lastRun ("<<last<< ")"<<std::endl; 
  Iterator->setRange(first,last);
  
  return true;
}
void HDQMInspector::setSeparator ( std::string const  in) [inline]

Definition at line 94 of file HDQMInspector.h.

References fSep, and recoMuon::in.

                                        {
    fSep = in;
    return;
  }
void HDQMInspector::setSkip0s ( bool const  in) [inline]

Definition at line 82 of file HDQMInspector.h.

References fSkip0s, and recoMuon::in.

                                 {
    fSkip0s = in;
    return;
  }
void HDQMInspector::setSkip99s ( bool const  in) [inline]

Definition at line 78 of file HDQMInspector.h.

References fSkip99s, and recoMuon::in.

                                  {
    fSkip99s = in;
    return;
  }
void HDQMInspector::setWhiteList ( std::string const &  ListItems)

Definition at line 154 of file HDQMInspector.cc.

References dtNoiseDBValidation_cfg::cerr, cmsRelvalreport::exit, i, python::multivaluedict::sort(), and whiteList.

{

  // Run over entire input string
  for (std::string::const_iterator Pos = ListItems.begin(); Pos != ListItems.end(); ) {

    // The rest of the string
    std::string Remainder(Pos, ListItems.end());

    // This entry will be from the beginning of the remainder to either a ","
    // or the end of the string
    std::string Entry = Remainder.substr(0, Remainder.find(","));

    // If we find a "-" we know it's a whitelist range
    if ( Entry.find("-") ) {

      // Get the first and last runs from this range
      int const FirstRun = atoi( Entry.substr(0, Entry.find("-")).c_str() );
      int const LastRun  = atoi( Entry.substr(Entry.find("-")+1).c_str() );

      // If you entered it stupidly we're going to stop here.
      if (FirstRun > LastRun) {
        std::cerr << "ERROR: FirstRun > LastRun in WhiteList" << std::endl;
        exit(1);
      }

      // For now the simplest thing to do is fill in gaps including each end
      for (int i = FirstRun; i <= LastRun; ++i) {
        whiteList.push_back(i);
      }

    } else {
      // If we didn't see a "-" just add it to the list
      whiteList.push_back( atoi(Entry.c_str()) );
    }

    // This is to make sure we are in the correct position as we go on.
    Pos += Entry.size();
    if (Pos != ListItems.end()) {
      Pos += 1;
    }

  }

  // sort the list for faster searching later
  std::sort(whiteList.begin(), whiteList.end());

  return;
}
void HDQMInspector::style ( ) [private]

Definition at line 30 of file HDQMInspector.cc.

Referenced by plot().

{
  TStyle* theStyle= new TStyle();
  theStyle->SetOptStat(0);
  //gROOT->SetStyle("Plain");
  theStyle->SetOptStat(0);
  theStyle->SetOptFit(111);
  theStyle->SetStatFont(12);
  theStyle->SetStatBorderSize(1);
  theStyle->SetCanvasColor(0);
  theStyle->SetCanvasBorderMode(0);
  theStyle->SetPadBorderMode(0);
  theStyle->SetPadColor(0);
  theStyle->SetLineWidth(1);
  theStyle->SetLineStyle(2);
  theStyle->SetPalette(1);
  theStyle->SetMarkerStyle(20);
  theStyle->SetMarkerColor(2);
  theStyle->SetLabelSize(0.05,"y");
  theStyle->SetLabelSize(0.04,"x");
  theStyle->SetTitleFontSize(0.2);
  theStyle->SetTitleW(0.9);
  theStyle->SetTitleH(0.06);
  theStyle->SetPadLeftMargin(0.12);   
  theStyle->SetPadTopMargin(0.13);
  theStyle->cd();
}
void HDQMInspector::unpackConditions ( std::string &  Conditions,
std::vector< DetIdItemList > &  vdetIdItemList 
) [private]

Definition at line 707 of file HDQMInspector.cc.

References gather_cfg::cout, HDQMInspector::DetIdItemList::detid, fSep, iDebug, HDQMInspector::DetIdItemList::items, and NULL.

Referenced by createTrend().

{
  char * pch;
  char delimiters[128]="><=+-*/&|() ";
  char copyConditions[1024];
  sprintf(copyConditions,"%s",Conditions.c_str());
  pch = strtok (copyConditions,delimiters);
  while (pch != NULL){
    if(strstr(pch,fSep.c_str())!=NULL){
      DetIdItemList detiditemlist;
      std::string itemD(pch);
      detiditemlist.detid=atol(itemD.substr(0,itemD.find(fSep)).c_str());
      detiditemlist.items.push_back(itemD.substr(itemD.find(fSep)+fSep.size())); // dhidas update +.size instead of "1"
      if (iDebug) {
        std::cout << "Found a Condition " << detiditemlist.items.back() << " for detId " << detiditemlist.detid << std::endl;
      }
      
      if(vdetIdItemList.size())
        if(vdetIdItemList.back().detid==detiditemlist.detid)
          vdetIdItemList.back().items.insert(vdetIdItemList.back().items.end(),detiditemlist.items.begin(),detiditemlist.items.end());
        else
          vdetIdItemList.push_back(detiditemlist);
      else
        vdetIdItemList.push_back(detiditemlist); 
    }
    pch = strtok (NULL,delimiters);
  }
}
size_t HDQMInspector::unpackItems ( std::string &  ListItems) [private]

Definition at line 690 of file HDQMInspector.cc.

References prof2calltree::count, gather_cfg::cout, and setItems().

Referenced by createTrend().

{
  std::string::size_type oldloc=0; 
  std::string::size_type loc = ListItems.find( ",", oldloc );
  size_t count=1;
  while( loc != std::string::npos ) {
    setItems(ListItems.substr(oldloc,loc-oldloc));
    oldloc=loc+1;
    loc=ListItems.find( ",", oldloc );
    count++; 
  } 
  //there is a single item
  setItems(ListItems.substr(oldloc,loc-oldloc));
  std::cout << std::endl;
  return count;
}

Member Data Documentation

std::vector<unsigned int> HDQMInspector::blackList [private]

Definition at line 125 of file HDQMInspector.h.

Referenced by createTrend(), and setBlackList().

std::string HDQMInspector::DBauth_ [private]

Definition at line 120 of file HDQMInspector.h.

Referenced by accessDB(), and setDB().

std::string HDQMInspector::DBName_ [private]

Definition at line 120 of file HDQMInspector.h.

Referenced by accessDB(), and setDB().

std::string HDQMInspector::DBTag_ [private]

Definition at line 120 of file HDQMInspector.h.

Referenced by accessDB(), and setDB().

Definition at line 139 of file HDQMInspector.h.

Referenced by plot().

Definition at line 144 of file HDQMInspector.h.

Referenced by closeFile(), and setDB().

std::string HDQMInspector::fSep [private]

Definition at line 141 of file HDQMInspector.h.

Referenced by ApplyConditions(), plot(), setItems(), setSeparator(), and unpackConditions().

bool HDQMInspector::fSkip0s [private]

Definition at line 137 of file HDQMInspector.h.

Referenced by plot(), and setSkip0s().

bool HDQMInspector::fSkip99s [private]

Definition at line 136 of file HDQMInspector.h.

Referenced by plot(), and setSkip99s().

int HDQMInspector::iDebug [private]
int HDQMInspector::iDoStat [private]

Definition at line 135 of file HDQMInspector.h.

Referenced by setDoStat().

std::vector<unsigned int> HDQMInspector::iovList [private]

Definition at line 124 of file HDQMInspector.h.

Referenced by createTrendLastRuns(), InitializeIOVList(), and setRange().

Definition at line 122 of file HDQMInspector.h.

Referenced by accessDB(), and ~HDQMInspector().

std::vector<unsigned int> HDQMInspector::vdetId_ [private]

Definition at line 132 of file HDQMInspector.h.

Referenced by getvDetId(), and plot().

Definition at line 130 of file HDQMInspector.h.

Referenced by createTrend(), plot(), and setItems().

std::vector<std::string> HDQMInspector::vlistItems_ [private]

Definition at line 131 of file HDQMInspector.h.

Referenced by getListItems(), and plot().

std::vector<unsigned int> HDQMInspector::vRun_ [private]

Definition at line 128 of file HDQMInspector.h.

Referenced by createTrend(), getRuns(), and plot().

std::vector<float> HDQMInspector::vSummary_ [private]

Definition at line 129 of file HDQMInspector.h.

Referenced by createTrend(), getSummary(), and plot().

std::vector<unsigned int> HDQMInspector::whiteList [private]

Definition at line 126 of file HDQMInspector.h.

Referenced by createTrend(), and setWhiteList().