CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/src/GeneratorInterface/RivetInterface/interface/RivetAnalyzer.h

Go to the documentation of this file.
00001 #ifndef GeneratorInterface_RivetInterface_RivetAnalyzer
00002 #define GeneratorInterface_RivetInterface_RivetAnalyzer
00003 
00004 #include "FWCore/Framework/interface/EDAnalyzer.h"
00005 #include "FWCore/Utilities/interface/InputTag.h"
00006 #include "Rivet/RivetAIDA.fhh"
00007 #include "Rivet/AnalysisHandler.hh"
00008 
00009 //DQM services
00010 #include "DQMServices/Core/interface/DQMStore.h"
00011 #include "FWCore/ServiceRegistry/interface/Service.h"
00012 #include "DQMServices/Core/interface/MonitorElement.h"
00013 
00014 #include "Rivet/RivetAIDA.hh"
00015 #include "LWH/AIHistogramFactory.h"
00016 #include "LWH/VariAxis.h"
00017 
00018 #include <vector>
00019 #include <string>
00020 
00021 namespace edm{
00022   class ParameterSet;
00023   class Event;
00024   class EventSetup;
00025   class InputTag;
00026 }
00027 
00028 class RivetAnalyzer : public edm::EDAnalyzer
00029 {
00030   public:
00031   RivetAnalyzer(const edm::ParameterSet&);
00032 
00033   virtual ~RivetAnalyzer();
00034 
00035   virtual void beginJob();
00036 
00037   virtual void endJob();  
00038 
00039   virtual void analyze(const edm::Event&, const edm::EventSetup&);
00040 
00041   virtual void beginRun(const edm::Run&, const edm::EventSetup&);
00042 
00043   virtual void endRun(const edm::Run&, const edm::EventSetup&);
00044   
00045   private:
00046 
00047   void normalizeTree(AIDA::ITree& tree);
00048   template<class AIDATYPE, class ROOTTYPE> ROOTTYPE* prebook(const AIDATYPE*, const std::string&);
00049   template<class AIDATYPE, class ROOTTYPE> ROOTTYPE* aida2root(const AIDATYPE*, const std::string&); 
00050   
00051 
00052   edm::InputTag            _hepmcCollection;
00053   bool                     _useExternalWeight;
00054   edm::InputTag            _genEventInfoCollection;
00055   Rivet::AnalysisHandler   _analysisHandler;   
00056   bool                     _isFirstEvent;
00057   std::string              _outFileName;
00058   bool                     _doFinalize;
00059   bool                     _produceDQM;
00060 
00061   DQMStore *dbe;
00062   std::vector<MonitorElement *> _mes;
00063 };
00064 
00065 template<class AIDATYPE, class ROOTTYPE> 
00066 ROOTTYPE* 
00067 RivetAnalyzer::prebook(const AIDATYPE* aidah, const std::string& name){
00068   ROOTTYPE* h = 0;
00069   if (aidah->axis().isFixedBinning() ) {//equidistant binning (easier case)
00070     int nbins = aidah->axis().bins();
00071     h = new ROOTTYPE(name.c_str(), name.c_str(), nbins, aidah->axis().lowerEdge(), aidah->axis().upperEdge());
00072   } else {
00073     int nbins = aidah->axis().bins();
00074     //need to dyn cast, IAxis lacks polymorfism
00075     const LWH::VariAxis* vax = dynamic_cast<const LWH::VariAxis*>(&aidah->axis());
00076     if (! vax ){
00077       throw cms::Exception("RivetAnalyzer") << "Cannot dynamix cast an AIDA axis to VariAxis ";
00078     }
00079     double* bins = new double[nbins+1];
00080     for (int i=0; i<nbins; ++i) {
00081       bins[i] = vax->binEdges(i).first;
00082     }
00083     bins[nbins] = vax->binEdges(nbins-1).second; //take last bin right border
00084     h = new ROOTTYPE(name.c_str(), name.c_str(), nbins, bins);
00085     delete bins;
00086   }
00087   return h; 
00088 }
00089 
00090 template<> 
00091 TH1F* RivetAnalyzer::aida2root(const AIDA::IHistogram1D* aidah, const std::string& name){
00092   /*TH1F* h = 0;
00093   if (aidah->axis().isFixedBinning() ) {//equidistant binning (easier case)
00094     int nbins = aidah->axis().bins();
00095     h = new TH1F(name.c_str(), name.c_str(), nbins, aidah->axis().lowerEdge(), aidah->axis().upperEdge());
00096   } else {
00097     int nbins = aidah->axis().bins();
00098     //need to dyn cast, IAxis lacks polymorfism
00099     const LWH::VariAxis* vax = dynamic_cast<const LWH::VariAxis*>(&aidah->axis());
00100     if (! vax ){
00101       throw cms::Exception("RivetAnalyzer") << "Cannot dynamix cast an AIDA axis to VariAxis ";
00102     }
00103     double* bins = new double[nbins+1];
00104     for (int i=0; i<nbins; ++i) {
00105       bins[i] = vax->binEdges(i).first;
00106     }
00107     bins[nbins] = vax->binEdges(nbins-1).second; //take last bin right border
00108     h = new TH1F(name.c_str(), name.c_str(), nbins, bins);
00109     delete bins;
00110   }
00111   */
00112   TH1F* h = prebook<AIDA::IHistogram1D, TH1F>(aidah, name);
00113   for (int i = 0; i < aidah->axis().bins(); ++i){
00114     h->SetBinContent(i+1, aidah->binHeight(i));
00115     h->SetBinError(i+1, aidah->binError(i));
00116   }  
00117   return h;
00118 }
00119 
00120 template<>
00121 TProfile* RivetAnalyzer::aida2root(const AIDA::IProfile1D* aidah, const std::string& name){
00122   TProfile* h = prebook<AIDA::IProfile1D, TProfile>(aidah, name);
00123   for (int i = 0; i < aidah->axis().bins(); ++i){
00124     h->SetBinContent(i+1, aidah->binMean(i));
00125     h->SetBinError(i+1, aidah->binRms(i));
00126   }
00127   return h;
00128 }
00129 
00130 
00131 #endif