CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/HLTriggerOffline/Higgs/src/HLTHiggsValidator.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:     HLTHiggsValidator
00004 // Class:       HLTHiggsValidator
00005 // 
00006 
00007 //
00008 // Jordi Duarte Campderros (based on the Jason Slaunwhite 
00009 // and Jeff Klukas coded from the HLTriggerOffline/Muon package
00010 //
00011 // $Id: HLTHiggsValidator.cc,v 1.7 2012/03/23 11:50:56 duarte Exp $
00012 //
00013 //
00014 
00015 // system include files
00016 //#include<memory>
00017 
00018 //#include "FWCore/Framework/interface/MakerMacros.h"
00019 
00020 #include "HLTriggerOffline/Higgs/interface/HLTHiggsValidator.h"
00021 #include "HLTriggerOffline/Higgs/src/EVTColContainer.cc"
00022 
00024 // Constructor
00025 HLTHiggsValidator::HLTHiggsValidator(const edm::ParameterSet& pset) :
00026         _pset(pset),
00027         _analysisnames(pset.getParameter<std::vector<std::string> >("analysis")),
00028         _collections(0),
00029         _dbe(0)
00030 {
00031         _collections = new EVTColContainer;
00032 }
00033 
00034 HLTHiggsValidator::~HLTHiggsValidator()
00035 {
00036         if( _collections != 0 )
00037         {
00038                 delete _collections;
00039                 _collections = 0;
00040         }
00041 }
00042 
00043 
00044 void HLTHiggsValidator::beginRun(const edm::Run & iRun, const edm::EventSetup & iSetup) 
00045 {
00046         for(size_t i = 0; i < _analysisnames.size() ; ++i)
00047         {
00048                 HLTHiggsSubAnalysis analyzer(_pset, _analysisnames.at(i));
00049                 _analyzers.push_back(analyzer);
00050         }
00051         // Call the Plotter beginRun (which stores the triggers paths..:)
00052         for(std::vector<HLTHiggsSubAnalysis>::iterator iter = _analyzers.begin(); 
00053                         iter != _analyzers.end(); ++iter) 
00054         {
00055                 iter->beginRun(iRun, iSetup);
00056         }
00057 }
00058         
00059 
00060 void HLTHiggsValidator::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00061 {
00062         static int eventNumber = 0;
00063         eventNumber++;
00064         LogTrace("HiggsValidation") << "In HLTHiggsSubAnalysis::analyze,  " 
00065                 << "Event: " << eventNumber;
00066         
00067         // Initialize the event collections
00068         this->_collections->reset();
00069 
00070         for(std::vector<HLTHiggsSubAnalysis>::iterator iter = _analyzers.begin(); 
00071                         iter != _analyzers.end(); ++iter) 
00072         {
00073                 iter->analyze(iEvent, iSetup, this->_collections);
00074         }
00075 }
00076 
00077 
00078 
00079 void HLTHiggsValidator::beginJob()
00080 {
00081 }
00082 
00083 void HLTHiggsValidator::endRun(const edm::Run & iRun, const edm::EventSetup& iSetup)
00084 {
00085         // vector<HLTMuonPlotter>::iterator iter;
00086         // for(std::vector<HLTHiggsPlotter>::iterator iter = _analyzers.begin(); 
00087         //                 iter != analyzers_.end(); ++iter) 
00088         // {
00089         //         iter->endRun(iRun, iSetup);
00090         // }
00091 }
00092 
00093 
00094 void HLTHiggsValidator::endJob()
00095 {
00096 }
00097 
00098 
00099 
00100 //define this as a plug-in
00101 //DEFINE_FWK_MODULE(HLTHiggsValidator);