Belle II Software development
TrackingHLTDQMModule.cc
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#include <tracking/modules/trackingDQM/TrackingHLTDQMModule.h>
10#include <tracking/modules/trackingDQM/TrackDQMEventProcessor.h>
11#include <tracking/dqmUtils/HistogramFactory.h>
12
13#include <TDirectory.h>
14
15using namespace Belle2;
16using namespace Belle2::HistogramFactory;
17using boost::format;
18
19//-----------------------------------------------------------------
20// Register the Module
21//-----------------------------------------------------------------
22
23REG_MODULE(TrackingHLTDQM);
24
25//-----------------------------------------------------------------
26// Implementation
27//-----------------------------------------------------------------
28
30{
32 setDescription("Data Quality Monitoring of the tracking run on HLT. ");
33}
34
35//------------------------------------------------------------------
36// Function to define histograms
37//-----------------------------------------------------------------
38
40{
42
43 // eventLevelTrackingInfo is currently only set by VXDTF2, if VXDTF2 is not in path the StoreObject is not there
44 m_eventLevelTrackingInfo.isOptional();
45
46 m_rawTTD.isOptional();
47}
48
50{
53
54 if (VXD::GeoCache::getInstance().getGeoTools()->getNumberOfLayers() == 0)
55 B2FATAL("Missing geometry for VXD.");
56
57 // Create a separate histogram directories and cd into it.
58 TDirectory* originalDirectory = gDirectory;
59
60 // There might be problems with nullptr if the directory with the same name already exists (but I am not sure because there isn't anything like that in AlignmentDQM)
61 TDirectory* TracksDQM = originalDirectory->GetDirectory("TrackingHLTDQM");
62 if (!TracksDQM)
63 TracksDQM = originalDirectory->mkdir("TrackingHLTDQM");
64
65 TracksDQM->cd();
67 DefineHits();
72
74
75 originalDirectory->cd();
76
77 for (auto change : m_histogramParameterChanges)
78 ProcessHistogramParameterChange(std::get<0>(change), std::get<1>(change), std::get<2>(change));
79}
80
82{
85 return;
86
87 bool runningOnHLT = true;
88
91
92 eventProcessor.Run();
93
94 if (m_eventLevelTrackingInfo.isValid()) {
95 m_trackingErrorFlags->Fill((double)m_eventLevelTrackingInfo->hasAnErrorFlag());
96
97 if (m_rawTTD.isValid()) {
98
99 for (auto& it : m_rawTTD) {
100
101 // get last injection time
102 auto difference = it.GetTimeSinceLastInjection(0);
103 // check time overflow, too long ago
104 if (difference != 0x7FFFFFFF) {
105
106 double timeSinceInj = it.GetTimeSinceLastInjection(0) / c_globalClock;
107 double timeInCycle = timeSinceInj - (int)(timeSinceInj / c_revolutionTime) * c_revolutionTime;
108
109 if (it.GetIsHER(0))
110 m_allVStimeHER->Fill(timeSinceInj, timeInCycle);
111 else
112 m_allVStimeLER->Fill(timeSinceInj, timeInCycle);
113
114 if (m_eventLevelTrackingInfo->hasAnErrorFlag()) {
115 if (it.GetIsHER(0))
116 m_abortVStimeHER->Fill(timeSinceInj, timeInCycle);
117 else
118 m_abortVStimeLER->Fill(timeSinceInj, timeInCycle);
119 } //has error flag
120 } //time overflow
121 }// loop on RawFTSW
122 } //RawFTSW is valid
123 } else
124 m_trackingErrorFlags->Fill(0.0);
125
126}
127
129{
130 // only monitor if any flag was set so only 2 bins needed
131 const char* flagTitle =
132 "Tracking error summary. Mean = errors/event (should be 0 or very close to 0); Error occured yes or no; Number of events";
133
134 m_trackingErrorFlags = Create("NumberTrackingErrorFlags", flagTitle, 2, -0.5, 1.5, "", "");
135
136 m_trackingErrorFlags->GetXaxis()->SetBinLabel(1, "No Error");
137 m_trackingErrorFlags->GetXaxis()->SetBinLabel(2, "Error occured");
138
139 //tracking abort VS time after HER/LER injection and time within a beam cycle
140 m_abortVStimeHER = new TH2F(
141 "TrkAbortVsTimeHER",
142 "Number of Events with Tracking Aborts vs HER injection;Time since last injection [#mus];Time in beam cycle [#mus];Events / bin",
143 100, 0, c_noInjectionTime, 50, 0, c_revolutionTime);
144
145 m_abortVStimeLER = new TH2F(
146 "TrkAbortVsTimeLER",
147 "Number of Events with Tracking Aborts vs LER injection;Time since last injection [#mus];Time in beam cycle [#mus];Events / bin",
148 100, 0, c_noInjectionTime, 50, 0, c_revolutionTime);
149
150 //tracking all VS time after HER/LER injection and time within a beam cycle
151 m_allVStimeHER = new TH2F(
152 "allEvtsVsTimeHER",
153 "Number Of Events vs HER injection;Time since last injection [#mus];Time in beam cycle [#mus];Events / bin",
154 100, 0, c_noInjectionTime, 50, 0, c_revolutionTime);
155
156 m_allVStimeLER = new TH2F(
157 "allEvtsVsTimeLER",
158 "Number of Events vs LER injection;Time since last injection [#mus];Time in beam cycle [#mus];Events / bin",
159 100, 0, c_noInjectionTime, 50, 0, c_revolutionTime);
160
161}
virtual void Run()
Call this to start processing the event data and filling histograms.
This class serves as a base for the TrackDQMModule and AlignDQMModule (and possibly other DQM histogr...
bool histogramsDefined
True if the defineHisto() was called.
std::vector< std::tuple< std::string, std::string, std::string > > m_histogramParameterChanges
Used for changing parameters of histograms via the ProcessHistogramParameterChange function.
virtual void initialize() override
Initializer.
void ProcessHistogramParameterChange(const std::string &name, const std::string &parameter, const std::string &value)
Process one change in histogram parameters.
virtual void event() override
This method is called for each event.
virtual void DefineTrackFitStatus()
Define histograms which require FitStatus.
std::string m_tracksStoreArrayName
StoreArray name where Tracks are written.
virtual TH1F * Create(std::string name, std::string title, int nbinsx, double xlow, double xup, std::string xTitle, std::string yTitle)
Function to create TH1F and add it to the vector of histograms (m_histograms).
virtual void DefineHits()
Define histograms with numbers of hits.
virtual void DefineMomentumCoordinates()
Define histograms with track momentum Pt.
virtual void DefineTracks()
All the following Define- functions should be used in the defineHisto() function to define histograms...
virtual void DefineMomentumAngles()
Define histograms with track momentum Pt.
void runningOnHLT()
function called when the module is run on HLT
std::string m_recoTracksStoreArrayName
StoreArray name where RecoTracks are written.
virtual void DefineHelixParametersAndCorrelations()
Define histograms with helix parameters and their correlations.
virtual void defineHisto() override
Histogram definitions such as TH1(), TH2(), TNtuple(), TTree()....
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
The purpose of this class is to process one event() in TrackDQMModule.
void initialize() override
Module functions.
StoreObjPtr< EventLevelTrackingInfo > m_eventLevelTrackingInfo
Acccess to the EventLevelTrackingInfo object in the datastore.
void event() override
fill of the histograms happens here
StoreArray< RawFTSW > m_rawTTD
Input array for DAQ Status.
TH1F * m_trackingErrorFlags
Monitors the Error flags set by the tracking code.
TH2F * m_allVStimeLER
number of events as a function of time after injection and time within a beam cycle - LER
static constexpr double c_revolutionTime
Beam revolution time in microseconds (approximated).
static constexpr double c_globalClock
Approximated global clock frequency in MHz.
TH2F * m_abortVStimeLER
abort rate as a function of time after injection and time within a beam cycle - LER
TH2F * m_allVStimeHER
number of events as a function of time after injection and time within a beam cycle - HER
virtual void DefineAbortFlagsHistograms()
All the following Define- functions should be used in the defineHisto() function to define histograms...
TH2F * m_abortVStimeHER
abort rate as a function of time after injection and time within a beam cycle - HER
static constexpr double c_noInjectionTime
Defines the range of the x axis of the 2D time histogram.
virtual void defineHisto() override
Histogram definitions such as TH1(), TH2(), TNtuple(), TTree()....
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.