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