Belle II Software  release-08-01-10
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 <TVector3.h>
16 //#include <limits>
17 #include <math.h>
18 
19 namespace Belle2 {
25  namespace PXD {
26 
30  const int Z_Si = 14;
31  const double A_Si = 28.085;
32  const double rho_Si = 2.3290 * Unit::g_cm3;
41  inline double xiBeta2_L(const int Z = Z_Si,
42  const double A = A_Si,
43  const double rho = rho_Si,
44  const int z = 1)
45  {
46  const double K = 0.307075 * Unit::MeV * pow(Unit::cm, 2);
47  return K / 2 * Z / A * z * z * rho;
48  }
49 
56  inline double hbarWp(const int Z = Z_Si,
57  const double A = A_Si,
58  const double rho = rho_Si)
59  {
60  return std::sqrt(rho * Z / A) * 28.816 * Unit::eV;
61  }
62 
69  inline double getDeltaP(const double mom, const double length, const double mass = Const::electronMass)
70  {
71  double betaGamma = mom / mass;
72  if (betaGamma <= 100) return 0.0; // requirement of the formula.
73  double beta2 = 1. / (1. + 1. / betaGamma / betaGamma);
74  double xi = xiBeta2_L() * length / beta2;
75  return xi * log(2 * mass * xi / pow(hbarWp(), 2) + 0.2);
76  }
77 
79  inline unsigned short getPXDModuleID(const VxdID& sensorID)
80  {
81  return sensorID.getLayerNumber() * 1000 +
82  sensorID.getLadderNumber() * 10 +
83  sensorID.getSensorNumber();
84  }
85 
87  inline VxdID getVxdIDFromPXDModuleID(const unsigned short& id)
88  {
89  return VxdID(id / 1000, (id % 1000) / 10, id % 10);
90  }
91 
94 
101  std::shared_ptr<TrackState> getTrackStateOnModule(const VXD::SensorInfoBase& pxdSensorInfo,
102  RecoTrack& recoTrack, double lambda = 0.0);
103 
110  bool isCloseToBorder(int u, int v, int checkDistance);
118  bool isDefectivePixelClose(int u, int v, int checkDistance, const VxdID& moduleID);
119 
126  inline bool isClusterAtUEdge(VxdID id, unsigned int umin, unsigned int umax)
127  {
128  unsigned int uedgemax = Belle2::VXD::GeoCache::getInstance().get(id).getUCells();
129  return (umin == 0 || umax == (uedgemax - 1));
130  }
137  inline bool isClusterAtVEdge(VxdID id, unsigned int vmin, unsigned int vmax)
138  {
139  unsigned int vedgemax = Belle2::VXD::GeoCache::getInstance().get(id).getVCells();
140  return ((id.getSensorNumber() == 1 && vmax == (vedgemax - 1))
141  || (id.getSensorNumber() == 2 && vmin == 0));
142  }
149  inline bool isClusterAtLadderJoint(VxdID id, unsigned int vmin, unsigned int vmax)
150  {
151  unsigned int vedgemax = Belle2::VXD::GeoCache::getInstance().get(id).getVCells();
152  return ((id.getSensorNumber() == 2 && vmax == (vedgemax - 1))
153  || (id.getSensorNumber() == 1 && vmin == 0));
154  }
155 
162  inline unsigned short getBinU(VxdID id, unsigned int uid, unsigned int vid, unsigned short nBinsU)
163  {
164  unsigned int drainsPerBin = 4 * Belle2::VXD::GeoCache::getInstance().get(id).getUCells() / nBinsU;
165  return (uid * 4 + vid % 4) / drainsPerBin;
166  }
172  inline unsigned short getBinV(VxdID id, unsigned int vid, unsigned short nBinsV)
173  {
174  unsigned int rowsPerBin = Belle2::VXD::GeoCache::getInstance().get(id).getVCells() / nBinsV;
175  return vid / rowsPerBin;
176  }
177 
178  } // end namespace PXD
180 } // end namespace Belle2
#define K(x)
macro autogenerated by FFTW
static const double electronMass
electron mass
Definition: Const.h:676
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
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
static const SensorInfoBase & get(Belle2::VxdID id)
Return a reference to the SensorInfo of a given SensorID.
Definition: GeoCache.h:139
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
#StateOnPlane with additional covariance matrix.
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
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:162
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:56
bool isCloseToBorder(int u, int v, int checkDistance)
Helper function to check if a pixel is close to the border.
Definition: PXDUtilities.cc:51
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:31
const double rho_Si
Silicon density in g cm^-3.
Definition: PXDUtilities.h:32
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:172
unsigned short getPXDModuleID(const VxdID &sensorID)
Helper function to get DHE id like module id from VxdID.
Definition: PXDUtilities.h:79
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:137
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:41
genfit::MeasuredStateOnPlane TrackState
Typedef TrackState (genfit::MeasuredStateOnPlane)
Definition: PXDUtilities.h:93
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:69
const int Z_Si
Const and Const expressions Only valid when g_mol is the default unit.
Definition: PXDUtilities.h:30
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:126
VxdID getVxdIDFromPXDModuleID(const unsigned short &id)
Helper function to get VxdID from DHE id like module iid.
Definition: PXDUtilities.h:87
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
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:149
Abstract base class for different kinds of events.