CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/DQMServices/ClientConfig/plugins/DQMGenericClient.cc

Go to the documentation of this file.
00001 /*
00002  *  Class:DQMGenericClient 
00003  *
00004  *
00005  *  $Date: 2011/07/15 08:33:35 $
00006  *  $Revision: 1.28 $
00007  * 
00008  *  \author Junghwan Goh - SungKyunKwan University
00009  */
00010 
00011 #include "DQMServices/ClientConfig/interface/DQMGenericClient.h"
00012 
00013 #include "DQMServices/ClientConfig/interface/FitSlicesYTool.h"
00014 #include "DQMServices/Core/interface/DQMStore.h"
00015 #include "DQMServices/Core/interface/MonitorElement.h"
00016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00017 #include "FWCore/ServiceRegistry/interface/Service.h"
00018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00019 
00020 #include <TH1F.h>
00021 #include <TClass.h>
00022 #include <TString.h>
00023 #include <TPRegexp.h>
00024 
00025 #include <cmath>
00026 #include <boost/tokenizer.hpp>
00027 
00028 using namespace std;
00029 using namespace edm;
00030 
00031 typedef MonitorElement ME;
00032 
00033 TPRegexp metacharacters("[\\^\\$\\.\\*\\+\\?\\|\\(\\)\\{\\}\\[\\]]");
00034 TPRegexp nonPerlWildcard("\\w\\*|^\\*");
00035 
00036 DQMGenericClient::DQMGenericClient(const ParameterSet& pset)
00037 {
00038   typedef std::vector<edm::ParameterSet> VPSet;
00039   typedef std::vector<std::string> vstring;
00040   typedef boost::escaped_list_separator<char> elsc;
00041 
00042   elsc commonEscapes("\\", " \t", "\'");
00043 
00044   verbose_ = pset.getUntrackedParameter<unsigned int>("verbose", 0);
00045 
00046   // Parse efficiency commands
00047   vstring effCmds = pset.getParameter<vstring>("efficiency");
00048   for ( vstring::const_iterator effCmd = effCmds.begin();
00049         effCmd != effCmds.end(); ++effCmd )
00050   {
00051     if ( effCmd->empty() ) continue;
00052 
00053     boost::tokenizer<elsc> tokens(*effCmd, commonEscapes);
00054 
00055     vector<string> args;
00056     for(boost::tokenizer<elsc>::const_iterator iToken = tokens.begin();
00057         iToken != tokens.end(); ++iToken) {
00058       if ( iToken->empty() ) continue;
00059       args.push_back(*iToken);
00060     }
00061 
00062     if ( args.size() < 4 ) {
00063       LogInfo("DQMGenericClient") << "Wrong input to effCmds\n";
00064       continue;
00065     }
00066 
00067     EfficOption opt;
00068     opt.name = args[0];
00069     opt.title = args[1];
00070     opt.numerator = args[2];
00071     opt.denominator = args[3];
00072     opt.isProfile = false;
00073 
00074     const string typeName = args.size() == 4 ? "eff" : args[4];
00075     if ( typeName == "eff" ) opt.type = 1;
00076     else if ( typeName == "fake" ) opt.type = 2;
00077     else opt.type = 0;
00078  
00079     efficOptions_.push_back(opt);
00080   }
00081   
00082   VPSet efficSets = pset.getUntrackedParameter<VPSet>("efficiencySets", VPSet());
00083   for ( VPSet::const_iterator efficSet = efficSets.begin();
00084         efficSet != efficSets.end(); ++efficSet )
00085   {
00086     EfficOption opt;
00087     opt.name = efficSet->getUntrackedParameter<string>("name");
00088     opt.title = efficSet->getUntrackedParameter<string>("title");
00089     opt.numerator = efficSet->getUntrackedParameter<string>("numerator");
00090     opt.denominator = efficSet->getUntrackedParameter<string>("denominator");
00091     opt.isProfile = false;
00092 
00093     const string typeName = efficSet->getUntrackedParameter<string>("typeName", "eff");
00094     if ( typeName == "eff" ) opt.type = 1;
00095     else if ( typeName == "fake" ) opt.type = 2;
00096     else opt.type = 0;
00097 
00098     efficOptions_.push_back(opt);
00099   }
00100 
00101   // Parse profiles
00102   vstring profileCmds = pset.getUntrackedParameter<vstring>("efficiencyProfile", vstring());
00103   for ( vstring::const_iterator profileCmd = profileCmds.begin();
00104         profileCmd != profileCmds.end(); ++profileCmd )
00105   {
00106     if ( profileCmd->empty() ) continue;
00107 
00108     boost::tokenizer<elsc> tokens(*profileCmd, commonEscapes);
00109 
00110     vector<string> args;
00111     for(boost::tokenizer<elsc>::const_iterator iToken = tokens.begin();
00112         iToken != tokens.end(); ++iToken) {
00113       if ( iToken->empty() ) continue;
00114       args.push_back(*iToken);
00115     }
00116 
00117     if ( args.size() < 4 ) {
00118       LogInfo("DQMGenericClient") << "Wrong input to profileCmds\n";
00119       continue;
00120     }
00121 
00122     EfficOption opt;
00123     opt.name = args[0];
00124     opt.title = args[1];
00125     opt.numerator = args[2];
00126     opt.denominator = args[3];
00127     opt.isProfile = true;
00128 
00129     const string typeName = args.size() == 4 ? "eff" : args[4];
00130     if ( typeName == "eff" ) opt.type = 1;
00131     else if ( typeName == "fake" ) opt.type = 2;
00132     else opt.type = 0;
00133  
00134     efficOptions_.push_back(opt);
00135   }
00136 
00137   VPSet profileSets = pset.getUntrackedParameter<VPSet>("efficiencyProfileSets", VPSet());
00138   for ( VPSet::const_iterator profileSet = profileSets.begin();
00139         profileSet != profileSets.end(); ++profileSet )
00140   {
00141     EfficOption opt;
00142     opt.name = profileSet->getUntrackedParameter<string>("name");
00143     opt.title = profileSet->getUntrackedParameter<string>("title");
00144     opt.numerator = profileSet->getUntrackedParameter<string>("numerator");
00145     opt.denominator = profileSet->getUntrackedParameter<string>("denominator");
00146     opt.isProfile = true;
00147 
00148     const string typeName = profileSet->getUntrackedParameter<string>("typeName", "eff");
00149     if ( typeName == "eff" ) opt.type = 1;
00150     else if ( typeName == "fake" ) opt.type = 2;
00151     else opt.type = 0;
00152 
00153     efficOptions_.push_back(opt);
00154   }
00155 
00156   // Parse resolution commands
00157   vstring resCmds = pset.getParameter<vstring>("resolution");
00158   for ( vstring::const_iterator resCmd = resCmds.begin();
00159         resCmd != resCmds.end(); ++resCmd )
00160   {
00161     if ( resCmd->empty() ) continue;
00162     boost::tokenizer<elsc> tokens(*resCmd, commonEscapes);
00163 
00164     vector<string> args;
00165     for(boost::tokenizer<elsc>::const_iterator iToken = tokens.begin();
00166         iToken != tokens.end(); ++iToken) {
00167       if ( iToken->empty() ) continue;
00168       args.push_back(*iToken);
00169     }
00170 
00171     if ( args.size() != 3 ) {
00172       LogInfo("DQMGenericClient") << "Wrong input to resCmds\n";
00173       continue;
00174     }
00175 
00176     ResolOption opt;
00177     opt.namePrefix = args[0];
00178     opt.titlePrefix = args[1];
00179     opt.srcName = args[2];
00180 
00181     resolOptions_.push_back(opt);
00182   }
00183 
00184   VPSet resolSets = pset.getUntrackedParameter<VPSet>("resolutionSets", VPSet());
00185   for ( VPSet::const_iterator resolSet = resolSets.begin();
00186         resolSet != resolSets.end(); ++resolSet )
00187   {
00188     ResolOption opt;
00189     opt.namePrefix = resolSet->getUntrackedParameter<string>("namePrefix");
00190     opt.titlePrefix = resolSet->getUntrackedParameter<string>("titlePrefix");
00191     opt.srcName = resolSet->getUntrackedParameter<string>("srcName");
00192 
00193     resolOptions_.push_back(opt);
00194   }
00195 
00196   // Parse Normalization commands
00197   vstring normCmds = pset.getUntrackedParameter<vstring>("normalization", vstring());
00198   for ( vstring::const_iterator normCmd = normCmds.begin();
00199         normCmd != normCmds.end(); ++normCmd )
00200   {
00201     if ( normCmd->empty() ) continue;
00202     boost::tokenizer<elsc> tokens(*normCmd, commonEscapes);
00203 
00204     vector<string> args;
00205     for(boost::tokenizer<elsc>::const_iterator iToken = tokens.begin();
00206         iToken != tokens.end(); ++iToken) {
00207       if ( iToken->empty() ) continue;
00208       args.push_back(*iToken);
00209     }
00210 
00211     if ( args.empty() or args.size() > 2 ) {
00212       LogInfo("DQMGenericClient") << "Wrong input to normCmds\n";
00213       continue;
00214     }
00215 
00216     NormOption opt;
00217     opt.name = args[0];
00218     opt.normHistName = args.size() == 2 ? args[1] : args[0];
00219 
00220     normOptions_.push_back(opt);
00221   }
00222 
00223   VPSet normSets = pset.getUntrackedParameter<VPSet>("normalizationSets", VPSet());
00224   for ( VPSet::const_iterator normSet = normSets.begin();
00225         normSet != normSets.end(); ++normSet )
00226   {
00227     NormOption opt;
00228     opt.name = normSet->getUntrackedParameter<string>("name");
00229     opt.normHistName = normSet->getUntrackedParameter<string>("normalizedTo", opt.name);
00230 
00231     normOptions_.push_back(opt);
00232   }
00233 
00234   // Cumulative distributions
00235   vstring cdCmds = pset.getUntrackedParameter<vstring>("cumulativeDists", vstring());
00236   for ( vstring::const_iterator cdCmd = cdCmds.begin();
00237         cdCmd != cdCmds.end(); ++cdCmd )
00238   {
00239     if ( cdCmd->empty() ) continue;
00240     boost::tokenizer<elsc> tokens(*cdCmd, commonEscapes);
00241 
00242     vector<string> args;
00243     for(boost::tokenizer<elsc>::const_iterator iToken = tokens.begin();
00244         iToken != tokens.end(); ++iToken) {
00245       if ( iToken->empty() ) continue;
00246       args.push_back(*iToken);
00247     }
00248 
00249     if ( args.size() != 1 ) {
00250       LogInfo("DQMGenericClient") << "Wrong input to cdCmds\n";
00251       continue;
00252     }
00253 
00254     CDOption opt;
00255     opt.name = args[0];
00256 
00257     cdOptions_.push_back(opt);
00258   }
00259 
00260   VPSet cdSets = pset.getUntrackedParameter<VPSet>("cumulativeDistSets", VPSet());
00261   for ( VPSet::const_iterator cdSet = cdSets.begin();
00262         cdSet != cdSets.end(); ++cdSet )
00263   {
00264     CDOption opt;
00265     opt.name = cdSet->getUntrackedParameter<string>("name");
00266 
00267     cdOptions_.push_back(opt);
00268   }
00269 
00270   outputFileName_ = pset.getUntrackedParameter<string>("outputFileName", "");
00271   subDirs_ = pset.getUntrackedParameter<vstring>("subDirs");
00272 
00273   resLimitedFit_ = pset.getUntrackedParameter<bool>("resolutionLimitedFit",false);
00274   isWildcardUsed_ = false;
00275 }
00276 
00277 void DQMGenericClient::endRun(const edm::Run& r, const edm::EventSetup& c) {
00278 
00279   typedef vector<string> vstring;
00280 
00281   // Update 2009-09-23
00282   // Migrated all code from endJob to this function
00283   // endJob is not necessarily called in the proper sequence
00284   // and does not necessarily book histograms produced in
00285   // that step.
00286   // It more robust to do the histogram manipulation in
00287   // this endRun function
00288 
00289 
00290   
00291   theDQM = 0;
00292   theDQM = Service<DQMStore>().operator->();
00293 
00294   if ( ! theDQM ) {
00295     LogInfo("DQMGenericClient") << "Cannot create DQMStore instance\n";
00296     return;
00297   }
00298 
00299   // Process wildcard in the sub-directory
00300   set<string> subDirSet;
00301 
00302   for(vstring::const_iterator iSubDir = subDirs_.begin();
00303       iSubDir != subDirs_.end(); ++iSubDir) {
00304     string subDir = *iSubDir;
00305 
00306     if ( subDir[subDir.size()-1] == '/' ) subDir.erase(subDir.size()-1);
00307 
00308     if ( TString(subDir).Contains(metacharacters) ) {
00309       isWildcardUsed_ = true;
00310 
00311       const string::size_type shiftPos = subDir.rfind('/');
00312       const string searchPath = subDir.substr(0, shiftPos);
00313       const string pattern    = subDir.substr(shiftPos + 1, subDir.length());
00314       //std::cout << "\n\n\n\nLooking for all subdirs of " << subDir << std::endl;
00315       
00316       findAllSubdirectories (searchPath, &subDirSet, pattern);
00317 
00318     }
00319     else {
00320       subDirSet.insert(subDir);
00321     }
00322   }
00323 
00324   for(set<string>::const_iterator iSubDir = subDirSet.begin();
00325       iSubDir != subDirSet.end(); ++iSubDir) {
00326     const string& dirName = *iSubDir;
00327 
00328     for ( vector<EfficOption>::const_iterator efficOption = efficOptions_.begin();
00329           efficOption != efficOptions_.end(); ++efficOption )
00330     {
00331       computeEfficiency(dirName, efficOption->name, efficOption->title, 
00332                         efficOption->numerator, efficOption->denominator,
00333                         efficOption->type, efficOption->isProfile);
00334     }
00335 
00336     for ( vector<ResolOption>::const_iterator resolOption = resolOptions_.begin();
00337           resolOption != resolOptions_.end(); ++resolOption )
00338     {
00339       computeResolution(dirName, resolOption->namePrefix, resolOption->titlePrefix, resolOption->srcName);
00340     }
00341 
00342     for ( vector<NormOption>::const_iterator normOption = normOptions_.begin();
00343           normOption != normOptions_.end(); ++normOption )
00344     {
00345       normalizeToEntries(dirName, normOption->name, normOption->normHistName);
00346     }
00347 
00348     for ( vector<CDOption>::const_iterator cdOption = cdOptions_.begin();
00349           cdOption != cdOptions_.end(); ++cdOption )
00350     {
00351       makeCumulativeDist(dirName, cdOption->name);
00352     }
00353   }
00354 
00355   //if ( verbose_ > 0 ) theDQM->showDirStructure();
00356 
00357   if ( ! outputFileName_.empty() ) theDQM->save(outputFileName_);
00358   
00359 }
00360 
00361 void DQMGenericClient::endJob()
00362 {
00363 
00364   // Update 2009-09-23
00365   // Migrated all code from here to endRun
00366 
00367   LogTrace ("DQMGenericClient") << "inside of DQMGenericClient::endJob()"
00368                                 << endl;
00369 
00370 }
00371 
00372 void DQMGenericClient::computeEfficiency(const string& startDir, const string& efficMEName, const string& efficMETitle, 
00373                                          const string& recoMEName, const string& simMEName, const int type, const bool makeProfile)
00374 {
00375   if ( ! theDQM->dirExists(startDir) ) {
00376     if ( verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_) ) {
00377       LogInfo("DQMGenericClient") << "computeEfficiency() : "
00378                                   << "Cannot find sub-directory " << startDir << endl; 
00379     }
00380     return;
00381   }
00382 
00383   theDQM->cd();
00384 
00385   ME* simME  = theDQM->get(startDir+"/"+simMEName);
00386   ME* recoME = theDQM->get(startDir+"/"+recoMEName);
00387 
00388   if ( !simME ) {
00389     if ( verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_) ) {
00390       LogInfo("DQMGenericClient") << "computeEfficiency() : "
00391                                   << "No sim-ME '" << simMEName << "' found\n";
00392     }
00393     return;
00394   }
00395 
00396   if ( !recoME ) {
00397     if ( verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_) ) {
00398       LogInfo("DQMGenericClient") << "computeEfficiency() : " 
00399                                   << "No reco-ME '" << recoMEName << "' found\n";
00400     }
00401     return;
00402   }
00403 
00404   // Treat everything as the base class, TH1
00405   
00406   TH1* hSim  = simME ->getTH1();
00407   TH1* hReco = recoME->getTH1();
00408   
00409   if ( !hSim || !hReco ) {
00410     if ( verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_) ) {
00411       LogInfo("DQMGenericClient") << "computeEfficiency() : "
00412                                   << "Cannot create TH1 from ME\n";
00413     }
00414     return;
00415   }
00416 
00417   string efficDir = startDir;
00418   string newEfficMEName = efficMEName;
00419   string::size_type shiftPos;
00420   if ( string::npos != (shiftPos = efficMEName.rfind('/')) ) {
00421     efficDir += "/"+efficMEName.substr(0, shiftPos);
00422     newEfficMEName.erase(0, shiftPos+1);
00423   }
00424   theDQM->setCurrentFolder(efficDir);
00425 
00426   if (makeProfile) {
00427     TProfile * efficHist = (hReco->GetXaxis()->GetXbins()->GetSize()==0) ?
00428       new TProfile(newEfficMEName.c_str(), efficMETitle.c_str(),
00429                    hReco->GetXaxis()->GetNbins(),
00430                    hReco->GetXaxis()->GetXmin(),
00431                    hReco->GetXaxis()->GetXmax()) :
00432       new TProfile(newEfficMEName.c_str(), efficMETitle.c_str(),
00433                    hReco->GetXaxis()->GetNbins(),
00434                    hReco->GetXaxis()->GetXbins()->GetArray());
00435 #if ROOT_VERSION_CODE >= ROOT_VERSION(5,27,0)
00436     for (int i=1; i <= hReco->GetNbinsX(); i++) {
00437       const double nReco = hReco->GetBinContent(i);
00438       const double nSim = hSim->GetBinContent(i);
00439     
00440       if ( nSim == 0 || nReco > nSim ) continue;
00441       const double effVal = nReco/nSim;
00442 
00443       const double errLo = TEfficiency::ClopperPearson((int)hSim->GetBinContent(i), 
00444                                                        (int)hReco->GetBinContent(i),
00445                                                        0.683,false);
00446       const double errUp = TEfficiency::ClopperPearson((int)hSim->GetBinContent(i), 
00447                                                        (int)hReco->GetBinContent(i),
00448                                                        0.683,true);
00449       const double errVal = (effVal - errLo > errUp - effVal) ? effVal - errLo : errUp - effVal;
00450       efficHist->SetBinContent(i, effVal);
00451       efficHist->SetBinEntries(i, 1);
00452       efficHist->SetBinError(i, sqrt(effVal * effVal + errVal * errVal));
00453     }
00454 #else
00455     for (int i=1; i <= hReco->GetNbinsX(); i++) {
00456       TGraphAsymmErrorsWrapper asymm;
00457       std::pair<double, double> efficiencyWithError;
00458       efficiencyWithError = asymm.efficiency((int)hReco->GetBinContent(i), 
00459                                              (int)hSim->GetBinContent(i));
00460       double effVal = efficiencyWithError.first;
00461       double errVal = efficiencyWithError.second;
00462       if ((int)hSim->GetBinContent(i) > 0) {
00463         efficHist->SetBinContent(i, effVal);
00464         efficHist->SetBinEntries(i, 1);
00465         efficHist->SetBinError(i, sqrt(effVal * effVal + errVal * errVal));
00466       }
00467     }
00468 #endif
00469     theDQM->bookProfile(newEfficMEName.c_str(),efficHist);
00470     delete efficHist;  
00471   }
00472 
00473   else {
00474 
00475     TH1* efficHist = (TH1*)hSim->Clone(newEfficMEName.c_str());
00476     efficHist->SetTitle(efficMETitle.c_str());
00477 
00478     // Here is where you have trouble --- you need
00479     // to understand what type of hist you have.
00480 
00481     ME* efficME = 0;
00482 
00483     // Parse the class name
00484     // This works, but there might be a better way
00485     TClass * myHistClass = efficHist->IsA();
00486     TString histClassName = myHistClass->GetName();
00487   
00488     if (histClassName == "TH1F"){
00489       efficME = theDQM->book1D(newEfficMEName, (TH1F*)efficHist);
00490     } else if (histClassName == "TH2F"){
00491       efficME = theDQM->book2D(newEfficMEName, (TH2F*)efficHist);    
00492     } else if (histClassName == "TH3F"){
00493       efficME = theDQM->book3D(newEfficMEName, (TH3F*)efficHist);    
00494     } 
00495   
00496 
00497     if ( !efficME ) {
00498       LogInfo("DQMGenericClient") << "computeEfficiency() : "
00499                                   << "Cannot book effic-ME from the DQM\n";
00500       return;
00501     }
00502 
00503     // Update: 2009-9-16 slaunwhj
00504     // call the most generic efficiency function
00505     // works up to 3-d histograms
00506 
00507     generic_eff (hSim, hReco, efficME, type);
00508   
00509     //   const int nBin = efficME->getNbinsX();
00510     //   for(int bin = 0; bin <= nBin; ++bin) {
00511     //     const float nSim  = simME ->getBinContent(bin);
00512     //     const float nReco = recoME->getBinContent(bin);
00513     //     float eff =0;
00514     //     if (type=="fake")eff = nSim ? 1-nReco/nSim : 0.;
00515     //     else eff= nSim ? nReco/nSim : 0.;
00516     //     const float err = nSim && eff <= 1 ? sqrt(eff*(1-eff)/nSim) : 0.;
00517     //     efficME->setBinContent(bin, eff);
00518     //     efficME->setBinError(bin, err);
00519     //   }
00520     efficME->setEntries(simME->getEntries());
00521 
00522   }
00523 
00524   // Global efficiency
00525   ME* globalEfficME = theDQM->get(efficDir+"/globalEfficiencies");
00526   if ( !globalEfficME ) globalEfficME = theDQM->book1D("globalEfficiencies", "Global efficiencies", 1, 0, 1);
00527   if ( !globalEfficME ) {
00528     LogInfo("DQMGenericClient") << "computeEfficiency() : "
00529                               << "Cannot book globalEffic-ME from the DQM\n";
00530     return;
00531   }
00532   TH1F* hGlobalEffic = globalEfficME->getTH1F();
00533   if ( !hGlobalEffic ) {
00534     LogInfo("DQMGenericClient") << "computeEfficiency() : "
00535                               << "Cannot create TH1F from ME, globalEfficME\n";
00536     return;
00537   }
00538 
00539   const float nSimAll = hSim->GetEntries();
00540   const float nRecoAll = hReco->GetEntries();
00541   float efficAll=0; 
00542   if ( type == 1 ) efficAll = nSimAll ? nRecoAll/nSimAll : 0;
00543   else if ( type == 2 ) efficAll = nSimAll ? 1-nRecoAll/nSimAll : 0;
00544   const float errorAll = nSimAll && efficAll < 1 ? sqrt(efficAll*(1-efficAll)/nSimAll) : 0;
00545 
00546   const int iBin = hGlobalEffic->Fill(newEfficMEName.c_str(), 0);
00547   hGlobalEffic->SetBinContent(iBin, efficAll);
00548   hGlobalEffic->SetBinError(iBin, errorAll);
00549 }
00550 
00551 void DQMGenericClient::computeResolution(const string& startDir, const string& namePrefix, const string& titlePrefix,
00552                                          const std::string& srcName)
00553 {
00554   if ( ! theDQM->dirExists(startDir) ) {
00555     if ( verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_) ) {
00556       LogInfo("DQMGenericClient") << "computeResolution() : "
00557                                   << "Cannot find sub-directory " << startDir << endl;
00558     }
00559     return;
00560   }
00561 
00562   theDQM->cd();
00563 
00564   ME* srcME = theDQM->get(startDir+"/"+srcName);
00565   if ( !srcME ) {
00566     if ( verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_) ) {
00567       LogInfo("DQMGenericClient") << "computeResolution() : "
00568                                   << "No source ME '" << srcName << "' found\n";
00569     }
00570     return;
00571   }
00572 
00573   TH2F* hSrc = srcME->getTH2F();
00574   if ( !hSrc ) {
00575     if ( verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_) ) {
00576       LogInfo("DQMGenericClient") << "computeResolution() : "
00577                                   << "Cannot create TH2F from source-ME\n";
00578     }
00579     return;
00580   }
00581 
00582   const int nBin = hSrc->GetNbinsX();
00583   const double xMin = hSrc->GetXaxis()->GetXmin();
00584   const double xMax = hSrc->GetXaxis()->GetXmax();
00585 
00586   string newDir = startDir;
00587   string newPrefix = namePrefix;
00588   string::size_type shiftPos;
00589   if ( string::npos != (shiftPos = namePrefix.rfind('/')) ) {
00590     newDir += "/"+namePrefix.substr(0, shiftPos);
00591     newPrefix.erase(0, shiftPos+1);
00592   }
00593 
00594   theDQM->setCurrentFolder(newDir);
00595 
00596   ME* meanME = theDQM->book1D(newPrefix+"_Mean", titlePrefix+" Mean", nBin, xMin, xMax);
00597   ME* sigmaME = theDQM->book1D(newPrefix+"_Sigma", titlePrefix+" Sigma", nBin, xMin, xMax);
00598 //  ME* chi2ME  = theDQM->book1D(namePrefix+"_Chi2" , titlePrefix+" #Chi^{2}", nBin, xMin, xMax); // N/A
00599 
00600   if (! resLimitedFit_ ) {
00601     FitSlicesYTool fitTool(srcME);
00602     fitTool.getFittedMeanWithError(meanME);
00603     fitTool.getFittedSigmaWithError(sigmaME);
00605   } else {
00606     limitedFit(srcME,meanME,sigmaME);
00607   }
00608 }
00609 
00610 void DQMGenericClient::normalizeToEntries(const std::string& startDir, const std::string& histName, const std::string& normHistName) 
00611 {
00612   if ( ! theDQM->dirExists(startDir) ) {
00613     if ( verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_) ) {
00614       LogInfo("DQMGenericClient") << "normalizeToEntries() : "
00615                                      << "Cannot find sub-directory " << startDir << endl;
00616     }
00617     return;
00618   }
00619 
00620   theDQM->cd();
00621 
00622   ME* element = theDQM->get(startDir+"/"+histName);
00623   ME* normME = theDQM->get(startDir+"/"+normHistName);
00624 
00625   if ( !element ) {
00626     if ( verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_) ) {
00627       LogInfo("DQMGenericClient") << "normalizeToEntries() : "
00628                                      << "No such element '" << histName << "' found\n";
00629     }
00630     return;
00631   }
00632 
00633   if ( !normME ) {
00634     if ( verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_) ) {
00635       LogInfo("DQMGenericClient") << "normalizeToEntries() : "
00636                                      << "No such element '" << normHistName << "' found\n";
00637     }
00638     return;
00639   }
00640 
00641   TH1F* hist  = element->getTH1F();
00642   if ( !hist) {
00643     if ( verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_) ) {
00644       LogInfo("DQMGenericClient") << "normalizeToEntries() : "
00645                                      << "Cannot create TH1F from ME\n";
00646     }
00647     return;
00648   }
00649 
00650   TH1F* normHist = normME->getTH1F();
00651   if ( !normHist ) {
00652     if ( verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_) ) {
00653       LogInfo("DQMGenericClient") << "normalizeToEntries() : "
00654                                      << "Cannot create TH1F from ME\n";
00655     }
00656     return;
00657   }
00658 
00659   const double entries = normHist->GetEntries();
00660   if ( entries != 0 ) {
00661     hist->Scale(1./entries);
00662   }
00663   else {
00664     LogInfo("DQMGenericClient") << "normalizeToEntries() : " 
00665                                    << "Zero entries in histogram\n";
00666   }
00667 
00668   return;
00669 }
00670 
00671 void DQMGenericClient::makeCumulativeDist(const std::string& startDir, const std::string& cdName) 
00672 {
00673   if ( ! theDQM->dirExists(startDir) ) {
00674     if ( verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_) ) {
00675       LogInfo("DQMGenericClient") << "makeCumulativeDist() : "
00676                                      << "Cannot find sub-directory " << startDir << endl;
00677     }
00678     return;
00679   }
00680 
00681   theDQM->cd();
00682 
00683   ME* element_cd = theDQM->get(startDir+"/"+cdName);
00684 
00685   if ( !element_cd ) {
00686     if ( verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_) ) {
00687       LogInfo("DQMGenericClient") << "makeCumulativeDist() : "
00688                                      << "No such element '" << cdName << "' found\n";
00689     }
00690     return;
00691   }
00692 
00693   TH1F* cd  = element_cd->getTH1F();
00694 
00695   if ( !cd ) {
00696     if ( verbose_ >= 2 || (verbose_ == 1 && !isWildcardUsed_) ) {
00697       LogInfo("DQMGenericClient") << "makeCumulativeDist() : "
00698                                      << "Cannot create TH1F from ME\n";
00699     }
00700     return;
00701   }
00702 
00703   int n_bins = cd->GetNbinsX() + 1;
00704 
00705   for (int i = 1; i <= n_bins; i++) {
00706     cd->SetBinContent(i,cd->GetBinContent(i) + cd->GetBinContent(i-1));
00707   }
00708 
00709   return;
00710 }
00711 
00712 void DQMGenericClient::limitedFit(MonitorElement * srcME, MonitorElement * meanME, MonitorElement * sigmaME)
00713 {
00714   TH2F * histo = srcME->getTH2F();
00715 
00716   static int i = 0;
00717   i++;
00718 
00719   // Fit slices projected along Y from bins in X 
00720   double cont_min = 100;    //Minimum number of entries
00721   Int_t binx =  histo->GetXaxis()->GetNbins();
00722 
00723   for (int i = 1; i <= binx ; i++) {
00724     TString iString(i);
00725     TH1 *histoY =  histo->ProjectionY(" ", i, i);
00726     double cont = histoY->GetEntries();
00727 
00728     if (cont >= cont_min) {
00729       float minfit = histoY->GetMean() - histoY->GetRMS();
00730       float maxfit = histoY->GetMean() + histoY->GetRMS();
00731       
00732       TF1 *fitFcn = new TF1(TString("g")+histo->GetName()+iString,"gaus",minfit,maxfit);
00733       double x1,x2;
00734       fitFcn->GetRange(x1,x2);
00735 
00736       histoY->Fit(fitFcn,"QR0","",x1,x2);
00737 
00738 //      histoY->Fit(fitFcn->GetName(),"RME");
00739       double *par = fitFcn->GetParameters();
00740       double *err = fitFcn->GetParErrors();
00741 
00742       meanME->setBinContent(i,par[1]);
00743       meanME->setBinError(i,err[1]);
00744 
00745       sigmaME->setBinContent(i,par[2]);
00746       sigmaME->setBinError(i,err[2]);
00747       if(fitFcn) delete fitFcn;
00748       if(histoY) delete histoY;
00749     }
00750     else {
00751       if(histoY) delete histoY;
00752       continue;
00753     }
00754   }
00755 }
00756 
00757 //=================================
00758 
00759 void DQMGenericClient::findAllSubdirectories (std::string dir, std::set<std::string> * myList, TString pattern = "") {
00760   if (!theDQM->dirExists(dir)) {
00761     LogError("DQMGenericClient") << " DQMGenericClient::findAllSubdirectories ==> Missing folder " << dir << " !!!"; 
00762     return;
00763   }
00764   if (pattern != "") {
00765     if (pattern.Contains(nonPerlWildcard)) pattern.ReplaceAll("*",".*");
00766     TPRegexp regexp(pattern);
00767     theDQM->cd(dir);
00768     vector <string> foundDirs = theDQM->getSubdirs();
00769     for(vector<string>::const_iterator iDir = foundDirs.begin();
00770         iDir != foundDirs.end(); ++iDir) {
00771       TString dirName = iDir->substr(iDir->rfind('/') + 1, iDir->length());
00772       if (dirName.Contains(regexp))
00773         findAllSubdirectories ( *iDir, myList);
00774     }
00775   }
00776   //std::cout << "Looking for directory " << dir ;
00777   else if (theDQM->dirExists(dir)){
00778     //std::cout << "... it exists! Inserting it into the list ";
00779     myList->insert(dir);
00780     //std::cout << "... now list has size " << myList->size() << std::endl;
00781     theDQM->cd(dir);
00782     findAllSubdirectories (dir, myList, "*");
00783   } else {
00784     //std::cout << "... DOES NOT EXIST!!! Skip bogus dir" << std::endl;
00785     
00786     LogInfo ("DQMGenericClient") << "Trying to find sub-directories of " << dir
00787                                  << " failed because " << dir  << " does not exist";
00788                                  
00789   }
00790   return;
00791 }
00792 
00793 
00794 void DQMGenericClient::generic_eff (TH1* denom, TH1* numer, MonitorElement* efficiencyHist, const int type) {
00795   for (int iBinX = 1; iBinX < denom->GetNbinsX()+1; iBinX++){
00796     for (int iBinY = 1; iBinY < denom->GetNbinsY()+1; iBinY++){
00797       for (int iBinZ = 1; iBinZ < denom->GetNbinsZ()+1; iBinZ++){
00798 
00799         int globalBinNum = denom->GetBin(iBinX, iBinY, iBinZ);
00800         
00801         float numerVal = numer->GetBinContent(globalBinNum);
00802         float denomVal = denom->GetBinContent(globalBinNum);
00803 
00804         float effVal = 0;
00805 
00806         // fake eff is in use
00807         if (type == 2 ) {          
00808           effVal = denomVal ? (1 - numerVal / denomVal) : 0;
00809         } else {
00810           effVal = denomVal ? numerVal / denomVal : 0;
00811         }
00812 
00813         float errVal = (denomVal && (effVal <=1)) ? sqrt(effVal*(1-effVal)/denomVal) : 0;
00814 
00815         LogDebug ("DQMGenericClient") << "(iBinX, iBinY, iBinZ)  = "
00816              << iBinX << ", "
00817              << iBinY << ", "
00818              << iBinZ << "), global bin =  "  << globalBinNum
00819              << "eff = " << numerVal << "  /  " << denomVal
00820              << " =  " << effVal 
00821              << " ... setting the error for that bin ... " << endl
00822              << endl;
00823 
00824         
00825         efficiencyHist->setBinContent(globalBinNum, effVal);
00826         efficiencyHist->setBinError(globalBinNum, errVal);
00827       }
00828     }
00829   }
00830 
00831   //efficiencyHist->setMinimum(0.0);
00832   //efficiencyHist->setMaximum(1.0);
00833 }
00834 
00835 
00836 /* vim:set ts=2 sts=2 sw=2 expandtab: */