Belle II Software development
CDCDedxCorrectionModule.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 <cdc/modules/CDCDedxCorrection/CDCDedxCorrectionModule.h>
10
11#include <cdc/geometry/CDCGeometryPar.h>
12
13#include <TMath.h>
14
15using namespace Belle2;
16using namespace CDC;
17
18REG_MODULE(CDCDedxCorrection);
19
21{
22
23 setDescription("Apply hit level corrections to the dE/dx measurements.");
24 addParam("relativeCorrections", m_relative, "If true, apply corrections relative", true);
25 addParam("momentumCor", m_momCor, "Boolean to apply momentum correction", false);
26 addParam("momentumCorFromDB", m_useDBMomCor, "Boolean to apply momentum correction from DB", false);
27 addParam("scaleCor", m_scaleCor, "Boolean to apply scale correction", false);
28 addParam("cosineCor", m_cosineCor, "Boolean to apply cosine correction", false);
29 addParam("wireGain", m_wireGain, "Boolean to apply wire gains", false);
30 addParam("runGain", m_runGain, "Boolean to apply run gain", false);
31 addParam("timeGain", m_timeGain, "Boolean to apply time gain", false);
32 addParam("twoDCell", m_twoDCell, "Boolean to apply 2D correction", false);
33 addParam("oneDCell", m_oneDCell, "Boolean to apply 1D correction", false);
34 addParam("cosineEdge", m_cosineEdge, "Boolean to apply cosine edge correction", true);
35 addParam("nonlADC", m_nonlADC, "Boolean to apply non-linear adc correction", true);
36 addParam("removeLowest", m_removeLowest, "portion of events with low dE/dx", double(0.05));
37 addParam("removeHighest", m_removeHighest, "portion of events with high dE/dx", double(0.25));
38}
39
41
43{
44 // register in datastore
45 m_cdcDedxTracks.isRequired();
46
47 // make sure the calibration constants are reasonable
48 // run gains
49 if (m_DBRunGain->getRunGain() == 0) {
50 B2WARNING("Run gain is zero for this run");
51 }
52
53 // wire gains
54 for (unsigned int i = 0; i < c_nSenseWires; ++i) {
55 if (m_DBWireGains->getWireGain(i) == 0)
56 B2WARNING("Wire gain is zero for this wire: " << i);
57 }
58
59 // cosine correction (store the bin edges for extrapolation)
60 for (unsigned int il = 0; il < c_maxNSenseLayers; ++il) {
61 unsigned int ncosbins = m_DBCosineCor->getSize(il);
62
63 for (unsigned int i = 0; i < ncosbins; ++i) {
64 double gain = m_DBCosineCor->getMean(il, i);
65 if (gain == 0)
66 B2ERROR("Cosine gain is zero for this bin " << i);
67 }
68 }
69
70 // get the hadron correction parameters
71 if (!m_DBHadronCor || m_DBHadronCor->getSize() == 0) {
72 B2WARNING("No hadron correction parameters!");
73 for (int i = 0; i < 4; ++i)
74 m_hadronpars.push_back(0.0);
75 m_hadronpars.push_back(1.0);
76 } else m_hadronpars = m_DBHadronCor->getHadronPars();
77
78 int jwire = -1;
79 B2INFO("Creating CDCGeometryPar object");
81
82 for (unsigned int il = 0; il < c_maxNSenseLayers; il++) {
83 int activewires = 0;
84 m_lgainavg[il] = 0.0;
85
86 for (unsigned int iw = 0; iw < cdcgeo.nWiresInLayer(il); ++iw) {
87 jwire++;
88 if (m_DBWireGains->getWireGain(jwire) > 0.) {
89 //active wire only
90 m_lgainavg[il] += m_DBWireGains->getWireGain(jwire);
91 activewires++;
92 }
93 }
94 if (activewires > 0) m_lgainavg[il] /= activewires;
95 else m_lgainavg[il] = 1.0;
96 }
97}
98
99
101{
102
103 // **************************************************
104 // LOOP OVER EACH DEDX MEASUREMENT (TRACK LEVEL)
105 // **************************************************
106 for (auto& dedxTrack : m_cdcDedxTracks) {
107 if (dedxTrack.size() == 0) {
108 B2WARNING("No good hits on this track...");
109 continue;
110 }
111
112 // **************************************************
113 // LOOP OVER EACH DEDX MEASUREMENT (HIT LEVEL)
114 // **************************************************
115 // hit level
116 int nhits = dedxTrack.size();
117 double costh = dedxTrack.getCosTheta();
118 std::vector<double> newLayerHits;
119 double newLayerDe = 0, newLayerDx = 0;
120
121 if (costh < TMath::Cos(150.0 * TMath::DegToRad()))continue; //-0.866
122 if (costh > TMath::Cos(17.0 * TMath::DegToRad())) continue; //0.95
123
124
125 double injtime = dedxTrack.getInjectionTime();
126 double injring = dedxTrack.getInjectionRing();
127
128 for (int i = 0; i < nhits; ++i) {
129
130 //pay attention to deadwire or gain uses
131 //getADCount is already corrected w/ non linear ADC payload
132 //getADCbasecount is now incorrect ADC
133 int jadcbase = dedxTrack.getADCBaseCount(i);
134 double jLayer = dedxTrack.getHitLayer(i);
135 double jWire = dedxTrack.getWire(i);
136 double jNDocaRS = dedxTrack.getDocaRS(i) / dedxTrack.getCellHalfWidth(i);
137 double jEntaRS = dedxTrack.getEntaRS(i);
138 double jPath = dedxTrack.getPath(i);
139
140 double correction = dedxTrack.m_scale * dedxTrack.m_runGain * dedxTrack.m_timeGain * dedxTrack.getCosineCorrection(
141 i) * dedxTrack.m_cosEdgeCor *
142 dedxTrack.getTwoDCorrection(i) * dedxTrack.getOneDCorrection(i) * dedxTrack.getNonLADCCorrection(i);
143 if (dedxTrack.getWireGain(i) > 0)correction *= dedxTrack.getWireGain(i); //also keep dead wire
144
145 //Modify hit level dedx
146 double newhitdedx = (m_relative) ? 1 / correction : 1.0;
147 StandardCorrection(jadcbase, jLayer, jWire, jNDocaRS, jEntaRS, jPath, costh, injring, injtime, newhitdedx);
148 dedxTrack.setDedx(i, newhitdedx);
149
150 // do track level dedx and modify after loop over hits
151 // rel const -> upto 6 from calibrated GT and 2 are direct from dedx track (no rel cal for them now)
152 // abs const -> upto 6 from calibrated GT and 2 are direct from default GT
153 if (m_relative) {
154 //prewire gains need old tracks (old bad wire) and post need track new wg (new dead wire)
155 //get same base adc + rel correction factor
156 correction *= GetCorrection(jadcbase, jLayer, jWire, jNDocaRS, jEntaRS, costh, injring, injtime);
157 if (!m_DBWireGains && dedxTrack.getWireGain(i) == 0)correction = 0;
158 } else {
159 //get modified adc + abs correction factor
160 correction = GetCorrection(jadcbase, jLayer, jWire, jNDocaRS, jEntaRS, costh, injring, injtime);
161 }
162
163 // combine hits across layers
164 if (correction != 0) {
165 newLayerDe += jadcbase / correction;
166 newLayerDx += jPath;
167 }
168
169 if (i + 1 < nhits && dedxTrack.getHitLayer(i + 1) == jLayer) {
170 continue;
171 } else {
172 if (newLayerDx != 0)newLayerHits.push_back(newLayerDe / newLayerDx * std::sqrt(1 - costh * costh));
173 newLayerDe = 0;
174 newLayerDx = 0;
175 }
176 }
177
178 // recalculate the truncated means
179 calculateMeans(&(dedxTrack.m_dedxAvg),
180 &(dedxTrack.m_dedxAvgTruncatedNoSat),
181 &(dedxTrack.m_dedxAvgTruncatedErr),
182 newLayerHits);
183
184 dedxTrack.m_dedxAvgTruncated = dedxTrack.m_dedxAvgTruncatedNoSat;
185 HadronCorrection(costh, dedxTrack.m_dedxAvgTruncated);
186 } // end loop over tracks
187}
188
190{
191
192 B2INFO("CDCDedxCorrectionModule exiting...");
193}
194
196{
197 double gain = m_DBRunGain->getRunGain();
198 if (gain != 0) {
199 dedx = dedx / gain;
200 } else dedx = 0;
201}
202
203void CDCDedxCorrectionModule::TimeGainCorrection(double& dedx, double ring, double time) const
204{
205 double gain = m_DBInjectTime->getCorrection("mean", ring, time);
206 if (gain != 0) {
207 dedx = dedx / gain;
208 } else dedx = 0;
209}
210
211void CDCDedxCorrectionModule::WireGainCorrection(int wireID, double& dedx, int layer) const
212{
213 double gain = m_DBWireGains->getWireGain(wireID);
214 if (gain != 0) dedx = dedx / gain;
215 else {
216 //rel-abs method needs all wire for cal but w/ this method post calis (e.g.final RG)
217 //will also see all hitdedx but that is not an issue for track level calibration
218 if (m_relative)dedx = 0;
219 else {
220 if (m_lgainavg.at(layer) > 0)dedx = dedx / m_lgainavg.at(layer);
221 }
222 }
223}
224
225void CDCDedxCorrectionModule::TwoDCorrection(int layer, double doca, double enta, double& dedx) const
226{
227
228 double gain = (m_DB2DCell) ? m_DB2DCell->getMean(layer, doca, enta) : 1.0;
229 if (gain != 0) dedx = dedx / gain;
230 else dedx = 0;
231}
232
233
234void CDCDedxCorrectionModule::OneDCorrection(int layer, double enta, double& dedx) const
235{
236 double gain = (m_DB1DCell) ? m_DB1DCell->getMean(layer, enta) : 1.0;
237 if (gain != 0) dedx = dedx / gain;
238 else dedx = 0;
239}
240
241void CDCDedxCorrectionModule::CosineCorrection(unsigned int layer, double costh, double& dedx) const
242{
243 double coscor = m_DBCosineCor->getMean(layer, costh);
244 if (coscor != 0) dedx = dedx / coscor;
245 else dedx = 0;
246}
247
248void CDCDedxCorrectionModule::CosineEdgeCorrection(double costh, double& dedx) const
249{
250
251 double cosedgecor = m_DBCosEdgeCor->getMean(costh);
252 if (cosedgecor != 0) dedx = dedx / cosedgecor;
253 else dedx = 0;
254}
255
256
257void CDCDedxCorrectionModule::HadronCorrection(double costheta, double& dedx) const
258{
259 dedx = D2I(costheta, I2D(costheta, 1.00) / 1.00 * dedx);
260}
261
262void CDCDedxCorrectionModule::StandardCorrection(int adc, int layer, int wireID, double doca, double enta, double length,
263 double costheta, double ring, double time, double& dedx) const
264{
265
266 if (!m_relative && m_nonlADC)
267 adc = m_DBNonlADC->getCorrectedADC(adc, layer);
268
269 dedx *= adc * std::sqrt(1 - costheta * costheta) / length;
270
271 if (m_scaleCor) {
272 double scale = m_DBScaleFactor->getScaleFactor();
273 if (scale != 0) dedx = dedx / scale;
274 else dedx = 0;
275 }
276
277 if (m_runGain)
278 RunGainCorrection(dedx);
279
280 if (m_timeGain)
281 TimeGainCorrection(dedx, ring, time);
282
283 if (m_wireGain)
284 WireGainCorrection(wireID, dedx, layer);
285
286 if (m_cosineCor)
287 CosineCorrection(layer, costheta, dedx);
288
289 //only if constants are abs are for specific costh
290 if (!m_relative && m_cosineEdge && (costheta <= -0.850 || costheta >= 0.950))
291 CosineEdgeCorrection(costheta, dedx);
292
293 if (m_twoDCell)
294 TwoDCorrection(layer, doca, enta, dedx);
295
296 if (m_oneDCell)
297 OneDCorrection(layer, enta, dedx);
298}
299
300
301double CDCDedxCorrectionModule::GetCorrection(int& adc, int layer, int wireID, double doca, double enta, double costheta,
302 double ring, double time) const
303{
304 double correction = 1.0;
305 if (m_scaleCor) correction *= m_DBScaleFactor->getScaleFactor();
306 if (m_runGain) correction *= m_DBRunGain->getRunGain();
307 if (m_timeGain) correction *= m_DBInjectTime->getCorrection("mean", ring, time);
308 if (m_wireGain) correction *= m_DBWireGains->getWireGain(wireID);
309 if (m_twoDCell) correction *= m_DB2DCell->getMean(layer, doca, enta);
310 if (m_oneDCell) correction *= m_DB1DCell->getMean(layer, enta);
311 if (m_cosineCor) correction *= m_DBCosineCor->getMean(layer, costheta);
312
313 //these last two are only for abs constant
314 if (!m_relative) {
315 if (m_nonlADC)adc = m_DBNonlADC->getCorrectedADC(adc, layer);
316 if (m_cosineEdge && (costheta <= -0.850 || costheta >= 0.950)) {
317 correction *= m_DBCosEdgeCor->getMean(costheta);
318 }
319 }
320 return correction;
321}
322
323
324double CDCDedxCorrectionModule::D2I(const double cosTheta, const double D) const
325{
326 double absCosTheta = fabs(cosTheta);
327 double projection = pow(absCosTheta, m_hadronpars[3]) + m_hadronpars[2];
328 if (projection == 0) {
329 B2WARNING("Something wrong with dE/dx hadron constants!");
330 return D;
331 }
332
333 double chargeDensity = D / projection;
334 double numerator = 1 + m_hadronpars[0] * chargeDensity;
335 double denominator = 1 + m_hadronpars[1] * chargeDensity;
336
337 if (denominator == 0) {
338 B2WARNING("Something wrong with dE/dx hadron constants!");
339 return D;
340 }
341
342 double I = D * m_hadronpars[4] * numerator / denominator;
343 return I;
344}
345
346double CDCDedxCorrectionModule::I2D(const double cosTheta, const double I) const
347{
348 double absCosTheta = fabs(cosTheta);
349 double projection = pow(absCosTheta, m_hadronpars[3]) + m_hadronpars[2];
350
351 if (projection == 0 || m_hadronpars[4] == 0) {
352 B2WARNING("Something wrong with dE/dx hadron constants!");
353 return I;
354 }
355
356 double a = m_hadronpars[0] / projection;
357 double b = 1 - m_hadronpars[1] / projection * (I / m_hadronpars[4]);
358 double c = -1.0 * I / m_hadronpars[4];
359
360 if (b == 0 && a == 0) {
361 B2WARNING("both a and b coefficiants for hadron correction are 0");
362 return I;
363 }
364
365 double discr = b * b - 4.0 * a * c;
366 if (discr < 0) {
367 B2WARNING("negative discriminant; return uncorrectecd value");
368 return I;
369 }
370
371 double D = (a != 0) ? (-b + sqrt(discr)) / (2.0 * a) : -c / b;
372 if (D < 0) {
373 B2WARNING("D is less 0! will try another solution");
374 D = (a != 0) ? (-b - sqrt(discr)) / (2.0 * a) : -c / b;
375 if (D < 0) {
376 B2WARNING("D is still less 0! just return uncorrectecd value");
377 return I;
378 }
379 }
380
381 return D;
382}
383
384
385void CDCDedxCorrectionModule::calculateMeans(double* mean, double* truncatedMean, double* truncatedMeanErr,
386 const std::vector<double>& dedx) const
387{
388 // Calculate the truncated average by skipping the lowest & highest
389 // events in the array of dE/dx values
390 std::vector<double> sortedDedx = dedx;
391 std::sort(sortedDedx.begin(), sortedDedx.end());
392 sortedDedx.erase(std::remove(sortedDedx.begin(), sortedDedx.end(), 0), sortedDedx.end());
393 sortedDedx.shrink_to_fit();
394
395 double truncatedMeanTmp = 0.0;
396 double meanTmp = 0.0;
397 double sumOfSquares = 0.0;
398 int numValuesTrunc = 0;
399 const int numDedx = sortedDedx.size();
400
401 // add a factor of 0.51 here to make sure we are rounding appropriately...
402 const int lowEdgeTrunc = int(numDedx * m_removeLowest + 0.51);
403 const int highEdgeTrunc = int(numDedx * (1 - m_removeHighest) + 0.51);
404 for (int i = 0; i < numDedx; i++) {
405 meanTmp += sortedDedx[i];
406 if (i >= lowEdgeTrunc and i < highEdgeTrunc) {
407 truncatedMeanTmp += sortedDedx[i];
408 sumOfSquares += sortedDedx[i] * sortedDedx[i];
409 numValuesTrunc++;
410 }
411 }
412
413 if (numDedx != 0) {
414 meanTmp /= numDedx;
415 }
416 if (numValuesTrunc != 0) {
417 truncatedMeanTmp /= numValuesTrunc;
418 } else {
419 truncatedMeanTmp = meanTmp;
420 }
421
422 *mean = meanTmp;
423 *truncatedMean = truncatedMeanTmp;
424
425 if (numValuesTrunc > 1) {
426 *truncatedMeanErr = sqrt(sumOfSquares / double(numValuesTrunc) - truncatedMeanTmp * truncatedMeanTmp) / double(
427 numValuesTrunc - 1);
428 } else {
429 *truncatedMeanErr = 0;
430 }
431}
bool m_twoDCell
boolean to apply 2D correction
double D2I(const double cosTheta, const double D) const
Saturation correction: convert the measured ionization (D) to actual ionization (I)
bool m_runGain
boolean to apply run gains
std::vector< double > m_hadronpars
hadron saturation parameters
bool m_cosineEdge
boolean to apply cosine edge
DBObjPtr< CDCDedxRunGain > m_DBRunGain
Run gain DB object.
bool m_useDBMomCor
boolean to apply momentum correction from DB
void TwoDCorrection(int layer, double doca, double enta, double &dedx) const
Perform a 2D correction.
bool m_momCor
boolean to apply momentum correction
virtual void initialize() override
Initialize the Module.
bool m_scaleCor
boolean to apply scale factor
DBObjPtr< CDCDedxHadronCor > m_DBHadronCor
hadron saturation parameters
virtual void event() override
This method is called for each event.
double m_removeHighest
upper bound for truncated mean
DBObjPtr< CDCDedxCosineEdge > m_DBCosEdgeCor
cosine edge calibration
double I2D(const double cosTheta, const double I) const
Saturation correction: convert the actual ionization (I) to measured ionization (D)
std::array< double, 56 > m_lgainavg
average calibration factor for the layer
void OneDCorrection(int layer, double enta, double &dedx) const
Perform a wire gain correction.
DBObjPtr< CDCDedxADCNonLinearity > m_DBNonlADC
hadron saturation non linearity
DBObjPtr< CDCDedx1DCell > m_DB1DCell
1D correction DB object
virtual void terminate() override
End of the event processing.
bool m_oneDCell
boolean to apply 1D correction
void HadronCorrection(double costheta, double &dedx) const
Perform a hadron saturation correction.
void calculateMeans(double *mean, double *truncatedMean, double *truncatedMeanErr, const std::vector< double > &dedx) const
Recalculate the dE/dx mean values after corrections.
void StandardCorrection(int adc, int layer, int wireID, double doca, double enta, double length, double costheta, double ring, double time, double &dedx) const
Perform a standard set of corrections.
DBObjPtr< CDCDedx2DCell > m_DB2DCell
2D correction DB object
bool m_nonlADC
boolean to apply non linear ADC
bool m_relative
boolean to apply relative or absolute correction
void CosineCorrection(unsigned int superlayer, double costheta, double &dedx) const
Perform the cosine correction.
DBObjPtr< CDCDedxCosineCor > m_DBCosineCor
Electron saturation correction DB object.
StoreArray< CDCDedxTrack > m_cdcDedxTracks
Store array: CDCDedxTrack.
bool m_timeGain
boolean to apply injection time gains
void CosineEdgeCorrection(double costh, double &dedx) const
Perform the cosine edge correction.
DBObjPtr< CDCDedxWireGain > m_DBWireGains
Wire gain DB object.
void WireGainCorrection(int wireID, double &dedx, int layer) const
Perform a wire gain correction.
bool m_cosineCor
boolean to apply cosine correction
void RunGainCorrection(double &dedx) const
Perform a run gain correction.
double GetCorrection(int &adc, int layer, int wireID, double doca, double enta, double costheta, double ring, double time) const
Get the standard set of corrections.
double m_removeLowest
lower bound for truncated mean
DBObjPtr< CDCDedxInjectionTime > m_DBInjectTime
time gain/reso DB object
DBObjPtr< CDCDedxScaleFactor > m_DBScaleFactor
Scale factor to make electrons ~1.
void TimeGainCorrection(double &dedx, double ring, double time) const
Perform a injection time gain correction.
bool m_wireGain
boolean to apply wire gains
The Class for CDC Geometry Parameters.
unsigned nWiresInLayer(int layerId) const
Returns wire numbers in a layer.
static CDCGeometryPar & Instance(const CDCGeometry *=nullptr)
Static method to get a reference to the CDCGeometryPar instance.
void setDescription(const std::string &description)
Sets the description of the module.
Definition Module.cc:214
Module()
Constructor.
Definition Module.cc:30
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:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
double sqrt(double a)
sqrt for double
Definition beamHelpers.h:28
Abstract base class for different kinds of events.