CMS 3D CMS Logo

TotemDAQTriggerDQMSource.cc
Go to the documentation of this file.
1 /****************************************************************************
2 *
3 * This is a part of TOTEM offline software.
4 * Authors:
5 * Jan Kašpar (jan.kaspar@gmail.com)
6 * Rafał Leszko (rafal.leszko@gmail.com)
7 *
8 ****************************************************************************/
9 
16 
19 
21 
22 #include <string>
23 
24 //----------------------------------------------------------------------------------------------------
25 
27 public:
29  ~TotemDAQTriggerDQMSource() override;
30 
31 protected:
32  void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
33  void analyze(edm::Event const &e, edm::EventSetup const &eSetup) override;
34 
35 private:
36  unsigned int verbosity;
37 
39 
43 };
44 
45 //----------------------------------------------------------------------------------------------------
46 //----------------------------------------------------------------------------------------------------
47 
48 using namespace std;
49 using namespace edm;
50 
51 //----------------------------------------------------------------------------------------------------
52 
54  : verbosity(ps.getUntrackedParameter<unsigned int>("verbosity", 0)) {
55  tokenFEDInfo = consumes<vector<TotemFEDInfo>>(ps.getParameter<edm::InputTag>("tagFEDInfo"));
56 }
57 
58 //----------------------------------------------------------------------------------------------------
59 
61 
62 //----------------------------------------------------------------------------------------------------
63 
65  ibooker.cd();
66 
67  ibooker.setCurrentFolder("CTPPS/DAQ/");
68 
69  daq_bx_diff = ibooker.book1D("bx_diff", ";OptoRx_{i}.BX - OptoRx_{j}.BX", 21, -10.5, +10.5);
70  daq_event_bx_diff = ibooker.book1D("daq_event_bx_diff", ";OptoRx_{i}.BX - Event.BX", 21, -10.5, +10.5);
72  ibooker.book2D("daq_event_bx_diff_vs_fed", ";OptoRx.ID;OptoRx.BX - Event.BX", 10, 575.5, 585.5, 21, -10.5, +10.5);
73 }
74 
75 //----------------------------------------------------------------------------------------------------
76 
78  // get input
80  event.getByToken(tokenFEDInfo, fedInfo);
81 
82  // check validity
83  bool daqValid = fedInfo.isValid();
84 
85  if (!daqValid) {
86  if (verbosity) {
87  LogPrint("TotemDAQTriggerDQMSource")
88  << "WARNING in TotemDAQTriggerDQMSource::analyze > some of the inputs are not valid.\n"
89  << " fedInfo.isValid = " << fedInfo.isValid();
90  }
91  }
92 
93  // DAQ plots
94  if (daqValid) {
95  for (auto &it1 : *fedInfo) {
96  daq_event_bx_diff->Fill(it1.bx() - event.bunchCrossing());
97  daq_event_bx_diff_vs_fed->Fill(it1.fedId(), it1.bx() - event.bunchCrossing());
98 
99  for (auto &it2 : *fedInfo) {
100  if (it2.fedId() <= it1.fedId())
101  continue;
102 
103  daq_bx_diff->Fill(it2.bx() - it1.bx());
104  }
105  }
106  }
107 }
108 
109 //----------------------------------------------------------------------------------------------------
110 
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
TotemDAQTriggerDQMSource(const edm::ParameterSet &ps)
void analyze(edm::Event const &e, edm::EventSetup const &eSetup) override
edm::EDGetTokenT< std::vector< TotemFEDInfo > > tokenFEDInfo
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void Fill(long long x)
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Log< level::Warning, true > LogPrint
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:212
bool isValid() const
Definition: HandleBase.h:70
HLT enums.
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
Definition: event.py:1
Definition: Run.h:45