Belle II Software  release-05-02-19
ECLDigitizerModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2010 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Poyuan Chen *
7  * Guglielmo De Nardo (denardo@na.infn.it) *
8  * Alex Bobrov *
9  * *
10  * This software is provided "as is" without any warranty. *
11  **************************************************************************/
12 //This module
13 #include <ecl/modules/eclDigitizer/ECLDigitizerModule.h>
14 
15 // Root
16 #include <TRandom.h>
17 #include <TFile.h>
18 #include <TTree.h>
19 
20 //Framework
21 #include <framework/logging/Logger.h>
22 #include <framework/utilities/FileSystem.h>
23 #include <framework/gearbox/Unit.h>
24 
25 //ECL
26 #include <ecl/digitization/algorithms.h>
27 #include <ecl/digitization/shaperdsp.h>
28 #include <ecl/digitization/ECLCompress.h>
29 #include <ecl/geometry/ECLGeometryPar.h>
30 #include <ecl/dbobjects/ECLCrystalCalib.h>
31 #include <ecl/dbobjects/ECLWaveformData.h>
32 #include <ecl/dataobjects/ECLHit.h>
33 #include <ecl/dataobjects/ECLSimHit.h>
34 #include <ecl/dataobjects/ECLDigit.h>
35 #include <ecl/dataobjects/ECLDsp.h>
36 #include <ecl/dataobjects/ECLDspWithExtraMCInfo.h>
37 #include <ecl/dataobjects/ECLTrig.h>
38 #include <ecl/dataobjects/ECLWaveforms.h>
39 #include <ecl/utility/ECLDspUtilities.h>
40 
41 using namespace std;
42 using namespace Belle2;
43 using namespace ECL;
44 
45 //-----------------------------------------------------------------
46 // Register the Module
47 //-----------------------------------------------------------------
48 REG_MODULE(ECLDigitizer)
49 
50 //-----------------------------------------------------------------
51 // Implementation
52 //-----------------------------------------------------------------
53 
55 {
56  //Set module properties
57  setDescription("Creates ECLDigiHits from ECLHits.");
58  setPropertyFlags(c_ParallelProcessingCertified);
59  addParam("Background", m_background, "Flag to use the Digitizer configuration with backgrounds (default: false)", false);
60  addParam("Calibration", m_calibration, "Flag to use the Digitizer for Waveform fit Covariance Matrix calibration (default: false)",
61  false);
62  addParam("DiodeDeposition", m_inter,
63  "Flag to take into account energy deposition in photodiodes; Default diode is sensitive detector (default: true)", true);
64  addParam("WaveformMaker", m_waveformMaker, "Flag to produce background waveform digits (default: false)", false);
65  addParam("CompressionAlgorithm", m_compAlgo, "Waveform compression algorithm (default: 0u)", 0u);
66  addParam("eclWaveformsName", m_eclWaveformsName, "Name of the output/input collection (digitized waveforms)", string(""));
67  addParam("HadronPulseShapes", m_HadronPulseShape, "Flag to include hadron component in pulse shape construction (default: true)",
68  true);
69  addParam("ADCThreshold", m_ADCThreshold, "ADC threshold for waveform fits (default: 25)", 25);
70  addParam("WaveformThresholdOverride", m_WaveformThresholdOverride,
71  "If gt 0 value is applied to all crystals for waveform saving threshold. If lt 0 dbobject is used. (GeV)", -1.0);
72  addParam("StoreDspWithExtraMCInfo", m_storeDspWithExtraMCInfo,
73  "Flag to store Dsp with extra MC information in addition to normal Dsp (default: false)", false);
74  addParam("DspWithExtraMCInfoThreshold", m_DspWithExtraMCInfoThreshold,
75  "Threshold above with to store Dsp with extra MC information [GeV]",
76  0.02);
77  addParam("DspDataTest", m_dspDataTest,
78  "Use DSP coefficients from the database for the processing. This "
79  "will significantly reduce performance so this option is to be "
80  "used for testing only.", false);
81 
82 }
83 
84 ECLDigitizerModule::~ECLDigitizerModule()
85 {
86 }
87 
88 void ECLDigitizerModule::initialize()
89 {
90  m_eclDsps.registerInDataStore();
91  if (m_storeDspWithExtraMCInfo) m_eclDspsWithExtraMCInfo.registerInDataStore();
92  m_eclDigits.registerInDataStore();
93  m_eclTrigs.registerInDataStore();
94 
95  if (m_HadronPulseShape) {
96  B2DEBUG(20,
97  "Hadron pulse shapes for ECL simulations are enabled. Pulse shape simulations use techniques validated with test beam data documented in: S. Longo and J. M. Roney 2018 JINST 13 P03018");
98  } else {
99  B2DEBUG(20, "Hadron pulse shapes for ECL simulations are disabled.");
100  }
101 
102  if (m_inter) {
103  B2DEBUG(20, "Diode-crossing pulse shapes for ECL simulations are enabled.");
104  } else {
105  B2DEBUG(20, "Diode-crossing pulse shapes for ECL simulations are disabled.");
106  }
107 
108  m_eclDiodeHits.registerInDataStore("ECLDiodeHits");
109 
110  m_eclDsps.registerRelationTo(m_eclDigits);
111  if (m_storeDspWithExtraMCInfo)
112  m_eclDspsWithExtraMCInfo.registerRelationTo(m_eclDigits);
113  m_eclDigits.registerRelationTo(m_eclHits);
114  if (m_waveformMaker)
115  m_eclWaveforms.registerInDataStore(m_eclWaveformsName);
116 
117  readDSPDB();
118 
119  m_adc.resize(EclConfiguration::m_nch);
120 
121  EclConfiguration::get().setBackground(m_background);
122 }
123 
124 void ECLDigitizerModule::beginRun()
125 {
126  const EclConfiguration& ec = EclConfiguration::get();
128  Ael("ECLCrystalElectronics"),
129  Aen("ECLCrystalEnergy"),
130  Tel("ECLCrystalElectronicsTime"),
131  Ten("ECLCrystalTimeOffset"),
132  Tct("ECLCrateTimeOffset"),
133  Tmc("ECLMCTimeOffset"),
134  Awave("ECL_FPGA_StoreWaveform");
135  double ns_per_tick = 1.0 / (4.0 * ec.m_rf) * 1e3;// ~0.49126819903043308239 ns/tick
136 
137  calibration_t def = {1, 0};
138  m_calib.assign(8736, def);
139 
140  if (Ael) for (int i = 0; i < 8736; i++) m_calib[i].ascale /= Ael->getCalibVector()[i];
141  if (Aen) for (int i = 0; i < 8736; i++) m_calib[i].ascale /= Aen->getCalibVector()[i] * 20000.0;
142 
143  if (Tel) for (int i = 0; i < 8736; i++) m_calib[i].tshift += Tel->getCalibVector()[i] * ns_per_tick;
144  if (Ten) for (int i = 0; i < 8736; i++) m_calib[i].tshift += Ten->getCalibVector()[i] * ns_per_tick;
145  if (Tct) for (int i = 0; i < 8736; i++) m_calib[i].tshift += Tct->getCalibVector()[i] * ns_per_tick;
146  if (Tmc) for (int i = 0; i < 8736; i++) m_calib[i].tshift += Tmc->getCalibVector()[i] * ns_per_tick;
147 
148  m_Awave.assign(8736, -1);
149  if (m_WaveformThresholdOverride < 0) {
150  if (Awave) for (int i = 0; i < 8736; i++) m_Awave[i] = Awave->getCalibVector()[i];
151  } else {
152  //If m_WaveformThresholdOverride > 0 override ECL_FPGA_StoreWaveform;
153  for (int i = 0; i < 8736; i++) m_Awave[i] = m_WaveformThresholdOverride * m_calib[i].ascale; // convert GeV to ADC
154  }
155 
156  if (m_HadronPulseShape) callbackHadronSignalShapes();
157 
158 }
159 
160 // interface to C shape fitting function function
161 void ECLDigitizerModule::shapeFitterWrapper(const int j, const int* FitA, const int ttrig,
162  int& m_lar, int& m_ltr, int& m_lq, int& m_chi) const
163 {
164  const crystallinks_t& t = m_tbl[j]; //lookup table [0,8735]
165  const fitparams_t& r = m_fitparams[t.ifunc];
166 
167  ECLShapeFit result;
168 
169  if (!m_dspDataTest) {
170  short int* id = (short int*)m_idn[t.idn].id;
171 
172  int A0 = (int) * (id + 0) - 128;
173  int Askip = (int) * (id + 1) - 128;
174 
175  int Ahard = (int) * (id + 2);
176  int k_a = (int) * ((unsigned char*)id + 26);
177  int k_b = (int) * ((unsigned char*)id + 27);
178  int k_c = (int) * ((unsigned char*)id + 28);
179  int k_16 = (int) * ((unsigned char*)id + 29);
180  int k1_chi = (int) * ((unsigned char*)id + 24);
181  int k2_chi = (int) * ((unsigned char*)id + 25);
182 
183  int chi_thres = (int) * (id + 15);
184 
185  int trg_time = ttrig;
186 
187  result = lftda_((int*)r.f, (int*)r.f1, (int*)r.fg41, (int*)r.fg43,
188  (int*)r.fg31, (int*)r.fg32, (int*)r.fg33, (int*)FitA,
189  trg_time, A0, Ahard, Askip, k_a, k_b, k_c, k_16, k1_chi,
190  k2_chi, chi_thres);
191  } else {
192  std::vector<int> adc(31);
193  for (int i = 0; i < 31; i++) adc[i] = FitA[i];
194  // NOTE: ttrig_packed is any even number in set {8*Z, 8*Z + 2, 8*Z + 4}
195  // where Z is an integer number in 0..23 range
196  // ttrig = ttrig_packed - 2 * (ttrig_packed / 8);
197  // ttrig = 6*Z + q
198  int ttrig_packed = ttrig / 6 * 8 + ttrig % 6;
199  result = ECLDspUtilities::shapeFitter(j + 1, adc, ttrig_packed);
200  }
201 
202  m_lar = result.amp;
203  m_ltr = result.time;
204  m_lq = result.quality;
205  m_chi = result.chi2;
206 
207  //== Set precision of chi^2 to be the same as in the raw data.
208  int discarded_bits = 0;
209  if ((m_chi & 0x7800000) != 0) {
210  m_chi = 0x7800000;
211  } else if ((m_chi & 0x0600000) != 0) {
212  discarded_bits = 14;
213  } else if ((m_chi & 0x0180000) != 0) {
214  discarded_bits = 12;
215  } else if ((m_chi & 0x0060000) != 0) {
216  discarded_bits = 10;
217  } else if ((m_chi & 0x0018000) != 0) {
218  discarded_bits = 8;
219  } else if ((m_chi & 0x0006000) != 0) {
220  discarded_bits = 6;
221  } else if ((m_chi & 0x0001800) != 0) {
222  discarded_bits = 4;
223  } else if ((m_chi & 0x0000600) != 0) {
224  discarded_bits = 2;
225  }
226  if (discarded_bits > 0) {
227  m_chi >>= discarded_bits;
228  m_chi <<= discarded_bits;
229  }
230 }
231 
232 int ECLDigitizerModule::shapeSignals()
233 {
234  const EclConfiguration& ec = EclConfiguration::get();
235  ECLGeometryPar* eclp = ECLGeometryPar::Instance();
236 
237  const double trgtick = ec.s_clock / ec.m_rf / ec.m_ntrg; // trigger tick
238  const int DeltaT = gRandom->Uniform(0, double(ec.m_ntrg) / 2.); // trigger decision can come in any time ???
239  const double timeInt = 2.* (double) DeltaT * trgtick;
240  const int ttrig = 2 * DeltaT;
241  const double timeOffset = timeInt - ec.s_clock / (2 * ec.m_rf);
242 
243  const auto eclTrig = m_eclTrigs.appendNew();
244  eclTrig->setTimeTrig(timeInt); //t0 (us)= (1520 - m_ltr)*24.*
245 
246  // clear the storage for the event
247  memset(m_adc.data(), 0, m_adc.size()*sizeof(adccounts_t));
248 
249  const double E2GeV = 1 / Unit::GeV; // convert Geant energy units to GeV
250  const double T2us = 1 / Unit::us; // convert Geant time units to microseconds
251 
252  // emulate response for ECL hits after ADC measurements
253  for (const auto& hit : m_eclSimHits) {
254  int j = hit.getCellId() - 1; //0~8735
255  double hitE = hit.getEnergyDep() * m_calib[j].ascale * E2GeV;
256  double hitTimeAve = (hit.getFlightTime() + m_calib[j].tshift + eclp->time2sensor(j, hit.getPosition())) * T2us;
257 
258  m_adc[j].energyConversion = m_calib[j].ascale * E2GeV * 20000;
259  m_adc[j].flighttime += hit.getFlightTime() * hit.getEnergyDep(); //true time weighted by energy
260  m_adc[j].timeshift += m_calib[j].tshift * hit.getEnergyDep();
261  m_adc[j].timetosensor += eclp->time2sensor(j, hit.getPosition()) * hit.getEnergyDep();
262  m_adc[j].totalHadronDep += hit.getHadronEnergyDep(); // true deposited energy hadron component (in GeV)
263  m_adc[j].totalDep += hit.getEnergyDep(); //true deposited energy (in GeV)
264 
265  if (m_HadronPulseShape == true) {
266  double hitHadronE = hit.getHadronEnergyDep() * m_calib[j].ascale * E2GeV;
267  m_adc[j].AddHit(hitE - hitHadronE, hitTimeAve + timeOffset, m_ss_HadronShapeSimulations[0]);//Gamma Component
268  m_adc[j].AddHit(hitHadronE, hitTimeAve + timeOffset, m_ss_HadronShapeSimulations[1]); //Hadron Component
269  } else {
270  m_adc[j].AddHit(hitE, hitTimeAve + timeOffset, m_ss[m_tbl[j].iss]);
271  }
272  }
273 
274  // add only background hits
275  for (const auto& hit : m_eclHits) {
276  if (hit.getBackgroundTag() == BackgroundMetaData::bg_none) continue;
277  int j = hit.getCellId() - 1; //0~8735
278  double hitE = hit.getEnergyDep() * m_calib[j].ascale * E2GeV;
279  double hitTimeAve = (hit.getTimeAve() + m_calib[j].tshift) * T2us;
280  m_adc[j].AddHit(hitE, hitTimeAve + timeOffset, m_ss[m_tbl[j].iss]);
281  }
282 
283  // internuclear counter effect -- charged particle crosses diode and produces signal
284  if (m_inter) {
285  // ionisation energy in Si is I = 3.6x10^-6 MeV for electron-hole pair
286  // 5000 pairs in the diode per 1 MeV deposited in the crystal attached to the diode
287  // conversion factor to get equvalent energy deposition in the crystal to sum up it with deposition in crystal
288  const double diodeEdep2crystalEdep = E2GeV * (1 / (5000 * 3.6e-6));
289  for (const auto& hit : m_eclDiodeHits) {
290  int j = hit.getCellId() - 1; //0~8735
291  double hitE = hit.getEnergyDep() * m_calib[j].ascale * diodeEdep2crystalEdep;
292  double hitTimeAve = (hit.getTimeAve() + m_calib[j].tshift) * T2us;
293 
294  adccounts_t& a = m_adc[j];
295  // cout << "internuclearcountereffect " << j << " " << hit.getEnergyDep() << " " << hit.getTimeAve() << " " << a.total << endl;
296  // for (int i = 0; i < ec.m_nsmp; i++) cout << i << " " << a.c[i] << endl;
297  if (m_HadronPulseShape) {
298  a.AddHit(hitE, hitTimeAve + timeOffset, m_ss_HadronShapeSimulations[2]); // diode component
299  } else {
300  a.AddHit(hitE, hitTimeAve + timeOffset, m_ss[1]); // m_ss[1] is the sampled diode response
301  }
302  // for (int i = 0; i < ec.m_nsmp; i++) cout << i << " " << a.c[i] << endl;
303  }
304  }
305 
306  if (m_calibration) {
307  // This has been added by Alex Bobrov for calibration
308  // of covariance matrix artificially generate 100 MeV in time for each crystal
309  double hitE = 0.1, hitTimeAve = 0.0;
310  for (int j = 0; j < ec.m_nch; j++)
311  m_adc[j].AddHit(hitE, hitTimeAve + timeOffset, m_ss[m_tbl[j].iss]);
312  }
313 
314  return ttrig;
315 }
316 
317 void ECLDigitizerModule::makeElectronicNoiseAndPedestal(int J, int* FitA)
318 {
319  const EclConfiguration& ec = EclConfiguration::get();
320  float z[ec.m_nsmp], AdcNoise[ec.m_nsmp]; // buffers with electronic noise
321  // Noise generation
322  for (int i = 0; i < ec.m_nsmp; i++) z[i] = gRandom->Gaus(0, 1);
323  m_noise[m_tbl[J].inoise].generateCorrelatedNoise(z, AdcNoise);
324  for (int i = 0; i < ec.m_nsmp; i++) FitA[i] = 20 * AdcNoise[i] + 3000;
325 }
326 
327 void ECLDigitizerModule::makeWaveforms()
328 {
329  const EclConfiguration& ec = EclConfiguration::get();
330  BitStream out(ec.m_nch * ec.m_nsmp * 18 / 32);
331  out.putNBits(m_compAlgo & 0xff, 8);
332  ECLCompress* comp = selectAlgo(m_compAlgo & 0xff);
333  if (comp == NULL)
334  B2FATAL("Unknown compression algorithm: " << m_compAlgo);
335 
336  int FitA[ec.m_nsmp]; // buffer for the waveform fitter
337  // loop over entire calorimeter
338  for (int j = 0; j < ec.m_nch; j++) {
339  adccounts_t& a = m_adc[j];
340  makeElectronicNoiseAndPedestal(j, FitA);
341  for (int i = 0; i < ec.m_nsmp; i++) {
342  int A = 20000 * a.c[i] + FitA[i];
343  FitA[i] = max(0, min(A, (1 << 18) - 1));
344  }
345  comp->compress(out, FitA);
346  }
347  out.resize();
348 
349  ECLWaveforms* wf = new ECLWaveforms;
350  m_eclWaveforms.assign(wf);
351 
352  std::swap(out.getStore(), wf->getStore());
353  if (comp) delete comp;
354 }
355 
356 void ECLDigitizerModule::event()
357 {
358  const EclConfiguration& ec = EclConfiguration::get();
359  const int ttrig = shapeSignals();
360 
361  // We want to produce background waveforms in simulation first than
362  // dump to a disk, read from the disk to test before real data
363  if (m_waveformMaker) { makeWaveforms(); return; }
364 
365  // make relation between cellid and eclhits
366  struct ch_t {int cell, id;};
367  vector<ch_t> hitmap;
368  for (const auto& hit : m_eclHits) {
369  int j = hit.getCellId() - 1; //0~8735
370  if (hit.getBackgroundTag() == BackgroundMetaData::bg_none) hitmap.push_back({j, hit.getArrayIndex()});
371  // cout<<"C:"<<hit.getBackgroundTag()<<" "<<hit.getCellId()<<" "<<hit.getEnergyDep()<<" "<<hit.getTimeAve()<<endl;
372  }
373 
374  bool isBGOverlay = m_eclWaveforms.isValid();
375  BitStream out;
376  ECLCompress* comp = NULL;
377 
378  // check background overlay
379  if (isBGOverlay) {
380  std::swap(out.getStore(), m_eclWaveforms->getStore());
381  out.setPos(0);
382  unsigned int compAlgo = out.getNBits(8);
383  comp = selectAlgo(compAlgo);
384  if (comp == NULL)
385  B2FATAL("Unknown compression algorithm: " << compAlgo);
386  }
387 
388  int FitA[ec.m_nsmp]; // buffer for the waveform fitter
389 
390  // loop over entire calorimeter
391  for (int j = 0; j < ec.m_nch; j++) {
392  adccounts_t& a = m_adc[j];
393 
394  //normalize the MC true arrival times
395  if (m_adc[j].totalDep > 0) {
396  m_adc[j].flighttime /= m_adc[j].totalDep;
397  m_adc[j].timeshift /= m_adc[j].totalDep;
398  m_adc[j].timetosensor /= m_adc[j].totalDep;
399  }
400 
401  // if background waveform is here there is no need to generate
402  // electronic noise since it is already in the waveform
403  if (isBGOverlay) {
404  comp->uncompress(out, FitA);
405  } else {
406  // Signal amplitude should be above 100 keV
407  if (a.total < 0.0001) continue;
408  makeElectronicNoiseAndPedestal(j, FitA);
409  }
410 
411  for (int i = 0; i < ec.m_nsmp; i++) {
412  int A = 20000 * a.c[i] + FitA[i];
413  FitA[i] = max(0, min(A, (1 << 18) - 1));
414  }
415 
416  int energyFit = 0; // fit output : Amplitude 18 bits
417  int tFit = 0; // fit output : T_ave 12 bits
418  int qualityFit = 0; // fit output : quality 2 bits
419  int chi = 0; // fit output : chi square it is not available in the experiment
420 
421  shapeFitterWrapper(j, FitA, ttrig, energyFit, tFit, qualityFit, chi);
422 
423  if (energyFit > m_ADCThreshold) {
424  int CellId = j + 1;
425 
426  //note energyFit and m_Awave is in ADC units
427  if (energyFit > m_Awave[CellId - 1]) {
428  //only save waveforms above ADC threshold
429  const auto eclDsp = m_eclDsps.appendNew();
430  eclDsp->setCellId(CellId);
431  eclDsp->setDspA(FitA);
432  }
433 
434  // only store extra MC info if requested and above threshold
435  if (m_storeDspWithExtraMCInfo and a.totalDep >= m_DspWithExtraMCInfoThreshold) {
436  const auto eclDspWithExtraMCInfo = m_eclDspsWithExtraMCInfo.appendNew();
437  eclDspWithExtraMCInfo->setCellId(CellId);
438  eclDspWithExtraMCInfo->setDspA(FitA);
439  eclDspWithExtraMCInfo->setEnergyDep(a.totalDep);
440  eclDspWithExtraMCInfo->setHadronEnergyDep(a.totalHadronDep);
441  eclDspWithExtraMCInfo->setFlightTime(a.flighttime);
442  eclDspWithExtraMCInfo->setTimeShift(a.timeshift);
443  eclDspWithExtraMCInfo->setTimeToSensor(a.timetosensor);
444  eclDspWithExtraMCInfo->setEnergyConversion(a.energyConversion * 20000);
445  }
446 
447  const auto eclDigit = m_eclDigits.appendNew();
448  eclDigit->setCellId(CellId); // cellId in range from 1 to 8736
449  eclDigit->setAmp(energyFit); // E (GeV) = energyFit/20000;
450  eclDigit->setTimeFit(tFit); // t0 (us)= (1520 - m_ltr)*24.*12/508/(3072/2) ;
451  eclDigit->setQuality(qualityFit);
452  if (qualityFit == 2)
453  eclDigit->setChi(chi);
454  else eclDigit->setChi(0);
455  for (const auto& hit : hitmap)
456  if (hit.cell == j) eclDigit->addRelationTo(m_eclHits[hit.id]);
457 
458  // set relation to DspWithExtraInfo
459  for (auto& DspWithExtraMCInfo : m_eclDspsWithExtraMCInfo) {
460  if (eclDigit->getCellId() == DspWithExtraMCInfo.getCellId()) DspWithExtraMCInfo.addRelationTo(eclDigit);
461  }
462  }
463  } //store each crystal hit
464  if (comp) delete comp;
465 }
466 
467 void ECLDigitizerModule::endRun()
468 {
469 }
470 
471 void ECLDigitizerModule::terminate()
472 {
473 }
474 
475 void ECLDigitizerModule::callbackHadronSignalShapes()
476 {
477 
478  if (m_waveformParametersMC.hasChanged()) {
479 
480  m_ss_HadronShapeSimulations.resize(3);
481 
482  //read MC templates from database
483  float photonParsPSD[10];
484  float hadronParsPSD[10];
485  float diodeParsPSD[10];
486  for (int i = 0; i < 10; i++) {
487  photonParsPSD[i] = m_waveformParametersMC->getPhotonParameters()[i + 1];
488  hadronParsPSD[i] = m_waveformParametersMC->getHadronParameters()[i + 1];
489  diodeParsPSD[i] = m_waveformParametersMC->getDiodeParameters()[i + 1];
490  }
491 
492  //Initialize templates for hadron shape simulations
493  m_ss_HadronShapeSimulations[0].InitSample(photonParsPSD, m_waveformParametersMC->getPhotonParameters()[0]);
494  m_ss_HadronShapeSimulations[1].InitSample(hadronParsPSD, m_waveformParametersMC->getHadronParameters()[0]);
495  m_ss_HadronShapeSimulations[2].InitSample(diodeParsPSD, m_waveformParametersMC->getDiodeParameters()[0]);
496 
497  }
498 
499 }
500 
501 void ECLDigitizerModule::readDSPDB()
502 {
503  const EclConfiguration& ec = EclConfiguration::get();
504 
505  string dataFileName;
506  if (m_background) {
507  dataFileName = FileSystem::findFile("/data/ecl/ECL-WF-BG.root");
508  B2DEBUG(150, "ECLDigitizer: Reading configuration data with background from: " << dataFileName);
509  } else {
510  dataFileName = FileSystem::findFile("/data/ecl/ECL-WF.root");
511  B2DEBUG(150, "ECLDigitizer: Reading configuration data without background from: " << dataFileName);
512  }
513  assert(! dataFileName.empty());
514 
515  TFile rootfile(dataFileName.c_str());
516  TTree* tree = (TTree*) rootfile.Get("EclWF");
517  TTree* tree2 = (TTree*) rootfile.Get("EclAlgo");
518  TTree* tree3 = (TTree*) rootfile.Get("EclNoise");
519 
520  if (tree == 0 || tree2 == 0 || tree3 == 0) B2FATAL("Data not found");
521 
522  m_tbl.resize(ec.m_nch);
523 
524  const int maxncellid = 512;
525  int ncellId;
526  vector<int> cellId(maxncellid);//[ncellId] buffer for crystal identification number
527 
528  tree->SetBranchAddress("ncellId", &ncellId);
529  tree->SetBranchAddress("cellId", cellId.data());
530 
531  vector<int> eclWaveformDataTable(ec.m_nch);
532  for (Long64_t j = 0, jmax = tree->GetEntries(); j < jmax; j++) {
533  tree->GetEntry(j);
534  assert(ncellId <= maxncellid);
535  for (int i = 0; i < ncellId; ++i)
536  eclWaveformDataTable[cellId[i] - 1] = j;
537  }
538  B2DEBUG(150, "ECLDigitizer: " << tree->GetEntries() << " sets of wave form covariance matricies will be used.");
539 
540  ECLWFAlgoParams* algo = new ECLWFAlgoParams;
541  tree2->SetBranchAddress("Algopars", &algo);
542  tree2->SetBranchAddress("ncellId", &ncellId);
543  tree2->SetBranchAddress("cellId", cellId.data());
544  Long64_t jmax2 = tree2->GetEntries();
545  vector<ECLWFAlgoParams> eclWFAlgoParams;
546  eclWFAlgoParams.reserve(jmax2);
547  for (Long64_t j = 0; j < jmax2; j++) {
548  tree2->GetEntry(j);
549  assert(ncellId <= maxncellid);
550  eclWFAlgoParams.push_back(*algo);
551  for (int i = 0; i < ncellId; ++i)
552  m_tbl[cellId[i] - 1].idn = j;
553  }
554  if (algo) delete algo;
555  B2DEBUG(150, "ECLDigitizer: " << eclWFAlgoParams.size() << " parameter sets of fitting algorithm were read.");
556 
557  ECLNoiseData* noise = new ECLNoiseData;
558  tree3->SetBranchAddress("NoiseM", &noise);
559  tree3->SetBranchAddress("ncellId", &ncellId);
560  tree3->SetBranchAddress("cellId", cellId.data());
561 
562  Long64_t jmax3 = tree3->GetEntries();
563  m_noise.reserve(jmax3);
564  for (Long64_t j = 0; j < jmax3; j++) {
565  tree3->GetEntry(j);
566  assert(ncellId <= maxncellid);
567  m_noise.push_back(*noise);
568  if (ncellId == 0) {
569  for (int i = 0; i < ec.m_nch; i++)
570  m_tbl[i].inoise = 0;
571  break;
572  } else {
573  for (int i = 0; i < ncellId; ++i)
574  m_tbl[cellId[i] - 1].inoise = j;
575  }
576  }
577  if (noise) delete noise;
578  B2DEBUG(150, "ECLDigitizer: " << eclWFAlgoParams.size() << " noise matricies were loaded.");
579 
580  // repack fitting algorithm parameters
581  m_idn.resize(eclWFAlgoParams.size());
582  for (int i = 0, imax = eclWFAlgoParams.size(); i < imax; i++)
583  repack(eclWFAlgoParams[i], m_idn[i]);
584 
585  vector<uint_pair_t> pairIdx;
586  for (int i = 0; i < ec.m_nch; i++) {
587  unsigned int wfIdx = eclWaveformDataTable[i];
588  unsigned int algoIdx = m_tbl[i].idn;
589  uint_pair_t p(wfIdx, algoIdx);
590  vector<uint_pair_t>::iterator ip = find(pairIdx.begin(), pairIdx.end(), p);
591  if (ip != pairIdx.end()) { // crystal i already have the same parameters as before
592  m_tbl[i].ifunc = ip - pairIdx.begin();
593  } else { // new combination of parameters
594  m_tbl[i].ifunc = pairIdx.size();
595  pairIdx.push_back(p);
596  }
597  }
598 
599  // now we know how many distinct elements of the (Waveform x AlgoParams) matrix should be
600  m_fitparams.resize(pairIdx.size());
601 
602  ECLWaveformData* eclWFData = new ECLWaveformData;
603  tree->SetBranchAddress("CovarianceM", &eclWFData);
604  tree->SetBranchStatus("ncellId", 0); // do not read it at the moment
605  tree->SetBranchStatus("cellId", 0); // do not read it at the moment
606 
607 // fill parameters for each (Waveform x AlgoParams) parameter set
608  for (unsigned int ip = 0; ip < pairIdx.size(); ip++) {
609  const uint_pair_t& p = pairIdx[ip];
610  tree->GetEntry(p.first); // retrieve data to eclWFData pointer
611  getfitparams(*eclWFData, eclWFAlgoParams[p.second], m_fitparams[ip]);
612  }
613  B2DEBUG(150, "ECLDigitizer: " << m_fitparams.size() << " fitting crystals groups were created.");
614 
615  // at the moment there is only one sampled signal shape in the pool
616  // since all shaper parameters are the same for all crystals
617  m_ss.resize(2);
618  float MP[10]; eclWFData->getWaveformParArray(MP);
619  m_ss[0].InitSample(MP, 27.7221);
620  // parameters vector from ps.dat file, time offset 0.5 usec added to
621  // have peak position with parameters from ps.dat roughly in the
622  // same place as in current MC
623  // double crystal_params[10] = {0.5, 0.301298, 0.631401, 0.470911, 0.903988, -0.11734200/19.5216, 2.26567, 0.675393, 0.683995, 0.0498786};
624  // m_ss[0].InitSample(crystal_params, 0.9999272*19.5216);
625  for (int i = 0; i < ec.m_nch; i++) m_tbl[i].iss = 0;
626  // one sampled diode response in the pool, parameters vector from
627  // pg.dat file, time offset 0.5 usec added to have peak position with
628  // parameters from ps.dat roughly in the same place as in current MC
629  double diode_params[] = {0 + 0.5, 0.100002, 0.756483, 0.456153, 0.0729031, 0.3906 / 9.98822, 2.85128, 0.842469, 0.854184, 0.110284};
630  m_ss[1].InitSample(diode_params, 0.9569100 * 9.98822);
631 
632  B2DEBUG(150, "ECLDigitizer: " << m_ss.size() << " sampled signal templates were created.");
633 
634  if (eclWFData) delete eclWFData;
635 
636  rootfile.Close();
637 }
638 
639 void ECLDigitizerModule::repack(const ECLWFAlgoParams& eclWFAlgo, algoparams_t& t)
640 {
641  // filling short int array
642  t.id[ 0] = eclWFAlgo.getlAT() + 128;
643  t.id[ 1] = eclWFAlgo.getsT() + 128;
644  t.id[ 2] = eclWFAlgo.gethT();
645  t.id[ 3] = 17952;
646  t.id[ 4] = 19529;
647  t.id[ 5] = 69;
648  t.id[ 6] = 0;
649  t.id[ 7] = 0;
650  t.id[ 8] = 257;
651  t.id[ 9] = -1;
652  t.id[10] = 0;
653  t.id[11] = 0;
654 
655  // filling unsigned char array
656  t.ic[12 * 2 + 0] = eclWFAlgo.getk1();
657  t.ic[12 * 2 + 1] = eclWFAlgo.getk2();
658  t.ic[13 * 2 + 0] = eclWFAlgo.getka();
659  t.ic[13 * 2 + 1] = eclWFAlgo.getkb();
660  t.ic[14 * 2 + 0] = eclWFAlgo.getkc();
661  t.ic[14 * 2 + 1] = eclWFAlgo.gety0s() - 16;
662 
663  // again, filling short int array
664  t.id[15] = eclWFAlgo.getcT();
665 }
666 
667 void ECLDigitizerModule::getfitparams(const ECLWaveformData& eclWFData, const ECLWFAlgoParams& eclWFAlgo, fitparams_t& p)
668 {
669  const EclConfiguration& ec = EclConfiguration::get();
670 
671  double ssd[16][16];
672  eclWFData.getMatrix(ssd);
673  vector<double> MP(10);
674  eclWFData.getWaveformParArray(MP.data());
675 
676  // shape function parameters
677  int iff = 1 << 14; //Alex Bobrov for correct chi square
678 
679  int ia = 1 << eclWFAlgo.getka();
680  int ib = 1 << eclWFAlgo.getkb();
681  int ic = 1 << eclWFAlgo.getkc();
682 
683  int_array_192x16_t& ref_f = p.f;
684  int_array_192x16_t& ref_f1 = p.f1;
685  int_array_192x16_t& ref_fg31 = p.fg31;
686  int_array_192x16_t& ref_fg32 = p.fg32;
687  int_array_192x16_t& ref_fg33 = p.fg33;
688  int_array_24x16_t& ref_fg41 = p.fg41;
689  int_array_24x16_t& ref_fg43 = p.fg43;
690 
691  ShaperDSP_t dsp(MP);
692  dsp.settimestride(ec.m_step);
693  dsp.setseedoffset(ec.m_step / ec.m_ndt);
694  dsp.settimeseed(-(ec.m_step - (ec.m_step / ec.m_ndt)));
695  vector<dd_t> f(16);
696  for (int k = 0; k < 2 * ec.m_ndt; k++, dsp.nextseed()) { // calculate fit parameters around 0 +- 1 ADC tick
697  dsp.fillvector(f);
698 
699  double g0g0 = 0, g0g1 = 0, g1g1 = 0, g0g2 = 0, g1g2 = 0, g2g2 = 0;
700  double sg0[16], sg1[16], sg2[16];
701  for (int j = 0; j < 16; j++) {
702  double g0 = 0, g1 = 0, g2 = 0;
703  for (int i = 0; i < 16; i++) {
704  g0 += ssd[j][i] * f[i].first;
705  g1 += ssd[j][i] * f[i].second;
706  g2 += ssd[j][i];
707  }
708  g0g0 += g0 * f[j].first;
709  g0g1 += g1 * f[j].first;
710  g1g1 += g1 * f[j].second;
711  g0g2 += g0;
712  g1g2 += g1;
713  g2g2 += g2;
714  sg0[j] = g0;
715  sg1[j] = g1;
716  sg2[j] = g2;
717  }
718 
719  double a00 = g1g1 * g2g2 - g1g2 * g1g2;
720  double a11 = g0g0 * g2g2 - g0g2 * g0g2;
721  double a22 = g0g0 * g1g1 - g0g1 * g0g1;
722  double a01 = g1g2 * g0g2 - g0g1 * g2g2;
723  double a02 = g0g1 * g1g2 - g1g1 * g0g2;
724  double a12 = g0g1 * g0g2 - g0g0 * g1g2;
725 
726  double igg2 = 1 / (a11 * g1g1 + g0g1 * a01 + g1g2 * a12);
727 
728  // to fixed point representation
729  const double isd = 3. / 4., sd = 1 / isd ; // conversion factor (???)
730  for (int i = 0; i < 16; i++) {
731  double w = i ? 1.0 : 1. / 16.;
732 
733  ref_f [k][i] = lrint(f[i].first * iff * w);
734  ref_f1 [k][i] = lrint(f[i].second * iff * w * sd);
735 
736  double fg31 = (a00 * sg0[i] + a01 * sg1[i] + a02 * sg2[i]) * igg2;
737  double fg32 = (a01 * sg0[i] + a11 * sg1[i] + a12 * sg2[i]) * igg2;
738  double fg33 = (a02 * sg0[i] + a12 * sg1[i] + a22 * sg2[i]) * igg2;
739 
740  ref_fg31[k][i] = lrint(fg31 * ia * w);
741  ref_fg32[k][i] = lrint(fg32 * ib * w * isd);
742  ref_fg33[k][i] = lrint(fg33 * ic * w);
743  }
744 
745  //first approximation without time correction
746  int jk = 23 + ((48 - k) >> 2);
747  if (jk >= 0 && jk < 24 && (48 - k) % 4 == 0) {
748 
749  double igg1 = 1 / a11;
750  // to fixed point
751  for (int i = 0; i < 16; i++) {
752  double w = i ? 1.0 : 1. / 16.;
753 
754  double fg41 = (g2g2 * sg0[i] - g0g2 * sg2[i]) * igg1;
755  double fg43 = (g0g0 * sg2[i] - g0g2 * sg0[i]) * igg1;
756  ref_fg41[jk][i] = lrint(fg41 * ia * w);
757  ref_fg43[jk][i] = lrint(fg43 * ic * w);
758  }
759  }
760  }
761 }
Belle2::ECL::ShaperDSP_t::settimeseed
void settimeseed(double)
set initial time
Definition: shaperdsp.cc:363
Belle2::ECL::EclConfiguration::m_nch
static constexpr int m_nch
total number of electronic channels (crystals) in calorimeter
Definition: EclConfiguration.h:44
Belle2::ECLNoiseData
Container for constant matrix used to generate electronic noise.
Definition: ECLWaveformData.h:178
Belle2::ECL::EclConfiguration::m_ntrg
static constexpr int m_ntrg
number of trigger counts per ADC clock tick
Definition: EclConfiguration.h:49
Belle2::ECL::ShaperDSP_t
Class include function that calculate electronic response from energy deposit
Definition: shaperdsp.h:30
Belle2::ECL::EclConfiguration::m_nsmp
static constexpr int m_nsmp
number of ADC measurements for signal fitting
Definition: EclConfiguration.h:51
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::ECLDigitizerModule::uint_pair_t
std::pair< unsigned int, unsigned int > uint_pair_t
a pair of unsigned ints
Definition: ECLDigitizerModule.h:118
Belle2::ECL::BitStream
Bit stream struct.
Definition: ECLCompress.h:39
Belle2::ECL::ShaperDSP_t::fillvector
void fillvector(std::vector< double > &) const
fill vector with response function values and its derivative
Definition: shaperdsp.cc:409
Belle2::ECLDigitizerModule::int_array_24x16_t
fitparams_t::int_array_24x16_t int_array_24x16_t
weighting coefficients amplitude calculation.
Definition: ECLDigitizerModule.h:117
Belle2::ECL::EclConfiguration::algoparams_t
a struct for the parameters of the algorithm
Definition: EclConfiguration.h:97
Belle2::ECL::ShaperDSP_t::setseedoffset
void setseedoffset(double)
set timeoffset
Definition: shaperdsp.cc:358
Belle2::ECL::ShaperDSP_t::settimestride
void settimestride(double)
set grid step for function calculation
Definition: shaperdsp.cc:353
Belle2::ECLWaveforms::getStore
std::vector< unsigned int > & getStore()
get data
Definition: ECLWaveforms.h:38
Belle2::ECL::ECLShapeFit
ShaperDSP fit results from _lftda function.
Definition: ECLDspEmulator.h:30
Belle2::ECL::ECLCompress::uncompress
virtual void uncompress(BitStream &in, int *adc)=0
Decompress the ECL waveform.
Belle2::ECLWaveformData::getWaveformParArray
void getWaveformParArray(float P[10]) const
Getter method for waveform shape parameters as one dimentional array of floats.
Definition: ECLWaveformData.h:122
Belle2::DBObjPtr
Class for accessing objects in the database.
Definition: DBObjPtr.h:31
Belle2::Module
Base class for Modules.
Definition: Module.h:74
Belle2::ECL::EclConfiguration
singleton class to hold the ECL configuration
Definition: EclConfiguration.h:31
Belle2::ECL::EclConfiguration::fitparams_t
a struct for the fit parameters
Definition: EclConfiguration.h:89
Belle2::ECLDigitWaveformParametersForMC
DB object to store photon, hadron and diode shape parameters used in simulations.
Definition: ECLDigitWaveformParametersForMC.h:41
Belle2::ECL::ShaperDSP_t::nextseed
void nextseed()
substruct toffset to tzero
Definition: shaperdsp.cc:368
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::ECLDigitizerModule::int_array_192x16_t
fitparams_t::int_array_192x16_t int_array_192x16_t
weighting coefficients for time and amplitude calculation
Definition: ECLDigitizerModule.h:115
Belle2::ECL::ECLCompress
Abstract class (interface) for ECL waveform compression/decompression to/from the BitStream storage.
Definition: ECLCompress.h:118
Belle2::ECL::EclConfiguration::s_clock
static constexpr double s_clock
digitization clock in RF units
Definition: EclConfiguration.h:48
Belle2::ECL::EclConfiguration::m_rf
static constexpr double m_rf
accelerating RF, http://ptep.oxfordjournals.org/content/2013/3/03A006.full.pdf
Definition: EclConfiguration.h:45
Belle2::ECL::EclConfiguration::m_ndt
static constexpr int m_ndt
number of points per ADC tick where signal fit procedure parameters are evaluated
Definition: EclConfiguration.h:57
Belle2::ECLWaveformData
ECLWaveformData - container for inverse covariant matrix and shape parameters for time and amplitude ...
Definition: ECLWaveformData.h:39
Belle2::ECLDigitizerModule
The ECLDigitizer module.
Definition: ECLDigitizerModule.h:79
Belle2::ECL::ECLGeometryPar
The Class for ECL Geometry Parameters.
Definition: ECLGeometryPar.h:45
Belle2::ECLWaveforms
Class to store ECL waveforms for entire calorimeter.
Definition: ECLWaveforms.h:32
Belle2::ECL::EclConfiguration::m_step
static constexpr double m_step
time between points in internal units t_{asrto}*m_rf/2.
Definition: EclConfiguration.h:47
Belle2::ECLWaveformData::getMatrix
void getMatrix(float M[16][16]) const
Getter method for all matrix as two dimentional array (floats)
Definition: ECLWaveformData.h:66
Belle2::ECLDigitizerModule::calibration_t
calibration constants per channel
Definition: ECLDigitizerModule.h:141
Belle2::ECL::EclConfiguration::adccounts_t
a struct for the ADC count
Definition: EclConfiguration.h:75
Belle2::ECL::ECLCompress::compress
virtual void compress(BitStream &out, const int *adc)=0
Compress the ECL waveform.
Belle2::ECLWFAlgoParams
Container for constant parameters used in waveform fits.
Definition: ECLWaveformData.h:150