Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:06:35

0001 // -*- C++ -*-
0002 //
0003 // Package:    EventWithHistoryProducerFromL1ABC
0004 // Class:      EventWithHistoryProducerFromL1ABC
0005 //
0006 /**\class EventWithHistoryProducerFromL1ABC EventWithHistoryProducerFromL1ABC.cc DPGAnalysis/SiStripTools/src/EventWithHistoryProducerFromL1ABC.cc
0007 
0008  Description: <one line class summary>
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Andrea Venturi
0015 //         Created:  Tue Jun 30 15:26:20 CET 2009
0016 //
0017 //
0018 
0019 // system include files
0020 #include <memory>
0021 #include <map>
0022 
0023 // user include files
0024 #include "FWCore/Framework/interface/Frameworkfwd.h"
0025 #include "FWCore/Framework/interface/stream/EDProducer.h"
0026 
0027 #include "FWCore/Framework/interface/Event.h"
0028 #include "FWCore/Framework/interface/Run.h"
0029 #include "FWCore/Framework/interface/MakerMacros.h"
0030 
0031 #include "FWCore/Utilities/interface/InputTag.h"
0032 
0033 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0034 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0035 
0036 #include "DataFormats/Scalers/interface/L1AcceptBunchCrossing.h"
0037 #include "DataFormats/TCDS/interface/TCDSRecord.h"
0038 #include "DPGAnalysis/SiStripTools/interface/EventWithHistory.h"
0039 //
0040 // class declarations
0041 //
0042 
0043 class EventWithHistoryProducerFromL1ABC : public edm::stream::EDProducer<> {
0044 public:
0045   explicit EventWithHistoryProducerFromL1ABC(const edm::ParameterSet&);
0046   ~EventWithHistoryProducerFromL1ABC() override;
0047 
0048 private:
0049   void beginRun(const edm::Run&, const edm::EventSetup&) override;
0050   void produce(edm::Event&, const edm::EventSetup&) override;
0051   void endRun(const edm::Run&, const edm::EventSetup&) override;
0052 
0053   // ----------member data ---------------------------
0054 
0055   edm::EDGetTokenT<L1AcceptBunchCrossingCollection> _l1abccollectionToken;
0056   edm::EDGetTokenT<TCDSRecord> _tcdsRecordToken;
0057   const bool _forceNoOffset;
0058   const bool _forceSCAL;
0059   std::map<edm::EventNumber_t, long long> _offsets;
0060   long long _curroffset;
0061   edm::EventNumber_t _curroffevent;
0062 };
0063 
0064 //
0065 // constants, enums and typedefs
0066 //
0067 
0068 //
0069 // static data member definitions
0070 //
0071 
0072 //
0073 // constructors and destructor
0074 //
0075 EventWithHistoryProducerFromL1ABC::EventWithHistoryProducerFromL1ABC(const edm::ParameterSet& iConfig)
0076     : _l1abccollectionToken(
0077           mayConsume<L1AcceptBunchCrossingCollection>(iConfig.getParameter<edm::InputTag>("l1ABCCollection"))),
0078       _tcdsRecordToken(mayConsume<TCDSRecord>(iConfig.getParameter<edm::InputTag>("tcdsRecordLabel"))),
0079       _forceNoOffset(iConfig.getUntrackedParameter<bool>("forceNoOffset", false)),
0080       _forceSCAL(iConfig.getParameter<bool>("forceSCAL")),
0081       _offsets(),
0082       _curroffset(0),
0083       _curroffevent(0) {
0084   if (_forceNoOffset)
0085     edm::LogWarning("NoOffsetComputation") << "Orbit and BX offset will NOT be computed: Be careful!";
0086 
0087   produces<EventWithHistory>();
0088 
0089   //now do what ever other initialization is needed
0090 }
0091 
0092 EventWithHistoryProducerFromL1ABC::~EventWithHistoryProducerFromL1ABC() {
0093   // do anything here that needs to be done at desctruction time
0094   // (e.g. close files, deallocate resources etc.)
0095 }
0096 
0097 //
0098 // member functions
0099 //
0100 
0101 // ------------ method called to produce the data  ------------
0102 void EventWithHistoryProducerFromL1ABC::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0103   using namespace edm;
0104 
0105   if (iEvent.run() < 110878) {
0106     std::unique_ptr<EventWithHistory> pOut(new EventWithHistory(iEvent));
0107     iEvent.put(std::move(pOut));
0108 
0109   } else {
0110     Handle<L1AcceptBunchCrossingCollection> pIn;
0111     iEvent.getByToken(_l1abccollectionToken, pIn);
0112     Handle<TCDSRecord> tcds_pIn;
0113     iEvent.getByToken(_tcdsRecordToken, tcds_pIn);
0114     bool useTCDS(tcds_pIn.isValid() && !_forceSCAL);
0115 
0116     const auto* tcdsRecord = useTCDS ? tcds_pIn.product() : nullptr;
0117     // offset computation
0118 
0119     long long orbitoffset = 0;
0120     int bxoffset = 0;
0121 
0122     if (useTCDS) {
0123       if (!_forceNoOffset) {
0124         if (tcdsRecord->getL1aHistoryEntry(0).getIndex() == 0) {
0125           orbitoffset = (long long)tcdsRecord->getOrbitNr() - (long long)iEvent.orbitNumber();
0126           bxoffset = tcdsRecord->getBXID() - iEvent.bunchCrossing();
0127         }
0128       }
0129     } else {
0130       if (!_forceNoOffset) {
0131         for (L1AcceptBunchCrossingCollection::const_iterator l1abc = pIn->begin(); l1abc != pIn->end(); ++l1abc) {
0132           if (l1abc->l1AcceptOffset() == 0) {
0133             orbitoffset = (long long)l1abc->orbitNumber() - (long long)iEvent.orbitNumber();
0134             bxoffset = l1abc->bunchCrossing() - iEvent.bunchCrossing();
0135           }
0136         }
0137       }
0138     }
0139 
0140     std::unique_ptr<EventWithHistory> pOut(useTCDS ? new EventWithHistory(iEvent, *tcdsRecord, orbitoffset, bxoffset)
0141                                                    : new EventWithHistory(iEvent, *pIn, orbitoffset, bxoffset));
0142     iEvent.put(std::move(pOut));
0143 
0144     // monitor offset
0145 
0146     long long absbxoffset = orbitoffset * 3564 + bxoffset;
0147 
0148     if (_offsets.empty()) {
0149       _curroffset = absbxoffset;
0150       _curroffevent = iEvent.id().event();
0151       _offsets[iEvent.id().event()] = absbxoffset;
0152     } else {
0153       if (_curroffset != absbxoffset || iEvent.id().event() < _curroffevent) {
0154         if (_curroffset != absbxoffset) {
0155           edm::LogInfo("AbsoluteBXOffsetChanged")
0156               << "Absolute BX offset changed from " << _curroffset << " to " << absbxoffset << " at orbit "
0157               << iEvent.orbitNumber() << " and BX " << iEvent.bunchCrossing();
0158           if (useTCDS) {
0159             edm::LogVerbatim("AbsoluteBXOffsetChanged") << *tcdsRecord;  // Not sure about this
0160           } else {
0161             for (L1AcceptBunchCrossingCollection::const_iterator l1abc = pIn->begin(); l1abc != pIn->end(); ++l1abc) {
0162               edm::LogVerbatim("AbsoluteBXOffsetChanged") << *l1abc;
0163             }
0164           }
0165         }
0166 
0167         _curroffset = absbxoffset;
0168         _curroffevent = iEvent.id().event();
0169         _offsets[iEvent.id().event()] = absbxoffset;
0170       }
0171     }
0172   }
0173 }
0174 
0175 void EventWithHistoryProducerFromL1ABC::beginRun(const edm::Run&, const edm::EventSetup&) {
0176   // reset offset vector
0177 
0178   _offsets.clear();
0179   edm::LogInfo("AbsoluteBXOffsetReset") << "Absolute BX offset map reset";
0180 }
0181 
0182 void EventWithHistoryProducerFromL1ABC::endRun(const edm::Run&, const edm::EventSetup&) {
0183   // summary of absolute bx offset vector
0184 
0185   edm::LogInfo("AbsoluteBXOffsetSummary") << "Absolute BX offset summary:";
0186   for (std::map<edm::EventNumber_t, long long>::const_iterator offset = _offsets.begin(); offset != _offsets.end();
0187        ++offset) {
0188     edm::LogVerbatim("AbsoluteBXOffsetSummary") << offset->first << " " << offset->second;
0189   }
0190 }
0191 
0192 //define this as a plug-in
0193 DEFINE_FWK_MODULE(EventWithHistoryProducerFromL1ABC);