CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/FWCore/Modules/src/EventIDChecker.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    Modules
00004 // Class:      EventIDChecker
00005 //
00013 //
00014 // Original Author:  Chris Jones
00015 //         Created:  Tue Jun 16 15:42:17 CDT 2009
00016 //
00017 //
00018 
00019 // user include files
00020 #include "DataFormats/Provenance/interface/EventID.h"
00021 #include "FWCore/Framework/interface/EDAnalyzer.h"
00022 #include "FWCore/Framework/interface/Event.h"
00023 #include "FWCore/Framework/interface/MakerMacros.h"
00024 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
00025 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00026 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
00027 
00028 // system include files
00029 #include <algorithm>
00030 #include <memory>
00031 #include <vector>
00032 
00033 //
00034 // class decleration
00035 //
00036 
00037 class EventIDChecker : public edm::EDAnalyzer {
00038 public:
00039    explicit EventIDChecker(edm::ParameterSet const&);
00040    ~EventIDChecker();
00041     static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
00042 
00043 
00044 private:
00045    virtual void beginJob();
00046    virtual void analyze(edm::Event const&, edm::EventSetup const&);
00047    virtual void endJob();
00048    virtual void postForkReacquireResources(unsigned int iChildIndex, unsigned int iNumberOfChildren);
00049 
00050    // ----------member data ---------------------------
00051    std::vector<edm::EventID> ids_;
00052    unsigned int index_;
00053 
00054    unsigned int multiProcessSequentialEvents_;
00055    unsigned int numberOfEventsLeftBeforeSearch_;
00056    bool mustSearch_;
00057 };
00058 
00059 //
00060 // constants, enums and typedefs
00061 //
00062 
00063 //
00064 // static data member definitions
00065 //
00066 
00067 //
00068 // constructors and destructor
00069 //
00070 EventIDChecker::EventIDChecker(edm::ParameterSet const& iConfig) :
00071   ids_(iConfig.getUntrackedParameter<std::vector<edm::EventID> >("eventSequence")),
00072   index_(0),
00073   multiProcessSequentialEvents_(iConfig.getUntrackedParameter<unsigned int>("multiProcessSequentialEvents")),
00074   numberOfEventsLeftBeforeSearch_(0),
00075   mustSearch_(false)
00076 {
00077    //now do what ever initialization is needed
00078 
00079 }
00080 
00081 
00082 EventIDChecker::~EventIDChecker() {
00083 
00084    // do anything here that needs to be done at desctruction time
00085    // (e.g. close files, deallocate resources etc.)
00086 
00087 }
00088 
00089 
00090 //
00091 // member functions
00092 //
00093 
00094 namespace {
00095    struct CompareWithoutLumi {
00096       CompareWithoutLumi(edm::EventID const& iThis) : m_this(iThis) {
00097       }
00098       bool operator()(edm::EventID const& iOther) {
00099          return m_this.run() == iOther.run() && m_this.event() == iOther.event();
00100       }
00101       edm::EventID m_this;
00102    };
00103 }
00104 
00105 // ------------ method called to for each event  ------------
00106 void
00107 EventIDChecker::analyze(edm::Event const& iEvent, edm::EventSetup const&) {
00108    if(mustSearch_) {
00109       if(0 == numberOfEventsLeftBeforeSearch_) {
00110          numberOfEventsLeftBeforeSearch_ = multiProcessSequentialEvents_;
00111          //the event must be after the last event in our list since multicore doesn't go backwards
00112          std::vector<edm::EventID>::iterator itFind= std::find_if(ids_.begin()+index_, ids_.end(), CompareWithoutLumi(iEvent.id()));
00113          if(itFind == ids_.end()) {
00114             throw cms::Exception("MissedEvent") << "The event " << iEvent.id() << "is not in the list.\n";
00115          }
00116          index_ = itFind-ids_.begin();
00117       }
00118       --numberOfEventsLeftBeforeSearch_;
00119    }
00120 
00121    if(index_ >= ids_.size()) {
00122       throw cms::Exception("TooManyEvents") << "Was passes " << ids_.size() << " EventIDs but have processed more events than that\n";
00123    }
00124    if(iEvent.id().run() != ids_[index_].run() || iEvent.id().event() != ids_[index_].event()) {
00125       throw cms::Exception("WrongEvent") << "Was expecting event " << ids_[index_] << " but was given " << iEvent.id() << "\n";
00126    }
00127    ++index_;
00128 }
00129 
00130 // ------------ method called once each job just before starting event loop  ------------
00131 void
00132 EventIDChecker::beginJob() {
00133 }
00134 
00135 // ------------ method called once each job just after ending the event loop  ------------
00136 void
00137 EventIDChecker::endJob() {
00138 }
00139 
00140 // ------------ method called once each job for validation
00141 void
00142 EventIDChecker::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
00143   edm::ParameterSetDescription desc;
00144   desc.addUntracked<std::vector<edm::EventID> >("eventSequence");
00145   desc.addUntracked<unsigned int>("multiProcessSequentialEvents", 0U);
00146   descriptions.add("eventIDChecker", desc);
00147 }
00148 
00149 void
00150 EventIDChecker::postForkReacquireResources(unsigned int /*iChildIndex*/, unsigned int /*iNumberOfChildren*/) {
00151    mustSearch_ = true;
00152 }
00153 
00154 //define this as a plug-in
00155 DEFINE_FWK_MODULE(EventIDChecker);