CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch1/src/DPGAnalysis/SiStripTools/plugins/EventWithHistoryProducerFromL1ABC.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    EventWithHistoryProducerFromL1ABC
00004 // Class:      EventWithHistoryProducerFromL1ABC
00005 // 
00013 //
00014 // Original Author:  Andrea Venturi
00015 //         Created:  Tue Jun 30 15:26:20 CET 2009
00016 //
00017 //
00018 
00019 
00020 // system include files
00021 #include <memory>
00022 #include <map>
00023 
00024 // user include files
00025 #include "FWCore/Framework/interface/Frameworkfwd.h"
00026 #include "FWCore/Framework/interface/EDProducer.h"
00027 
00028 #include "FWCore/Framework/interface/Event.h"
00029 #include "FWCore/Framework/interface/Run.h"
00030 #include "FWCore/Framework/interface/MakerMacros.h"
00031 
00032 #include "FWCore/Utilities/interface/InputTag.h"
00033 
00034 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00035 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00036 
00037 #include "DataFormats/Scalers/interface/L1AcceptBunchCrossing.h"
00038 #include "DPGAnalysis/SiStripTools/interface/EventWithHistory.h"
00039 //
00040 // class decleration
00041 //
00042 
00043 class EventWithHistoryProducerFromL1ABC : public edm::EDProducer {
00044    public:
00045       explicit EventWithHistoryProducerFromL1ABC(const edm::ParameterSet&);
00046       ~EventWithHistoryProducerFromL1ABC();
00047 
00048    private:
00049       virtual void beginJob() ;
00050       virtual void beginRun(edm::Run&, const edm::EventSetup&) ;
00051       virtual void produce(edm::Event&, const edm::EventSetup&);
00052       virtual void endRun(edm::Run&, const edm::EventSetup&) ;
00053       virtual void endJob() ;
00054       
00055       // ----------member data ---------------------------
00056 
00057   edm::InputTag _l1abccollection;
00058   const bool _forceNoOffset;
00059   std::map<unsigned int, long long> _offsets;
00060   long long _curroffset;
00061   unsigned int _curroffevent;
00062 };
00063 
00064 //
00065 // constants, enums and typedefs
00066 //
00067 
00068 
00069 //
00070 // static data member definitions
00071 //
00072 
00073 //
00074 // constructors and destructor
00075 //
00076 EventWithHistoryProducerFromL1ABC::EventWithHistoryProducerFromL1ABC(const edm::ParameterSet& iConfig):
00077   _l1abccollection(iConfig.getParameter<edm::InputTag>("l1ABCCollection")),
00078   _forceNoOffset(iConfig.getUntrackedParameter<bool>("forceNoOffset",false)),
00079   _offsets(), _curroffset(0), _curroffevent(0)
00080 {
00081 
00082   if(_forceNoOffset) edm::LogWarning("NoOffsetComputation") << "Orbit and BX offset will NOT be computed: Be careful!";
00083 
00084   produces<EventWithHistory>();
00085    
00086    //now do what ever other initialization is needed
00087   
00088 }
00089 
00090 
00091 EventWithHistoryProducerFromL1ABC::~EventWithHistoryProducerFromL1ABC()
00092 {
00093  
00094    // do anything here that needs to be done at desctruction time
00095    // (e.g. close files, deallocate resources etc.)
00096 
00097 }
00098 
00099 
00100 //
00101 // member functions
00102 //
00103 
00104 // ------------ method called to produce the data  ------------
00105 void
00106 EventWithHistoryProducerFromL1ABC::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00107 {
00108    using namespace edm;
00109 
00110    if(iEvent.run() < 110878 ) {
00111 
00112      std::auto_ptr<EventWithHistory> pOut(new EventWithHistory(iEvent));
00113      iEvent.put(pOut);
00114 
00115    }
00116    else {
00117 
00118      Handle<L1AcceptBunchCrossingCollection > pIn;
00119      iEvent.getByLabel(_l1abccollection,pIn);
00120      
00121      // offset computation
00122      
00123      long long orbitoffset = 0;
00124      int bxoffset = 0;
00125      if(!_forceNoOffset) {
00126        for(L1AcceptBunchCrossingCollection::const_iterator l1abc=pIn->begin();l1abc!=pIn->end();++l1abc) {
00127          if(l1abc->l1AcceptOffset()==0) {
00128            orbitoffset = (long long)l1abc->orbitNumber() - (long long)iEvent.orbitNumber();
00129            bxoffset = l1abc->bunchCrossing() - iEvent.bunchCrossing();
00130          }
00131        }
00132      }
00133      
00134      
00135      std::auto_ptr<EventWithHistory> pOut(new EventWithHistory(iEvent,*pIn,orbitoffset,bxoffset));
00136      iEvent.put(pOut);
00137      
00138      // monitor offset
00139      
00140      long long absbxoffset = orbitoffset*3564 + bxoffset;
00141      
00142      if(_offsets.size()==0) {
00143        _curroffset = absbxoffset;
00144        _curroffevent = iEvent.id().event();
00145        _offsets[iEvent.id().event()] = absbxoffset;
00146      }
00147      else {
00148        if(_curroffset != absbxoffset || iEvent.id().event() < _curroffevent ) {
00149 
00150          if( _curroffset != absbxoffset) {
00151            edm::LogInfo("AbsoluteBXOffsetChanged") << "Absolute BX offset changed from " 
00152                                                    << _curroffset << " to "
00153                                                    << absbxoffset << " at orbit "
00154                                                    << iEvent.orbitNumber() << " and BX "
00155                                                    << iEvent.bunchCrossing();
00156            for(L1AcceptBunchCrossingCollection::const_iterator l1abc=pIn->begin();l1abc!=pIn->end();++l1abc) {
00157              edm::LogVerbatim("AbsoluteBXOffsetChanged") << *l1abc;
00158            }
00159          }
00160 
00161          _curroffset = absbxoffset;
00162          _curroffevent = iEvent.id().event();
00163          _offsets[iEvent.id().event()] = absbxoffset;
00164        }
00165      }
00166    }
00167 }
00168 
00169 // ------------ method called once each job just before starting event loop  ------------
00170 void 
00171 EventWithHistoryProducerFromL1ABC::beginJob()
00172 {
00173 }
00174 
00175 void 
00176 EventWithHistoryProducerFromL1ABC::beginRun(edm::Run&, const edm::EventSetup&)
00177 {
00178   // reset offset vector
00179 
00180   _offsets.clear();
00181   edm::LogInfo("AbsoluteBXOffsetReset") << "Absolute BX offset map reset";
00182 
00183 }
00184 
00185 void 
00186 EventWithHistoryProducerFromL1ABC::endRun(edm::Run&, const edm::EventSetup&)
00187 {
00188   // summary of absolute bx offset vector
00189 
00190   edm::LogInfo("AbsoluteBXOffsetSummary") << "Absolute BX offset summary:";
00191   for(std::map<unsigned int, long long>::const_iterator offset=_offsets.begin();offset!=_offsets.end();++offset) {
00192     edm::LogVerbatim("AbsoluteBXOffsetSummary") << offset->first << " " << offset->second;
00193   }
00194 
00195 }
00196 
00197 // ------------ method called once each job just after ending the event loop  ------------
00198 void 
00199 EventWithHistoryProducerFromL1ABC::endJob() {
00200 }
00201 
00202 //define this as a plug-in
00203 DEFINE_FWK_MODULE(EventWithHistoryProducerFromL1ABC);