Belle II Software  release-05-02-19
DQMHistAnalysisSVDDose.h
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2021 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Ludovico Massaccesi *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #pragma once
12 
13 #ifdef _BELLE2_EPICS
14 #include "cadef.h"
15 #endif
16 
17 #include <dqm/analysis/modules/DQMHistAnalysis.h>
18 #include <TCanvas.h>
19 #include <TString.h>
20 #include <TPaveText.h>
21 #include <TH2.h>
22 #include <TH1.h>
23 #include <string>
24 #include <vector>
25 
26 namespace Belle2 {
48  class DQMHistAnalysisSVDDoseModule : public DQMHistAnalysisModule {
49  public:
50  DQMHistAnalysisSVDDoseModule();
51  virtual ~DQMHistAnalysisSVDDoseModule();
52 
53  private:
57  typedef struct SensorGroup {
58  TString nameSuffix;
59  TString titleSuffix;
60  const char* pvMiddle;
61  int nStrips;
62  } SensorGroup;
63 
64 #ifdef _BELLE2_EPICS
65 
66  typedef struct MyPV {
67  chid mychid;
68  double lastNHits = 0.0;
69  double lastNEvts = 0.0;
70  } MyPV;
71 #endif
72 
73  void initialize() override final;
74  void beginRun() override final;
75  void event() override final;
76  void endRun() override final;
77 
78  void updateCanvases();
79 
81  template<typename T>
82  inline T* findHistT(TString name) { return dynamic_cast<T*>(findHist(name.Data())); }
83 
94  template<typename T>
95  static T* divide(T* num, T* den, float scale = 1.0f)
96  {
97  TString name = TString("occu_from_") + num->GetName();
98  T* res = (T*)num->Clone(name);
99  for (int i = 0; i < num->GetNcells(); i++) {
100  float n = num->GetBinContent(i), d = den->GetBinContent(i), e = num->GetBinError(i);
101  res->SetBinContent(i, d ? scale * n / d : 0.0f);
102  res->SetBinError(i, d ? scale * e / d : 0.0f);
103  }
104  return res;
105  }
106 
108  static void carryOverflowOver(TH1F* h);
109 
110  // Steerable data members (parameters)
111  std::string m_pvPrefix;
112  bool m_useEpics;
113  double m_epicsUpdateSeconds;
114  std::string m_pvSuffix;
115  std::string m_deltaTPVSuffix;
116  std::string m_statePVSuffix;
118  // Data members for outputs
120  TPaveText* m_legend = nullptr;
121  // Canvases & output histos for Poisson trigger (TTYP_POIS) events
122  std::vector<TCanvas*> m_c_instOccu;
123  std::vector<TCanvas*> m_c_occuLER;
124  std::vector<TCanvas*> m_c_occuHER;
125  std::vector<TCanvas*> m_c_occuLER1;
126  std::vector<TCanvas*> m_c_occuHER1;
127  std::vector<TH2F*> m_h_occuLER;
128  std::vector<TH2F*> m_h_occuHER;
129  std::vector<TH1F*> m_h_occuLER1;
130  std::vector<TH1F*> m_h_occuHER1;
131  // Canvases & output histos for all events
132  std::vector<TCanvas*> m_c_instOccuAll;
133  std::vector<TCanvas*> m_c_occuLERAll;
134  std::vector<TCanvas*> m_c_occuHERAll;
135  std::vector<TCanvas*> m_c_occuLER1All;
136  std::vector<TCanvas*> m_c_occuHER1All;
137  std::vector<TH2F*> m_h_occuLERAll;
138  std::vector<TH2F*> m_h_occuHERAll;
139  std::vector<TH1F*> m_h_occuLER1All;
140  std::vector<TH1F*> m_h_occuHER1All;
142 #ifdef _BELLE2_EPICS
143  std::vector<MyPV> m_myPVs;
144  double m_lastPVUpdate = -1.0;
145  chid m_timeSinceLastPVUpdateChan;
146  struct dbr_ctrl_enum m_stateCtrl;
147  chid m_stateChan;
148 #endif
149 
153  static const std::vector<SensorGroup> c_sensorGroups;
154  };
156 }
Belle2::DQMHistAnalysisSVDDoseModule::m_c_instOccuAll
std::vector< TCanvas * > m_c_instOccuAll
Canvases for the instantaneous occupancy.
Definition: DQMHistAnalysisSVDDose.h:140
Belle2::DQMHistAnalysisSVDDoseModule::m_monObj
MonitoringObject * m_monObj
Monitoring object for MiraBelle.
Definition: DQMHistAnalysisSVDDose.h:127
Belle2::DQMHistAnalysisSVDDoseModule::m_h_occuHERAll
std::vector< TH2F * > m_h_occuHERAll
Histograms for the occu.
Definition: DQMHistAnalysisSVDDose.h:146
Belle2::DQMHistAnalysisSVDDoseModule::m_c_occuLER1All
std::vector< TCanvas * > m_c_occuLER1All
Canvases for the 1D occu.
Definition: DQMHistAnalysisSVDDose.h:143
Belle2::DQMHistAnalysisSVDDoseModule::m_pvPrefix
std::string m_pvPrefix
Prefix for EPICS PVs.
Definition: DQMHistAnalysisSVDDose.h:119
Belle2::DQMHistAnalysisSVDDoseModule::SensorGroup
struct Belle2::DQMHistAnalysisSVDDoseModule::SensorGroup SensorGroup
A struct to define the sensors group we average over.
Belle2::DQMHistAnalysisSVDDoseModule::endRun
void endRun() override final
This method is called if the current run ends.
Definition: DQMHistAnalysisSVDDose.cc:213
Belle2::DQMHistAnalysisSVDDoseModule::m_statePVSuffix
std::string m_statePVSuffix
Suffix of the state PV.
Definition: DQMHistAnalysisSVDDose.h:124
Belle2::DQMHistAnalysisSVDDoseModule::m_c_occuLERAll
std::vector< TCanvas * > m_c_occuLERAll
Canvases for the occu.
Definition: DQMHistAnalysisSVDDose.h:141
Belle2::DQMHistAnalysisModule::findHist
static TH1 * findHist(const std::string &histname)
Find histogram.
Definition: DQMHistAnalysis.cc:83
Belle2::DQMHistAnalysisSVDDoseModule::m_useEpics
bool m_useEpics
Whether to update EPICS PVs.
Definition: DQMHistAnalysisSVDDose.h:120
Belle2::DQMHistAnalysisSVDDoseModule::SensorGroup::pvMiddle
const char * pvMiddle
Middle part of the PV name.
Definition: DQMHistAnalysisSVDDose.h:68
Belle2::DQMHistAnalysisSVDDoseModule::m_h_occuLER1All
std::vector< TH1F * > m_h_occuLER1All
Histograms for the 1D occu.
Definition: DQMHistAnalysisSVDDose.h:147
Belle2::DQMHistAnalysisSVDDoseModule::m_c_occuLER1
std::vector< TCanvas * > m_c_occuLER1
Canvases for the 1D occu.
Definition: DQMHistAnalysisSVDDose.h:133
Belle2::DQMHistAnalysisSVDDoseModule::m_h_occuHER1
std::vector< TH1F * > m_h_occuHER1
Histograms for the 1D occu.
Definition: DQMHistAnalysisSVDDose.h:138
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::DQMHistAnalysisSVDDoseModule::m_h_occuHER1All
std::vector< TH1F * > m_h_occuHER1All
Histograms for the 1D occu.
Definition: DQMHistAnalysisSVDDose.h:148
Belle2::DQMHistAnalysisSVDDoseModule::m_h_occuLER
std::vector< TH2F * > m_h_occuLER
Histograms for the occu.
Definition: DQMHistAnalysisSVDDose.h:135
Belle2::DQMHistAnalysisSVDDoseModule::m_epicsUpdateSeconds
double m_epicsUpdateSeconds
Minimum interval between successive PV updates.
Definition: DQMHistAnalysisSVDDose.h:121
Belle2::DQMHistAnalysisSVDDoseModule::m_h_occuLERAll
std::vector< TH2F * > m_h_occuLERAll
Histograms for the occu.
Definition: DQMHistAnalysisSVDDose.h:145
Belle2::DQMHistAnalysisSVDDoseModule::initialize
void initialize() override final
Initialize the Module.
Definition: DQMHistAnalysisSVDDose.cc:46
Belle2::DQMHistAnalysisSVDDoseModule::m_legend
TPaveText * m_legend
Legend of the inst.
Definition: DQMHistAnalysisSVDDose.h:128
Belle2::DQMHistAnalysisSVDDoseModule::m_h_occuHER
std::vector< TH2F * > m_h_occuHER
Histograms for the occu.
Definition: DQMHistAnalysisSVDDose.h:136
Belle2::DQMHistAnalysisSVDDoseModule::SensorGroup::titleSuffix
TString titleSuffix
Suffix for the title of the canvases.
Definition: DQMHistAnalysisSVDDose.h:67
Belle2::DQMHistAnalysisSVDDoseModule::m_pvSuffix
std::string m_pvSuffix
Suffix for EPICS PVs.
Definition: DQMHistAnalysisSVDDose.h:122
Belle2::DQMHistAnalysisSVDDoseModule::carryOverflowOver
static void carryOverflowOver(TH1F *h)
Carries the content of the overflow bin into the last bin.
Definition: DQMHistAnalysisSVDDose.cc:450
Belle2::DQMHistAnalysisSVDDoseModule::SensorGroup::nameSuffix
TString nameSuffix
Suffix of the name of the histograms.
Definition: DQMHistAnalysisSVDDose.h:66
Belle2::DQMHistAnalysisSVDDoseModule::m_c_occuHER1
std::vector< TCanvas * > m_c_occuHER1
Canvases for the 1D occu.
Definition: DQMHistAnalysisSVDDose.h:134
Belle2::DQMHistAnalysisSVDDoseModule::m_c_occuLER
std::vector< TCanvas * > m_c_occuLER
Canvases for the occu.
Definition: DQMHistAnalysisSVDDose.h:131
Belle2::DQMHistAnalysisSVDDoseModule::findHistT
T * findHistT(TString name)
Utility method.
Definition: DQMHistAnalysisSVDDose.h:90
Belle2::DQMHistAnalysisSVDDoseModule::beginRun
void beginRun() override final
Called when entering a new run.
Definition: DQMHistAnalysisSVDDose.cc:153
Belle2::DQMHistAnalysisSVDDoseModule::m_deltaTPVSuffix
std::string m_deltaTPVSuffix
Suffix of the update-time monitoring PV.
Definition: DQMHistAnalysisSVDDose.h:123
Belle2::DQMHistAnalysisSVDDoseModule::SensorGroup::nStrips
int nStrips
Total number of strips in the sensor group.
Definition: DQMHistAnalysisSVDDose.h:69
Belle2::DQMHistAnalysisSVDDoseModule::m_c_instOccu
std::vector< TCanvas * > m_c_instOccu
Canvases for the instantaneous occupancy.
Definition: DQMHistAnalysisSVDDose.h:130
Belle2::DQMHistAnalysisSVDDoseModule::m_c_occuHER
std::vector< TCanvas * > m_c_occuHER
Canvases for the occu.
Definition: DQMHistAnalysisSVDDose.h:132
Belle2::DQMHistAnalysisSVDDoseModule::m_h_occuLER1
std::vector< TH1F * > m_h_occuLER1
Histograms for the 1D occu.
Definition: DQMHistAnalysisSVDDose.h:137
Belle2::DQMHistAnalysisSVDDoseModule::m_c_occuHERAll
std::vector< TCanvas * > m_c_occuHERAll
Canvases for the occu.
Definition: DQMHistAnalysisSVDDose.h:142
Belle2::DQMHistAnalysisSVDDoseModule::event
void event() override final
This method is the core of the module.
Definition: DQMHistAnalysisSVDDose.cc:168
Belle2::DQMHistAnalysisSVDDoseModule::c_sensorGroups
static const std::vector< SensorGroup > c_sensorGroups
List of sensors groups.
Definition: DQMHistAnalysisSVDDose.h:161
Belle2::MonitoringObject
MonitoringObject is a basic object to hold data for the run-dependency monitoring Run summary TCanvas...
Definition: MonitoringObject.h:41
Belle2::DQMHistAnalysisSVDDoseModule::divide
static T * divide(T *num, T *den, float scale=1.0f)
Divide two histograms, ignoring errors on the second histogram.
Definition: DQMHistAnalysisSVDDose.h:103
Belle2::DQMHistAnalysisSVDDoseModule::m_c_occuHER1All
std::vector< TCanvas * > m_c_occuHER1All
Canvases for the 1D occu.
Definition: DQMHistAnalysisSVDDose.h:144