Belle II Software  release-05-01-25
DQMHistAnalysisPXDCM.h
1 //+
2 // File : DQMHistAnalysisPXDCM.h
3 // Description : DAQM Analysis for PXD Common Modes
4 //
5 // Author : Bjoern Spruck, University Mainz
6 // Date : 2018
7 //-
8 
9 #pragma once
10 
11 #ifdef _BELLE2_EPICS
12 // EPICS
13 #include "cadef.h"
14 #endif
15 
16 #include <dqm/analysis/modules/DQMHistAnalysis.h>
17 #include <vxd/dataobjects/VxdID.h>
18 
19 #include <vector>
20 #include <map>
21 #include <TH2.h>
22 #include <TCanvas.h>
23 #include <TLine.h>
24 
25 namespace Belle2 {
33 
34  // Public functions
35  public:
36 
41  private:
42 
44  void initialize(void) override final;
45 
47  void beginRun(void) override final;
48  void event(void) override final;
49  void terminate(void) override final;
50 
51  // Data members
55  std::string m_pvPrefix;
57  int m_minEntries = 10000;
59  double m_warnMeanAdhoc{};
61  double m_errorMeanAdhoc{};
62 
64  std::vector<VxdID> m_PXDModules;
65 
67  TH2D* m_hCommonMode = nullptr;
69  TH2D* m_hCommonModeDelta = nullptr;
71  TH2D* m_hCommonModeOld = nullptr;
73  TCanvas* m_cCommonMode = nullptr;
75  TCanvas* m_cCommonModeDelta = nullptr;
77  TLine* m_line1 = nullptr;
79  TLine* m_line2 = nullptr;
81 // TLine* m_line3 = nullptr;
82 
85 
87  bool m_useEpics;
88 
89 #ifdef _BELLE2_EPICS
90  std::vector <chid> mychid;
93  std::map <VxdID, chid> mychid_mean;
94 #endif
95  };
97 } // end namespace Belle2
98 
Belle2::DQMHistAnalysisPXDCMModule::DQMHistAnalysisPXDCMModule
DQMHistAnalysisPXDCMModule()
Constructor.
Definition: DQMHistAnalysisPXDCM.cc:27
Belle2::DQMHistAnalysisPXDCMModule
DQM Histogram Analysis for PXD Common Modes.
Definition: DQMHistAnalysisPXDCM.h:32
Belle2::DQMHistAnalysisPXDCMModule::m_warnMeanAdhoc
double m_warnMeanAdhoc
warn level for mean adhoc plot
Definition: DQMHistAnalysisPXDCM.h:59
Belle2::DQMHistAnalysisPXDCMModule::m_minEntries
int m_minEntries
Update entry intervall.
Definition: DQMHistAnalysisPXDCM.h:57
Belle2::DQMHistAnalysisPXDCMModule::m_hCommonModeOld
TH2D * m_hCommonModeOld
histogram covering all modules
Definition: DQMHistAnalysisPXDCM.h:71
Belle2::DQMHistAnalysisPXDCMModule::m_hCommonMode
TH2D * m_hCommonMode
histogram covering all modules
Definition: DQMHistAnalysisPXDCM.h:67
Belle2::DQMHistAnalysisPXDCMModule::beginRun
void beginRun(void) override final
Module functions to be called from event process.
Definition: DQMHistAnalysisPXDCM.cc:131
Belle2::DQMHistAnalysisPXDCMModule::m_hCommonModeDelta
TH2D * m_hCommonModeDelta
histogram covering all modules
Definition: DQMHistAnalysisPXDCM.h:69
Belle2::DQMHistAnalysisPXDCMModule::m_cCommonMode
TCanvas * m_cCommonMode
Final Canvas.
Definition: DQMHistAnalysisPXDCM.h:73
Belle2::DQMHistAnalysisPXDCMModule::m_pvPrefix
std::string m_pvPrefix
prefix for EPICS PVs
Definition: DQMHistAnalysisPXDCM.h:55
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::DQMHistAnalysisPXDCMModule::m_errorMeanAdhoc
double m_errorMeanAdhoc
error level for mean adhoc plot
Definition: DQMHistAnalysisPXDCM.h:61
Belle2::DQMHistAnalysisPXDCMModule::m_monObj
MonitoringObject * m_monObj
Line in the Canvas to guide the eye.
Definition: DQMHistAnalysisPXDCM.h:84
Belle2::DQMHistAnalysisPXDCMModule::m_histogramDirectoryName
std::string m_histogramDirectoryName
name of histogram directory
Definition: DQMHistAnalysisPXDCM.h:53
Belle2::DQMHistAnalysisPXDCMModule::m_useEpics
bool m_useEpics
flag if to export to EPICS
Definition: DQMHistAnalysisPXDCM.h:87
Belle2::DQMHistAnalysisPXDCMModule::event
void event(void) override final
This method is the core of the module.
Definition: DQMHistAnalysisPXDCM.cc:146
Belle2::DQMHistAnalysisPXDCMModule::initialize
void initialize(void) override final
Module functions to be called from main process.
Definition: DQMHistAnalysisPXDCM.cc:51
Belle2::DQMHistAnalysisPXDCMModule::terminate
void terminate(void) override final
This method is called at the end of the event processing.
Definition: DQMHistAnalysisPXDCM.cc:332
Belle2::DQMHistAnalysisPXDCMModule::~DQMHistAnalysisPXDCMModule
~DQMHistAnalysisPXDCMModule()
Destructor.
Definition: DQMHistAnalysisPXDCM.cc:42
Belle2::DQMHistAnalysisPXDCMModule::m_PXDModules
std::vector< VxdID > m_PXDModules
IDs of all PXD Modules to iterate over.
Definition: DQMHistAnalysisPXDCM.h:64
Belle2::DQMHistAnalysisPXDCMModule::m_line1
TLine * m_line1
Line in the Canvas to guide the eye.
Definition: DQMHistAnalysisPXDCM.h:77
Belle2::DQMHistAnalysisPXDCMModule::m_line2
TLine * m_line2
Line in the Canvas to guide the eye.
Definition: DQMHistAnalysisPXDCM.h:79
Belle2::MonitoringObject
MonitoringObject is a basic object to hold data for the run-dependency monitoring Run summary TCanvas...
Definition: MonitoringObject.h:41
Belle2::DQMHistAnalysisPXDCMModule::m_cCommonModeDelta
TCanvas * m_cCommonModeDelta
Final Canvas.
Definition: DQMHistAnalysisPXDCM.h:75
Belle2::DQMHistAnalysisModule
The base class for the histogram analysis module.
Definition: DQMHistAnalysis.h:27