CMS 3D CMS Logo

EventIDChecker.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Modules
4 // Class: EventIDChecker
5 //
13 //
14 // Original Author: Chris Jones
15 // Created: Tue Jun 16 15:42:17 CDT 2009
16 //
17 //
18 
19 // user include files
27 
28 // system include files
29 #include <algorithm>
30 #include <memory>
31 #include <vector>
32 
33 //
34 // class decleration
35 //
36 
38 public:
39  explicit EventIDChecker(edm::ParameterSet const&);
40  ~EventIDChecker() override;
41  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
42 
43 private:
44  void beginJob() override;
45  void analyze(edm::Event const&, edm::EventSetup const&) override;
46  void endJob() override;
47 
48  // ----------member data ---------------------------
49  std::vector<edm::EventID> ids_;
50  unsigned int index_;
51 
55 };
56 
57 //
58 // constants, enums and typedefs
59 //
60 
61 //
62 // static data member definitions
63 //
64 
65 //
66 // constructors and destructor
67 //
69  : ids_(iConfig.getUntrackedParameter<std::vector<edm::EventID> >("eventSequence")),
70  index_(0),
71  multiProcessSequentialEvents_(iConfig.getUntrackedParameter<unsigned int>("multiProcessSequentialEvents")),
74  //now do what ever initialization is needed
75 }
76 
78  // do anything here that needs to be done at desctruction time
79  // (e.g. close files, deallocate resources etc.)
80 }
81 
82 //
83 // member functions
84 //
85 
86 namespace {
87  struct CompareWithoutLumi {
88  CompareWithoutLumi(edm::EventID const& iThis) : m_this(iThis) {}
89  bool operator()(edm::EventID const& iOther) {
90  return m_this.run() == iOther.run() && m_this.event() == iOther.event();
91  }
92  edm::EventID m_this;
93  };
94 } // namespace
95 
96 // ------------ method called to for each event ------------
98  if (mustSearch_) {
101  //the event must be after the last event in our list since multicore doesn't go backwards
102  std::vector<edm::EventID>::iterator itFind =
103  std::find_if(ids_.begin() + index_, ids_.end(), CompareWithoutLumi(iEvent.id()));
104  if (itFind == ids_.end()) {
105  throw cms::Exception("MissedEvent") << "The event " << iEvent.id() << "is not in the list.\n";
106  }
107  index_ = itFind - ids_.begin();
108  }
110  }
111 
112  if (index_ >= ids_.size()) {
113  throw cms::Exception("TooManyEvents")
114  << "Was passes " << ids_.size() << " EventIDs but have processed more events than that\n";
115  }
116  if (iEvent.id().run() != ids_[index_].run() || iEvent.id().event() != ids_[index_].event()) {
117  throw cms::Exception("WrongEvent") << "Was expecting event " << ids_[index_] << " but was given " << iEvent.id()
118  << "\n";
119  }
120  ++index_;
121 }
122 
123 // ------------ method called once each job just before starting event loop ------------
125 
126 // ------------ method called once each job just after ending the event loop ------------
128 
129 // ------------ method called once each job for validation
132  desc.addUntracked<std::vector<edm::EventID> >("eventSequence");
133  desc.addUntracked<unsigned int>("multiProcessSequentialEvents", 0U);
134  descriptions.add("eventIDChecker", desc);
135 }
136 
137 //define this as a plug-in
RunNumber_t run() const
Definition: EventID.h:39
EventNumber_t event() const
Definition: EventID.h:41
~EventIDChecker() override
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
unsigned int multiProcessSequentialEvents_
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void analyze(edm::Event const &, edm::EventSetup const &) override
unsigned int index_
std::vector< edm::EventID > ids_
unsigned int numberOfEventsLeftBeforeSearch_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
EventIDChecker(edm::ParameterSet const &)
edm::EventID id() const
Definition: EventBase.h:59
HLT enums.
void endJob() override
void beginJob() override