Belle II Software  release-05-02-19
TOPWaveformQualityPlotterModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2017 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Jan Strube *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 // basf2
12 #include <framework/core/HistoModule.h>
13 #include <top/modules/TOPWaveformQualityPlotter/TOPWaveformQualityPlotterModule.h>
14 
15 // stl
16 #include <utility>
17 
18 //ROOT
19 #include "TDirectory.h"
20 #include "TCanvas.h"
21 #include "TGraphErrors.h"
22 
23 using namespace std;
24 
25 namespace Belle2 {
31  REG_MODULE(TOPWaveformQualityPlotter)
32 
34  : HistoModule()
35  {
36  setDescription("TOP DQM histogram module");
37  addParam("histogramDirectoryName", m_histogramDirectoryName,
38  "histogram directory in ROOT file", string("TOP"));
39  addParam("drawWaves", m_DRAWWAVES, "option to draw waveforms", true);
40  addParam("debugHistos", m_DEBUGGING, "option to draw debug histograms", true);
41  addParam("noisemaps", m_NOISE, "option to draw noisemaps", false);
42  }
43 
44 
45  void TOPWaveformQualityPlotterModule::defineHisto()
46  {
47  TDirectory* oldDir = gDirectory;
48  m_directory = oldDir->mkdir(m_histogramDirectoryName.c_str());
49  m_directory->cd();
50  m_samples = new TH1F("ADCvalues", "ADC values ", 100, -50, 50);
51  m_samples->GetXaxis()->SetTitle("ADC Value");
52  m_samples->GetYaxis()->SetTitle("Number of Samples");
53  m_scrod_id = new TH1F("scrodID", "scrodID", 100, 0, 100);
54  m_asic = new TH1F("IRSX", "IRSX", 4, 0, 4);
55  m_carrier = new TH1F("carrier", "asic col", 4, 0, 4);
56  m_asic_ch = new TH1F("asicCh", "channel", 8, 0, 8);
57  m_errorFlag = new TH1F("errorFlag", "errorFlag", 1000, 0, 1000);
58  m_asic_win = new TH1F("window", "window", 4, 0, 4);
59  m_entries = new TH1F("entries", "entries", 100, 0, 2600);
60  m_moduleID = new TH1F("moduleID", "moduleID", 16, 1, 17);
61  m_pixelID = new TH1F("pixelID", "pixelID", 512, 1, 513);
62  oldDir->cd();
63  }
64 
65 
66  void TOPWaveformQualityPlotterModule::initialize()
67  {
68  // Register histograms (calls back defineHisto)
69  REG_HISTOGRAM;
70  //Get Waveform from datastore
71  m_waveform.isRequired();
72  }
73 
74 
75  void TOPWaveformQualityPlotterModule::basicDebuggingPlots(const TOPRawWaveform& v)
76  {
77  int scrodid = v.getScrodID();
78  int asicid = v.getASICNumber();
79  int channelid = v.getASICChannel();
80  int carrierid = v.getCarrierNumber();
81  m_scrod_id->Fill(scrodid);
82  m_asic_ch->Fill(channelid);
83  m_asic->Fill(asicid);
84  m_carrier->Fill(carrierid);
85  m_asic_win->Fill(v.getStorageWindow());
86  m_entries->Fill(v.getWaveform().size());
87  m_moduleID->Fill(v.getModuleID());
88  m_pixelID->Fill(v.getPixelID());
89 
90  if (m_hitmap.find(scrodid) == m_hitmap.end()) { // cppcheck-suppress stlFindInsert
91  m_hitmap[scrodid] = new TH2F((string("scrod ") + to_string(scrodid) + string("Hitmap")).c_str(),
92  (string("scrod ") + to_string(scrodid) + string("carrier vs. asic;asic;carrier")).c_str(), 4, 0, 4, 4, 0, 4);
93  }
94  m_hitmap[scrodid]->Fill(asicid, carrierid);
95  const vector<short>& waveform = v.getWaveform();
96  if (waveform.empty()) {
97  return;
98  }
99  for (short adc : waveform) {
100  m_samples->Fill(adc);
101  }
102  }
103 
104  void
105  TOPWaveformQualityPlotterModule::drawWaveforms(const TOPRawWaveform& v)
106  {
107  vector<short> waveform = v.getWaveform();
108  if (waveform.empty()) {
109  return;
110  }
111  int scrodid = v.getScrodID();
112  // skip broken events
113  if (scrodid == 0) {
114  return;
115  }
116  int asicNumber = v.getASICNumber();
117  int carrierNumber = v.getCarrierNumber();
118  int iChannel = v.getASICChannel();
119  string gname = string("scrod_") + to_string(scrodid) + string("_carrier") + to_string(carrierNumber) + string("_asic") + to_string(
120  asicNumber) + to_string(iChannel);
121  if (m_waveformHists[scrodid][carrierNumber][asicNumber].find(iChannel) ==
122  m_waveformHists[scrodid][carrierNumber][asicNumber].end()) {
123  // FIXME ? This assumes exactly 256 samples in a waveform (FE data)
124  auto h = new TProfile(gname.c_str(), gname.c_str(), 256, 0, 256);
125  h->Sumw2(false); // unweighted filling.
126  m_waveformHists[scrodid][carrierNumber][asicNumber][iChannel] = h;
127  }
128  // FIXME assumes 256 samples in a waveform
129  for (size_t i = 0; i < 256; ++i) {
130  // exit broken waveforms early
131  if (i >= waveform.size()) {
132  break;
133  }
134  m_waveformHists[scrodid][carrierNumber][asicNumber][iChannel]->Fill(i + 0.5, iChannel * 1500 + waveform[i]);
135  }
136  }
137 
138 
139  void TOPWaveformQualityPlotterModule::event()
140  {
141  if (not m_waveform) {
142  return;
143  }
144  for (auto evtwave : m_waveform) {
145  if (m_DRAWWAVES) {
146  drawWaveforms(evtwave);
147  }
148  if (m_DEBUGGING) {
149  basicDebuggingPlots(evtwave);
150  }
151  if (m_NOISE) {
152  auto channelID = evtwave.getChannel();
153  const vector<short> v_samples = evtwave.getWaveform();
154  size_t nsamples = v_samples.size();
155  if (m_channelNoiseMap.find(channelID) == m_channelNoiseMap.end()) {
156  string idName = string("noise_") + to_string(channelID);
157  m_channelNoiseMap.insert(make_pair(channelID, new TH1F(idName.c_str(), idName.c_str(), 200, -100, 100)));
158  m_channelEventMap.insert(make_pair(channelID, m_iEvent));
159  }
160  TH1F* noise = m_channelNoiseMap[channelID];
161  // Plot all samples in common histogram for quick sanity check
162  for (size_t s = 0; s < nsamples; s++) {
163  double adc = v_samples.at(s);
164  m_samples->Fill(adc);
165  if (s < nsamples - 1) {
166  noise->Fill(v_samples[s + 1] - adc);
167  }
168  }
169  }
170  }
171  m_iEvent += 1;
172  return;
173  }
174 
175  void TOPWaveformQualityPlotterModule::endRun()
176  {
177  if (m_DRAWWAVES) {
178  // Each waveform was stored in a TH1F
179  // This now gets transformed to a TGraph
180  // All 8 channels for a given ASIC are put in the same TMultiGraph
181  // Then we make one canvas of several TMultigraphs for each scrod
182  // FIXME: This is going to get called at the end of the run; the mem leak I am introducing here should be fixed by somebody more patient with ROOT memory management than me.
183  for (auto scrod_it : m_waveformHists) {
184  int scrodid = scrod_it.first;
185  string name = string("scrod_") + to_string(scrodid);
186  TCanvas* c = new TCanvas(name.c_str(), name.c_str());
187  c->Divide(4, 4);
188  int canvasPad = 0;
189  for (auto carrier_it : scrod_it.second) {
190  int carrierNumber = carrier_it.first;
191  for (auto asic_it : carrier_it.second) {
192  int asicNumber = asic_it.first;
193  canvasPad += 1;
194  string gname = string("scrod_") + to_string(scrodid) + string("_carrier") + to_string(carrierNumber) + string("_asic") + to_string(
195  asicNumber);
196  TMultiGraph* mg = new TMultiGraph(gname.c_str(), gname.c_str());
197  for (auto channel_it : asic_it.second) {
198  // read the data from the TProfile
199  TGraphErrors* g = new TGraphErrors(channel_it.second);
200  g->SetMarkerStyle(7);
201  mg->Add(g);
202  }
203  // FIXME nSamples is hard-coded to 256
204  TH2F* h = new TH2F(gname.c_str(), gname.c_str(), 256, 0, 256, 8, -500, -500 + 8 * 1500);
205  for (int ibin = 0; ibin < 8; ibin++) {
206  h->GetYaxis()->SetBinLabel(ibin + 1, to_string(ibin).c_str());
207  }
208  h->GetYaxis()->SetTickSize(0);
209  h->GetYaxis()->SetTickLength(0);
210  h->SetStats(0);
211  c->cd(canvasPad);
212  h->Draw();
213  mg->Draw("P");
214  }
215  }
216  m_directory->WriteTObject(c);
217  }
218  }
219  }
220 
222 } // end Belle2 namespace
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::TOPWaveformQualityPlotterModule
Plots and histograms of waveforms and feature extracted parameters.
Definition: TOPWaveformQualityPlotterModule.h:48
Belle2::TOPRawWaveform
Class to store raw data waveforms.
Definition: TOPRawWaveform.h:34
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::HistoModule
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29