CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/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 using namespace edm;
00015 using namespace dqmTnP;
00016 
00017 class DQMGenericTnPClient : public edm::EDAnalyzer{
00018   public:
00019     DQMGenericTnPClient(const edm::ParameterSet& pset);
00020     virtual ~DQMGenericTnPClient();
00021     virtual void analyze(const edm::Event& event, const edm::EventSetup& eventSetup) {};
00022     virtual void endRun(const edm::Run &run, const edm::EventSetup &setup);
00023     void calculateEfficiency(const ParameterSet& pset);
00024   private:
00025     DQMStore * dqmStore;
00026     TFile * plots;
00027     string myDQMrootFolder;
00028     bool verbose;
00029     const VParameterSet efficiencies;
00030     GaussianPlusLinearFitter *GPLfitter;
00031     VoigtianPlusExponentialFitter *VPEfitter;
00032 };
00033 
00034 DQMGenericTnPClient::DQMGenericTnPClient(const edm::ParameterSet& pset):
00035   myDQMrootFolder( pset.getUntrackedParameter<string>("MyDQMrootFolder") ),
00036   verbose( pset.getUntrackedParameter<bool>("Verbose",false) ),
00037   efficiencies( pset.getUntrackedParameter<VParameterSet>("Efficiencies") )
00038 {
00039   TString savePlotsInRootFileName = pset.getUntrackedParameter<string>("SavePlotsInRootFileName","");
00040   plots = savePlotsInRootFileName!="" ? new TFile(savePlotsInRootFileName,"recreate") : 0;
00041   GPLfitter = new GaussianPlusLinearFitter(verbose);
00042   VPEfitter = new VoigtianPlusExponentialFitter(verbose);
00043 }
00044 
00045 void DQMGenericTnPClient::endRun(const edm::Run &run, const edm::EventSetup &setup){
00046   dqmStore = Service<DQMStore>().operator->();
00047   if( !dqmStore ){
00048     LogError("DQMGenericTnPClient")<<"Could not find DQMStore service\n";
00049     return;
00050   }
00051   dqmStore->setCurrentFolder(myDQMrootFolder);
00052   //loop over all efficiency tasks
00053   for(VParameterSet::const_iterator pset = efficiencies.begin(); pset!=efficiencies.end(); pset++){
00054     calculateEfficiency(*pset);
00055   }
00056 }
00057   
00058 void DQMGenericTnPClient::calculateEfficiency(const ParameterSet& pset){
00059   //get hold of numerator and denominator histograms
00060   string allMEname = myDQMrootFolder+"/"+pset.getUntrackedParameter<string>("DenominatorMEname");
00061   string passMEname = myDQMrootFolder+"/"+pset.getUntrackedParameter<string>("NumeratorMEname");
00062   MonitorElement *allME = dqmStore->get(allMEname);
00063   MonitorElement *passME = dqmStore->get(passMEname);
00064   if(allME==0 || passME==0){
00065     LogDebug("DQMGenericTnPClient")<<"Could not find MEs: "<<allMEname<<" or "<<passMEname<<endl;
00066     return;
00067   }
00068   TH1 *all = allME->getTH1();
00069   TH1 *pass = passME->getTH1();
00070   //setup the fitter  
00071   string fitFunction = pset.getUntrackedParameter<string>("FitFunction");
00072   AbstractFitter *fitter = 0;
00073   if(fitFunction=="GaussianPlusLinear"){
00074     GPLfitter->setup(
00075       pset.getUntrackedParameter<double>("ExpectedMean"),
00076       pset.getUntrackedParameter<double>("FitRangeLow"),
00077       pset.getUntrackedParameter<double>("FitRangeHigh"),
00078       pset.getUntrackedParameter<double>("ExpectedSigma")
00079     );
00080     fitter = GPLfitter;
00081   }else if(fitFunction=="VoigtianPlusExponential"){
00082     VPEfitter->setup(
00083       pset.getUntrackedParameter<double>("ExpectedMean"),
00084       pset.getUntrackedParameter<double>("FitRangeLow"),
00085       pset.getUntrackedParameter<double>("FitRangeHigh"),
00086       pset.getUntrackedParameter<double>("ExpectedSigma"),
00087       pset.getUntrackedParameter<double>("FixedWidth")
00088     );
00089     fitter = VPEfitter;
00090   }else{
00091     LogError("DQMGenericTnPClient")<<"Fit function: "<<fitFunction<<" does not exist"<<endl;
00092     return;
00093   }
00094   //check dimensions
00095   int dimensions = all->GetDimension();
00096   int massDimension = pset.getUntrackedParameter<int>("MassDimension");
00097   if(massDimension>dimensions){
00098     LogError("DQMGenericTnPClient")<<"Monitoring Element "<<allMEname<<" has smaller dimension than "<<massDimension<<endl;
00099     return;
00100   }
00101   //figure out directory and efficiency names  
00102   string effName = pset.getUntrackedParameter<string>("EfficiencyMEname");
00103   string effDir = myDQMrootFolder;
00104   string::size_type slashPos = effName.rfind('/');
00105   if ( string::npos != slashPos ) {
00106     effDir += "/"+effName.substr(0, slashPos);
00107     effName.erase(0, slashPos+1);
00108   }
00109   dqmStore->setCurrentFolder(effDir);
00110   TString prefix(effDir.c_str());
00111   prefix.ReplaceAll('/','_');
00112   //calculate and book efficiency
00113   if(dimensions==2){
00114     TProfile* eff = 0;
00115     TProfile* effChi2 = 0;
00116     TString error = fitter->calculateEfficiency((TH2*)pass, (TH2*)all, massDimension, eff, effChi2, plots?prefix+effName.c_str():"");
00117     if(error!=""){
00118       LogError("DQMGenericTnPClient")<<error<<endl;
00119       return;
00120     }
00121     dqmStore->bookProfile(effName,eff);
00122     dqmStore->bookProfile(effName+"Chi2",effChi2);
00123     delete eff;
00124     delete effChi2;
00125   }else if(dimensions==3){
00126     TProfile2D* eff = 0;
00127     TProfile2D* effChi2 = 0;
00128     TString error = fitter->calculateEfficiency((TH3*)pass, (TH3*)all, massDimension, eff, effChi2, plots?prefix+effName.c_str():"");
00129     if(error!=""){
00130       LogError("DQMGenericTnPClient")<<error<<endl;
00131       return;
00132     }
00133     dqmStore->bookProfile2D(effName,eff);
00134     dqmStore->bookProfile2D(effName+"Chi2",effChi2);
00135     delete eff;
00136     delete effChi2;
00137   }
00138 }
00139 
00140 DQMGenericTnPClient::~DQMGenericTnPClient(){
00141   delete GPLfitter;
00142   if(plots){
00143     plots->Close();
00144   }
00145 }
00146 
00147 //define this as a plug-in
00148 DEFINE_FWK_MODULE(DQMGenericTnPClient);