File indexing completed on 2024-04-06 12:07:00
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010 #include "FWCore/Framework/interface/ESHandle.h"
0011 #include "FWCore/Framework/interface/MakerMacros.h"
0012 #include "FWCore/Framework/interface/Event.h"
0013 #include "FWCore/Framework/interface/EventSetup.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 #include "FWCore/Utilities/interface/InputTag.h"
0016
0017 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0018 #include "DQMServices/Core/interface/DQMStore.h"
0019
0020 #include "DataFormats/CTPPSDigi/interface/TotemFEDInfo.h"
0021
0022 #include <string>
0023
0024
0025
0026 class TotemDAQTriggerDQMSource : public DQMEDAnalyzer {
0027 public:
0028 TotemDAQTriggerDQMSource(const edm::ParameterSet &ps);
0029 ~TotemDAQTriggerDQMSource() override;
0030
0031 protected:
0032 void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
0033 void analyze(edm::Event const &e, edm::EventSetup const &eSetup) override;
0034
0035 private:
0036 unsigned int verbosity;
0037
0038 edm::EDGetTokenT<std::vector<TotemFEDInfo>> tokenFEDInfo;
0039
0040 MonitorElement *daq_bx_diff;
0041 MonitorElement *daq_event_bx_diff;
0042 MonitorElement *daq_event_bx_diff_vs_fed;
0043 };
0044
0045
0046
0047
0048 using namespace std;
0049 using namespace edm;
0050
0051
0052
0053 TotemDAQTriggerDQMSource::TotemDAQTriggerDQMSource(const edm::ParameterSet &ps)
0054 : verbosity(ps.getUntrackedParameter<unsigned int>("verbosity", 0)) {
0055 tokenFEDInfo = consumes<vector<TotemFEDInfo>>(ps.getUntrackedParameter<edm::InputTag>("tagFEDInfo"));
0056 }
0057
0058
0059
0060 TotemDAQTriggerDQMSource::~TotemDAQTriggerDQMSource() {}
0061
0062
0063
0064 void TotemDAQTriggerDQMSource::bookHistograms(DQMStore::IBooker &ibooker, edm::Run const &, edm::EventSetup const &) {
0065 ibooker.cd();
0066
0067 ibooker.setCurrentFolder("CTPPS/DAQ/");
0068
0069 daq_bx_diff = ibooker.book1D("bx_diff", ";OptoRx_{i}.BX - OptoRx_{j}.BX", 21, -10.5, +10.5);
0070 daq_event_bx_diff = ibooker.book1D("daq_event_bx_diff", ";OptoRx_{i}.BX - Event.BX", 21, -10.5, +10.5);
0071 daq_event_bx_diff_vs_fed =
0072 ibooker.book2D("daq_event_bx_diff_vs_fed", ";OptoRx.ID;OptoRx.BX - Event.BX", 10, 575.5, 585.5, 21, -10.5, +10.5);
0073 }
0074
0075
0076
0077 void TotemDAQTriggerDQMSource::analyze(edm::Event const &event, edm::EventSetup const &eventSetup) {
0078
0079 Handle<vector<TotemFEDInfo>> fedInfo;
0080 event.getByToken(tokenFEDInfo, fedInfo);
0081
0082
0083 bool daqValid = fedInfo.isValid();
0084
0085 if (!daqValid) {
0086 if (verbosity) {
0087 LogPrint("TotemDAQTriggerDQMSource")
0088 << "WARNING in TotemDAQTriggerDQMSource::analyze > some of the inputs are not valid.\n"
0089 << " fedInfo.isValid = " << fedInfo.isValid();
0090 }
0091 }
0092
0093
0094 if (daqValid) {
0095 for (auto &it1 : *fedInfo) {
0096 daq_event_bx_diff->Fill(it1.bx() - event.bunchCrossing());
0097 daq_event_bx_diff_vs_fed->Fill(it1.fedId(), it1.bx() - event.bunchCrossing());
0098
0099 for (auto &it2 : *fedInfo) {
0100 if (it2.fedId() <= it1.fedId())
0101 continue;
0102
0103 daq_bx_diff->Fill(it2.bx() - it1.bx());
0104 }
0105 }
0106 }
0107 }
0108
0109
0110
0111 DEFINE_FWK_MODULE(TotemDAQTriggerDQMSource);