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