CMS 3D CMS Logo

EventWithHistory.cc
Go to the documentation of this file.
1 #include <map>
7 
9 
11 
13  : TinyEvent(event, orbit, bx), _prevse() {}
14 
16  : TinyEvent(event, orbit, bx), _prevse() {}
17 
19 
20 EventWithHistory::EventWithHistory(const std::vector<edm::EventAuxiliary>& veaux)
21  : TinyEvent((!veaux.empty()) ? veaux[veaux.size() - 1] : TinyEvent()), _prevse() {
22  for (std::vector<edm::EventAuxiliary>::const_reverse_iterator eaux = veaux.rbegin(); eaux != veaux.rend(); eaux++) {
23  if (eaux != veaux.rbegin()) {
24  _prevse.push_back(*eaux);
25  }
26  }
27 }
28 
30  const L1AcceptBunchCrossingCollection& l1abcc,
31  const long long orbitoffset,
32  const int bxoffset)
33  : TinyEvent(), _prevse() {
34  std::map<int, TinyEvent> tmpmap;
35 
36  for (L1AcceptBunchCrossingCollection::const_iterator l1abc = l1abcc.begin(); l1abc != l1abcc.end(); ++l1abc) {
37  if (event.id().event() > (edm::EventNumber_t)(-1 * l1abc->l1AcceptOffset())) {
38  edm::EventNumber_t evnumb = event.id().event() + l1abc->l1AcceptOffset();
39  if (orbitoffset < (long long)l1abc->orbitNumber()) {
40  unsigned int neworbit = l1abc->orbitNumber() - orbitoffset;
41  int newbx = l1abc->bunchCrossing() - bxoffset;
42 
43  /*
44  the lines below assumes that the BX number is between 0 and 3563. If this is not the case it will jump to 0 and to the next orbit in case of
45  evets with BX=3564
46  */
47  while (newbx > 3563) {
48  ++neworbit;
49  newbx -= 3564;
50  }
51  while (newbx < 0) {
52  --neworbit;
53  newbx += 3564;
54  }
55 
56  if (l1abc->eventType() != 0) {
57  TinyEvent tmpse(evnumb, neworbit, newbx);
58  tmpmap[l1abc->l1AcceptOffset()] = tmpse;
59  } else {
60  edm::LogWarning("L1AcceptBunchCrossingNoType") << "L1AcceptBunchCrossing with no type found: ";
61  for (L1AcceptBunchCrossingCollection::const_iterator debu = l1abcc.begin(); debu != l1abcc.end(); ++debu) {
62  edm::LogPrint("L1AcceptBunchCrossingNoType") << *debu;
63  }
64  }
65  } else {
66  edm::LogError("L1AcceptBunchCrossingOffsetTooLarge")
67  << " Too large orbit offset " << orbitoffset << " " << l1abc->orbitNumber();
68  }
69  } else {
70  edm::LogInfo("L1AcceptBunchCrossingNegativeEvent") << "L1AcceptBunchCrossing with negative event: ";
71  for (L1AcceptBunchCrossingCollection::const_iterator debu = l1abcc.begin(); debu != l1abcc.end(); ++debu) {
72  edm::LogVerbatim("L1AcceptBunchCrossingNegativeEvent") << *debu;
73  }
74  }
75  }
76  // look for the event itself
77  if (tmpmap.find(0) != tmpmap.end()) {
78  TinyEvent::operator=(tmpmap[0]);
79 
80  // loop on the rest of the map and stop when it is missing
81  // check that the events are in the right order and break if not
82 
83  int counter = -1;
84  while (tmpmap.find(counter) != tmpmap.end()) {
85  if (tmpmap[counter + 1].deltaBX(tmpmap[counter]) > 0) {
86  _prevse.push_back(tmpmap[counter]);
87  --counter;
88  } else {
89  edm::LogWarning("L1AcceptBunchCrossingNotInOrder")
90  << "L1AcceptBunchCrossing not in order: orbit " << event.orbitNumber() << " BX " << event.bunchCrossing()
91  << " orbit offset " << orbitoffset << " bx offset " << bxoffset << " :";
92  for (L1AcceptBunchCrossingCollection::const_iterator debu = l1abcc.begin(); debu != l1abcc.end(); ++debu) {
93  edm::LogPrint("L1AcceptBunchCrossingNotInOrder") << *debu;
94  }
95  break;
96  }
97  }
98  } else {
99  TinyEvent::operator=(event);
100  edm::LogWarning("L1AcceptBunchCrossingNoCollection") << " L1AcceptBunchCrossing with offset=0 not found "
101  << " likely L1ABCCollection is empty ";
102  }
103 }
104 
106 
108  if (this != &he) {
110  _prevse = he._prevse;
111  }
112  return *this;
113 }
114 
115 // int EventWithHistory::operator<(const EventWithHistory& other) const { return TinyEvent::operator<(other); }
116 
118  int equal = TinyEvent::operator==(other);
119 
120  // depth is not checked anymore
121 
122  // equal = equal && (depth() == other.depth());
123 
124  if (equal) {
125  for (unsigned int i = 0; i < ((depth() < other.depth()) ? depth() : other.depth()); i++) {
126  equal = equal && (_prevse[i] == other._prevse[i]);
127  }
128  }
129 
130  return equal;
131 }
132 
133 int EventWithHistory::add(const EventWithHistory& he, const int idepth) {
134  if (!add((const TinyEvent&)he, idepth))
135  return 0;
136 
137  for (std::vector<TinyEvent>::const_iterator ev = he._prevse.begin(); ev != he._prevse.end(); ev++) {
138  if (!add(*ev, idepth))
139  return 0;
140  }
141  return 1;
142 }
143 
144 int EventWithHistory::add(const TinyEvent& se, const int idepth) {
145  bool isfuture = (idepth < 0);
146  const unsigned int adepth = abs(idepth);
147 
148  // protect against the possibility of filling with past and future history
149 
150  if (depth() > 0 && ((isfuture && !isFutureHistory()) || (!isfuture && isFutureHistory())))
151  return 0;
152 
153  if (adepth == 0)
154  return 0;
155  if (_prevse.size() >= adepth)
156  return 0;
157 
158  if (_prevse.empty()) {
159  if ((!isfuture && isNextOf(se)) || (isfuture && se.isNextOf(*this))) {
160  _prevse.push_back(se);
161  return 1;
162  } else {
163  return 0;
164  }
165  } else {
166  if ((!isfuture && _prevse[_prevse.size() - 1].isNextOf(se)) ||
167  (isfuture && se.isNextOf(_prevse[_prevse.size() - 1]))) {
168  _prevse.push_back(se);
169  return 1;
170  } else {
171  return 0;
172  }
173  }
174  return 0;
175 }
176 
178 const unsigned int EventWithHistory::orbit() const { return TinyEvent::_orbit; }
179 const int EventWithHistory::bx() const { return TinyEvent::_bx; }
180 
181 const TinyEvent* EventWithHistory::get(const unsigned int ev) const {
182  if (ev == 0)
183  return this;
184  if (ev <= _prevse.size())
185  return &_prevse[ev - 1];
186  return nullptr;
187 }
188 
189 unsigned int EventWithHistory::depth() const { return _prevse.size(); }
190 
191 bool EventWithHistory::isFutureHistory() const { return (depth() > 0 && _prevse[0].isNextOf(*this)); }
192 
193 long long EventWithHistory::deltaBX(const unsigned int ev2, const unsigned int ev1) const {
194  if (ev2 == ev1)
195  return 0;
196 
197  if (ev2 < ev1 && ev1 <= _prevse.size()) {
198  if (ev2 == 0)
199  return TinyEvent::deltaBX(_prevse[ev1 - 1]);
200  return _prevse[ev2 - 1].deltaBX(_prevse[ev1 - 1]);
201  }
202 
203  return -1;
204 }
205 
206 long long EventWithHistory::deltaBX(const unsigned int ev1) const { return deltaBX(0, ev1); }
207 
208 long long EventWithHistory::deltaBX() const { return deltaBX(0, 1); }
209 
210 long long EventWithHistory::deltaBX(const TinyEvent& se) const { return TinyEvent::deltaBX(se); }
211 
212 long long EventWithHistory::absoluteBX(const unsigned int ev1) const {
213  if (ev1 == 0)
214  return TinyEvent::absoluteBX();
215  if (ev1 <= _prevse.size())
216  return _prevse[ev1 - 1].absoluteBX();
217 
218  return -1;
219 }
220 
222 
223 long long EventWithHistory::absoluteBXinCycle(const unsigned int ev1, const int bx0) const {
224  if (ev1 == 0)
225  return TinyEvent::absoluteBXinCycle(bx0);
226  if (ev1 <= _prevse.size())
227  return _prevse[ev1 - 1].absoluteBXinCycle(bx0);
228 
229  return -1;
230 }
231 
232 long long EventWithHistory::absoluteBXinCycle(const int bx0) const { return TinyEvent::absoluteBXinCycle(bx0); }
233 
234 long long EventWithHistory::deltaBXinCycle(const unsigned int ev2, const unsigned int ev1, const int bx0) const {
235  if (ev2 == ev1 && ev1 <= _prevse.size()) {
236  if (ev2 == 0)
237  return TinyEvent::deltaBXinCycle(*this, bx0);
238  return _prevse[ev2 - 1].deltaBXinCycle(_prevse[ev1 - 1], bx0);
239  }
240 
241  if (ev2 < ev1 && ev1 <= _prevse.size()) {
242  if (ev2 == 0)
243  return TinyEvent::deltaBXinCycle(_prevse[ev1 - 1], bx0);
244  return _prevse[ev2 - 1].deltaBXinCycle(_prevse[ev1 - 1], bx0);
245  }
246 
247  return -1;
248 }
249 
250 long long EventWithHistory::deltaBXinCycle(const unsigned int ev1, const int bx0) const {
251  return deltaBXinCycle(0, ev1, bx0);
252 }
253 
254 long long EventWithHistory::deltaBXinCycle(const int bx0) const { return deltaBXinCycle(0, 1, bx0); }
255 
256 long long EventWithHistory::deltaBXinCycle(const TinyEvent& se, const int bx0) const {
257  return TinyEvent::deltaBXinCycle(se, bx0);
258 }
unsigned int depth() const
size
Write out results.
long long absoluteBXinCycle(const int bx0) const
Definition: TinyEvent.h:47
EventNumber_t event() const
Definition: EventID.h:40
long long deltaBX() const
const unsigned int orbit() const
const TinyEvent * get(const unsigned int ev) const
unsigned int _bx
Definition: TinyEvent.h:66
unsigned long long EventNumber_t
long long absoluteBXinCycle(const unsigned int ev1, const int bx0) const
bool ev
bool equal(const T &first, const T &second)
Definition: Equal.h:32
std::vector< TinyEvent > _prevse
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< L1AcceptBunchCrossing > L1AcceptBunchCrossingCollection
int add(const EventWithHistory &he, const int idepth)
TinyEvent & operator=(const TinyEvent &se)
Definition: TinyEvent.h:27
long long absoluteBX() const
Definition: TinyEvent.h:45
long long deltaBX(const TinyEvent &se) const
Definition: TinyEvent.h:49
const int bx() const
const edm::EventNumber_t event() const
int isNextOf(const TinyEvent &se) const
Definition: TinyEvent.h:43
long long deltaBXinCycle(const unsigned int ev2, const unsigned int ev1, const int bx0) const
bool isFutureHistory() const
edm::EventID id() const
Definition: EventBase.h:59
long long absoluteBX() const
static std::atomic< unsigned int > counter
edm::EventNumber_t _event
Definition: TinyEvent.h:64
EventWithHistory & operator=(const EventWithHistory &he)
unsigned int _orbit
Definition: TinyEvent.h:65
int operator==(const TinyEvent &other) const
Definition: TinyEvent.h:39
int operator==(const EventWithHistory &other) const
Definition: event.py:1
long long deltaBXinCycle(const TinyEvent &se, const int bx0) const
Definition: TinyEvent.h:58