File indexing completed on 2024-04-06 12:06:35
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include <memory>
0021
0022
0023 #include "TH2F.h"
0024 #include "TProfile.h"
0025 #include "DPGAnalysis/SiStripTools/interface/RunHistogramManager.h"
0026
0027 #include <vector>
0028 #include <string>
0029
0030 #include "FWCore/Framework/interface/Frameworkfwd.h"
0031 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0032
0033 #include "FWCore/Framework/interface/Event.h"
0034 #include "FWCore/Framework/interface/Run.h"
0035 #include "FWCore/Framework/interface/MakerMacros.h"
0036
0037 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0038 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0039
0040 #include "FWCore/Utilities/interface/InputTag.h"
0041
0042 #include "DataFormats/Scalers/interface/L1AcceptBunchCrossing.h"
0043
0044
0045
0046
0047 class L1ABCDebugger : public edm::one::EDAnalyzer<edm::one::WatchRuns> {
0048 public:
0049 explicit L1ABCDebugger(const edm::ParameterSet&);
0050 ~L1ABCDebugger() override;
0051
0052 private:
0053 void beginJob() override;
0054 void analyze(const edm::Event&, const edm::EventSetup&) override;
0055 void beginRun(const edm::Run&, const edm::EventSetup&) override;
0056 void endRun(const edm::Run&, const edm::EventSetup&) override {}
0057 void endJob() override;
0058
0059
0060
0061 edm::EDGetTokenT<L1AcceptBunchCrossingCollection> m_l1abccollectionToken;
0062 const unsigned int m_maxLS;
0063 const unsigned int m_LSfrac;
0064
0065 RunHistogramManager m_rhm;
0066
0067 TH2F** m_hoffsets;
0068 TProfile** m_horboffvsorb;
0069 TProfile** m_hbxoffvsorb;
0070 };
0071
0072
0073
0074
0075
0076
0077
0078
0079
0080
0081
0082
0083 L1ABCDebugger::L1ABCDebugger(const edm::ParameterSet& iConfig)
0084 : m_l1abccollectionToken(
0085 consumes<L1AcceptBunchCrossingCollection>(iConfig.getParameter<edm::InputTag>("l1ABCCollection"))),
0086 m_maxLS(iConfig.getUntrackedParameter<unsigned int>("maxLSBeforeRebin", 250)),
0087 m_LSfrac(iConfig.getUntrackedParameter<unsigned int>("startingLSFraction", 16)),
0088 m_rhm(consumesCollector()) {
0089
0090
0091 m_hoffsets = m_rhm.makeTH2F(
0092 "offsets", "Orbit vs BX offsets between SCAL and Event", 2 * 3564 + 1, -3564.5, 3564.5, 201, -100.5, 100.5);
0093 m_horboffvsorb =
0094 m_rhm.makeTProfile("orboffvsorb", "SCAL Orbit offset vs orbit number", m_LSfrac * m_maxLS, 0, m_maxLS * 262144);
0095 m_hbxoffvsorb =
0096 m_rhm.makeTProfile("bxoffvsorb", "SCAL BX offset vs orbit number", m_LSfrac * m_maxLS, 0, m_maxLS * 262144);
0097 }
0098
0099 L1ABCDebugger::~L1ABCDebugger() {
0100
0101
0102 }
0103
0104
0105
0106
0107
0108
0109 void L1ABCDebugger::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0110 using namespace edm;
0111
0112 Handle<L1AcceptBunchCrossingCollection> pIn;
0113 iEvent.getByToken(m_l1abccollectionToken, pIn);
0114
0115
0116 for (L1AcceptBunchCrossingCollection::const_iterator l1abc = pIn->begin(); l1abc != pIn->end(); ++l1abc) {
0117 if (l1abc->l1AcceptOffset() == 0) {
0118 if (m_hoffsets && *m_hoffsets)
0119 (*m_hoffsets)
0120 ->Fill((int)l1abc->bunchCrossing() - (int)iEvent.bunchCrossing(),
0121 (long long)l1abc->orbitNumber() - (long long)iEvent.orbitNumber());
0122 if (m_horboffvsorb && *m_horboffvsorb)
0123 (*m_horboffvsorb)->Fill(iEvent.orbitNumber(), (long long)l1abc->orbitNumber() - (long long)iEvent.orbitNumber());
0124 if (m_hbxoffvsorb && *m_hbxoffvsorb)
0125 (*m_hbxoffvsorb)->Fill(iEvent.orbitNumber(), (int)l1abc->bunchCrossing() - (int)iEvent.bunchCrossing());
0126 }
0127 }
0128
0129
0130
0131 edm::LogInfo("L1ABCDebug") << "Dump of L1AcceptBunchCrossing Collection for event in orbit " << iEvent.orbitNumber()
0132 << " and BX " << iEvent.bunchCrossing();
0133
0134 for (L1AcceptBunchCrossingCollection::const_iterator l1abc = pIn->begin(); l1abc != pIn->end(); ++l1abc) {
0135 edm::LogVerbatim("L1ABCDebug") << *l1abc;
0136 }
0137 }
0138
0139 void L1ABCDebugger::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
0140 m_rhm.beginRun(iRun);
0141
0142 if (m_hoffsets && *m_hoffsets) {
0143 (*m_hoffsets)->GetXaxis()->SetTitle("#Delta BX (SCAL-Event)");
0144 (*m_hoffsets)->GetYaxis()->SetTitle("#Delta orbit (SCAL-Event)");
0145 }
0146 if (m_horboffvsorb && *m_horboffvsorb) {
0147 (*m_horboffvsorb)->GetXaxis()->SetTitle("Orbit");
0148 (*m_horboffvsorb)->GetYaxis()->SetTitle("#Delta orbit (SCAL-Event)");
0149 (*m_horboffvsorb)->SetCanExtend(TH1::kXaxis);
0150 }
0151 if (m_hbxoffvsorb && *m_hbxoffvsorb) {
0152 (*m_hbxoffvsorb)->GetXaxis()->SetTitle("Orbit");
0153 (*m_hbxoffvsorb)->GetYaxis()->SetTitle("#Delta BX (SCAL-Event)");
0154 (*m_hbxoffvsorb)->SetCanExtend(TH1::kXaxis);
0155 }
0156 }
0157
0158 void L1ABCDebugger::beginJob() {}
0159
0160
0161 void L1ABCDebugger::endJob() {}
0162
0163
0164 DEFINE_FWK_MODULE(L1ABCDebugger);