Belle II Software  release-06-02-00
SVDClusterPosition.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 <framework/logging/Logger.h>
10 #include <svd/reconstruction/SVDClusterPosition.h>
11 #include <svd/reconstruction/SVDCoG3Time.h>
12 #include <svd/reconstruction/SVDELS3Time.h>
13 #include <svd/reconstruction/SVDELS3Charge.h>
14 #include <svd/reconstruction/SVDMaxSampleCharge.h>
15 #include <svd/reconstruction/SVDSumSamplesCharge.h>
16 #include <svd/reconstruction/SVDMaxSumAlgorithm.h>
17 #include <vxd/geometry/GeoCache.h>
18 #include <svd/geometry/SensorInfo.h>
19 
20 #include <TMath.h>
21 
22 using namespace std;
23 
24 namespace Belle2 {
30  namespace SVD {
31 
32 
33  void SVDClusterPosition::applyCoGPosition(const Belle2::SVD::RawCluster& rawCluster, double& position, double& positionError)
34  {
35 
36  const VXD::GeoCache& geo = VXD::GeoCache::getInstance();
37  const VXD::SensorInfoBase& info = geo.getSensorInfo(rawCluster.getSensorID());
38 
39  position = 0;
40  double charge = 0;
41 
42  //take the strips in the rawCluster
43  std::vector<Belle2::SVD::StripInRawCluster> strips = rawCluster.getStripsInRawCluster();
44 
45  for (auto aStrip : strips) {
46 
47  double stripPos = rawCluster.isUSide() ? info.getUCellPosition(aStrip.cellID) : info.getVCellPosition(aStrip.cellID);
48 
49  //getting the charge of the strip
50  double stripCharge = aStrip.charge;
51 
52  position += stripPos * stripCharge;
53  charge += stripCharge;
54  }
55 
56  position /= charge;
57 
58  //now compute position error
59  positionError = 0;
60  double pitch = rawCluster.isUSide() ? info.getUPitch() : info.getVPitch();
61  double sumStripCharge = getSumOfStripCharges(rawCluster);
62 
63  positionError = m_CoGOnlyErr.getPositionError(rawCluster.getSensorID(), rawCluster.isUSide(), 0,
64  sumStripCharge / getClusterNoise(rawCluster), rawCluster.getSize(), pitch);
65 
66  }
67 
68 
69 
70  void SVDClusterPosition::applyAHTPosition(const Belle2::SVD::RawCluster& rawCluster, double& position, double& positionError)
71  {
72 
73  // NOTE:
74  // tail and head are the two strips at the edge of the cluster
75 
76  const VXD::GeoCache& geo = VXD::GeoCache::getInstance();
77  const VXD::SensorInfoBase& info = geo.getSensorInfo(rawCluster.getSensorID());
78  double pitch = rawCluster.isUSide() ? info.getUPitch() : info.getVPitch();
79 
80  //take the strips in the rawCluster
81  std::vector<Belle2::SVD::StripInRawCluster> strips = rawCluster.getStripsInRawCluster();
82 
83  //informations about the head strip
84  int headStripCellID = strips.at(strips.size() - 1).cellID;
85  double headStripCharge = strips.at(strips.size() - 1).charge;
86  //informations about the tail strip
87  int tailStripCellID = strips.at(0).cellID;
88  double tailStripCharge = strips.at(0).charge;
89 
90  // average strip charge of the center of the cluster
91 
92  double centreCharge = (getSumOfStripCharges(rawCluster) - tailStripCharge - headStripCharge) / (strips.size() - 2);
93 
94  tailStripCharge = (tailStripCharge < centreCharge) ? tailStripCharge : centreCharge;
95  headStripCharge = (headStripCharge < centreCharge) ? headStripCharge : centreCharge;
96  double tailPos = rawCluster.isUSide() ? info.getUCellPosition(tailStripCellID) : info.getVCellPosition(tailStripCellID);
97  double headPos = rawCluster.isUSide() ? info.getUCellPosition(headStripCellID) : info.getVCellPosition(headStripCellID);
98  position = 0.5 * (tailPos + headPos + (headStripCharge - tailStripCharge) / centreCharge * pitch);
99 
100  //now compute position error
101  double cutAdjacent = m_ClusterCal.getMinAdjSNR(rawCluster.getSensorID(), rawCluster.isUSide());
102  double sn = centreCharge / cutAdjacent / getClusterNoise(rawCluster);
103 
104  // Rough estimates of Landau noise
105  double landauHead = tailStripCharge / centreCharge;
106  double landauTail = headStripCharge / centreCharge;
107  positionError = 0.5 * pitch * sqrt(1.0 / sn / sn +
108  0.5 * landauHead * landauHead +
109  0.5 * landauTail * landauTail);
110 
111  }
112 
113  double SVDClusterPosition::getSumOfStripCharges(const Belle2::SVD::RawCluster& rawCluster)
114  {
115 
116  double sumStripCharge = 0;
117 
118  //take the strips in the rawCluster
119  std::vector<Belle2::SVD::StripInRawCluster> strips = rawCluster.getStripsInRawCluster();
120 
121  // compute the sum of strip charges
122  for (auto aStrip : strips) {
123 
124  double stripCharge = aStrip.charge;
125  sumStripCharge += stripCharge;
126  }
127  return sumStripCharge;
128  }
129 
130  double SVDClusterPosition::getClusterNoise(const Belle2::SVD::RawCluster& rawCluster)
131  {
132 
133  double clusterNoise = 0;
134 
135  //take the strips in the rawCluster
136  std::vector<Belle2::SVD::StripInRawCluster> strips = rawCluster.getStripsInRawCluster();
137 
138  // compute the cluster noise as sum in quadrature of strip noise
139  for (auto aStrip : strips) {
140 
141  float averageNoiseInElectrons = m_NoiseCal.getNoiseInElectrons(rawCluster.getSensorID(), rawCluster.isUSide(), aStrip.cellID);
142  clusterNoise += averageNoiseInElectrons * averageNoiseInElectrons;
143  }
144  return sqrt(clusterNoise);
145  }
146 
147  double SVDClusterPosition::getAverageStripNoise(const Belle2::SVD::RawCluster& rawCluster)
148  {
149 
150  double averageNoise = 0;
151 
152  //take the strips in the rawCluster
153  std::vector<Belle2::SVD::StripInRawCluster> strips = rawCluster.getStripsInRawCluster();
154 
155  // compute the average strip noise
156  for (auto aStrip : strips) {
157 
158  float averageNoiseInElectrons = m_NoiseCal.getNoiseInElectrons(rawCluster.getSensorID(), rawCluster.isUSide(), aStrip.cellID);
159  averageNoise += averageNoiseInElectrons;
160  }
161  return averageNoise / strips.size();
162  }
163 
164  void SVDClusterPosition::reconstructStrips(Belle2::SVD::RawCluster& rawCluster)
165  {
166 
167 
168  std::vector<Belle2::SVD::StripInRawCluster> strips = rawCluster.getStripsInRawCluster();
169 
170  //loop on strips
171  for (int i = 0; i < (int)strips.size(); i++) {
172 
173  Belle2::SVD::StripInRawCluster strip = strips.at(i);
174 
175  RawCluster tmp(rawCluster.getSensorID(), rawCluster.isUSide(), 0, 0);
176  if (tmp.add(rawCluster.getSensorID(), rawCluster.isUSide(), strip)) {
177 
178  // time computation
179  if (m_stripTimeAlgo.compare("dontdo") != 0) {
180 
181  double time = 0;
182  double timeError = 0;
183  int firstFrame = 0;
184 
185  if (m_stripTimeAlgo == "ELS3") {
186  SVDELS3Time ct;
187  ct.computeClusterTime(tmp, time, timeError, firstFrame);
188 
189  } else if (m_stripTimeAlgo == "CoG3") {
190  SVDCoG3Time ct;
191  ct.computeClusterTime(tmp, time, timeError, firstFrame);
192  }
193  rawCluster.setStripTime(i, time);
194  }
195 
196  // charge computation
197  // may be not needed in case the cluster position is computed using APV samples
198  // and not using reconstructed strips
199  if (m_stripChargeAlgo.compare("dontdo") != 0) {
200  double charge = 0;
201  double SNR;
202  double seedCharge;
203 
204  if (m_stripChargeAlgo == "ELS3") {
205  // ELS3 can return non-sense values for off-time or low-charge strips)
206  // if returned charge is negative or more than 30% different than MaxSample, we use MaxSample
207  // without notice to the user!
208 
209  SVDELS3Charge cc;
210  cc.computeClusterCharge(tmp, charge, SNR, seedCharge);
211 
212  double maxSample_charge = 0;
213  SVDMaxSampleCharge maxSample_cc;
214  maxSample_cc.computeClusterCharge(tmp, maxSample_charge, SNR, seedCharge);
215 
216  if ((abs(charge - maxSample_charge) / maxSample_charge > 0.3) || charge < 0)
217  rawCluster.setStripCharge(i, maxSample_charge);
218  else
219  rawCluster.setStripCharge(i, charge);
220 
221  } else if (m_stripChargeAlgo == "SumSamples") {
223  cc.computeClusterCharge(tmp, charge, SNR, seedCharge);
224  rawCluster.setStripCharge(i, charge);
225  } else {
226  // MaxSample is used when the algorithm is not recognized
228  cc.computeClusterCharge(tmp, charge, SNR, seedCharge);
229  rawCluster.setStripCharge(i, charge);
230  }
231  }
232  } else
233  B2ERROR("this should not happen...");
234 
235  }
236 
237  }
238 
239  } //SVD namespace
241 } //Belle2 namespace
Class representing a raw cluster candidate during clustering of the SVD.
Definition: RawCluster.h:33
void setStripCharge(int index, double charge)
set the strip charge
Definition: RawCluster.h:127
bool isUSide() const
Definition: RawCluster.h:85
VxdID getSensorID() const
Definition: RawCluster.h:80
const std::vector< StripInRawCluster > getStripsInRawCluster() const
Definition: RawCluster.h:110
void setStripTime(int index, double time)
set the strip time
Definition: RawCluster.h:133
Derived Class representing the SVD cluster time computed with the CoG3 algorithm.
Definition: SVDCoG3Time.h:22
void computeClusterTime(Belle2::SVD::RawCluster &rawCluster, double &time, double &timeError, int &firstFrame) override
computes the cluster time, timeError and FirstFrame with the CoG3 algorithm
Definition: SVDCoG3Time.cc:21
Derived Class representing the SVD cluster charge computed with the ELS3 algorithm.
Definition: SVDELS3Charge.h:26
Derived Class representing the SVD cluster time computed with the ELS3 algorithm.
Definition: SVDELS3Time.h:24
void computeClusterTime(Belle2::SVD::RawCluster &rawCluster, double &time, double &timeError, int &firstFrame) override
computes the cluster time, timeError and FirstFrame with the ELS3 algorithm
Definition: SVDELS3Time.cc:22
Derived Class representing the SVD cluster charge computed summing the max sample of each strip.
void computeClusterCharge(Belle2::SVD::RawCluster &rawCluster, double &charge, double &SNR, double &seedCharge) override
compute the cluster charge, charge error and SNR with MaxSample
Derived Class representing the SVD cluster charge computed summing the samples of each strip.
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:39
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:66
Base class to provide Sensor Information for PXD and SVD.
Abstract base class for different kinds of events.
structure containing the relevant informations of each strip of the raw cluster
Definition: RawCluster.h:20