File indexing completed on 2024-04-06 12:06:42
0001 #include <map>
0002 #include "FWCore/Framework/interface/Event.h"
0003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0004 #include "DataFormats/Provenance/interface/EventAuxiliary.h"
0005 #include "DataFormats/Scalers/interface/L1AcceptBunchCrossing.h"
0006 #include "DPGAnalysis/SiStripTools/interface/EventWithHistory.h"
0007 #include "DataFormats/TCDS/interface/TCDSRecord.h"
0008
0009 EventWithHistory::EventWithHistory() : TinyEvent(), _prevse() {}
0010
0011 EventWithHistory::EventWithHistory(const TinyEvent& se) : TinyEvent(se), _prevse() {}
0012
0013 EventWithHistory::EventWithHistory(const edm::EventNumber_t event, const int orbit, const int bx)
0014 : TinyEvent(event, orbit, bx), _prevse() {}
0015
0016 EventWithHistory::EventWithHistory(const edm::EventNumber_t event, const unsigned int orbit, const int bx)
0017 : TinyEvent(event, orbit, bx), _prevse() {}
0018
0019 EventWithHistory::EventWithHistory(const edm::Event& event) : TinyEvent(event), _prevse() {}
0020
0021 EventWithHistory::EventWithHistory(const std::vector<edm::EventAuxiliary>& veaux)
0022 : TinyEvent((!veaux.empty()) ? veaux[veaux.size() - 1] : TinyEvent()), _prevse() {
0023 for (std::vector<edm::EventAuxiliary>::const_reverse_iterator eaux = veaux.rbegin(); eaux != veaux.rend(); eaux++) {
0024 if (eaux != veaux.rbegin()) {
0025 _prevse.push_back(*eaux);
0026 }
0027 }
0028 }
0029
0030 EventWithHistory::EventWithHistory(const edm::Event& event,
0031 const L1AcceptBunchCrossingCollection& l1abcc,
0032 const long long orbitoffset,
0033 const int bxoffset)
0034 : TinyEvent(), _prevse() {
0035 std::map<int, TinyEvent> tmpmap;
0036
0037 for (L1AcceptBunchCrossingCollection::const_iterator l1abc = l1abcc.begin(); l1abc != l1abcc.end(); ++l1abc) {
0038 if (event.id().event() > (edm::EventNumber_t)(-1 * l1abc->l1AcceptOffset())) {
0039 edm::EventNumber_t evnumb = event.id().event() + l1abc->l1AcceptOffset();
0040 if (orbitoffset < (long long)l1abc->orbitNumber()) {
0041 unsigned int neworbit = l1abc->orbitNumber() - orbitoffset;
0042 int newbx = l1abc->bunchCrossing() - bxoffset;
0043
0044
0045
0046
0047
0048 while (newbx > 3563) {
0049 ++neworbit;
0050 newbx -= 3564;
0051 }
0052 while (newbx < 0) {
0053 --neworbit;
0054 newbx += 3564;
0055 }
0056
0057 if (l1abc->eventType() != 0) {
0058 TinyEvent tmpse(evnumb, neworbit, newbx);
0059 tmpmap[l1abc->l1AcceptOffset()] = tmpse;
0060 } else {
0061 edm::LogWarning("L1AcceptBunchCrossingNoType") << "L1AcceptBunchCrossing with no type found: ";
0062 for (L1AcceptBunchCrossingCollection::const_iterator debu = l1abcc.begin(); debu != l1abcc.end(); ++debu) {
0063 edm::LogPrint("L1AcceptBunchCrossingNoType") << *debu;
0064 }
0065 }
0066 } else {
0067 edm::LogError("L1AcceptBunchCrossingOffsetTooLarge")
0068 << " Too large orbit offset " << orbitoffset << " " << l1abc->orbitNumber();
0069 }
0070 } else {
0071 edm::LogInfo("L1AcceptBunchCrossingNegativeEvent") << "L1AcceptBunchCrossing with negative event: ";
0072 for (L1AcceptBunchCrossingCollection::const_iterator debu = l1abcc.begin(); debu != l1abcc.end(); ++debu) {
0073 edm::LogVerbatim("L1AcceptBunchCrossingNegativeEvent") << *debu;
0074 }
0075 }
0076 }
0077
0078 if (tmpmap.find(0) != tmpmap.end()) {
0079 TinyEvent::operator=(tmpmap[0]);
0080
0081
0082
0083
0084 int counter = -1;
0085 while (tmpmap.find(counter) != tmpmap.end()) {
0086 if (tmpmap[counter + 1].deltaBX(tmpmap[counter]) > 0) {
0087 _prevse.push_back(tmpmap[counter]);
0088 --counter;
0089 } else {
0090 edm::LogWarning("L1AcceptBunchCrossingNotInOrder")
0091 << "L1AcceptBunchCrossing not in order: orbit " << event.orbitNumber() << " BX " << event.bunchCrossing()
0092 << " orbit offset " << orbitoffset << " bx offset " << bxoffset << " :";
0093 for (L1AcceptBunchCrossingCollection::const_iterator debu = l1abcc.begin(); debu != l1abcc.end(); ++debu) {
0094 edm::LogPrint("L1AcceptBunchCrossingNotInOrder") << *debu;
0095 }
0096 break;
0097 }
0098 }
0099 } else {
0100 TinyEvent::operator=(event);
0101 edm::LogWarning("L1AcceptBunchCrossingNoCollection") << " L1AcceptBunchCrossing with offset=0 not found "
0102 << " likely L1ABCCollection is empty ";
0103 }
0104 }
0105
0106 EventWithHistory::EventWithHistory(const edm::Event& event,
0107 const TCDSRecord& tcdsRecord,
0108 const long long orbitoffset,
0109 const int bxoffset)
0110 : TinyEvent(), _prevse() {
0111 std::map<int, TinyEvent> tmpmap;
0112
0113
0114 auto l1aHistory = tcdsRecord.getFullL1aHistory();
0115 const L1aInfo infoForIndex0 = L1aInfo(0, tcdsRecord.getOrbitNr(), tcdsRecord.getBXID(), tcdsRecord.getEventType());
0116 l1aHistory.push_back(infoForIndex0);
0117
0118 for (auto l1a : l1aHistory) {
0119 int l1AcceptOffset(l1a.getIndex());
0120
0121 if (event.id().event() > (edm::EventNumber_t)(-1 * l1AcceptOffset)) {
0122 edm::EventNumber_t evnumb = event.id().event() + l1AcceptOffset;
0123 if (orbitoffset < (long long)l1a.getOrbitNr()) {
0124 unsigned int neworbit = l1a.getOrbitNr() - orbitoffset;
0125 int newbx = l1a.getBXID() - bxoffset;
0126
0127
0128
0129
0130
0131
0132 while (newbx > 3563) {
0133 ++neworbit;
0134 newbx -= 3564;
0135 }
0136 while (newbx < 0) {
0137 --neworbit;
0138 newbx += 3564;
0139 }
0140
0141 if (tcdsRecord.getEventType() != 0) {
0142 TinyEvent tmpse(evnumb, neworbit, newbx);
0143 tmpmap[l1AcceptOffset] = tmpse;
0144 } else {
0145 edm::LogWarning("L1AcceptBunchCrossingNoType") << "L1AcceptBunchCrossing with no type found: ";
0146 edm::LogPrint("L1AcceptBunchCrossingNoType") << &tcdsRecord;
0147 }
0148 } else {
0149 edm::LogError("L1AcceptBunchCrossingOffsetTooLarge")
0150 << " Too large orbit offset " << orbitoffset << " " << tcdsRecord.getOrbitNr();
0151 }
0152 } else {
0153 edm::LogInfo("L1AcceptBunchCrossingNegativeEvent") << "L1AcceptBunchCrossing with negative event: ";
0154 edm::LogVerbatim("L1AcceptBunchCrossingNegativeEvent") << &tcdsRecord;
0155 }
0156 }
0157
0158
0159 if (tmpmap.find(0) != tmpmap.end()) {
0160 TinyEvent::operator=(tmpmap[0]);
0161
0162
0163
0164
0165 int counter = -1;
0166 while (tmpmap.find(counter) != tmpmap.end()) {
0167 if (tmpmap[counter + 1].deltaBX(tmpmap[counter]) > 0) {
0168 _prevse.push_back(tmpmap[counter]);
0169 --counter;
0170 } else {
0171 edm::LogWarning("L1AcceptBunchCrossingNotInOrder")
0172 << "L1AcceptBunchCrossing not in order: orbit " << event.orbitNumber() << " BX " << event.bunchCrossing()
0173 << " orbit offset " << orbitoffset << " bx offset " << bxoffset << " :";
0174 edm::LogPrint("L1AcceptBunchCrossingNotInOrder") << &tcdsRecord;
0175 break;
0176 }
0177 }
0178 } else {
0179 TinyEvent::operator=(event);
0180 edm::LogWarning("L1AcceptBunchCrossingNoCollection") << " L1AcceptBunchCrossing with offset=0 not found "
0181 << " likely TCDSRecord is empty ";
0182 }
0183 }
0184
0185 EventWithHistory::EventWithHistory(const EventWithHistory& he) : TinyEvent(he), _prevse(he._prevse) {}
0186
0187 EventWithHistory& EventWithHistory::operator=(const EventWithHistory& he) {
0188 if (this != &he) {
0189 TinyEvent::operator=(he);
0190 _prevse = he._prevse;
0191 }
0192 return *this;
0193 }
0194
0195
0196
0197 int EventWithHistory::operator==(const EventWithHistory& other) const {
0198 int equal = TinyEvent::operator==(other);
0199
0200
0201
0202
0203
0204 if (equal) {
0205 for (unsigned int i = 0; i < ((depth() < other.depth()) ? depth() : other.depth()); i++) {
0206 equal = equal && (_prevse[i] == other._prevse[i]);
0207 }
0208 }
0209
0210 return equal;
0211 }
0212
0213 int EventWithHistory::add(const EventWithHistory& he, const int idepth) {
0214 if (!add((const TinyEvent&)he, idepth))
0215 return 0;
0216
0217 for (std::vector<TinyEvent>::const_iterator ev = he._prevse.begin(); ev != he._prevse.end(); ev++) {
0218 if (!add(*ev, idepth))
0219 return 0;
0220 }
0221 return 1;
0222 }
0223
0224 int EventWithHistory::add(const TinyEvent& se, const int idepth) {
0225 bool isfuture = (idepth < 0);
0226 const unsigned int adepth = abs(idepth);
0227
0228
0229
0230 if (depth() > 0 && ((isfuture && !isFutureHistory()) || (!isfuture && isFutureHistory())))
0231 return 0;
0232
0233 if (adepth == 0)
0234 return 0;
0235 if (_prevse.size() >= adepth)
0236 return 0;
0237
0238 if (_prevse.empty()) {
0239 if ((!isfuture && isNextOf(se)) || (isfuture && se.isNextOf(*this))) {
0240 _prevse.push_back(se);
0241 return 1;
0242 } else {
0243 return 0;
0244 }
0245 } else {
0246 if ((!isfuture && _prevse[_prevse.size() - 1].isNextOf(se)) ||
0247 (isfuture && se.isNextOf(_prevse[_prevse.size() - 1]))) {
0248 _prevse.push_back(se);
0249 return 1;
0250 } else {
0251 return 0;
0252 }
0253 }
0254 return 0;
0255 }
0256
0257 const edm::EventNumber_t EventWithHistory::event() const { return TinyEvent::_event; }
0258 const unsigned int EventWithHistory::orbit() const { return TinyEvent::_orbit; }
0259 const int EventWithHistory::bx() const { return TinyEvent::_bx; }
0260
0261 const TinyEvent* EventWithHistory::get(const unsigned int ev) const {
0262 if (ev == 0)
0263 return this;
0264 if (ev <= _prevse.size())
0265 return &_prevse[ev - 1];
0266 return nullptr;
0267 }
0268
0269 unsigned int EventWithHistory::depth() const { return _prevse.size(); }
0270
0271 bool EventWithHistory::isFutureHistory() const { return (depth() > 0 && _prevse[0].isNextOf(*this)); }
0272
0273 long long EventWithHistory::deltaBX(const unsigned int ev2, const unsigned int ev1) const {
0274 if (ev2 == ev1)
0275 return 0;
0276
0277 if (ev2 < ev1 && ev1 <= _prevse.size()) {
0278 if (ev2 == 0)
0279 return TinyEvent::deltaBX(_prevse[ev1 - 1]);
0280 return _prevse[ev2 - 1].deltaBX(_prevse[ev1 - 1]);
0281 }
0282
0283 return -1;
0284 }
0285
0286 long long EventWithHistory::deltaBX(const unsigned int ev1) const { return deltaBX(0, ev1); }
0287
0288 long long EventWithHistory::deltaBX() const { return deltaBX(0, 1); }
0289
0290 long long EventWithHistory::deltaBX(const TinyEvent& se) const { return TinyEvent::deltaBX(se); }
0291
0292 long long EventWithHistory::absoluteBX(const unsigned int ev1) const {
0293 if (ev1 == 0)
0294 return TinyEvent::absoluteBX();
0295 if (ev1 <= _prevse.size())
0296 return _prevse[ev1 - 1].absoluteBX();
0297
0298 return -1;
0299 }
0300
0301 long long EventWithHistory::absoluteBX() const { return TinyEvent::absoluteBX(); }
0302
0303 long long EventWithHistory::absoluteBXinCycle(const unsigned int ev1, const int bx0) const {
0304 if (ev1 == 0)
0305 return TinyEvent::absoluteBXinCycle(bx0);
0306 if (ev1 <= _prevse.size())
0307 return _prevse[ev1 - 1].absoluteBXinCycle(bx0);
0308
0309 return -1;
0310 }
0311
0312 long long EventWithHistory::absoluteBXinCycle(const int bx0) const { return TinyEvent::absoluteBXinCycle(bx0); }
0313
0314 long long EventWithHistory::deltaBXinCycle(const unsigned int ev2, const unsigned int ev1, const int bx0) const {
0315 if (ev2 == ev1 && ev1 <= _prevse.size()) {
0316 if (ev2 == 0)
0317 return TinyEvent::deltaBXinCycle(*this, bx0);
0318 return _prevse[ev2 - 1].deltaBXinCycle(_prevse[ev1 - 1], bx0);
0319 }
0320
0321 if (ev2 < ev1 && ev1 <= _prevse.size()) {
0322 if (ev2 == 0)
0323 return TinyEvent::deltaBXinCycle(_prevse[ev1 - 1], bx0);
0324 return _prevse[ev2 - 1].deltaBXinCycle(_prevse[ev1 - 1], bx0);
0325 }
0326
0327 return -1;
0328 }
0329
0330 long long EventWithHistory::deltaBXinCycle(const unsigned int ev1, const int bx0) const {
0331 return deltaBXinCycle(0, ev1, bx0);
0332 }
0333
0334 long long EventWithHistory::deltaBXinCycle(const int bx0) const { return deltaBXinCycle(0, 1, bx0); }
0335
0336 long long EventWithHistory::deltaBXinCycle(const TinyEvent& se, const int bx0) const {
0337 return TinyEvent::deltaBXinCycle(se, bx0);
0338 }