Belle II Software release-09-00-00
EventT0Validation.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 <reconstruction/modules/EventT0Validation/EventT0Validation.h>
10
11#include <TEfficiency.h>
12
13using namespace Belle2;
14
15REG_MODULE(EventT0Validation);
16
17//---------------------------------
19{
20 setPropertyFlags(c_ParallelProcessingCertified); // parallel processing
21 setDescription("Make data quality monitoring plots for EventT0 for bhabha, mu mu, and hadron samples for different trigger (time) sources.");
22 addParam("RootFileName", m_RootFileName, "Name of the ROOT output file.", m_RootFileName);
23}
24
25//---------------------------------
27
28
29//---------------------------------
31{
32 m_outputFile = new TFile(m_RootFileName.c_str(), "RECREATE");
33 if (m_outputFile != nullptr) {
34 m_outputFile->cd();
35
36 int nBins = 400 ;
37 double minT0 = -100 ;
38 double maxT0 = 100 ;
39
40 m_histECLEventT0 = new TH1F("m_histECLEventT0", "ECL EventT0;EventT0 [ns];events / 0.5 ns", nBins, minT0, maxT0);
41 setPlotMetaData(m_histECLEventT0, "Distribution of ECL EventT0s.", "Should be centered around 0.", m_contact);
42 m_histSVDEventT0 = new TH1F("m_histSVDEventT0", "SVD EventT0;EventT0 [ns];events / 0.5 ns", nBins, minT0, maxT0);
43 setPlotMetaData(m_histSVDEventT0, "Distribution of SVD EventT0s.", "Should be centered around 0.", m_contact);
44 m_histTOPEventT0 = new TH1F("m_histTOPEventT0", "TOP EventT0;EventT0 [ns];events / 0.5 ns", nBins, minT0, maxT0);
45 setPlotMetaData(m_histTOPEventT0, "Distribution of TOP EventT0s.", "Should be centered around 0.", m_contact);
46 m_histCDCEventT0 = new TH1F("m_histCDCEventT0", "CDC EventT0;EventT0 [ns];events / 0.5 ns", nBins, minT0, maxT0);
47 setPlotMetaData(m_histCDCEventT0, "Distribution of CDC EventT0s.", "Should be centered around 0.", m_contact);
48 m_histCDCHitBasedEventT0 = new TH1F("m_histCDCHitBasedEventT0", "CDC hit based EventT0;EventT0 [ns];events / 0.5 ns",
49 nBins, minT0, maxT0);
50 setPlotMetaData(m_histCDCHitBasedEventT0, "Distribution of CDC hit based EventT0s.", "Should be centered around 0.", m_contact);
51 m_histCDCChi2EventT0 = new TH1F("m_histCDCChi2EventT0", "CDC FullGrid #chi^{2} EventT0;EventT0 [ns];events / 0.5 ns",
52 nBins, minT0, maxT0);
53 setPlotMetaData(m_histCDCChi2EventT0, "Distribution of CDC FullGrid #chi^{2} EventT0s.", "Should be centered around 0.", m_contact);
54 m_histCDCGridEventT0 = new TH1F("m_histCDCGridEventT0", "CDC Grid search EventT0;EventT0 [ns];events / 0.5 ns", nBins, minT0,
55 maxT0);
56 setPlotMetaData(m_histCDCGridEventT0, "Distribution of CDC Grid search EventT0s.", "Should be centered around 0.", m_contact);
57
59 new TH1D("m_histAlgorithmSourceCounts",
60 "Number of events with EventT0 from each algorithm;Algorithm;Count",
61 11, 0, 11);
63 "Number of events in which an EventT0 was found by each algorithm. "\
64 "Some of the CDC algorithms are only executed if no SVD EventT0 is found.",
65 "Values should be around 1000 (1000 events in validation).", m_contact);
67 new TH1D("m_histAlgorithmSourceCountsActive",
68 "Number of events with EventT0 from each algorithm where it was active;Algorithm;Count",
69 11, 0, 11);
71 "Number of events in which an EventT0 was found by each algorithm if they are executed. "\
72 "Some of the CDC algorithms are only executed if no SVD EventT0 is found.",
73 "Values should be around 1000 (1000 events in validation).", m_contact);
74
75 for (uint i = 0; i < 11; i++) {
76 m_histAlgorithmSourceCounts->GetXaxis()->SetBinLabel(i + 1, c_eventT0Algorithms[i]);
77 m_histAlgorithmSourceCountsActive->GetXaxis()->SetBinLabel(i + 1, c_eventT0Algorithms[i]);
78 }
79 m_histAlgorithmSourceCounts->GetXaxis()->CenterTitle(kTRUE);
80 m_histAlgorithmSourceCounts->GetXaxis()->SetTitleOffset(1.2);
81 m_histAlgorithmSourceCountsActive->GetXaxis()->CenterTitle(kTRUE);
82 m_histAlgorithmSourceCountsActive->GetXaxis()->SetTitleOffset(1.2);
83 }
84
85 m_eventT0.isRequired();
86}
87
88
89
90//---------------------------------
92{
93 if (!m_eventT0.isValid()) {
94 B2WARNING("Missing EventT0, EventT0Validation is skipped.");
95 return;
96 }
97
98 m_histECLEventT0->Reset();
99 m_histSVDEventT0->Reset();
100 m_histTOPEventT0->Reset();
101 m_histCDCEventT0->Reset();
103 m_histCDCChi2EventT0->Reset();
104 m_histCDCGridEventT0->Reset();
105
108
109}
110
111
112//---------------------------------
114{
115 // Determine if there is a valid event t0 to use and then extract the information about it
116 if (!m_eventT0.isValid()) {
117 B2WARNING("Missing EventT0, EventT0Validation is skipped.");
118 return ;
119 }
120
121
122 // Set the different EventT0 values, default is -1000 in case there are no information based on a given detector
123 const double eventT0ECL =
124 m_eventT0->hasTemporaryEventT0(Const::EDetector::ECL) ? m_eventT0->getBestECLTemporaryEventT0()->eventT0 : -1000;
125 const double eventT0CDC =
126 m_eventT0->hasTemporaryEventT0(Const::EDetector::CDC) ? m_eventT0->getBestCDCTemporaryEventT0()->eventT0 : -1000;
127 const double eventT0TOP =
128 m_eventT0->hasTemporaryEventT0(Const::EDetector::TOP) ? m_eventT0->getBestTOPTemporaryEventT0()->eventT0 : -1000;
129 const double eventT0SVD =
130 m_eventT0->hasTemporaryEventT0(Const::EDetector::SVD) ? m_eventT0->getBestSVDTemporaryEventT0()->eventT0 : -1000;
131
132 const auto getCDCEventT0sForAlgorithm = [cdcEventT0s = m_eventT0->getTemporaryEventT0s(Const::EDetector::CDC)](
133 const std::string & algorithm) {
134 std::vector<EventT0::EventT0Component> temporaries;
135 temporaries.reserve(cdcEventT0s.size());
136 for (const auto& evtt0 : cdcEventT0s) {
137 if (evtt0.algorithm == algorithm) {
138 temporaries.push_back(evtt0);
139 }
140 }
141 return temporaries;
142 };
143
144 m_histECLEventT0->Fill(eventT0ECL);
145 m_histSVDEventT0->Fill(eventT0SVD);
146 m_histTOPEventT0->Fill(eventT0TOP);
147 m_histCDCEventT0->Fill(eventT0CDC);
148
149 const auto hitBasedCDCT0 = getCDCEventT0sForAlgorithm("hit based");
150 const auto chi2CDCT0 = getCDCEventT0sForAlgorithm("chi2");
151 const auto gridCDCT0 = getCDCEventT0sForAlgorithm("grid");
152 m_histCDCHitBasedEventT0->Fill(hitBasedCDCT0.empty() ? -1000 : hitBasedCDCT0.back().eventT0);
153 m_histCDCChi2EventT0->Fill(chi2CDCT0.empty() ? -1000 : chi2CDCT0.back().eventT0);
154 m_histCDCGridEventT0->Fill(gridCDCT0.empty() ? -1000 : gridCDCT0.back().eventT0);
155
156 B2DEBUG(20, "eventT0ECL = " << eventT0ECL << " ns") ;
157 B2DEBUG(20, "eventT0CDC = " << eventT0CDC << " ns") ;
158 B2DEBUG(20, "eventT0TOP = " << eventT0TOP << " ns") ;
159 B2DEBUG(20, "eventT0SVD = " << eventT0SVD << " ns") ;
160
161 const bool hasECLEventT0 = m_eventT0->hasTemporaryEventT0(Const::EDetector::ECL);
162 const bool hasSVDEventT0 = m_eventT0->hasTemporaryEventT0(Const::EDetector::SVD);
163 const bool hasTOPEventT0 = m_eventT0->hasTemporaryEventT0(Const::EDetector::TOP);
164 const bool hasCDCEventT0 = m_eventT0->hasTemporaryEventT0(Const::EDetector::CDC);
165 const bool hasCDCHitBasedEventT0 = not hitBasedCDCT0.empty();
166 const bool hasCDCFullGridChi2EventT0 = not chi2CDCT0.empty();
167 const bool hasCDCGridEventT0 = not gridCDCT0.empty();
168
171 m_histAlgorithmSourceCounts->Fill(c_eventT0Algorithms[1], hasECLEventT0);
172 m_histAlgorithmSourceCounts->Fill(c_eventT0Algorithms[2], hasSVDEventT0);
173 m_histAlgorithmSourceCounts->Fill(c_eventT0Algorithms[3], hasTOPEventT0);
174 m_histAlgorithmSourceCounts->Fill(c_eventT0Algorithms[4], hasCDCEventT0);
175 m_histAlgorithmSourceCounts->Fill(c_eventT0Algorithms[5], hasCDCHitBasedEventT0);
176 m_histAlgorithmSourceCounts->Fill(c_eventT0Algorithms[6], hasCDCHitBasedEventT0);
177 m_histAlgorithmSourceCounts->Fill(c_eventT0Algorithms[7], hasCDCFullGridChi2EventT0);
178 m_histAlgorithmSourceCounts->Fill(c_eventT0Algorithms[8], hasCDCFullGridChi2EventT0);
179 m_histAlgorithmSourceCounts->Fill(c_eventT0Algorithms[9], hasCDCGridEventT0);
180 m_histAlgorithmSourceCounts->Fill(c_eventT0Algorithms[10], hasCDCGridEventT0);
181
188 m_histAlgorithmSourceCountsActive->Fill(c_eventT0Algorithms[5], 1); // We always execute hit based search
189 m_histAlgorithmSourceCountsActive->Fill(c_eventT0Algorithms[6], 1); // We always execute hit based search
190 // We only execute the chi2 algorithm if no SVD value is found, but this is the "all" column
192 // We only execute the chi2 algorithm if no SVD value is found
193 m_histAlgorithmSourceCountsActive->Fill(c_eventT0Algorithms[8], not hasSVDEventT0);
194 // We only execute the chi2 algorithm if no SVD value is found, but this is the "all" column
196 // We only execute the grid algorithm if no SVD value is found
197 m_histAlgorithmSourceCountsActive->Fill(c_eventT0Algorithms[10], not hasSVDEventT0);
198}
199
200
202{
203 if (m_outputFile != nullptr) {
204 m_outputFile->cd();
205
206 TEfficiency* algorithmEffi = new TEfficiency(*m_histAlgorithmSourceCounts, *m_histAlgorithmSourceCountsActive);
207 algorithmEffi->SetTitle("Efficiency of finding an EventT0 per Algorithm");
208 algorithmEffi->Write("EventT0AlgorithmEfficiency");
209 algorithmEffi->GetListOfFunctions()->Add(new TNamed("Contact", m_contact));
210 algorithmEffi->GetListOfFunctions()->Add(new TNamed("Description", "Efficiencies of the various EventT0 algorithms."));
211 algorithmEffi->GetListOfFunctions()->Add(new TNamed("Check", "Efficiency should be close to 1 for all active algorithms."));
212 algorithmEffi->GetListOfFunctions()->Add(new TNamed("MetaOptions", "shifter"));
213 algorithmEffi->Write();
214
215 m_outputFile->Write();
216 m_outputFile->Close();
217 }
218
219}
220
221
222void EventT0ValidationModule::setPlotMetaData(TH1* hist, const std::string& description, const std::string& check,
223 const std::string& contact, const std::string& shifter)
224{
225 hist->GetListOfFunctions()->Add(new TNamed("Contact", contact));
226 hist->GetListOfFunctions()->Add(new TNamed("Description", description));
227 hist->GetListOfFunctions()->Add(new TNamed("Check", check));
228 hist->GetListOfFunctions()->Add(new TNamed("MetaOptions", shifter));
229 hist->Write();
230}
TH1F * m_histTOPEventT0
TOP based EventT0 histograms.
TFile * m_outputFile
Output ROOT file.
virtual ~EventT0ValidationModule()
Destructor.
StoreObjPtr< EventT0 > m_eventT0
Store array for event t0.
TH1F * m_histCDCEventT0
CDC based EventT0 histograms.
std::string m_contact
For simplicity, just set the contact once.
virtual void initialize() override
Initialize the module.
TH1D * m_histAlgorithmSourceCounts
Count of events with EventT0 from a given algorithm, numerator for efficiency calculation.
virtual void event() override
This method is called for each event.
TH1F * m_histCDCHitBasedEventT0
CDC hit based EventT0 histograms.
TH1F * m_histCDCChi2EventT0
CDC chi2 based EventT0 histograms.
virtual void endRun() override
This method is called for each run.
TH1F * m_histCDCGridEventT0
CDC grid EventT0 histograms.
std::string m_RootFileName
Name of the output ROOT file.
TH1F * m_histECLEventT0
ECL based EventT0 histograms.
virtual void beginRun() override
This method is called for each run.
const char * c_eventT0Algorithms[11]
EventT0 algorithms for which to calculate fractions of abundance.
EventT0ValidationModule()
Default constructor.
TH1D * m_histAlgorithmSourceCountsActive
Counts of events with EventT0 from a given algorithm, denominator for efficiency calculation.
void setPlotMetaData(TH1 *hist, const std::string &description, const std::string &check, const std::string &contact, const std::string &shifter="shifter")
Set additional plot info: contact, description, shifter.
TH1F * m_histSVDEventT0
SVD based EventT0 histograms.
Base class for Modules.
Definition: Module.h:72
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
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#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.