Belle II Software development
SVDCoGOnlyPosition Class Reference

Derived Class representing the SVD cluster position computed with the CoGOnly algorithm. More...

#include <SVDCoGOnlyPosition.h>

Inheritance diagram for SVDCoGOnlyPosition:
SVDClusterPosition

Public Member Functions

void computeClusterPosition (Belle2::SVD::RawCluster &rawCluster, double &position, double &positionError) override
 computes the cluster position and position error with the CoG algorithm
 
virtual ~SVDCoGOnlyPosition ()
 virtual destructor
 
void applyCoGPosition (const Belle2::SVD::RawCluster &rawCluster, double &position, double &positionError)
 CoG Position Algorithm.
 
void applyAHTPosition (const Belle2::SVD::RawCluster &rawCluster, double &position, double &positionError)
 AHT Position Algorithm.
 
void reconstructStrips (Belle2::SVD::RawCluster &rawCluster)
 reconstruct strips
 
void applyUnfolding (Belle2::SVD::RawCluster &rawCluster)
 Apply cluster charges unfolding.
 
void set_stripChargeAlgo (const std::string &user_stripChargeAlgo)
 set which algorithm to use for strip charge in cluster position reconstruction
 
void set_stripTimeAlgo (const std::string &user_stripTimeAlgo)
 set which algorithm to use for strip time in cluster position reconstruction, 'dontdo' will skip it
 

Protected Member Functions

double getSumOfStripCharges (const Belle2::SVD::RawCluster &rawCluster)
 helper, returns the sum of the strip charges
 
double getClusterNoise (const Belle2::SVD::RawCluster &rawCluster)
 helper, returns the sum in quadrature of the strip noise
 
double getAverageStripNoise (const Belle2::SVD::RawCluster &rawCluster)
 helper, returns the average strip noise
 

Protected Attributes

SVDCoGOnlyPositionError m_CoGOnlyErr
 CoGOnly Position Error.
 
SVDCoGOnlyErrorScaleFactors m_CoGOnlyCal
 Scaling Factors for the CoGOnly algorithm.
 
SVDOldDefaultErrorScaleFactors m_OldDefaultCal
 Scaling Factors for the OldDefault algorithm.
 
SVDClustering m_ClusterCal
 SVD clustering parameters.
 
SVDNoiseCalibrations m_NoiseCal
 Noise calibrations for the position error.
 

Private Attributes

std::string m_stripChargeAlgo
 algorithm used to reconstruct strip charge for cluster position
 
std::string m_stripTimeAlgo
 algorithm used to reconstruct strip time for cluster position
 

Detailed Description

Derived Class representing the SVD cluster position computed with the CoGOnly algorithm.

Definition at line 22 of file SVDCoGOnlyPosition.h.

Constructor & Destructor Documentation

◆ ~SVDCoGOnlyPosition()

virtual ~SVDCoGOnlyPosition ( )
inlinevirtual

virtual destructor

Definition at line 35 of file SVDCoGOnlyPosition.h.

35{};

Member Function Documentation

◆ applyAHTPosition()

void applyAHTPosition ( const Belle2::SVD::RawCluster rawCluster,
double &  position,
double &  positionError 
)
inherited

AHT Position Algorithm.

Definition at line 71 of file SVDClusterPosition.cc.

72 {
73
74 // NOTE:
75 // tail and head are the two strips at the edge of the cluster
76
77 const VXD::GeoCache& geo = VXD::GeoCache::getInstance();
78 const VXD::SensorInfoBase& info = geo.getSensorInfo(rawCluster.getSensorID());
79 double pitch = rawCluster.isUSide() ? info.getUPitch() : info.getVPitch();
80
81 //take the strips in the rawCluster
82 std::vector<Belle2::SVD::StripInRawCluster> strips = rawCluster.getStripsInRawCluster();
83
84 //informations about the head strip
85 int headStripCellID = strips.at(strips.size() - 1).cellID;
86 double headStripCharge = strips.at(strips.size() - 1).charge;
87 //informations about the tail strip
88 int tailStripCellID = strips.at(0).cellID;
89 double tailStripCharge = strips.at(0).charge;
90
91 // average strip charge of the center of the cluster
92
93 double centreCharge = (getSumOfStripCharges(rawCluster) - tailStripCharge - headStripCharge) / (strips.size() - 2);
94
95 tailStripCharge = (tailStripCharge < centreCharge) ? tailStripCharge : centreCharge;
96 headStripCharge = (headStripCharge < centreCharge) ? headStripCharge : centreCharge;
97 double tailPos = rawCluster.isUSide() ? info.getUCellPosition(tailStripCellID) : info.getVCellPosition(tailStripCellID);
98 double headPos = rawCluster.isUSide() ? info.getUCellPosition(headStripCellID) : info.getVCellPosition(headStripCellID);
99 position = 0.5 * (tailPos + headPos + (headStripCharge - tailStripCharge) / centreCharge * pitch);
100
101 //now compute position error
102 double cutAdjacent = m_ClusterCal.getMinAdjSNR(rawCluster.getSensorID(), rawCluster.isUSide());
103 double sn = centreCharge / cutAdjacent / getClusterNoise(rawCluster);
104
105 // Rough estimates of Landau noise
106 double landauHead = tailStripCharge / centreCharge;
107 double landauTail = headStripCharge / centreCharge;
108 positionError = 0.5 * pitch * sqrt(1.0 / sn / sn +
109 0.5 * landauHead * landauHead +
110 0.5 * landauTail * landauTail);
111
112 }
double getMinAdjSNR(const Belle2::VxdID &sensorID, const bool &isU) const
Return the minimum SNR for the adjacent.
Definition: SVDClustering.h:75
const std::vector< StripInRawCluster > getStripsInRawCluster() const
Definition: RawCluster.h:110
bool isUSide() const
Definition: RawCluster.h:85
VxdID getSensorID() const
Definition: RawCluster.h:80
SVDClustering m_ClusterCal
SVD clustering parameters.
double getSumOfStripCharges(const Belle2::SVD::RawCluster &rawCluster)
helper, returns the sum of the strip charges
double getClusterNoise(const Belle2::SVD::RawCluster &rawCluster)
helper, returns the sum in quadrature of the strip noise
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28

◆ applyCoGPosition()

void applyCoGPosition ( const Belle2::SVD::RawCluster rawCluster,
double &  position,
double &  positionError 
)
inherited

CoG Position Algorithm.

Definition at line 34 of file SVDClusterPosition.cc.

35 {
36
37 const VXD::GeoCache& geo = VXD::GeoCache::getInstance();
38 const VXD::SensorInfoBase& info = geo.getSensorInfo(rawCluster.getSensorID());
39
40 position = 0;
41 double charge = 0;
42
43 //take the strips in the rawCluster
44 std::vector<Belle2::SVD::StripInRawCluster> strips = rawCluster.getStripsInRawCluster();
45
46 for (auto aStrip : strips) {
47
48 double stripPos = rawCluster.isUSide() ? info.getUCellPosition(aStrip.cellID) : info.getVCellPosition(aStrip.cellID);
49
50 //getting the charge of the strip
51 double stripCharge = aStrip.charge;
52
53 position += stripPos * stripCharge;
54 charge += stripCharge;
55 }
56
57 position /= charge;
58
59 //now compute position error
60 positionError = 0;
61 double pitch = rawCluster.isUSide() ? info.getUPitch() : info.getVPitch();
62 double sumStripCharge = getSumOfStripCharges(rawCluster);
63
64 positionError = m_CoGOnlyErr.getPositionError(rawCluster.getSensorID(), rawCluster.isUSide(), 0,
65 sumStripCharge / getClusterNoise(rawCluster), rawCluster.getSize(), pitch);
66
67 }
double getPositionError(const Belle2::VxdID &sensorID, const bool &isU, const unsigned short &strip, const double &clSNR, const int &clSize, const float &pitch) const
Return the position error.
SVDCoGOnlyPositionError m_CoGOnlyErr
CoGOnly Position Error.
double charge(int pdgCode)
Returns electric charge of a particle with given pdg code.
Definition: EvtPDLUtil.cc:44

◆ applyUnfolding()

void applyUnfolding ( Belle2::SVD::RawCluster rawCluster)
inherited

Apply cluster charges unfolding.

Definition at line 240 of file SVDClusterPosition.cc.

241 {
242
243
244 std::vector<Belle2::SVD::StripInRawCluster> strips = rawCluster.getStripsInRawCluster();
245 double unfoldingCoefficient = m_ClusterCal.getUnfoldingCoeff(rawCluster.getSensorID(), rawCluster.isUSide());
246 unsigned int Size = strips.size();
247 double threshold = 0;
248 Eigen::VectorXd Charges(Size);
249 Eigen::MatrixXd Couplings(Size, Size);
250 // Unfolding Matrix
251 for (unsigned int i = 0; i < Size; ++i) {
252 for (unsigned int j = 0; j < Size; ++j) {
253 if (i == j) {Couplings(i, j) = 1 - 2 * unfoldingCoefficient;}
254 else if (j == i + 1) {Couplings(i, j) = unfoldingCoefficient;}
255 else if (j == i - 1) {Couplings(i, j) = unfoldingCoefficient;}
256 else {Couplings(i, j) = 0;}
257 }
258
259 Belle2::SVD::StripInRawCluster strip = strips.at(i);
260
261 RawCluster tmp(rawCluster.getSensorID(), rawCluster.isUSide(), 0, 0);
262 if (tmp.add(rawCluster.getSensorID(), rawCluster.isUSide(), strip)) {
263 //Fill a vector with strip charges
264 Charges(i) = strip.charge;
265 }
266
267 }
268 //Apply the unfolding
269 Charges = Couplings.inverse() * Charges;
270 for (unsigned i = 0; i < Size; i++) {
271 if (Charges(i) < threshold) {Charges(i) = 0;}
272 rawCluster.setStripCharge(i, Charges(i));
273 }
274
275 }
double getUnfoldingCoeff(const Belle2::VxdID &sensorID, const bool &isU) const
Return the unfolding coefficient for the strip charge.
Definition: SVDClustering.h:97
void setStripCharge(int index, double charge)
set the strip charge
Definition: RawCluster.h:127
structure containing the relevant informations of each strip of the raw cluster
Definition: RawCluster.h:20
double charge
strip charge
Definition: RawCluster.h:26

◆ computeClusterPosition()

void computeClusterPosition ( Belle2::SVD::RawCluster rawCluster,
double &  position,
double &  positionError 
)
overridevirtual

computes the cluster position and position error with the CoG algorithm

Implements SVDClusterPosition.

Definition at line 23 of file SVDCoGOnlyPosition.cc.

24 {
25 reconstructStrips(rawCluster);
26 applyUnfolding(rawCluster);
27 applyCoGPosition(rawCluster, position, positionError);
28
29 //apply scale factors for the position errors
30 positionError = m_CoGOnlyCal.getCorrectedClusterPositionError(rawCluster.getSensorID(), rawCluster.isUSide(), rawCluster.getSize(),
31 positionError);
32
33 }
double getCorrectedClusterPositionError(const Belle2::VxdID &sensorID, const bool &isU, const int &size, const double &raw_error) const
Return the corrected cluster position error.
SVDCoGOnlyErrorScaleFactors m_CoGOnlyCal
Scaling Factors for the CoGOnly algorithm.
void applyUnfolding(Belle2::SVD::RawCluster &rawCluster)
Apply cluster charges unfolding.
void applyCoGPosition(const Belle2::SVD::RawCluster &rawCluster, double &position, double &positionError)
CoG Position Algorithm.
void reconstructStrips(Belle2::SVD::RawCluster &rawCluster)
reconstruct strips

◆ getAverageStripNoise()

double getAverageStripNoise ( const Belle2::SVD::RawCluster rawCluster)
protectedinherited

helper, returns the average strip noise

Definition at line 148 of file SVDClusterPosition.cc.

149 {
150
151 double averageNoise = 0;
152
153 //take the strips in the rawCluster
154 std::vector<Belle2::SVD::StripInRawCluster> strips = rawCluster.getStripsInRawCluster();
155
156 // compute the average strip noise
157 for (auto aStrip : strips) {
158
159 float averageNoiseInElectrons = m_NoiseCal.getNoiseInElectrons(rawCluster.getSensorID(), rawCluster.isUSide(), aStrip.cellID);
160 averageNoise += averageNoiseInElectrons;
161 }
162 return averageNoise / strips.size();
163 }
float getNoiseInElectrons(const VxdID &sensorID, const bool &isU, const unsigned short &strip) const
This method provides the correct noise conversion into electrons, taking into account that the noise ...
SVDNoiseCalibrations m_NoiseCal
Noise calibrations for the position error.

◆ getClusterNoise()

double getClusterNoise ( const Belle2::SVD::RawCluster rawCluster)
protectedinherited

helper, returns the sum in quadrature of the strip noise

Definition at line 131 of file SVDClusterPosition.cc.

132 {
133
134 double clusterNoise = 0;
135
136 //take the strips in the rawCluster
137 std::vector<Belle2::SVD::StripInRawCluster> strips = rawCluster.getStripsInRawCluster();
138
139 // compute the cluster noise as sum in quadrature of strip noise
140 for (auto aStrip : strips) {
141
142 float averageNoiseInElectrons = m_NoiseCal.getNoiseInElectrons(rawCluster.getSensorID(), rawCluster.isUSide(), aStrip.cellID);
143 clusterNoise += averageNoiseInElectrons * averageNoiseInElectrons;
144 }
145 return sqrt(clusterNoise);
146 }

◆ getSumOfStripCharges()

double getSumOfStripCharges ( const Belle2::SVD::RawCluster rawCluster)
protectedinherited

helper, returns the sum of the strip charges

Definition at line 114 of file SVDClusterPosition.cc.

115 {
116
117 double sumStripCharge = 0;
118
119 //take the strips in the rawCluster
120 std::vector<Belle2::SVD::StripInRawCluster> strips = rawCluster.getStripsInRawCluster();
121
122 // compute the sum of strip charges
123 for (auto aStrip : strips) {
124
125 double stripCharge = aStrip.charge;
126 sumStripCharge += stripCharge;
127 }
128 return sumStripCharge;
129 }

◆ reconstructStrips()

void reconstructStrips ( Belle2::SVD::RawCluster rawCluster)
inherited

reconstruct strips

Definition at line 165 of file SVDClusterPosition.cc.

166 {
167
168
169 std::vector<Belle2::SVD::StripInRawCluster> strips = rawCluster.getStripsInRawCluster();
170
171 //loop on strips
172 for (int i = 0; i < (int)strips.size(); i++) {
173
174 Belle2::SVD::StripInRawCluster strip = strips.at(i);
175
176 RawCluster tmp(rawCluster.getSensorID(), rawCluster.isUSide(), 0, 0);
177 if (tmp.add(rawCluster.getSensorID(), rawCluster.isUSide(), strip)) {
178
179 // time computation
180 if (m_stripTimeAlgo.compare("dontdo") != 0) {
181
182 double time = 0;
183 double timeError = 0;
184 int firstFrame = 0;
185
186 if (m_stripTimeAlgo == "ELS3") {
187 SVDELS3Time ct;
188 ct.computeClusterTime(tmp, time, timeError, firstFrame);
189
190 } else if (m_stripTimeAlgo == "CoG3") {
191 SVDCoG3Time ct;
192 ct.computeClusterTime(tmp, time, timeError, firstFrame);
193 }
194 rawCluster.setStripTime(i, time);
195 }
196
197 // charge computation
198 // may be not needed in case the cluster position is computed using APV samples
199 // and not using reconstructed strips
200 if (m_stripChargeAlgo.compare("dontdo") != 0) {
201 double charge = 0;
202 double SNR;
203 double seedCharge;
204
205 if (m_stripChargeAlgo == "ELS3") {
206 // ELS3 can return non-sense values for off-time or low-charge strips)
207 // if returned charge is negative or more than 30% different than MaxSample, we use MaxSample
208 // without notice to the user!
209
210 SVDELS3Charge cc;
211 cc.computeClusterCharge(tmp, charge, SNR, seedCharge);
212
213 double maxSample_charge = 0;
214 SVDMaxSampleCharge maxSample_cc;
215 maxSample_cc.computeClusterCharge(tmp, maxSample_charge, SNR, seedCharge);
216
217 if ((abs(charge - maxSample_charge) / maxSample_charge > 0.3) || charge < 0)
218 rawCluster.setStripCharge(i, maxSample_charge);
219 else
220 rawCluster.setStripCharge(i, charge);
221
222 } else if (m_stripChargeAlgo == "SumSamples") {
223 SVDSumSamplesCharge cc;
224 cc.computeClusterCharge(tmp, charge, SNR, seedCharge);
225 rawCluster.setStripCharge(i, charge);
226 } else {
227 // MaxSample is used when the algorithm is not recognized
228 SVDMaxSampleCharge cc;
229 cc.computeClusterCharge(tmp, charge, SNR, seedCharge);
230 rawCluster.setStripCharge(i, charge);
231 }
232 }
233 } else
234 B2ERROR("this should not happen...");
235
236 }
237
238 }
void setStripTime(int index, double time)
set the strip time
Definition: RawCluster.h:133
std::string m_stripTimeAlgo
algorithm used to reconstruct strip time for cluster position
std::string m_stripChargeAlgo
algorithm used to reconstruct strip charge for cluster position

◆ set_stripChargeAlgo()

void set_stripChargeAlgo ( const std::string &  user_stripChargeAlgo)
inlineinherited

set which algorithm to use for strip charge in cluster position reconstruction

Definition at line 59 of file SVDClusterPosition.h.

59{m_stripChargeAlgo = user_stripChargeAlgo;}

◆ set_stripTimeAlgo()

void set_stripTimeAlgo ( const std::string &  user_stripTimeAlgo)
inlineinherited

set which algorithm to use for strip time in cluster position reconstruction, 'dontdo' will skip it

Definition at line 62 of file SVDClusterPosition.h.

62{m_stripTimeAlgo = user_stripTimeAlgo;}

Member Data Documentation

◆ m_ClusterCal

SVDClustering m_ClusterCal
protectedinherited

SVD clustering parameters.

Definition at line 77 of file SVDClusterPosition.h.

◆ m_CoGOnlyCal

SVDCoGOnlyErrorScaleFactors m_CoGOnlyCal
protectedinherited

Scaling Factors for the CoGOnly algorithm.

Definition at line 75 of file SVDClusterPosition.h.

◆ m_CoGOnlyErr

SVDCoGOnlyPositionError m_CoGOnlyErr
protectedinherited

CoGOnly Position Error.

Definition at line 74 of file SVDClusterPosition.h.

◆ m_NoiseCal

SVDNoiseCalibrations m_NoiseCal
protectedinherited

Noise calibrations for the position error.

Definition at line 78 of file SVDClusterPosition.h.

◆ m_OldDefaultCal

SVDOldDefaultErrorScaleFactors m_OldDefaultCal
protectedinherited

Scaling Factors for the OldDefault algorithm.

Definition at line 76 of file SVDClusterPosition.h.

◆ m_stripChargeAlgo

std::string m_stripChargeAlgo
privateinherited

algorithm used to reconstruct strip charge for cluster position

Definition at line 82 of file SVDClusterPosition.h.

◆ m_stripTimeAlgo

std::string m_stripTimeAlgo
privateinherited

algorithm used to reconstruct strip time for cluster position

Definition at line 83 of file SVDClusterPosition.h.


The documentation for this class was generated from the following files: