Belle II Software development
SVDDQMDoseModule.h
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#pragma once
10
11#include <framework/core/HistoModule.h>
12#include <framework/datastore/StoreArray.h>
13#include <framework/datastore/StoreObjPtr.h>
14#include <TH1.h>
15#include <TH2.h>
16#include <string>
17#include <vector>
18#include <functional>
19
20namespace Belle2 {
25 // Forward declarations to avoid the need of #includes here in the header
26 class RawFTSW;
27 class SVDShaperDigit;
28 class TRGSummary;
29 class VxdID;
30 // class SVDFADCMaskedStrips;
31
32 namespace SVD {
42 public:
45 c_HERInjection = 1,
46 c_LERInjection = 2,
47 c_NoInjection = 4
48 };
49
52 void initialize() override final;
53 void beginRun() override final;
54 void event() override final;
55 void defineHisto() override final;
57 private:
59 typedef struct SensorGroup {
61 const char* nameSuffix;
63 const char* titleSuffix;
65 int nBins;
67 double xMin;
69 double xMax;
71 std::function<bool(const VxdID&)> contains;
73 mutable int nStripsU = 0;
74 } SensorGroup;
75
76 // Steerable data members (parameters)
77 unsigned int m_eventFilter;
81 std::vector<int> m_trgTypes;
83 // Inputs
88 // Outputs (histograms)
90 TH2F* h_nEvtsVsTime = nullptr;
92 std::vector<TH2F*> m_groupNHitsU;
94 std::vector<TH1F*> m_groupOccupanciesU;
96 TH1F* h_nEvtsVsTime1 = nullptr;
98 std::vector<TH1F*> m_groupNHits1U;
100 TH2F* h_bunchNumVsNHits = nullptr;
101
102 // Other stuff
111 static constexpr double c_revolutionTime = 5120.0 / 508.0;
124 static constexpr double c_globalClock = 127.0;
125
127 static constexpr double c_defaultOccuMin = 0.0;
132 static const int c_defaultNBins = 90;
144 static constexpr double c_defaultOccuMax = 100.0 / 1536.0 * c_defaultNBins;
145
147 static const std::vector<SensorGroup> c_sensorGroups;
148 };
149 }
151}
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
The SVD dose-monitoring DQM module.
TH1F * h_nEvtsVsTime1
Hist of the total evts in each time bin (1D, time since inj.
unsigned int m_eventFilter
Bitmask for event type selection, see EEventType .
void initialize() override final
Overrides HistoModule::initialize.
TH2F * h_bunchNumVsNHits
Hist of bunch number vs number of fired strips (copied from SVDDQMInjection).
static const std::vector< SensorGroup > c_sensorGroups
List of interesting groups of sensors to average over.
StoreArray< RawFTSW > m_rawTTD
Input: DAQ status.
std::vector< TH1F * > m_groupNHits1U
Hists of the number of hits in each time bin (1D) per sensor group, U-side.
void defineHisto() override final
Overrides HistoModule::defineHisto.
StoreObjPtr< TRGSummary > m_trgSummary
Input: trigger type.
SVDDQMDoseModule()
Default constructor, defines parameters.
static constexpr double c_defaultOccuMax
Default maximum of the instantaneous occupancy histograms.
static constexpr double c_revolutionTime
Beam revolution time in microseconds (approximated).
void event() override final
Overrides HistoModule::event.
static constexpr double c_globalClock
Approximated global clock frequency in MHz.
double m_noInjectionTime
After this time (in microseconds) from last injection the event falls in the "No Injection" category.
std::vector< int > m_trgTypes
Trigger types to accept (all if the vector is empty).
std::string m_histogramDirectoryName
Name of the histograms' directory in the ROOT file.
std::vector< TH2F * > m_groupNHitsU
Hists of the number of hits in each time bin per sensor group, U-side.
StoreArray< SVDShaperDigit > m_digits
Input: raw hits.
std::vector< TH1F * > m_groupOccupanciesU
Hists of the instantaneous occupancy per sensor group (see c_sensorGroups), U-side.
TH2F * h_nEvtsVsTime
Hist of the total evts in each time bin (time since inj.
EEventType
Bits definition for the bitmask that selects the events to put in the histograms.
void beginRun() override final
Overrides HistoModule::beginRun.
static constexpr double c_defaultOccuMin
Default minimum of the instantaneous occupancy histograms.
std::string m_SVDShaperDigitsName
Name of the StoreArray of SVDShaperDigit to use.
static const int c_defaultNBins
Default number of bins of the instantaneous occupancy histograms.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
Abstract base class for different kinds of events.
A struct to define non-trivial histograms in a human-readable way.
double xMin
The lower limit for the instantaneous occupancy histo.
double xMax
The upper limit for the instantaneous occupancy histo.
std::function< bool(const VxdID &)> contains
Function that says if a sensor is in this group.
const char * titleSuffix
The title will be "SVD Instantaneous Occupancy @titleSuffix @side;Occupancy [%];Count / bin".
const char * nameSuffix
The name will be "SVDInstOccupancy_@nameSuffix@side".
int nBins
The number of bins for the instantaneous occupancy histo.