CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/DQMOffline/Trigger/plugins/DQMGenericTnPClient.cc

Go to the documentation of this file.
00001 #include "FWCore/Framework/interface/Frameworkfwd.h"
00002 #include "FWCore/Framework/interface/EDAnalyzer.h"
00003 #include "FWCore/Framework/interface/Event.h"
00004 #include "FWCore/Framework/interface/MakerMacros.h"
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 #include "FWCore/ServiceRegistry/interface/Service.h"
00007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00008 #include "FWCore/Utilities/interface/InputTag.h"
00009 #include "DQMServices/Core/interface/DQMStore.h"
00010 #include "DQMServices/Core/interface/MonitorElement.h"
00011 
00012 #include "DQMOffline/Trigger/plugins/GenericTnPFitter.h"
00013 
00014 #include<TString.h>
00015 #include<TPRegexp.h>
00016 
00017 using namespace edm;
00018 using namespace dqmTnP;
00019 using namespace std;
00020 
00021 
00022 typedef std::vector<std::string> vstring;
00023 
00024 class DQMGenericTnPClient : public edm::EDAnalyzer{
00025   public:
00026     DQMGenericTnPClient(const edm::ParameterSet& pset);
00027     virtual ~DQMGenericTnPClient();
00028     virtual void analyze(const edm::Event& event, const edm::EventSetup& eventSetup) {};
00029     virtual void endRun(const edm::Run &run, const edm::EventSetup &setup);
00030   void calculateEfficiency(std::string dirName, const ParameterSet& pset);
00031     void findAllSubdirectories (std::string dir, std::set<std::string> * myList, TString pattern);
00032   private:
00033     DQMStore * dqmStore;
00034     TFile * plots;
00035     vstring subDirs;
00036     std::string myDQMrootFolder;
00037     bool verbose;
00038     const VParameterSet efficiencies;
00039     GaussianPlusLinearFitter *GPLfitter;
00040     VoigtianPlusExponentialFitter *VPEfitter;
00041 };
00042 
00043 DQMGenericTnPClient::DQMGenericTnPClient(const edm::ParameterSet& pset):
00044   subDirs( pset.getUntrackedParameter<vstring>("subDirs", vstring()) ),
00045   myDQMrootFolder( pset.getUntrackedParameter<std::string>("MyDQMrootFolder", "") ),
00046   verbose( pset.getUntrackedParameter<bool>("Verbose",false) ),
00047   efficiencies( pset.getUntrackedParameter<VParameterSet>("Efficiencies") )
00048 {
00049   TString savePlotsInRootFileName = pset.getUntrackedParameter<string>("SavePlotsInRootFileName","");
00050   plots = savePlotsInRootFileName!="" ? new TFile(savePlotsInRootFileName,"recreate") : 0;
00051   GPLfitter = new GaussianPlusLinearFitter(verbose);
00052   VPEfitter = new VoigtianPlusExponentialFitter(verbose);
00053 }
00054 
00055 void DQMGenericTnPClient::endRun(const edm::Run &run, const edm::EventSetup &setup){
00056 
00057   TPRegexp metacharacters("[\\^\\$\\.\\*\\+\\?\\|\\(\\)\\{\\}\\[\\]]");
00058 
00059   dqmStore = Service<DQMStore>().operator->();
00060   if( !dqmStore ){
00061     LogError("DQMGenericTnPClient")<<"Could not find DQMStore service\n";
00062     return;
00063   }
00064   dqmStore->setCurrentFolder(myDQMrootFolder);
00065 
00066   set<std::string> subDirSet;
00067   
00068   if (myDQMrootFolder != "")
00069     subDirSet.insert(myDQMrootFolder);
00070   else {
00071     for(vstring::const_iterator iSubDir = subDirs.begin(); 
00072         iSubDir != subDirs.end(); ++iSubDir) {
00073       std::string subDir = *iSubDir;
00074       if ( subDir[subDir.size()-1] == '/' ) subDir.erase(subDir.size()-1);
00075       if ( TString(subDir).Contains(metacharacters) ) {
00076         const string::size_type shiftPos = subDir.rfind('/');
00077         const string searchPath = subDir.substr(0, shiftPos);
00078         const string pattern    = subDir.substr(shiftPos + 1, subDir.length());
00079         findAllSubdirectories (searchPath, &subDirSet, pattern);
00080       }
00081       else {
00082         subDirSet.insert(subDir);
00083       }
00084     }
00085   }
00086 
00087   for(set<string>::const_iterator iSubDir = subDirSet.begin();
00088       iSubDir != subDirSet.end(); ++iSubDir) {
00089     const string& dirName = *iSubDir;
00090     for(VParameterSet::const_iterator pset = efficiencies.begin(); 
00091         pset != efficiencies.end(); ++pset) {
00092       calculateEfficiency(dirName, *pset);
00093     }
00094   }
00095 
00096 }
00097   
00098 void DQMGenericTnPClient::calculateEfficiency(std::string dirName, const ParameterSet& pset){
00099   //get hold of numerator and denominator histograms
00100   string allMEname = dirName+"/"+pset.getUntrackedParameter<string>("DenominatorMEname");
00101   string passMEname = dirName+"/"+pset.getUntrackedParameter<string>("NumeratorMEname");
00102   MonitorElement *allME = dqmStore->get(allMEname);
00103   MonitorElement *passME = dqmStore->get(passMEname);
00104   if(allME==0 || passME==0){
00105     LogDebug("DQMGenericTnPClient")<<"Could not find MEs: "<<allMEname<<" or "<<passMEname<<endl;
00106     return;
00107   }
00108   TH1 *all = allME->getTH1();
00109   TH1 *pass = passME->getTH1();
00110   //setup the fitter  
00111   string fitFunction = pset.getUntrackedParameter<string>("FitFunction");
00112   AbstractFitter *fitter = 0;
00113   if(fitFunction=="GaussianPlusLinear"){
00114     GPLfitter->setup(
00115       pset.getUntrackedParameter<double>("ExpectedMean"),
00116       pset.getUntrackedParameter<double>("FitRangeLow"),
00117       pset.getUntrackedParameter<double>("FitRangeHigh"),
00118       pset.getUntrackedParameter<double>("ExpectedSigma")
00119     );
00120     fitter = GPLfitter;
00121   }else if(fitFunction=="VoigtianPlusExponential"){
00122     VPEfitter->setup(
00123       pset.getUntrackedParameter<double>("ExpectedMean"),
00124       pset.getUntrackedParameter<double>("FitRangeLow"),
00125       pset.getUntrackedParameter<double>("FitRangeHigh"),
00126       pset.getUntrackedParameter<double>("ExpectedSigma"),
00127       pset.getUntrackedParameter<double>("FixedWidth")
00128     );
00129     fitter = VPEfitter;
00130   }else{
00131     LogError("DQMGenericTnPClient")<<"Fit function: "<<fitFunction<<" does not exist"<<endl;
00132     return;
00133   }
00134   //check dimensions
00135   int dimensions = all->GetDimension();
00136   int massDimension = pset.getUntrackedParameter<int>("MassDimension");
00137   if(massDimension>dimensions){
00138     LogError("DQMGenericTnPClient")<<"Monitoring Element "<<allMEname<<" has smaller dimension than "<<massDimension<<endl;
00139     return;
00140   }
00141   //figure out directory and efficiency names  
00142   string effName = pset.getUntrackedParameter<string>("EfficiencyMEname");
00143   string effDir = dirName;
00144   string::size_type slashPos = effName.rfind('/');
00145   if ( string::npos != slashPos ) {
00146     effDir += "/"+effName.substr(0, slashPos);
00147     effName.erase(0, slashPos+1);
00148   }
00149   dqmStore->setCurrentFolder(effDir);
00150   TString prefix(effDir.c_str());
00151   prefix.ReplaceAll('/','_');
00152   //calculate and book efficiency
00153   if(dimensions==2){
00154     TProfile* eff = 0;
00155     TProfile* effChi2 = 0;
00156     TString error = fitter->calculateEfficiency((TH2*)pass, (TH2*)all, massDimension, eff, effChi2, plots?prefix+effName.c_str():"");
00157     if(error!=""){
00158       LogError("DQMGenericTnPClient")<<error<<endl;
00159       return;
00160     }
00161     dqmStore->bookProfile(effName,eff);
00162     dqmStore->bookProfile(effName+"Chi2",effChi2);
00163     delete eff;
00164     delete effChi2;
00165   }else if(dimensions==3){
00166     TProfile2D* eff = 0;
00167     TProfile2D* effChi2 = 0;
00168     TString error = fitter->calculateEfficiency((TH3*)pass, (TH3*)all, massDimension, eff, effChi2, plots?prefix+effName.c_str():"");
00169     if(error!=""){
00170       LogError("DQMGenericTnPClient")<<error<<endl;
00171       return;
00172     }
00173     dqmStore->bookProfile2D(effName,eff);
00174     dqmStore->bookProfile2D(effName+"Chi2",effChi2);
00175     delete eff;
00176     delete effChi2;
00177   }
00178 }
00179 
00180 DQMGenericTnPClient::~DQMGenericTnPClient(){
00181   delete GPLfitter;
00182   if(plots){
00183     plots->Close();
00184   }
00185 }
00186 
00187 void DQMGenericTnPClient::findAllSubdirectories (std::string dir, std::set<std::string> * myList, TString pattern = "") {
00188 
00189   if (!dqmStore->dirExists(dir)) {
00190     LogError("DQMGenericTnPClient") << " DQMGenericTnPClient::findAllSubdirectories ==> Missing folder " << dir << " !!!";
00191     return;
00192   }
00193   TPRegexp nonPerlWildcard("\\w\\*|^\\*");
00194   if (pattern != "") {
00195     if (pattern.Contains(nonPerlWildcard)) pattern.ReplaceAll("*",".*");
00196     TPRegexp regexp(pattern);
00197     dqmStore->cd(dir);
00198     vector <string> foundDirs = dqmStore->getSubdirs();
00199     for(vector<string>::const_iterator iDir = foundDirs.begin();
00200         iDir != foundDirs.end(); ++iDir) {
00201       TString dirName = iDir->substr(iDir->rfind('/') + 1, iDir->length());
00202       if (dirName.Contains(regexp))
00203         findAllSubdirectories ( *iDir, myList);
00204     }
00205   }
00206   else if (dqmStore->dirExists(dir)){
00207     myList->insert(dir);
00208     dqmStore->cd(dir);
00209     findAllSubdirectories (dir, myList, "*");
00210   } else {
00211     
00212     LogInfo ("DQMGenericClient") << "Trying to find sub-directories of " << dir
00213                                  << " failed because " << dir  << " does not exist";
00214                                  
00215   }
00216   return;
00217 }
00218 
00219 
00220 //define this as a plug-in
00221 DEFINE_FWK_MODULE(DQMGenericTnPClient);