CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/CommonTools/UtilAlgos/interface/HistoAnalyzer.h

Go to the documentation of this file.
00001 #ifndef UtilAlgos_HistoAnalyzer_h
00002 #define UtilAlgos_HistoAnalyzer_h
00003 
00013 #include "FWCore/Framework/interface/EDAnalyzer.h"
00014 #include "FWCore/Framework/interface/Event.h"
00015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00016 #include "FWCore/Utilities/interface/InputTag.h"
00017 #include "FWCore/ServiceRegistry/interface/Service.h"
00018 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00019 #include "CommonTools/UtilAlgos/interface/ExpressionHisto.h"
00020 
00021 template<typename C>
00022 class HistoAnalyzer : public edm::EDAnalyzer {
00023 public:
00025   HistoAnalyzer( const edm::ParameterSet& );
00027   ~HistoAnalyzer();
00028   
00029 protected:
00031   virtual void analyze( const edm::Event&, const edm::EventSetup& );
00032 
00033 private:
00035   edm::InputTag src_;
00037   bool usingWeights_;
00039   edm::InputTag weights_;
00041   std::vector<ExpressionHisto<typename C::value_type>* > vhistograms;
00042 };
00043 
00044 template<typename C>
00045 HistoAnalyzer<C>::HistoAnalyzer( const edm::ParameterSet& par ) : 
00046   src_(par.template getParameter<edm::InputTag>("src")),
00047   usingWeights_(par.exists("weights")),
00048   weights_(par.template getUntrackedParameter<edm::InputTag>("weights", edm::InputTag("fake")))
00049 {
00050    edm::Service<TFileService> fs;
00051    std::vector<edm::ParameterSet> histograms = 
00052                                    par.template getParameter<std::vector<edm::ParameterSet> >("histograms");
00053    std::vector<edm::ParameterSet>::const_iterator it = histograms.begin();
00054    std::vector<edm::ParameterSet>::const_iterator end = histograms.end();
00055 
00056    // create the histograms from the given parameter sets 
00057    for (; it!=end; ++it)
00058    {
00059       ExpressionHisto<typename C::value_type>* hist = new ExpressionHisto<typename C::value_type>(*it);
00060       hist->initialize(*fs);
00061       vhistograms.push_back(hist);
00062    }   
00063 
00064 }
00065 
00066 template<typename C>
00067 HistoAnalyzer<C>::~HistoAnalyzer() 
00068 {
00069    // delete all histograms and clear the vector of pointers
00070    typename std::vector<ExpressionHisto<typename C::value_type>* >::iterator it = vhistograms.begin(); 
00071    typename std::vector<ExpressionHisto<typename C::value_type>* >::iterator end = vhistograms.end();
00072    for (;it!=end; ++it){
00073      (*it)->~ExpressionHisto<typename C::value_type>();
00074    }
00075    vhistograms.clear(); 
00076 }
00077 
00078 template<typename C>
00079 void HistoAnalyzer<C>::analyze( const edm::Event& iEvent, const edm::EventSetup& ) 
00080 {
00081    edm::Handle<C> coll;
00082    iEvent.getByLabel( src_, coll);
00083    double weight = 1.0;
00084    if (usingWeights_) { 
00085        edm::Handle<double> weightColl;
00086        iEvent.getByLabel( weights_, weightColl ); 
00087        weight = *weightColl;
00088    }
00089 
00090    typename std::vector<ExpressionHisto<typename C::value_type>* >::iterator it = vhistograms.begin();
00091    typename std::vector<ExpressionHisto<typename C::value_type>* >::iterator end = vhistograms.end(); 
00092 
00093    for (;it!=end; ++it){
00094       uint32_t i = 0;
00095       for( typename C::const_iterator elem=coll->begin(); elem!=coll->end(); ++elem, ++i ) {
00096          if (!(*it)->fill( *elem, weight, i )) {
00097             break;
00098          }
00099       }
00100    }
00101 }
00102 
00103 #endif