CMS 3D CMS Logo

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