CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/DPGAnalysis/SiStripTools/src/EventWithHistory.cc

Go to the documentation of this file.
00001 #include <map>
00002 #include "FWCore/Framework/interface/Event.h"
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 #include "DataFormats/Provenance/interface/EventAuxiliary.h"
00005 #include "DataFormats/Scalers/interface/L1AcceptBunchCrossing.h"
00006 #include "DPGAnalysis/SiStripTools/interface/EventWithHistory.h"
00007 
00008 EventWithHistory::EventWithHistory(): TinyEvent(), _prevse() { }
00009 
00010 EventWithHistory::EventWithHistory(const TinyEvent& se): TinyEvent(se), _prevse() { }
00011 
00012 EventWithHistory::EventWithHistory(const unsigned int event, const int orbit, const int bx): 
00013   TinyEvent(event,orbit,bx), _prevse() { }
00014 
00015 EventWithHistory::EventWithHistory(const unsigned int event, const unsigned int orbit, const int bx): 
00016   TinyEvent(event,orbit,bx), _prevse() { }
00017 
00018 EventWithHistory::EventWithHistory(const edm::Event& event): 
00019   TinyEvent(event), _prevse() { }
00020 
00021 EventWithHistory::EventWithHistory(const std::vector<edm::EventAuxiliary>& veaux):
00022   TinyEvent((veaux.size()>0) ? veaux[veaux.size()-1]: TinyEvent()), _prevse() 
00023 {
00024   for(std::vector<edm::EventAuxiliary>::const_reverse_iterator eaux=veaux.rbegin();eaux!=veaux.rend();eaux++) {
00025     if(eaux!=veaux.rbegin()) {
00026       _prevse.push_back(*eaux);
00027     }
00028   }
00029 }
00030 
00031 EventWithHistory::EventWithHistory(const edm::Event& event, const L1AcceptBunchCrossingCollection& l1abcc,
00032                                    const long long orbitoffset, const int bxoffset):
00033   TinyEvent(), _prevse()
00034 {
00035   
00036   std::map<int,TinyEvent> tmpmap;
00037   
00038   for(L1AcceptBunchCrossingCollection::const_iterator l1abc=l1abcc.begin();l1abc!=l1abcc.end();++l1abc) {
00039     if(event.id().event() > (unsigned int)(-1*l1abc->l1AcceptOffset()) ) {
00040       unsigned int evnumb = event.id().event()+l1abc->l1AcceptOffset();
00041       if(orbitoffset < (long long)l1abc->orbitNumber()) {
00042         unsigned int neworbit = l1abc->orbitNumber() - orbitoffset;
00043         int newbx = l1abc->bunchCrossing() - bxoffset;
00044         
00045         while(newbx > 3563) {
00046           ++neworbit;
00047           newbx -= 3564;
00048         }
00049         while(newbx < 0) {
00050           --neworbit;
00051           newbx += 3564;
00052         }
00053         
00054         if(l1abc->eventType()!=0) {
00055           TinyEvent tmpse(evnumb,neworbit,newbx);
00056           tmpmap[l1abc->l1AcceptOffset()]=tmpse;
00057         }
00058         else {
00059           edm::LogWarning("L1AcceptBunchCrossingNoType") << "L1AcceptBunchCrossing with no type found: ";
00060           for(L1AcceptBunchCrossingCollection::const_iterator debu=l1abcc.begin();debu!=l1abcc.end();++debu) {
00061             edm::LogPrint("L1AcceptBunchCrossingNoType") << *debu;
00062           }
00063         }
00064       }
00065       else {
00066         edm::LogError("L1AcceptBunchCrossingOffsetTooLarge") << " Too large orbit offset "
00067                                                              << orbitoffset << " " 
00068                                                              << l1abc->orbitNumber();
00069       }
00070     }
00071     else {
00072       edm::LogInfo("L1AcceptBunchCrossingNegativeEvent") << "L1AcceptBunchCrossing with negative event: ";
00073       for(L1AcceptBunchCrossingCollection::const_iterator debu=l1abcc.begin();debu!=l1abcc.end();++debu) {
00074         edm::LogVerbatim("L1AcceptBunchCrossingNegativeEvent") << *debu;
00075       }
00076     }
00077   }
00078   // look for the event itself
00079   if(tmpmap.find(0)!=tmpmap.end()) {
00080     
00081     TinyEvent::operator=(tmpmap[0]);
00082   
00083     // loop on the rest of the map and stop when it is missing
00084     // check that the events are in the right order and break if not
00085 
00086     int counter=-1;
00087     while(tmpmap.find(counter)!=tmpmap.end()) {
00088 
00089       if(tmpmap[counter+1].deltaBX(tmpmap[counter]) > 0) {
00090         _prevse.push_back(tmpmap[counter]);
00091         --counter;
00092       }
00093       else {
00094         edm::LogWarning("L1AcceptBunchCrossingNotInOrder") << "L1AcceptBunchCrossing not in order: ";
00095         for(L1AcceptBunchCrossingCollection::const_iterator debu=l1abcc.begin();debu!=l1abcc.end();++debu) {
00096           edm::LogPrint("L1AcceptBunchCrossingNotInOrder") << *debu;
00097         }
00098         break;
00099       }
00100     }
00101   }
00102   else {
00103     TinyEvent::operator=(event);
00104     edm::LogWarning("L1AcceptBunchCrossingNoCollection") << " L1AcceptBunchCrossing with offset=0 not found "
00105                                                          << " likely L1ABCCollection is empty ";
00106   }
00107   
00108 }
00109 
00110 EventWithHistory::EventWithHistory(const EventWithHistory& he): TinyEvent(he), _prevse(he._prevse) { }
00111 
00112 EventWithHistory& EventWithHistory::operator=(const EventWithHistory& he) {
00113 
00114   if(this != &he) {
00115     TinyEvent::operator=(he);
00116     _prevse = he._prevse;
00117   }
00118   return *this;
00119 }
00120 
00121 // int EventWithHistory::operator<(const EventWithHistory& other) const { return TinyEvent::operator<(other); }
00122 
00123 int EventWithHistory::operator==(const EventWithHistory& other) const {
00124 
00125   int equal = TinyEvent::operator==(other);
00126 
00127   // depth is not checked anymore
00128 
00129   //  equal = equal && (depth() == other.depth());
00130 
00131   if(equal) {
00132     for(unsigned int i=0;i<((depth()<other.depth())?depth():other.depth());i++) {
00133       equal = equal && (_prevse[i] == other._prevse[i]);
00134     }
00135   }
00136 
00137   return equal;
00138 }
00139 
00140 int EventWithHistory::add(const EventWithHistory& he, const int idepth) {
00141 
00142   if(!add((const TinyEvent&) he,idepth)) return 0;
00143 
00144   for(std::vector<TinyEvent>::const_iterator ev=he._prevse.begin();ev!=he._prevse.end();ev++) {
00145     if(!add(*ev,idepth)) return 0;
00146   }
00147   return 1;
00148 }
00149 
00150 int EventWithHistory::add(const TinyEvent& se, const int idepth) {
00151 
00152   bool isfuture = (idepth <0);
00153   const unsigned int adepth = abs(idepth);
00154 
00155   // protect against the possibility of filling with past and future history
00156 
00157   if( depth()>0 && ((isfuture && !isFutureHistory()) || (!isfuture && isFutureHistory()))) return 0; 
00158 
00159   if(adepth==0) return 0;
00160   if(_prevse.size()>= adepth) return 0;
00161 
00162   if(_prevse.size()==0) {
00163     if((!isfuture && isNextOf(se)) || (isfuture && se.isNextOf(*this))) {
00164       _prevse.push_back(se);
00165       return 1;
00166     }
00167     else {
00168       return 0;
00169     }
00170   }
00171   else {
00172     if((!isfuture && _prevse[_prevse.size()-1].isNextOf(se)) || (isfuture && se.isNextOf(_prevse[_prevse.size()-1]))) {
00173       _prevse.push_back(se);
00174       return 1;
00175     }
00176     else {
00177       return 0;
00178     }
00179   }    
00180   return 0;
00181 }
00182 
00183 const unsigned int EventWithHistory::event() const { return TinyEvent::_event; } 
00184 const unsigned int EventWithHistory::orbit() const { return TinyEvent::_orbit; } 
00185 const int EventWithHistory::bx() const { return TinyEvent::_bx; } 
00186 
00187 const TinyEvent* EventWithHistory::get(const unsigned int ev) const {
00188 
00189   if(ev==0) return this;
00190   if(ev<=_prevse.size()) return &_prevse[ev-1];
00191   return 0;
00192 
00193 }
00194 
00195 unsigned int EventWithHistory::depth() const { return _prevse.size(); }
00196 
00197 bool EventWithHistory::isFutureHistory() const { 
00198 
00199   return (depth()>0 && _prevse[0].isNextOf(*this));
00200 
00201 }
00202  
00203 long long EventWithHistory::deltaBX(const unsigned int ev2, const unsigned int ev1) const {
00204 
00205   if(ev2==ev1) return 0;
00206 
00207   if(ev2<ev1 && ev1<=_prevse.size()) {
00208     if(ev2==0) return  TinyEvent::deltaBX(_prevse[ev1-1]);
00209     return _prevse[ev2-1].deltaBX(_prevse[ev1-1]);
00210   }
00211 
00212   return -1;
00213 }
00214 
00215 long long EventWithHistory::deltaBX(const unsigned int ev1) const {  return deltaBX(0,ev1); }
00216 
00217 long long EventWithHistory::deltaBX() const {  return deltaBX(0,1); }
00218 
00219 long long EventWithHistory::deltaBX(const TinyEvent& se) const { 
00220 
00221  return TinyEvent::deltaBX(se); 
00222 
00223 }
00224 
00225 long long EventWithHistory::absoluteBX(const unsigned int ev1) const {
00226 
00227   if(ev1==0) return TinyEvent::absoluteBX();
00228   if(ev1<=_prevse.size()) return _prevse[ev1-1].absoluteBX();
00229 
00230   return -1;
00231 
00232 }
00233 
00234 long long EventWithHistory::absoluteBX() const {
00235 
00236   return TinyEvent::absoluteBX();
00237 
00238 }
00239 
00240 long long EventWithHistory::absoluteBXinCycle(const unsigned int ev1, const int bx0) const {
00241 
00242   if(ev1==0) return TinyEvent::absoluteBXinCycle(bx0);
00243   if(ev1<=_prevse.size()) return _prevse[ev1-1].absoluteBXinCycle(bx0);
00244 
00245   return -1;
00246 
00247 }
00248 
00249 long long EventWithHistory::absoluteBXinCycle(const int bx0) const {
00250 
00251   return TinyEvent::absoluteBXinCycle(bx0);
00252 
00253 }
00254 
00255 long long EventWithHistory::deltaBXinCycle(const unsigned int ev2, const unsigned int ev1, const int bx0) const {
00256 
00257   if(ev2==ev1 && ev1<=_prevse.size()) {
00258     if(ev2==0) return TinyEvent::deltaBXinCycle(*this,bx0);
00259     return _prevse[ev2-1].deltaBXinCycle(_prevse[ev1-1],bx0);
00260   }
00261 
00262   if(ev2<ev1 && ev1<=_prevse.size()) {
00263     if(ev2==0) return  TinyEvent::deltaBXinCycle(_prevse[ev1-1],bx0);
00264     return _prevse[ev2-1].deltaBXinCycle(_prevse[ev1-1],bx0);
00265   }
00266 
00267   return -1;
00268 }
00269 
00270 long long EventWithHistory::deltaBXinCycle(const unsigned int ev1, const int bx0) const {
00271   return deltaBXinCycle(0,ev1,bx0);
00272 }
00273 
00274 long long EventWithHistory::deltaBXinCycle(const int bx0) const {
00275   return deltaBXinCycle(0,1,bx0);
00276 }
00277 
00278 long long EventWithHistory::deltaBXinCycle(const TinyEvent& se, const int bx0) const {
00279 
00280   return TinyEvent::deltaBXinCycle(se,bx0);
00281 
00282 }