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