Belle II Software development
PXDPerformanceCollectorModule.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 <pxd/modules/pxdPerformanceCollector/PXDPerformanceCollectorModule.h>
10
11#include <framework/datastore/RelationArray.h>
12
13#include <vxd/geometry/GeoCache.h>
14#include <pxd/geometry/SensorInfo.h>
15#include <pxd/reconstruction/PXDGainCalibrator.h>
16#include <pxd/utilities/PXDUtilities.h>
17
18#include <TTree.h>
19#include <TH1I.h>
20#include <TH2F.h>
21//#include <TEfficiency.h>
22
23#include <boost/format.hpp>
24
25using namespace std;
26using boost::format;
27using namespace Belle2;
28using namespace Belle2::PXD;
29
30//-----------------------------------------------------------------
31// Register the Module
32//-----------------------------------------------------------------
33REG_MODULE(PXDPerformanceCollector);
34
35//-----------------------------------------------------------------
36// Implementation
37//-----------------------------------------------------------------
38
40 , m_selectedEff(true), m_selectedRes(true)
41 , m_pxd2TrackEvent()
42 , m_deltaD0oSqrt2(0.0), m_deltaZ0oSqrt2(0.0)
43 , m_signal(0), m_estimated(0.0)
44 , m_run(0), m_exp(0)
45{
46 // Set module properties
47 setDescription("Calibration collector module for CDST data.");
49
50 addParam("minPt", m_minPt, "Minimum pt cut", float(1.0));
51 addParam("minClusterCharge", m_minClusterCharge, "Minimum cluster charge cut", int(0));
52 addParam("minClusterSize", m_minClusterSize, "Minimum cluster size cut ", int(2));
53 addParam("maxClusterSize", m_maxClusterSize, "Maximum cluster size cut ", int(6));
54 addParam("nBinsU", m_nBinsU, "Number of gain corrections per sensor along u side", int(4));
55 addParam("nBinsV", m_nBinsV, "Number of gain corrections per sensor along v side", int(6));
56 addParam("gainPayloadName", m_gainName, "Payload name for Gain to be read from DB", string(""));
57 addParam("fillChargeRatioHistogram", m_fillChargeRatioHistogram,
58 "Flag to fill Ratio (cluster charge to the expected MPV) histograms", bool(true));
59 addParam("fillChargeTree", m_fillChargeTree, "Flag to fill cluster charge with the estimated MPV to TTree", bool(false));
60 addParam("fillEventTree", m_fillEventTree, "Flag to fill event tree for validation", bool(false));
61
62 // additional parameters for validation. Considering modularAnalysis for more flexible controls.
63 addParam("minPt4Eff", m_minPt4Eff, "Minimum pt cut for efficiency validation", float(2.0));
64 addParam("maxAbsVx", m_maxAbsVx, "Minimum abs(Vx) cut in cm for efficiency validation", float(0.03));
65 addParam("maxAbsVy", m_maxAbsVy, "Minimum abs(Vy) cut in cm for efficiency validation", float(0.03));
66 addParam("maxAbsVz", m_maxAbsVz, "Minimum abs(Vz) cut in cm for efficiency validation", float(0.155));
67
68 addParam("minPt4Res", m_minPt4Res, "Minimum pt cut for resolution validation", float(1.0));
69 addParam("minSVDHits", m_minSVDHits, "Minimum number of SVD hits foor resolution validation", int(8));
70 addParam("minCDCHits", m_minCDCHits, "Minimum number of CDC hits foor resolution validation", int(21));
71 addParam("maxAbsLambda", m_maxAbsLambda, "Maximum absolute dip angle (lambda)", float(0.5));
72 addParam("minPBetaSinTheta3o2", m_minPBetaSinTheta3o2, "Minimum p*Beta*sin(theta_0)^{3/2}", float(2));
73 addParam("maxAbsZ0", m_maxAbsZ0, "Maximum abs(z0)", float(1));
74 addParam("maxAbsD0", m_maxAbsD0, "Maximum abs(d0)", float(0.3));
75
76
77}
78
79void PXDPerformanceCollectorModule::prepare() // Do your initialise() stuff here
80{
82
83 if (m_nBinsU == 0) {
84 B2WARNING("Number of bins along u side incremented from 0->1");
85 m_nBinsU = 1;
86 }
87
88 if (m_nBinsV == 0) {
89 B2WARNING("Number of bins along v side incremented from 0->1");
90 m_nBinsV = 1;
91 }
92
94 int nPXDSensors = gTools->getNumberOfPXDSensors();
95
96 //-------------------------------------------------------------------------------------
97 // PXDClusterCounter: Count the number of PXDClusters and store charge for each uBin/vBin pair
98 //-------------------------------------------------------------------------------------
99
100 auto hPXDClusterCounter = new TH1I("hPXDClusterCounter", "Number of clusters found in data sample",
101 m_nBinsU * m_nBinsV * nPXDSensors, 0,
102 m_nBinsU * m_nBinsV * nPXDSensors);
103 hPXDClusterCounter->GetXaxis()->SetTitle("bin id");
104 hPXDClusterCounter->GetYaxis()->SetTitle("Number of clusters");
105 auto hPXDClusterChargeRatio = new TH2F("hPXDClusterChargeRatio", "Charge ratio of clusters found in data sample",
106 m_nBinsU * m_nBinsV * nPXDSensors, 0,
107 m_nBinsU * m_nBinsV * nPXDSensors,
108 400, 0., 4.);
109 hPXDClusterChargeRatio->GetXaxis()->SetTitle("bin id");
110 hPXDClusterChargeRatio->GetYaxis()->SetTitle("Cluster charge ratio (relative to expected MPV)");
111 //auto hPXDTrackClusterCounter = new TH1I("hPXDTrackClusterCounter", "Number of clusters found in data sample",
112 //m_nBinsU * m_nBinsV * nPXDSensors, 0,
113 //m_nBinsU * m_nBinsV * nPXDSensors);
114 //hPXDTrackClusterCounter->GetXaxis()->SetTitle("bin id");
115 //hPXDTrackClusterCounter->GetYaxis()->SetTitle("Number of clusters");
116 //auto hPXDTrackPointCounter = new TH1I("hPXDTrackPointCounter", "Number of clusters found in data sample",
117 //m_nBinsU * m_nBinsV * nPXDSensors, 0,
118 //m_nBinsU * m_nBinsV * nPXDSensors);
119 //hPXDTrackPointCounter->GetXaxis()->SetTitle("bin id");
120 //hPXDTrackPointCounter->GetYaxis()->SetTitle("Number of clusters");
121 for (int iSensor = 0; iSensor < nPXDSensors; iSensor++) {
122 for (int uBin = 0; uBin < m_nBinsU; uBin++) {
123 for (int vBin = 0; vBin < m_nBinsV; vBin++) {
124 VxdID id = gTools->getSensorIDFromPXDIndex(iSensor);
125 string sensorDescr = id;
126 hPXDClusterCounter->GetXaxis()->SetBinLabel(iSensor * m_nBinsU * m_nBinsV + uBin * m_nBinsV + vBin + 1,
127 str(format("%1%_%2%_%3%") % sensorDescr % uBin % vBin).c_str());
129 hPXDClusterChargeRatio->GetXaxis()->SetBinLabel(iSensor * m_nBinsU * m_nBinsV + uBin * m_nBinsV + vBin + 1,
130 str(format("%1%_%2%_%3%") % sensorDescr % uBin % vBin).c_str());
131 //hPXDTrackClusterCounter->GetXaxis()->SetBinLabel(iSensor * m_nBinsU * m_nBinsV + uBin * m_nBinsV + vBin + 1,
132 //str(format("%1%_%2%_%3%") % sensorDescr % uBin % vBin).c_str());
133 //hPXDTrackPointCounter->GetXaxis()->SetBinLabel(iSensor * m_nBinsU * m_nBinsV + uBin * m_nBinsV + vBin + 1,
134 //str(format("%1%_%2%_%3%") % sensorDescr % uBin % vBin).c_str());
135 }
136 }
137 }
138
139 registerObject<TH1I>("PXDClusterCounter", hPXDClusterCounter);
141 registerObject<TH2F>("PXDClusterChargeRatio", hPXDClusterChargeRatio);
142
143 //-------------------------------------------------------------------------------------
144 // PXDTrackClusterCounter: Count the number of PXDClustersFrom tracks (the same track selection as for track points)
145 //-------------------------------------------------------------------------------------
146 auto hPXDTrackClusterCounter = (TH1I*)hPXDClusterCounter->Clone("hPXDTrackClusterCounter");
147 hPXDTrackClusterCounter->SetTitle("Number of track clusters");
148 hPXDTrackClusterCounter->GetYaxis()->SetTitle("Number of track clusters");
149 registerObject<TH1I>("PXDTrackClusterCounter", hPXDTrackClusterCounter);
150
151 //-------------------------------------------------------------------------------------
152 // PXDTrackPointCounter: Count the number of PXDClustersFrom tracks (the same track selection as for track points)
153 //-------------------------------------------------------------------------------------
154 auto hPXDTrackPointCounter = (TH1I*)hPXDClusterCounter->Clone("hPXDTrackPointCounter");
155 hPXDTrackPointCounter->SetTitle("Number of track points");
156 hPXDTrackPointCounter->GetYaxis()->SetTitle("Number of track points");
157 registerObject<TH1I>("PXDTrackPointCounter", hPXDTrackPointCounter);
158
159 //----------------------------------------------------------------------
160 // PXDTrees for gain calibration: One tree to store the calibration data for each grid bin
161 //----------------------------------------------------------------------
162
163 if (m_fillChargeTree) // only fill the tree when required
164 for (int iSensor = 0; iSensor < nPXDSensors; iSensor++) {
165 for (int uBin = 0; uBin < m_nBinsU; uBin++) {
166 for (int vBin = 0; vBin < m_nBinsV; vBin++) {
167 VxdID id = gTools->getSensorIDFromPXDIndex(iSensor);
168 auto layerNumber = id.getLayerNumber();
169 auto ladderNumber = id.getLadderNumber();
170 auto sensorNumber = id.getSensorNumber();
171 string treename = str(format("tree_%1%_%2%_%3%_%4%_%5%") % layerNumber % ladderNumber % sensorNumber % uBin % vBin);
172 auto tree = new TTree(treename.c_str(), treename.c_str());
173 tree->Branch<int>("signal", &m_signal);
174 tree->Branch<float>("estimated", &m_estimated);
175 registerObject<TTree>(treename, tree);
176 }
177 }
178 }
179
180 // TEfficiency for validation, no Reset() and may be questionable for merging
181 //auto effPXDLayer1 = new TEfficiency("PXDLayer1Efficiency", "Efficiency of PXD innner layer;#phi;#z [cm];#epsilon", 730, -M_PI, M_PI, 400, -3.2, 6.2);
182 //auto effPXDLayer2 = new TEfficiency("PXDLayer2Efficiency", "Efficiency of PXD outer layer;#phi;#z [cm];#epsilon", 128, 1.4, 2.5, 400, -4.2, 8.2);
183 //registerObject<TEfficiency>("PXDLayer1Efficiency", effPXDLayer1);
184 //registerObject<TEfficiency>("PXDLayer2Efficiency", effPXDLayer2);
185 auto hTotalHitsLayer1 = new TH2F("hTotalHitsLayer1", "Total number of hits from layer 1;#phi;z [cm]", 730, -M_PI, M_PI, 400,
186 -3.2, 6.2);
187 auto hPassedHitsLayer1 = new TH2F("hPassedHitsLayer1", "Passed number of hits from layer 1;#phi;z [cm]", 730, -M_PI, M_PI, 400,
188 -3.2, 6.2);
189 auto hTotalHitsLayer2 = new TH2F("hTotalHitsLayer2", "Total number of hits from layer 2;#phi;z [cm]", 128, 1.4, 2.5, 400,
190 -4.2, 8.2);
191 auto hPassedHitsLayer2 = new TH2F("hPassedHitsLayer2", "Passed number of hits from layer 2;#phi;z [cm]", 128, 1.4, 2.5, 400,
192 -4.2, 8.2);
193 registerObject<TH2F>("hTotalHitsLayer1", hTotalHitsLayer1);
194 registerObject<TH2F>("hPassedHitsLayer1", hPassedHitsLayer1);
195 registerObject<TH2F>("hTotalHitsLayer2", hTotalHitsLayer2);
196 registerObject<TH2F>("hPassedHitsLayer2", hPassedHitsLayer2);
197
198 // trees for correctd d0 and z0 to the IP
199 auto treeD0Z0 = new TTree("tree_d0z0", "TTree of corrected d0 and z0");
200 treeD0Z0->Branch<float>("d0", &m_deltaD0oSqrt2);
201 treeD0Z0->Branch<float>("z0", &m_deltaZ0oSqrt2);
202 registerObject<TTree>("tree_d0z0", treeD0Z0);
203
204 // dbtree
205 auto dbtree = new TTree("dbtree", "dbtree");
206 dbtree->Branch<int>("run", &m_run);
207 dbtree->Branch<int>("exp", &m_exp);
208 dbtree->Branch<PXDGainMapPar>("gainMap", &m_gainMap);
209 registerObject<TTree>("dbtree", dbtree);
210
211 if (m_fillEventTree) {
212 auto tree = new TTree("pxd", "PXD 2-track events");
213 tree->Branch<PXD2TrackEvent>("PXD2TrackEvent", &m_pxd2TrackEvent, 8000, 1);
214 registerObject<TTree>("pxd", tree);
215 }
216}
217
218void PXDPerformanceCollectorModule::startRun() // Do your beginRun() stuff here
219{
220 m_run = m_evtMetaData->getRun();
221 m_exp = m_evtMetaData->getExperiment();
222 if (m_gainName.length()) {
224 m_gainMap = *gainMap;
225 } else {
227 }
228 getObjectPtr<TTree>("dbtree")->Fill();
229}
230
231void PXDPerformanceCollectorModule::collect() // Do your event() stuff here
232{
233 // If no input, nothing to do
234 if (!m_pxd2TrackEvents) return;
235
236 // Beam spot
237 DBObjPtr<BeamSpot> beamSpotDB;
238 auto ip = ROOT::Math::XYZVector(beamSpotDB->getIPPosition());
239
240 // Actually only one event holder / event
241 for (auto& pxd2TrackEvent : m_pxd2TrackEvents) {
242 m_selectedRes = true;
243 m_selectedEff = true;
244
245 auto vertex = pxd2TrackEvent.getVertex();
246 vertex -= ip; // correct vertex relative to ip
247 if (fabs(vertex.X()) > m_maxAbsVx ||
248 fabs(vertex.Y()) > m_maxAbsVy ||
249 fabs(vertex.Z()) > m_maxAbsVz)
250 m_selectedEff = false;
251
252 // track level selection and collection
253 collectFromTrack(pxd2TrackEvent.getTrackP());
254 collectFromTrack(pxd2TrackEvent.getTrackM());
255
256 // event level collection
257 collectDeltaIP(pxd2TrackEvent);
258
259 if (m_fillEventTree) {
260 m_pxd2TrackEvent = pxd2TrackEvent;
261 getObjectPtr<TTree>("pxd")->Fill();
262 }
263 }
264
265}
266
268{
269 if (!m_selectedRes) return;
270 auto d0p_1 = event.getTrackP().d0p;
271 auto d0p_2 = event.getTrackM().d0p;
272 auto z0p_1 = event.getTrackP().z0p;
273 auto z0p_2 = event.getTrackM().z0p;
274 m_deltaD0oSqrt2 = (d0p_1 + d0p_2) / sqrt(2.);
275 m_deltaZ0oSqrt2 = (z0p_1 - z0p_2) / sqrt(2.);
276
277 // Fill the tree of impact parameters
278 getObjectPtr<TTree>("tree_d0z0")->Fill();
279}
280
282{
283 auto gTools = VXD::GeoCache::getInstance().getGeoTools();
284 bool selected4Gain = true;
285 bool selected4Eff = true;
286
287 if (track.pt < m_minPt) selected4Gain = false;
288 if (track.pt < m_minPt4Eff) selected4Eff = false; // just applied on track level
289
290 // Track level filtering for resolution validation
291 if (track.pt < m_minPt4Res) m_selectedRes = false;
292 if (track.nPXDHits < 1 || track.nCDCHits < m_minCDCHits || track.nSVDHits < m_minSVDHits)
293 m_selectedRes = false;
294 if (fabs(track.d0p) > m_maxAbsD0 || fabs(track.z0p) > m_maxAbsZ0)
295 m_selectedRes = false;
296 auto lambda0 = atan(track.tanLambda);
297 if (fabs(lambda0) > m_maxAbsLambda)
298 m_selectedRes = false;
299 auto sinTheta0 = 1. / sqrt(1. + pow(track.tanLambda, 2));
300 auto pBetaSinTheta3o2 = track.pt * 1.0 * sqrt(sinTheta0);
301 if (pBetaSinTheta3o2 < m_minPBetaSinTheta3o2)
302 m_selectedRes = false;
303
304 for (auto& trackCluster : track.trackClusters) {
305 bool selectedCluster = true;
306 auto cluster = trackCluster.cluster;
307 auto intersection = trackCluster.intersection;
308 auto usedInTrack = trackCluster.usedInTrack;
309
310 // Check for valid cluster and intersection
311 if (!usedInTrack || intersection.chargeMPV <= 0)
312 selectedCluster = false;
313
314 // Apply cluster selection cuts
315 if (cluster.charge < m_minClusterCharge || cluster.size < m_minClusterSize || cluster.size > m_maxClusterSize)
316 selectedCluster = false;
317
318 if (cluster.pxdID <= 0) {
319 B2FATAL("Unexpected cluster module id : " << cluster.pxdID);
320
321 }
322
323 // Fill tree or histograms for gain calibration
324 if (selected4Gain && selectedCluster) {
325
326 // Compute variables from cluster needed for gain estimation
327 m_signal = cluster.charge;
328 m_estimated = intersection.chargeMPV;
329
330 VxdID sensorID = PXD::getVxdIDFromPXDModuleID(cluster.pxdID);
331 const PXD::SensorInfo& Info = dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::getInstance().getSensorInfo(sensorID));
332 auto uID = Info.getUCellID(cluster.posU);
333 auto vID = Info.getVCellID(cluster.posV);
334 auto iSensor = gTools->getPXDSensorIndex(sensorID);
335 auto layerNumber = sensorID.getLayerNumber();
336 auto ladderNumber = sensorID.getLadderNumber();
337 auto sensorNumber = sensorID.getSensorNumber();
338 auto uBin = PXD::PXDGainCalibrator::getInstance().getBinU(sensorID, uID, vID, m_nBinsU);
339 auto vBin = PXD::PXDGainCalibrator::getInstance().getBinV(sensorID, vID, m_nBinsV);
340 // Calculate bin ID based on iSensor, uBin, vBin and number of bins in u/v
341 int binID = iSensor * m_nBinsU * m_nBinsV + uBin * m_nBinsV + vBin;
342
343 // Increment the counter
344 getObjectPtr<TH1I>("PXDClusterCounter")->Fill(binID);
345
346 // Fill variabels into tree
347 if (m_fillChargeTree) {
348 string treename = str(format("tree_%1%_%2%_%3%_%4%_%5%") % layerNumber % ladderNumber % sensorNumber % uBin % vBin);
349 getObjectPtr<TTree>(treename)->Fill();
350 }
351
352 // Fill cluster charge ratio histogram if enabled
354 double ratio = m_signal / m_estimated;
355 auto axis = getObjectPtr<TH2F>("PXDClusterChargeRatio")->GetYaxis();
356 double maxY = axis->GetBinCenter(axis->GetNbins());
357 // Manipulate too large ratio for better estimation on median.
358 getObjectPtr<TH2F>("PXDClusterChargeRatio")->Fill(binID, TMath::Min(ratio, maxY));
359 }
360 }
361
362 // Fill effciency
363 if (m_selectedEff && selected4Eff) {
364 auto x = intersection.x;
365 auto y = intersection.y;
366 auto phi = atan2(y, x);
367 auto z = intersection.z;
368
369 // Get uBin and vBin from a global point.
370 VxdID sensorID = PXD::getVxdIDFromPXDModuleID(cluster.pxdID);
371 const PXD::SensorInfo& Info = dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::getInstance().getSensorInfo(sensorID));
372 auto localPoint = Info.pointToLocal(ROOT::Math::XYZVector(x, y, z));
373 auto uID = Info.getUCellID(localPoint.X());
374 auto vID = Info.getVCellID(localPoint.Y());
375 auto iSensor = gTools->getPXDSensorIndex(sensorID);
376 auto uBin = PXD::PXDGainCalibrator::getInstance().getBinU(sensorID, uID, vID, m_nBinsU);
377 auto vBin = PXD::PXDGainCalibrator::getInstance().getBinV(sensorID, vID, m_nBinsV);
378
379 // Filling counters
380 getObjectPtr<TH1I>("PXDTrackPointCounter")->Fill(iSensor * m_nBinsU * m_nBinsV + uBin * m_nBinsV + vBin);
381 if (usedInTrack)
382 getObjectPtr<TH1I>("PXDTrackClusterCounter")->Fill(iSensor * m_nBinsU * m_nBinsV + uBin * m_nBinsV + vBin);
383
384 // Filling 2D histograms
385 if (cluster.pxdID < 2000) {
386 //getObjectPtr<TEfficiency>("PXDLayer1Efficiency")->Fill(usedInTrack,phi,z);
387 getObjectPtr<TH2F>("hTotalHitsLayer1")->Fill(phi, z);
388 if (usedInTrack)
389 getObjectPtr<TH2F>("hPassedHitsLayer1")->Fill(phi, z);
390 } else {
391 //getObjectPtr<TEfficiency>("PXDLayer2Efficiency")->Fill(usedInTrack,phi,z);
392 getObjectPtr<TH2F>("hTotalHitsLayer2")->Fill(phi, z);
393 if (usedInTrack)
394 getObjectPtr<TH2F>("hPassedHitsLayer2")->Fill(phi, z);
395 }
396 }
397
398 } // end loop trackClusters
399
400}
Calibration collector module base class.
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
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
Class PXD2TrackEvent: Event data container for performance and calibration studies.
The payload class for PXD gain corrections.
Definition: PXDGainMapPar.h:43
StoreArray< PXD2TrackEvent > m_pxd2TrackEvents
Required input PXD2TrackEvent
bool m_selectedRes
Flag of selection for resolution validation.
float m_estimated
Estimated cluster charge in ADU.
float m_maxAbsVy
Maximum absolute value for y coordinate of vertex.
int m_nBinsV
Number of corrections per sensor along v side.
PXDPerformanceCollectorModule()
Constructor: Sets the description, the properties and the parameters of the module.
bool m_fillChargeRatioHistogram
Flag to fill cluster charge ratio (relative to expected MPV) histograms.
float m_minPBetaSinTheta3o2
Minimum p*Beta*sin(theta_0)^{3/2}.
bool m_selectedEff
Flag of selection for efficiency validation.
std::string m_gainName
Payload name for Gain to be read from DB.
PXDGainMapPar m_gainMap
GainMap to be stored in dbtree.
int m_nBinsU
Number of corrections per sensor along u side.
StoreObjPtr< EventMetaData > m_evtMetaData
Required input EventMetaData.
float m_maxAbsLambda
Maximum absolute dip angle (lambda)
bool m_fillEventTree
Flag to fill event tree for validation.
float m_minPt4Res
Minimum pt cut for resolution monitoring.
int m_exp
Experiment number to be stored in dbtree.
int m_minCDCHits
Minimum number of CDC hits for resolution.
float m_maxAbsVx
Maximum absolute value for x coordinate of vertex.
bool m_fillChargeTree
Flag to fill cluster charge and its estimated MPV in TTree.
void collectFromTrack(const PXD2TrackEvent::baseType &track)
Collect info on track level.
int m_minSVDHits
Minimum number of SVD hits for resolution.
float m_maxAbsVz
Maximum absolute value for z coordinate of vertex.
float m_minPt4Eff
Minimum pt cut for efficiency monitoring.
void collectDeltaIP(const PXD2TrackEvent &event)
Collect info for impact parameter study on event level.
std::string m_store2TrackEventsName
Name of the collection to use for PXD2TrackEvents.
int m_run
Run number to be stored in dbtree.
unsigned short getBinV(VxdID id, unsigned int vid) const
Get gain correction bin along sensor v side.
unsigned short getBinU(VxdID id, unsigned int uid, unsigned int vid) const
Get gain correction bin along sensor u side.
static PXDGainCalibrator & getInstance()
Main (and only) way to access the PXDGainCalibrator.
Specific implementation of SensorInfo for PXD Sensors which provides additional pixel specific inform...
Definition: SensorInfo.h:23
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:67
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
const GeoTools * getGeoTools()
Return a raw pointer to a GeoTools object.
Definition: GeoCache.h:142
unsigned short getNumberOfPXDSensors() const
Get number of PXD sensors.
Definition: GeoTools.h:133
ROOT::Math::XYZVector pointToLocal(const ROOT::Math::XYZVector &global, bool reco=false) const
Convert a point from global to local coordinates.
int getVCellID(double v, bool clamp=false) const
Return the corresponding pixel/strip ID of a given v coordinate.
int getUCellID(double u, double v=0, bool clamp=false) const
Return the corresponding pixel/strip ID of a given u coordinate.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
baseType getSensorNumber() const
Get the sensor id.
Definition: VxdID.h:100
baseType getLadderNumber() const
Get the ladder id.
Definition: VxdID.h:98
baseType getLayerNumber() const
Get the layer id.
Definition: VxdID.h:96
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 atan(double a)
atan for double
Definition: beamHelpers.h:34
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Namespace to encapsulate code needed for simulation and reconstrucion of the PXD.
VxdID getVxdIDFromPXDModuleID(const unsigned short &id)
Helper function to get VxdID from DHE id like module iid.
Definition: PXDUtilities.h:86
Abstract base class for different kinds of events.
STL namespace.
Struct to hold variables from a track which contains a vector of data type like TrackCluster.