CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EventWithHistory.cc
Go to the documentation of this file.
1 #include <map>
7 
9 
11 
12 EventWithHistory::EventWithHistory(const unsigned int event, const int orbit, const int bx):
13  TinyEvent(event,orbit,bx), _prevse() { }
14 
15 EventWithHistory::EventWithHistory(const unsigned int 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.size()>0) ? 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() > (unsigned int)(-1*l1abc->l1AcceptOffset()) ) {
40  unsigned int 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  while(newbx > 3563) {
46  ++neworbit;
47  newbx -= 3564;
48  }
49  while(newbx < 0) {
50  --neworbit;
51  newbx += 3564;
52  }
53 
54  if(l1abc->eventType()!=0) {
55  TinyEvent tmpse(evnumb,neworbit,newbx);
56  tmpmap[l1abc->l1AcceptOffset()]=tmpse;
57  }
58  else {
59  edm::LogWarning("L1AcceptBunchCrossingNoType") << "L1AcceptBunchCrossing with no type found: ";
60  for(L1AcceptBunchCrossingCollection::const_iterator debu=l1abcc.begin();debu!=l1abcc.end();++debu) {
61  edm::LogPrint("L1AcceptBunchCrossingNoType") << *debu;
62  }
63  }
64  }
65  else {
66  edm::LogError("L1AcceptBunchCrossingOffsetTooLarge") << " Too large orbit offset "
67  << orbitoffset << " "
68  << l1abc->orbitNumber();
69  }
70  }
71  else {
72  edm::LogInfo("L1AcceptBunchCrossingNegativeEvent") << "L1AcceptBunchCrossing with negative event: ";
73  for(L1AcceptBunchCrossingCollection::const_iterator debu=l1abcc.begin();debu!=l1abcc.end();++debu) {
74  edm::LogVerbatim("L1AcceptBunchCrossingNegativeEvent") << *debu;
75  }
76  }
77  }
78  // look for the event itself
79  if(tmpmap.find(0)!=tmpmap.end()) {
80 
81  TinyEvent::operator=(tmpmap[0]);
82 
83  // loop on the rest of the map and stop when it is missing
84  // check that the events are in the right order and break if not
85 
86  int counter=-1;
87  while(tmpmap.find(counter)!=tmpmap.end()) {
88 
89  if(tmpmap[counter+1].deltaBX(tmpmap[counter]) > 0) {
90  _prevse.push_back(tmpmap[counter]);
91  --counter;
92  }
93  else {
94  edm::LogWarning("L1AcceptBunchCrossingNotInOrder") << "L1AcceptBunchCrossing not in order: ";
95  for(L1AcceptBunchCrossingCollection::const_iterator debu=l1abcc.begin();debu!=l1abcc.end();++debu) {
96  edm::LogPrint("L1AcceptBunchCrossingNotInOrder") << *debu;
97  }
98  break;
99  }
100  }
101  }
102  else {
103  TinyEvent::operator=(event);
104  edm::LogWarning("L1AcceptBunchCrossingNoCollection") << " L1AcceptBunchCrossing with offset=0 not found "
105  << " likely L1ABCCollection is empty ";
106  }
107 
108 }
109 
110 EventWithHistory::EventWithHistory(const EventWithHistory& he): TinyEvent(he), _prevse(he._prevse) { }
111 
113 
114  if(this != &he) {
116  _prevse = he._prevse;
117  }
118  return *this;
119 }
120 
121 // int EventWithHistory::operator<(const EventWithHistory& other) const { return TinyEvent::operator<(other); }
122 
124 
125  int equal = TinyEvent::operator==(other);
126 
127  // depth is not checked anymore
128 
129  // equal = equal && (depth() == other.depth());
130 
131  if(equal) {
132  for(unsigned int i=0;i<((depth()<other.depth())?depth():other.depth());i++) {
133  equal = equal && (_prevse[i] == other._prevse[i]);
134  }
135  }
136 
137  return equal;
138 }
139 
140 int EventWithHistory::add(const EventWithHistory& he, const int idepth) {
141 
142  if(!add((const TinyEvent&) he,idepth)) return 0;
143 
144  for(std::vector<TinyEvent>::const_iterator ev=he._prevse.begin();ev!=he._prevse.end();ev++) {
145  if(!add(*ev,idepth)) return 0;
146  }
147  return 1;
148 }
149 
150 int EventWithHistory::add(const TinyEvent& se, const int idepth) {
151 
152  bool isfuture = (idepth <0);
153  const unsigned int adepth = abs(idepth);
154 
155  // protect against the possibility of filling with past and future history
156 
157  if( depth()>0 && ((isfuture && !isFutureHistory()) || (!isfuture && isFutureHistory()))) return 0;
158 
159  if(adepth==0) return 0;
160  if(_prevse.size()>= adepth) return 0;
161 
162  if(_prevse.size()==0) {
163  if((!isfuture && isNextOf(se)) || (isfuture && se.isNextOf(*this))) {
164  _prevse.push_back(se);
165  return 1;
166  }
167  else {
168  return 0;
169  }
170  }
171  else {
172  if((!isfuture && _prevse[_prevse.size()-1].isNextOf(se)) || (isfuture && se.isNextOf(_prevse[_prevse.size()-1]))) {
173  _prevse.push_back(se);
174  return 1;
175  }
176  else {
177  return 0;
178  }
179  }
180  return 0;
181 }
182 
183 const unsigned int EventWithHistory::event() const { return TinyEvent::_event; }
184 const unsigned int EventWithHistory::orbit() const { return TinyEvent::_orbit; }
185 const int EventWithHistory::bx() const { return TinyEvent::_bx; }
186 
187 const TinyEvent* EventWithHistory::get(const unsigned int ev) const {
188 
189  if(ev==0) return this;
190  if(ev<=_prevse.size()) return &_prevse[ev-1];
191  return 0;
192 
193 }
194 
195 unsigned int EventWithHistory::depth() const { return _prevse.size(); }
196 
198 
199  return (depth()>0 && _prevse[0].isNextOf(*this));
200 
201 }
202 
203 long long EventWithHistory::deltaBX(const unsigned int ev2, const unsigned int ev1) const {
204 
205  if(ev2==ev1) return 0;
206 
207  if(ev2<ev1 && ev1<=_prevse.size()) {
208  if(ev2==0) return TinyEvent::deltaBX(_prevse[ev1-1]);
209  return _prevse[ev2-1].deltaBX(_prevse[ev1-1]);
210  }
211 
212  return -1;
213 }
214 
215 long long EventWithHistory::deltaBX(const unsigned int ev1) const { return deltaBX(0,ev1); }
216 
217 long long EventWithHistory::deltaBX() const { return deltaBX(0,1); }
218 
219 long long EventWithHistory::deltaBX(const TinyEvent& se) const {
220 
221  return TinyEvent::deltaBX(se);
222 
223 }
224 
225 long long EventWithHistory::absoluteBX(const unsigned int ev1) const {
226 
227  if(ev1==0) return TinyEvent::absoluteBX();
228  if(ev1<=_prevse.size()) return _prevse[ev1-1].absoluteBX();
229 
230  return -1;
231 
232 }
233 
234 long long EventWithHistory::absoluteBX() const {
235 
236  return TinyEvent::absoluteBX();
237 
238 }
239 
240 long long EventWithHistory::absoluteBXinCycle(const unsigned int ev1, const int bx0) const {
241 
242  if(ev1==0) return TinyEvent::absoluteBXinCycle(bx0);
243  if(ev1<=_prevse.size()) return _prevse[ev1-1].absoluteBXinCycle(bx0);
244 
245  return -1;
246 
247 }
248 
249 long long EventWithHistory::absoluteBXinCycle(const int bx0) const {
250 
251  return TinyEvent::absoluteBXinCycle(bx0);
252 
253 }
254 
255 long long EventWithHistory::deltaBXinCycle(const unsigned int ev2, const unsigned int ev1, const int bx0) const {
256 
257  if(ev2==ev1 && ev1<=_prevse.size()) {
258  if(ev2==0) return TinyEvent::deltaBXinCycle(*this,bx0);
259  return _prevse[ev2-1].deltaBXinCycle(_prevse[ev1-1],bx0);
260  }
261 
262  if(ev2<ev1 && ev1<=_prevse.size()) {
263  if(ev2==0) return TinyEvent::deltaBXinCycle(_prevse[ev1-1],bx0);
264  return _prevse[ev2-1].deltaBXinCycle(_prevse[ev1-1],bx0);
265  }
266 
267  return -1;
268 }
269 
270 long long EventWithHistory::deltaBXinCycle(const unsigned int ev1, const int bx0) const {
271  return deltaBXinCycle(0,ev1,bx0);
272 }
273 
274 long long EventWithHistory::deltaBXinCycle(const int bx0) const {
275  return deltaBXinCycle(0,1,bx0);
276 }
277 
278 long long EventWithHistory::deltaBXinCycle(const TinyEvent& se, const int bx0) const {
279 
280  return TinyEvent::deltaBXinCycle(se,bx0);
281 
282 }
unsigned int depth() const
long long absoluteBXinCycle(const int bx0) const
Definition: TinyEvent.h:63
EventNumber_t event() const
Definition: EventID.h:44
int i
Definition: DBlmapReader.cc:9
long long deltaBX() const
const unsigned int orbit() const
const TinyEvent * get(const unsigned int ev) const
unsigned int _bx
Definition: TinyEvent.h:89
#define abs(x)
Definition: mlp_lapack.h:159
long long absoluteBXinCycle(const unsigned int ev1, const int bx0) const
std::vector< TinyEvent > _prevse
std::vector< L1AcceptBunchCrossing > L1AcceptBunchCrossingCollection
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
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
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:56
long long absoluteBX() const
unsigned int _event
Definition: TinyEvent.h:87
EventWithHistory & operator=(const EventWithHistory &he)
unsigned int _orbit
Definition: TinyEvent.h:88
const unsigned int event() const
int operator==(const TinyEvent &other) const
Definition: TinyEvent.h:45
int operator==(const EventWithHistory &other) const
tuple size
Write out results.
long long deltaBXinCycle(const TinyEvent &se, const int bx0) const
Definition: TinyEvent.h:79