Belle II Software development
SVDSpacePointQICalibrationModule.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 <svd/modules/svdSpacePointCreator/SpacePointHelperFunctions.h>
10#include <svd/modules/svdSpacePointCreator/SVDSpacePointQICalibrationModule.h>
11//#include <svd/modules/svdSpacePointCreator/SVDSpacePointCreatorModule.h>
12
13#include <math.h>
14
15using namespace std;
16using namespace Belle2;
17
18
19REG_MODULE(SVDSpacePointQICalibration);
20
22 Module()
23{
24
25 setDescription("Generate PDFs used in assigning quality to SpacePoints.");
27
28 // 1. Collections.
29 addParam("SVDClusters", m_svdClustersName,
30 "SVDCluster collection name", string(""));
31
32 addParam("RecoTracks", m_recoTracksName,
33 "RecoTracks collection name", string(""));
34
35 // 2.Modification parameters:
36 addParam("NameOfInstance", m_nameOfInstance,
37 "allows the user to set an identifier for this module. Usefull if one wants to use several instances of that module", string(""));
38 addParam("binSize", m_binSize, "Number of bins in charge distribution.",
39 int(50));
40
41 addParam("maxClusterSize", m_maxClusterSize, "Max number of strips the PDF are separated into.",
42 int(5));
43
44 addParam("useLegacyNaming", m_useLegacyNaming,
45 "Use legacy pdf naming", bool(true));
46
47 addParam("outputFileName", m_outputFileName,
48 "Name of output file containing pdfs", std::string("spacePointQICalibration.root"));
49}
50
51
52
54{
55 // prepare all store- and relationArrays:
58
59
60
62 for (auto& layers : geo.getLayers(VXD::SensorInfoBase::SVD)) {
63 for (auto& ladders : geo.getLadders(layers)) {
64 for (auto& sensors : geo.getSensors(ladders)) {
65
66 for (int uSize = 1; uSize <= m_maxClusterSize; uSize++) {
67 for (int vSize = 1; vSize <= m_maxClusterSize; vSize++) {
68 std::string sensorName = "";
69 std::string errorStringName = "";
70
71
72 spPDFName(sensors, uSize, vSize, m_maxClusterSize, sensorName, errorStringName, m_useLegacyNaming);
73
74 std::string signalHistName = sensorName + "signal";
75 std::string backgroundHistName = sensorName + "background";
76
77 if (signalHistMap.count(sensorName) == 0) {
78 TH2F* sigHist = new TH2F(signalHistName.c_str(), "", m_binSize, 0, 200000, m_binSize, 0, 200000);
79 TH2F* bkgHist = new TH2F(backgroundHistName.c_str(), "", m_binSize, 0, 200000, m_binSize, 0, 200000);
80
81 signalHistMap[sensorName] = sigHist;
82 backgroundHistMap[sensorName] = bkgHist;
83 }
84 }
85 }
86 }
87 }
88 }
89 TH2F* timeSignal = new TH2F("timeSignal", "", 40, -100, 100, 40, -100, 100);
90 TH2F* timeBackground = new TH2F("timeBackground", "", 40, -100, 100, 40, -100, 100);
91 TH2F* sizeSignal = new TH2F("sizeSignal", "", 20, 0, 20, 20, 0, 20);
92 TH2F* sizeBackground = new TH2F("sizeBackground", "", 20, 0, 20, 20, 0, 20);
93
94 signalHistMap["timeSignal"] = timeSignal;
95 signalHistMap["sizeSignal"] = sizeSignal;
96 backgroundHistMap["timeBackground"] = timeBackground;
97 backgroundHistMap["sizeBackground"] = sizeBackground;
98
99}
100
101
102
104{
105 for (auto& track : m_recoTracks) {
106 if (track.wasFitSuccessful() == 1) {
107 if (track.hasSVDHits() == 1) {
108 for (auto& svdHit : track.getSVDHitList()) {
109 if (svdHit->isUCluster()) {
110 for (auto& svdHit2 : track.getSVDHitList()) {
111 if (svdHit2->isUCluster() == 0) {
112 if (svdHit->getSensorID() == svdHit2->getSensorID()) {
113 std::string pdfName;
114 std::string errorStringName;
115 spPDFName(svdHit->getSensorID(), svdHit->getSize(), svdHit2->getSize(), m_maxClusterSize, pdfName, errorStringName,
117 auto sigHist = signalHistMap.at(pdfName);
118 sigHist->Fill(svdHit->getCharge(), svdHit2->getCharge());
119 auto timeSigHist = signalHistMap.at("timeSignal");
120 timeSigHist->Fill(svdHit->getClsTime(), svdHit2->getClsTime());
121 auto sizeSigHist = signalHistMap.at("sizeSignal");
122 sizeSigHist->Fill(svdHit->getSize(), svdHit2->getSize());
123 }
124 }
125 }
126 }
127 }
128 }
129 }
130 }
131
132
133 for (auto& uCluster : m_svdClusters) {
134 if (uCluster.isUCluster() == 1) {
135 for (auto& vCluster : m_svdClusters) {
136 if (vCluster.isUCluster() == 0) {
137 if (uCluster.getSensorID() == vCluster.getSensorID()) {
138 std::string pdfName;
139 std::string errorStringName;
140 spPDFName(uCluster.getSensorID(), uCluster.getSize(), vCluster.getSize(), m_maxClusterSize, pdfName, errorStringName,
142 auto bkgHist = backgroundHistMap.at(pdfName);
143 bkgHist->Fill(uCluster.getCharge(), vCluster.getCharge());
144 auto timeBkgHist = backgroundHistMap.at("timeBackground");
145 timeBkgHist->Fill(uCluster.getClsTime(), vCluster.getClsTime());
146 auto sizeBkgHist = backgroundHistMap.at("sizeBackground");
147 sizeBkgHist->Fill(uCluster.getSize(), vCluster.getSize());
148 }
149 }
150 }
151 }
152 }
153}
154
155
156
158{
159 //Open rootfile
160 TFile* f = new TFile(m_outputFileName.c_str(), "RECREATE");
161 f->cd();
162 std::vector<std::string> usedSensors;
164 for (auto& layers : geo.getLayers(VXD::SensorInfoBase::SVD)) {
165 for (auto& ladders : geo.getLadders(layers)) {
166 for (auto& sensors : geo.getSensors(ladders)) {
167 for (int uSize = 1; uSize <= m_maxClusterSize; uSize++) {
168 for (int vSize = 1; vSize <= m_maxClusterSize; vSize++) {
169 std::string sensorName;
170 std::string errorName;
171
172
173 spPDFName(sensors, uSize, vSize, m_maxClusterSize, sensorName, errorName, m_useLegacyNaming);
174
175 if (std::find(usedSensors.begin(), usedSensors.end(), sensorName.c_str()) == usedSensors.end()) {
176 usedSensors.push_back(sensorName);
177 TH2F* probHist = new TH2F(sensorName.c_str(), "", m_binSize, 0, 200000, m_binSize, 0, 200000);
178 TH2F* errorHist = new TH2F(errorName.c_str(), "", m_binSize, 0, 200000, m_binSize, 0, 200000);
179 calculateProb(signalHistMap[sensorName], backgroundHistMap[sensorName], probHist);
180 calculateError(signalHistMap[sensorName], backgroundHistMap[sensorName], errorHist);
181
182 probHist->Write();
183 errorHist->Write();
184 }
185 }
186 }
187 }
188 }
189 }
190 TH2F* timeProb = new TH2F("timeProb", "", 40, -100, 100, 40, -100, 100);
191 TH2F* timeError = new TH2F("timeError", "", 40, -100, 100, 40, -100, 100);
192 TH2F* sizeProb = new TH2F("sizeProb", "", 20, 0, 20, 20, 0, 20);
193 TH2F* sizeError = new TH2F("sizeError", "", 20, 0, 20, 20, 0, 20);
194 calculateProb(signalHistMap["timeSignal"], backgroundHistMap["timeBackground"], timeProb);
195 calculateError(signalHistMap["timeSignal"], backgroundHistMap["timeBackground"], timeError);
196 calculateProb(signalHistMap["sizeSignal"], backgroundHistMap["sizeBackground"], sizeProb);
197 calculateError(signalHistMap["sizeSignal"], backgroundHistMap["sizeBackground"], sizeError);
198 timeProb->Write();
199 timeError->Write();
200 sizeProb->Write();
201 sizeError->Write();
202
203 f->Close();
204}
205
206
207void SVDSpacePointQICalibrationModule::calculateProb(TH2F* signal, TH2F* background, TH2F* probability)
208{
209 signal->Divide(background);
210 probability->Add(signal);
211
212}
213void SVDSpacePointQICalibrationModule::calculateError(TH2F* signal, TH2F* background, TH2F* error)
214{
215 int imax = signal->GetXaxis()->GetNbins();
216 int jmax = signal->GetYaxis()->GetNbins();
217 for (int i = 1; i <= imax; i++) {
218 for (int j = 1; j <= jmax; j++) {
219 int bkg = background->GetBinContent(i, j);
220 int sig = signal->GetBinContent(i, j);
221 double var = ((sig + 1) * (sig + 2)) / ((bkg + 2) * (bkg + 3)) -
222 ((sig + 1) * (sig + 1)) / ((bkg + 2) * (bkg + 2));
223 double err = sqrt(var);
224 error->SetBinContent(i, j, err);
225 }
226 }
227}
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
std::string m_svdClustersName
SVDCluster collection name.
std::string m_nameOfInstance
allows the user to set an identifier for this module.
virtual void initialize() override
Init the module.
StoreArray< SVDCluster > m_svdClusters
the storeArray for svdClusters as member, is faster than recreating link for each event
void calculateError(TH2F *signal, TH2F *background, TH2F *error)
compute error
std::map< std::string, TH2F * > backgroundHistMap
map for background histograms
std::string m_recoTracksName
RecoTracks collection name.
std::map< std::string, TH2F * > signalHistMap
map of signal histograms
int m_maxClusterSize
max numnber of strips the PDF are separated into
StoreArray< RecoTrack > m_recoTracks
RecoTrack store array.
void calculateProb(TH2F *signal, TH2F *background, TH2F *probability)
compute probability
int m_binSize
number of bins in charge distribution
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
const std::set< Belle2::VxdID > getLayers(SensorInfoBase::SensorType sensortype=SensorInfoBase::VXD)
Return a set of all known Layers.
Definition: GeoCache.cc:176
const std::set< Belle2::VxdID > & getSensors(Belle2::VxdID ladder) const
Return a set of all sensor IDs belonging to a given ladder.
Definition: GeoCache.cc:204
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
const std::set< Belle2::VxdID > & getLadders(Belle2::VxdID layer) const
Return a set of all ladder IDs belonging to a given layer.
Definition: GeoCache.cc:193
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
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
void spPDFName(const VxdID &sensor, int uSize, int vSize, int maxClusterSize, std::string &PDFName, std::string &errorPDFName, bool useLegacyNaming)
Function to set name of PDF for spacePoint quality estimation.
Abstract base class for different kinds of events.
STL namespace.