Belle II Software development
PXDPerformanceVariablesCollectorModule.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/pxdPerformanceVariablesCollector/PXDPerformanceVariablesCollectorModule.h>
10
11#include <framework/database/DBObjPtr.h>
12#include <analysis/dataobjects/ParticleList.h>
13#include <mdst/dbobjects/BeamSpot.h>
14#include <mdst/dataobjects/Track.h>
15#include <tracking/dataobjects/RecoTrack.h>
16
17#include <TTree.h>
18#include <TH1I.h>
19#include <TH2F.h>
20
21#include <boost/format.hpp>
22
23using namespace std;
24using boost::format;
25using namespace Belle2;
26using namespace Belle2::PXD;
27
28//-----------------------------------------------------------------
29// Register the Module
30//-----------------------------------------------------------------
31REG_MODULE(PXDPerformanceVariablesCollector);
32
33//-----------------------------------------------------------------
34// Implementation
35//-----------------------------------------------------------------
36
38 , m_selected4Eff(false)
39 , m_deltaD0oSqrt2(0.0), m_deltaZ0oSqrt2(0.0)
40 , m_signal(0), m_estimated(0.0)
41 , m_run(0), m_exp(0)
42{
43 // Set module properties
44 setDescription("Calibration collector module for CDST data.");
46
47 addParam("minClusterCharge", m_minClusterCharge, "Minimum cluster charge cut", int(0));
48 addParam("minClusterSize", m_minClusterSize, "Minimum cluster size cut ", int(2));
49 addParam("maxClusterSize", m_maxClusterSize, "Maximum cluster size cut ", int(6));
50 addParam("nBinsU", m_nBinsU, "Number of gain corrections per sensor along u side", int(4));
51 addParam("nBinsV", m_nBinsV, "Number of gain corrections per sensor along v side", int(6));
52 addParam("gainPayloadName", m_gainName, "Payload name for Gain to be read from DB", string(""));
53 addParam("useClusterPosition", m_useClusterPosition,
54 "Flag to use cluster position rather than track point to group pixels for gain calibration.", bool(false));
55 addParam("fillChargeRatioHistogram", m_fillChargeRatioHistogram,
56 "Flag to fill Ratio (cluster charge to the expected MPV) histograms", bool(true));
57 addParam("fillChargeTree", m_fillChargeTree, "Flag to fill cluster charge with the estimated MPV to TTree", bool(false));
58 addParam("maskedDistance", m_maskedDistance, "Distance inside which no masked pixel or sensor border is allowed", int(10));
59
60 // Particle list names
61 addParam("PList4GainName", m_PList4GainName, "Name of the particle list for gain calibration and efficiency study",
62 std::string("e+:gain"));
63 addParam("PList4EffName", m_PList4EffName, "Name of the particle list for event selection in efficiency study",
64 std::string("vpho:eff"));
65 addParam("PList4ResName", m_PList4ResName, "Name of the particle list for resolution study", std::string("vpho:res"));
66
67}
68
69void PXDPerformanceVariablesCollectorModule::prepare() // Do your initialise() stuff here
70{
71
72 if (m_nBinsU == 0) {
73 B2WARNING("Number of bins along u side incremented from 0->1");
74 m_nBinsU = 1;
75 }
76
77 if (m_nBinsV == 0) {
78 B2WARNING("Number of bins along v side incremented from 0->1");
79 m_nBinsV = 1;
80 }
81
83 int nPXDSensors = gTools->getNumberOfPXDSensors();
84
85 //-------------------------------------------------------------------------------------
86 // PXDClusterCounter: Count the number of PXDClusters and store charge for each uBin/vBin pair
87 //-------------------------------------------------------------------------------------
88
89 auto hPXDClusterCounter = new TH1I("hPXDClusterCounter", "Number of clusters found in data sample",
90 m_nBinsU * m_nBinsV * nPXDSensors, 0,
91 m_nBinsU * m_nBinsV * nPXDSensors);
92 hPXDClusterCounter->GetXaxis()->SetTitle("bin id");
93 hPXDClusterCounter->GetYaxis()->SetTitle("Number of clusters");
94 auto hPXDClusterChargeRatio = new TH2F("hPXDClusterChargeRatio", "Charge ratio of clusters found in data sample",
95 m_nBinsU * m_nBinsV * nPXDSensors, 0,
96 m_nBinsU * m_nBinsV * nPXDSensors,
97 400, 0., 4.);
98 hPXDClusterChargeRatio->GetXaxis()->SetTitle("bin id");
99 hPXDClusterChargeRatio->GetYaxis()->SetTitle("Cluster charge ratio (relative to expected MPV)");
100 for (int iSensor = 0; iSensor < nPXDSensors; iSensor++) {
101 for (int uBin = 0; uBin < m_nBinsU; uBin++) {
102 for (int vBin = 0; vBin < m_nBinsV; vBin++) {
103 VxdID id = gTools->getSensorIDFromPXDIndex(iSensor);
104 string sensorDescr = id;
105 hPXDClusterCounter->GetXaxis()->SetBinLabel(iSensor * m_nBinsU * m_nBinsV + uBin * m_nBinsV + vBin + 1,
106 str(format("%1%_%2%_%3%") % sensorDescr % uBin % vBin).c_str());
108 hPXDClusterChargeRatio->GetXaxis()->SetBinLabel(iSensor * m_nBinsU * m_nBinsV + uBin * m_nBinsV + vBin + 1,
109 str(format("%1%_%2%_%3%") % sensorDescr % uBin % vBin).c_str());
110 }
111 }
112 }
113
114 registerObject<TH1I>("PXDClusterCounter", hPXDClusterCounter);
116 registerObject<TH2F>("PXDClusterChargeRatio", hPXDClusterChargeRatio);
117
118 //-------------------------------------------------------------------------------------
119 // PXDTrackClusterCounter: Count the number of PXD clusters from tracks (the same track selection as for PXDTrackPointCounter)
120 //-------------------------------------------------------------------------------------
121 auto hPXDTrackClusterCounter = (TH1I*)hPXDClusterCounter->Clone("hPXDTrackClusterCounter");
122 hPXDTrackClusterCounter->SetTitle("Number of track clusters");
123 hPXDTrackClusterCounter->GetYaxis()->SetTitle("Number of track clusters");
124 registerObject<TH1I>("PXDTrackClusterCounter", hPXDTrackClusterCounter);
125
126 //-------------------------------------------------------------------------------------
127 // PXDTrackPointCounter: Count the number of PXD track points
128 //-------------------------------------------------------------------------------------
129 auto hPXDTrackPointCounter = (TH1I*)hPXDClusterCounter->Clone("hPXDTrackPointCounter");
130 hPXDTrackPointCounter->SetTitle("Number of track points");
131 hPXDTrackPointCounter->GetYaxis()->SetTitle("Number of track points");
132 registerObject<TH1I>("PXDTrackPointCounter", hPXDTrackPointCounter);
133
134 //-------------------------------------------------------------------------------------
135 // PXDSelTrackClusterCounter: Count the number of PXD clusters from tracks (the same track selection as for PXDSelTrackPointCounter)
136 //-------------------------------------------------------------------------------------
137 auto hPXDSelTrackClusterCounter = (TH1I*)hPXDClusterCounter->Clone("hPXDSelTrackClusterCounter");
138 hPXDSelTrackClusterCounter->SetTitle("Number of selected track clusters (the same selectrion as for PXDSelTrackPointCounter)");
139 hPXDSelTrackClusterCounter->GetYaxis()->SetTitle("Number of track clusters");
140 registerObject<TH1I>("PXDSelTrackClusterCounter", hPXDSelTrackClusterCounter);
141
142 //-------------------------------------------------------------------------------------
143 // PXDSelTrackPointCounter: Count the number of PXD track points if they are away from hot/dead pixels
144 //-------------------------------------------------------------------------------------
145 auto hPXDSelTrackPointCounter = (TH1I*)hPXDClusterCounter->Clone("hPXDSelTrackPointCounter");
146 hPXDSelTrackPointCounter->SetTitle("Number of selected track points excluding hot/dead regions");
147 hPXDSelTrackPointCounter->GetYaxis()->SetTitle("Number of track points");
148 registerObject<TH1I>("PXDSelTrackPointCounter", hPXDSelTrackPointCounter);
149
150 //----------------------------------------------------------------------
151 // PXDTrees for gain calibration: One tree to store the calibration data for each grid bin
152 //----------------------------------------------------------------------
153
154 if (m_fillChargeTree) // only fill the tree when required
155 for (int iSensor = 0; iSensor < nPXDSensors; iSensor++) {
156 for (int uBin = 0; uBin < m_nBinsU; uBin++) {
157 for (int vBin = 0; vBin < m_nBinsV; vBin++) {
158 VxdID id = gTools->getSensorIDFromPXDIndex(iSensor);
159 auto layerNumber = id.getLayerNumber();
160 auto ladderNumber = id.getLadderNumber();
161 auto sensorNumber = id.getSensorNumber();
162 string treename = str(format("tree_%1%_%2%_%3%_%4%_%5%") % layerNumber % ladderNumber % sensorNumber % uBin % vBin);
163 auto tree = new TTree(treename.c_str(), treename.c_str());
164 tree->Branch<int>("signal", &m_signal);
165 tree->Branch<float>("estimated", &m_estimated);
166 registerObject<TTree>(treename, tree);
167 }
168 }
169 }
170
171 auto hTotalHitsLayer1 = new TH2F("hTotalHitsLayer1", "Total number of hits from layer 1;#phi;z [cm]", 720, -M_PI, M_PI, 400,
172 -3.2, 6.2);
173 auto hPassedHitsLayer1 = new TH2F("hPassedHitsLayer1", "Passed number of hits from layer 1;#phi;z [cm]", 720, -M_PI, M_PI, 400,
174 -3.2, 6.2);
175 auto hTotalHitsLayer2 = new TH2F("hTotalHitsLayer2", "Total number of hits from layer 2;#phi;z [cm]", 720, -M_PI, M_PI, 400,
176 -4.5, 8.5);
177 auto hPassedHitsLayer2 = new TH2F("hPassedHitsLayer2", "Passed number of hits from layer 2;#phi;z [cm]", 720, -M_PI, M_PI, 400,
178 -4.5, 8.5);
179 registerObject<TH2F>("hTotalHitsLayer1", hTotalHitsLayer1);
180 registerObject<TH2F>("hPassedHitsLayer1", hPassedHitsLayer1);
181 registerObject<TH2F>("hTotalHitsLayer2", hTotalHitsLayer2);
182 registerObject<TH2F>("hPassedHitsLayer2", hPassedHitsLayer2);
183
184 // trees for correctd d0 and z0 to the IP
185 auto treeD0Z0 = new TTree("tree_d0z0", "TTree of delta d0 (z0) over sqrt(2)");
186 treeD0Z0->Branch<float>("d0", &m_deltaD0oSqrt2);
187 treeD0Z0->Branch<float>("z0", &m_deltaZ0oSqrt2);
188 registerObject<TTree>("tree_d0z0", treeD0Z0);
189
190 // dbtree
191 auto dbtree = new TTree("dbtree", "dbtree");
192 dbtree->Branch<int>("run", &m_run);
193 dbtree->Branch<int>("exp", &m_exp);
194 dbtree->Branch<PXDGainMapPar>("gainMap", &m_gainMap);
195 registerObject<TTree>("dbtree", dbtree);
196
197}
198
199void PXDPerformanceVariablesCollectorModule::startRun() // Do your beginRun() stuff here
200{
201 m_run = m_evtMetaData->getRun();
202 m_exp = m_evtMetaData->getExperiment();
203 if (m_gainName.length()) {
205 m_gainMap = *gainMap;
206 } else {
208 }
209 getObjectPtr<TTree>("dbtree")->Fill();
210}
211
212void PXDPerformanceVariablesCollectorModule::collect() // Do your event() stuff here
213{
214 // Update booleans for event selection
215 m_selected4Eff = false;
217 if (particles4Eff.isValid() && particles4Eff->getListSize() == 1)
218 m_selected4Eff = true;
219 const Particle* vpho4eff = particles4Eff->getParticle(0);
220
221 collectDeltaIP(); // using ParticleList(m_PList4ResName)
222
224 // This list is used to extract info for both gain calibration and efficiency monitoring
225 // as it has the loosest selection and we can reuse the for loop.
226 // Do nothing if the list is empty
227 if (!particles.isValid() || particles->getListSize() < 1)
228 return;
229
230 for (auto const& particle : *particles) {
231 const Track* trackPtr = particle.getTrack();
232 auto mass = particle.getMass();
233 if (!trackPtr) return;
234 auto recoTrackPtr = trackPtr->getRelated<RecoTrack>("");
235 if (!recoTrackPtr) return;
236 auto pxdIntercepts = recoTrackPtr->getRelationsTo<PXDIntercept>("");
237 for (auto const& pxdIntercept : pxdIntercepts) {
238 TrackCluster_t trackCluster;
239 // Function setValues() also returns a recoTrack pointer
240 if (!trackCluster.setValues(pxdIntercept, "", "PXDClustersFromTracks", mass))
241 continue;
242
243 auto const& cluster = trackCluster.cluster;
244 auto const& intersection = trackCluster.intersection;
245 auto const& usedInTrack = trackCluster.usedInTrack;
246
247
248 // Collect info for efficiency study
249 if (m_selected4Eff && intersection.inside) {
250 // check if the particle is vpho4eff's daughter
251 for (auto const& daughter : vpho4eff->getDaughters()) {
252 if (particle.isCopyOf(daughter))
253 collectEfficiencyVariables(trackCluster);
254 } // end of vpho4eff daughter loop
255 }
256
257 // Collect info for gain calibration
258 // Check for valid cluster and intersection
259 if (!usedInTrack || !intersection.inside || intersection.chargeMPV <= 0)
260 continue;
261 // Apply cluster selection cuts
262 if (cluster.charge < m_minClusterCharge || cluster.size < m_minClusterSize || cluster.size > m_maxClusterSize)
263 continue;
264 collectGainVariables(trackCluster);
265
266 } // end of pxdIntercepts loop
267 } // end of particles loop
268}
269
271{
273 // Requiring vpho -> l+l-
274 if (!particles.isValid() || particles->getListSize() != 1)
275 return;
276
277 const Particle* mother = particles->getParticle(0);
278 // Use the vertex of the mother particle (vph0 -> l+l-) instead of IP for impact parameter correction
279 // This leads to very small sigma_z0
280 //TVector3 vertex = mother->getVertex();
281
282 // Use beam spot info for the interaction point
283 DBObjPtr<BeamSpot> beamSpotDB;
284 auto ip = ROOT::Math::XYZVector(beamSpotDB->getIPPosition());
285 auto vertex = ip;
286
287 const Particle* part0 = mother->getDaughter(0);
288 const Particle* part1 = mother->getDaughter(1);
289 //const TrackFitResult* tr0 = part0->getTrackFitResult();
290 //const TrackFitResult* tr1 = part1->getTrackFitResult();
293
295 auto d0p_0 = m_track_struct.d0p;
296 auto z0p_0 = m_track_struct.z0p;
298 auto d0p_1 = m_track_struct.d0p;
299 auto z0p_1 = m_track_struct.z0p;
300
301 m_deltaD0oSqrt2 = (d0p_0 + d0p_1) / sqrt(2.);
302 m_deltaZ0oSqrt2 = (z0p_0 - z0p_1) / sqrt(2.);
303
304 // Fill the tree of impact parameters
305 getObjectPtr<TTree>("tree_d0z0")->Fill();
306}
307
309{
310 auto cluster = trackCluster.cluster;
311 auto intersection = trackCluster.intersection;
312
313 // Compute variables from the cluster for gain estimation
314 m_signal = cluster.charge;
315 m_estimated = intersection.chargeMPV;
316
317 int uBin(-1), vBin(-1);
318 int binID = 0;
319 VxdID sensorID = PXD::getVxdIDFromPXDModuleID(cluster.pxdID);
320 auto layerNumber = sensorID.getLayerNumber();
321 auto ladderNumber = sensorID.getLadderNumber();
322 auto sensorNumber = sensorID.getSensorNumber();
323 try {
324 binID = getBinID(trackCluster, uBin, vBin, m_useClusterPosition);
325 } catch (...) {
326 // It happends very rarely (2 track clusters out of all in s-proc3 HLT bhabha skim.)
327 // One could check hit reconstruction bias and the calculation of uBin/vBin in PXDGainCalibrator
328 B2ERROR("On module " << sensorID
329 << ": Failed to get bin ID for the track cluster at (u,v) = (" << cluster.posU << "," << cluster.posV
330 << ") and (uBin, vBin) = (" << uBin << "," << vBin << ").");
331 return;
332 }
333
334// Increment the counter
335 getObjectPtr<TH1I>("PXDClusterCounter")->Fill(binID);
336
337 // Fill variabels into tree
338 if (m_fillChargeTree) {
339 string treename = str(format("tree_%1%_%2%_%3%_%4%_%5%") % layerNumber % ladderNumber % sensorNumber % uBin % vBin);
340 getObjectPtr<TTree>(treename)->Fill();
341 }
342
343 // Fill cluster charge ratio histogram if enabled
345 double ratio = m_signal / m_estimated;
346 auto axis = getObjectPtr<TH2F>("PXDClusterChargeRatio")->GetYaxis();
347 double maxY = axis->GetBinCenter(axis->GetNbins());
348 // Manipulate too large ratio for better estimation on median.
349 getObjectPtr<TH2F>("PXDClusterChargeRatio")->Fill(binID, TMath::Min(ratio, maxY));
350 }
351}
352
354{
355 auto const& cluster = trackCluster.cluster;
356 auto const& tPoint = trackCluster.intersection;
357 auto const& usedInTrack = trackCluster.usedInTrack;
358
359 auto phi = atan2(tPoint.y, tPoint.x);
360 auto z = tPoint.z;
361
362 VxdID sensorID = PXD::getVxdIDFromPXDModuleID(cluster.pxdID);
363 const PXD::SensorInfo& Info = dynamic_cast<const PXD::SensorInfo&>(VXD::GeoCache::getInstance().getSensorInfo(sensorID));
364 auto localPoint = Info.pointToLocal(ROOT::Math::XYZVector(tPoint.x, tPoint.y, tPoint.z), true);
365 auto uID = Info.getUCellID(localPoint.X());
366 auto vID = Info.getVCellID(localPoint.Y());
367 auto iSensor = VXD::GeoCache::getInstance().getGeoTools()->getPXDSensorIndex(sensorID);
368 auto uBin = PXD::PXDGainCalibrator::getInstance().getBinU(sensorID, uID, vID, m_nBinsU);
369 auto vBin = PXD::PXDGainCalibrator::getInstance().getBinV(sensorID, vID, m_nBinsV);
370 if (uBin >= m_nBinsU || vBin >= m_nBinsV)
371 return;
372 auto binID = iSensor * m_nBinsU * m_nBinsV + uBin * m_nBinsV + vBin;
373
374 // Filling counters
375 getObjectPtr<TH1I>("PXDTrackPointCounter")->Fill(binID);
376 if (usedInTrack)
377 getObjectPtr<TH1I>("PXDTrackClusterCounter")->Fill(binID);
378
379 if (isSelected(sensorID, uID, vID)) {
380 getObjectPtr<TH1I>("PXDSelTrackPointCounter")->Fill(binID);
381 if (usedInTrack)
382 getObjectPtr<TH1I>("PXDSelTrackClusterCounter")->Fill(binID);
383 }
384
385 // Filling 2D histograms
386 if (cluster.pxdID < 2000) {
387 //getObjectPtr<TEfficiency>("PXDLayer1Efficiency")->Fill(usedInTrack,phi,z);
388 getObjectPtr<TH2F>("hTotalHitsLayer1")->Fill(phi, z);
389 if (usedInTrack)
390 getObjectPtr<TH2F>("hPassedHitsLayer1")->Fill(phi, z);
391 } else {
392 //getObjectPtr<TEfficiency>("PXDLayer2Efficiency")->Fill(usedInTrack,phi,z);
393 getObjectPtr<TH2F>("hTotalHitsLayer2")->Fill(phi, z);
394 if (usedInTrack)
395 getObjectPtr<TH2F>("hPassedHitsLayer2")->Fill(phi, z);
396 }
397
398}
Calibration collector module base class.
static const ChargedStable pion
charged pion particle
Definition: Const.h:661
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
The payload class for PXD gain corrections.
Definition: PXDGainMapPar.h:43
PXDIntercept stores the U,V coordinates and uncertainties of the intersection of a track with an PXD ...
Definition: PXDIntercept.h:22
bool isSelected(const VxdID &sensorID, const int &uID, const int &vID)
Helper function to select a track point according to the status of the related pixels.
std::string m_PList4ResName
Name of the particle list for resolution study.
int m_nBinsV
Number of corrections per sensor along v side.
bool m_fillChargeRatioHistogram
Flag to fill cluster charge ratio (relative to expected MPV) histograms.
std::string m_PList4GainName
Name of the particle list for gain calibration.
PXD::Track_t m_track_struct
Track struct for holding required variables.
std::string m_gainName
Payload name for Gain to be read from DB.
void collectGainVariables(const PXD::TrackCluster_t &trackCluster)
Collect variables for gain calibration.
bool m_useClusterPosition
Flag to use cluster position rather than track point to group pixels into bins for gain calibration.
std::string m_PList4EffName
Name of the particle list for efficiency study.
int m_nBinsU
Number of corrections per sensor along u side.
StoreObjPtr< EventMetaData > m_evtMetaData
Required input EventMetaData.
void collectEfficiencyVariables(const PXD::TrackCluster_t &trackCluster)
Collect variables for efficiency monitoring.
PXDPerformanceVariablesCollectorModule()
Constructor: Sets the description, the properties and the parameters of the module.
bool m_selected4Eff
Event selection for efficiency validation.
int m_maskedDistance
Distance inside which no dead pixel or module border is allowed.
int getBinID(const PXD::TrackCluster_t &trackCluster, int &uBin, int &vBin, bool useCluster=false)
Helper function to get binID, uBin and vBin from Cluster_t struct The binning is derived from PXD::PX...
bool m_fillChargeTree
Flag to fill cluster charge and its estimated MPV in TTree.
void collectDeltaIP()
Collect info for impact parameter study on event level.
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
Class to store reconstructed particles.
Definition: Particle.h:75
const Track * getTrack() const
Returns the pointer to the Track object that was used to create this Particle (ParticleType == c_Trac...
Definition: Particle.cc:845
std::vector< Belle2::Particle * > getDaughters() const
Returns a vector of pointers to daughter particles.
Definition: Particle.cc:637
const Particle * getDaughter(unsigned i) const
Returns a pointer to the i-th daughter particle.
Definition: Particle.cc:631
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
T * getRelated(const std::string &name="", const std::string &namedRelation="") const
Get the object to or from which this object has a relation.
RelationVector< TO > getRelationsTo(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from this object to another store array.
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
bool isValid() const
Check whether the object was created.
Definition: StoreObjPtr.h:111
Values of the result of a track fit with a given particle hypothesis.
Class that bundles various TrackFitResults.
Definition: Track.h:25
const TrackFitResult * getTrackFitResultWithClosestMass(const Const::ChargedStable &requestedType) const
Return the track fit for a fit hypothesis with the closest mass.
Definition: Track.cc:104
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
int getPXDSensorIndex(VxdID sensorID) const
Return index of sensor in plots.
Definition: GeoTools.h:224
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 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.
float d0p
Corrected impact parameter in r-phi.
void setTrackVariables(const TrackFitResult *tfrPtr, const ROOT::Math::XYZVector &ip)
update track level variables from TrackFitResult
float z0p
Corrected impact parameter in z.
Struct to hold variables for track clusters.
bool usedInTrack
True if the cluster is used in tracking.
RecoTrack * setValues(const PXDIntercept &pxdIntercept, const std::string &recoTracksName="", const std::string &pxdTrackClustersName="PXDClustersFromTracks", const double &mass=Const::electronMass)
Update values from a PXDIntercept.
Cluster_t cluster
Cluster associated to the track.
TrackPoint_t intersection
The track-module intersection.