Belle II Software development
PXDUtilities.h
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#pragma once
10
11#include <framework/gearbox/Unit.h>
12#include <framework/gearbox/Const.h>
13#include <tracking/dataobjects/RecoTrack.h>
14#include <genfit/MeasuredStateOnPlane.h>
15//#include <limits>
16#include <math.h>
17
18namespace Belle2 {
24 namespace PXD {
25
29 const int Z_Si = 14;
30 const double A_Si = 28.085;
31 const double rho_Si = 2.3290 * Unit::g_cm3;
40 inline double xiBeta2_L(const int Z = Z_Si,
41 const double A = A_Si,
42 const double rho = rho_Si,
43 const int z = 1)
44 {
45 const double K = 0.307075 * Unit::MeV * pow(Unit::cm, 2);
46 return K / 2 * Z / A * z * z * rho;
47 }
48
55 inline double hbarWp(const int Z = Z_Si,
56 const double A = A_Si,
57 const double rho = rho_Si)
58 {
59 return std::sqrt(rho * Z / A) * 28.816 * Unit::eV;
60 }
61
68 inline double getDeltaP(const double mom, const double length, const double mass = Const::electronMass)
69 {
70 double betaGamma = mom / mass;
71 if (betaGamma <= 100) return 0.0; // requirement of the formula.
72 double beta2 = 1. / (1. + 1. / betaGamma / betaGamma);
73 double xi = xiBeta2_L() * length / beta2;
74 return xi * log(2 * mass * xi / pow(hbarWp(), 2) + 0.2);
75 }
76
78 inline unsigned short getPXDModuleID(const VxdID& sensorID)
79 {
80 return sensorID.getLayerNumber() * 1000 +
81 sensorID.getLadderNumber() * 10 +
82 sensorID.getSensorNumber();
83 }
84
86 inline VxdID getVxdIDFromPXDModuleID(const unsigned short& id)
87 {
88 return VxdID(id / 1000, (id % 1000) / 10, id % 10);
89 }
90
92 typedef genfit::MeasuredStateOnPlane TrackState;
93
100 std::shared_ptr<TrackState> getTrackStateOnModule(const VXD::SensorInfoBase& pxdSensorInfo,
101 RecoTrack& recoTrack, double lambda = 0.0);
102
109 bool isCloseToBorder(int u, int v, int checkDistance);
117 bool isDefectivePixelClose(int u, int v, int checkDistance, const VxdID& moduleID);
118
125 inline bool isClusterAtUEdge(VxdID id, unsigned int umin, unsigned int umax)
126 {
127 unsigned int uedgemax = Belle2::VXD::GeoCache::getInstance().getSensorInfo(id).getUCells();
128 return (umin == 0 || umax == (uedgemax - 1));
129 }
136 inline bool isClusterAtVEdge(VxdID id, unsigned int vmin, unsigned int vmax)
137 {
138 unsigned int vedgemax = Belle2::VXD::GeoCache::getInstance().getSensorInfo(id).getVCells();
139 return ((id.getSensorNumber() == 1 && vmax == (vedgemax - 1))
140 || (id.getSensorNumber() == 2 && vmin == 0));
141 }
148 inline bool isClusterAtLadderJoint(VxdID id, unsigned int vmin, unsigned int vmax)
149 {
150 unsigned int vedgemax = Belle2::VXD::GeoCache::getInstance().getSensorInfo(id).getVCells();
151 return ((id.getSensorNumber() == 2 && vmax == (vedgemax - 1))
152 || (id.getSensorNumber() == 1 && vmin == 0));
153 }
154
161 inline unsigned short getBinU(VxdID id, unsigned int uid, unsigned int vid, unsigned short nBinsU)
162 {
163 unsigned int drainsPerBin = 4 * Belle2::VXD::GeoCache::getInstance().getSensorInfo(id).getUCells() / nBinsU;
164 return (uid * 4 + vid % 4) / drainsPerBin;
165 }
171 inline unsigned short getBinV(VxdID id, unsigned int vid, unsigned short nBinsV)
172 {
173 unsigned int rowsPerBin = Belle2::VXD::GeoCache::getInstance().getSensorInfo(id).getVCells() / nBinsV;
174 return vid / rowsPerBin;
175 }
176
177 } // end namespace PXD
179} // end namespace Belle2
#define K(x)
macro autogenerated by FFTW
static const double electronMass
electron mass
Definition: Const.h:685
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
static const double eV
[electronvolt]
Definition: Unit.h:112
static const double MeV
[megaelectronvolt]
Definition: Unit.h:114
static const double cm
Standard units with the value = 1.
Definition: Unit.h:47
static const double g_cm3
Practical units with the value set at 1.
Definition: Unit.h:60
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
Base class to provide Sensor Information for PXD and SVD.
int getVCells() const
Return number of pixel/strips in v direction.
int getUCells() const
Return number of pixel/strips in u direction.
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
VxdID getVxdIDFromPXDModuleID(const unsigned short &id)
Helper function to get VxdID from DHE id like module iid.
Definition: PXDUtilities.h:86
bool isCloseToBorder(int u, int v, int checkDistance)
Helper function to check if a pixel is close to the border.
Definition: PXDUtilities.cc:51
unsigned short getBinU(VxdID id, unsigned int uid, unsigned int vid, unsigned short nBinsU)
Function to return a bin number for equal sized binning in U.
Definition: PXDUtilities.h:161
std::shared_ptr< TrackState > getTrackStateOnModule(const VXD::SensorInfoBase &pxdSensorInfo, RecoTrack &recoTrack, double lambda=0.0)
Helper function to get a track state on a module.
Definition: PXDUtilities.cc:21
const double A_Si
Atomic mass of silicon in g mol^-1.
Definition: PXDUtilities.h:30
const double rho_Si
Silicon density in g cm^-3.
Definition: PXDUtilities.h:31
unsigned short getPXDModuleID(const VxdID &sensorID)
Helper function to get DHE id like module id from VxdID.
Definition: PXDUtilities.h:78
double xiBeta2_L(const int Z=Z_Si, const double A=A_Si, const double rho=rho_Si, const int z=1)
xi = (K/2)*(Z/A)*z*z*(rho*L)/beta2 in MeV
Definition: PXDUtilities.h:40
bool isClusterAtUEdge(VxdID id, unsigned int umin, unsigned int umax)
Helper function to check if one of the end pixels are at the edge of the sensor.
Definition: PXDUtilities.h:125
genfit::MeasuredStateOnPlane TrackState
Typedef TrackState (genfit::MeasuredStateOnPlane)
Definition: PXDUtilities.h:92
bool isClusterAtVEdge(VxdID id, unsigned int vmin, unsigned int vmax)
Helper function to check if one of the end pixels are at the edge of the sensor.
Definition: PXDUtilities.h:136
double hbarWp(const int Z=Z_Si, const double A=A_Si, const double rho=rho_Si)
hbarWp = sqrt(rho*Z/A)*28.816 in eV
Definition: PXDUtilities.h:55
bool isClusterAtLadderJoint(VxdID id, unsigned int vmin, unsigned int vmax)
Helper function to check if one of the end pixels are at the ladder joint.
Definition: PXDUtilities.h:148
double getDeltaP(const double mom, const double length, const double mass=Const::electronMass)
helper function to estimate the most probable energy loss for a given track length.
Definition: PXDUtilities.h:68
const int Z_Si
Const and Const expressions Only valid when g_mol is the default unit.
Definition: PXDUtilities.h:29
bool isDefectivePixelClose(int u, int v, int checkDistance, const VxdID &moduleID)
Helper function to check if a defective (hot/dead) pixel is close.
Definition: PXDUtilities.cc:61
unsigned short getBinV(VxdID id, unsigned int vid, unsigned short nBinsV)
Function to return a bin number for equal sized binning in V.
Definition: PXDUtilities.h:171
Abstract base class for different kinds of events.