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