Belle II Software  release-08-00-10
ECLDigitizerModule.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 /* Own header. */
10 #include <ecl/modules/eclDigitizer/ECLDigitizerModule.h>
11 
12 /* ECL headers. */
13 #include <ecl/dataobjects/ECLDigit.h>
14 #include <ecl/dataobjects/ECLDsp.h>
15 #include <ecl/dataobjects/ECLDspWithExtraMCInfo.h>
16 #include <ecl/dataobjects/ECLHit.h>
17 #include <ecl/dataobjects/ECLSimHit.h>
18 #include <ecl/dataobjects/ECLTrig.h>
19 #include <ecl/dataobjects/ECLWaveforms.h>
20 #include <ecl/dbobjects/ECLWaveformData.h>
21 #include <ecl/digitization/algorithms.h>
22 #include <ecl/digitization/BitStream.h>
23 #include <ecl/digitization/ECLCompress.h>
24 #include <ecl/digitization/shaperdsp.h>
25 #include <ecl/geometry/ECLGeometryPar.h>
26 #include <ecl/utility/ECLDspUtilities.h>
27 
28 /* Basf2 headers. */
29 #include <framework/gearbox/Unit.h>
30 #include <framework/logging/Logger.h>
31 #include <framework/utilities/FileSystem.h>
32 
33 /* ROOT headers. */
34 #include <TFile.h>
35 #include <TRandom.h>
36 #include <TTree.h>
37 
38 using namespace std;
39 using namespace Belle2;
40 using namespace ECL;
41 
42 //-----------------------------------------------------------------
43 // Register the Module
44 //-----------------------------------------------------------------
45 REG_MODULE(ECLDigitizer);
46 
47 //-----------------------------------------------------------------
48 // Implementation
49 //-----------------------------------------------------------------
50 
51 ECLDigitizerModule::ECLDigitizerModule() : Module(), m_waveformParametersMC("ECLDigitWaveformParametersForMC"),
52  m_waveformParameters("ECLWFParameters"),
53  m_algoParameters("ECLWFAlgoParams"),
54  m_noiseParameters("ECLWFNoiseParams")
55 {
56  //Set module properties
57  setDescription("Creates ECLDigiHits from ECLHits.");
59  addParam("TriggerTime", m_trigTime,
60  "Flag to use crate trigger times from beam background overlay if there are any (default: false)", false);
61  addParam("Background", m_background, "Flag to use the Digitizer configuration with backgrounds (default: false)", false);
62  addParam("Calibration", m_calibration, "Flag to use the Digitizer for Waveform fit Covariance Matrix calibration (default: false)",
63  false);
64  addParam("DiodeDeposition", m_inter,
65  "Flag to take into account energy deposition in photodiodes; Default diode is sensitive detector (default: true)", true);
66  addParam("WaveformMaker", m_waveformMaker, "Flag to produce background waveform digits (default: false)", false);
67  addParam("CompressionAlgorithm", m_compAlgo, "Waveform compression algorithm (default: 0u)", 0u);
68  addParam("eclWaveformsName", m_eclWaveformsName, "Name of the output/input collection (digitized waveforms)", string(""));
69  addParam("HadronPulseShapes", m_HadronPulseShape, "Flag to include hadron component in pulse shape construction (default: true)",
70  true);
71  addParam("ADCThreshold", m_ADCThreshold, "ADC threshold for waveform fits (default: 25)", 25);
72  addParam("WaveformThresholdOverride", m_WaveformThresholdOverride,
73  "If gt 0 value is applied to all crystals for waveform saving threshold. If lt 0 dbobject is used. (GeV)", -1.0);
74  addParam("StoreDspWithExtraMCInfo", m_storeDspWithExtraMCInfo,
75  "Flag to store Dsp with extra MC information in addition to normal Dsp (default: false)", false);
76  addParam("DspWithExtraMCInfoThreshold", m_DspWithExtraMCInfoThreshold,
77  "Threshold above with to store Dsp with extra MC information [GeV]",
78  0.02);
79  addParam("DspDataTest", m_dspDataTest,
80  "Use DSP coefficients from the database for the processing. This "
81  "will significantly reduce performance so this option is to be "
82  "used for testing only.", false);
83  addParam("useWaveformParameters", m_useWaveformParameters,
84  "Use ECLWF{Parameters,AlgoParams,NoiseParams} payloads", true);
85  addParam("unitscale", m_unitscale,
86  "Normalization coefficient for ECL signal shape. "
87  "If positive, use same static value for all ECL channels. "
88  "If negative, calculate dynamically at beginRun().", -1.0);
89 
90 }
91 
93 {
94 }
95 
97 {
98  m_eclDsps.registerInDataStore();
99  if (m_storeDspWithExtraMCInfo) m_eclDspsWithExtraMCInfo.registerInDataStore();
100  m_eclDigits.registerInDataStore();
101  m_eclTrigs.registerInDataStore();
102 
103  if (m_HadronPulseShape) {
104  B2DEBUG(20,
105  "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");
106  } else {
107  B2DEBUG(20, "Hadron pulse shapes for ECL simulations are disabled.");
108  }
109 
110  if (m_inter) {
111  B2DEBUG(20, "Diode-crossing pulse shapes for ECL simulations are enabled.");
112  } else {
113  B2DEBUG(20, "Diode-crossing pulse shapes for ECL simulations are disabled.");
114  }
115 
116  m_eclDiodeHits.registerInDataStore("ECLDiodeHits");
117 
118  m_eclDsps.registerRelationTo(m_eclDigits);
120  m_eclDspsWithExtraMCInfo.registerRelationTo(m_eclDigits);
121  m_eclDigits.registerRelationTo(m_eclHits);
122  if (m_waveformMaker)
123  m_eclWaveforms.registerInDataStore(m_eclWaveformsName);
124 
126 
128 }
129 
131 {
133  double ns_per_tick = 1.0 / (4.0 * ec.getRF()) * 1e3;// ~0.49126819903043308239 ns/tick
134 
136  readDSPDB();
137  m_loadOnce = false;
138  }
139 
140  calibration_t def = {1, 0};
142 
143  if (m_CrystalElectronics.isValid()) {
144  for (int i = 0; i < ECLElementNumbers::c_NCrystals; i++)
145  m_calib[i].ascale /= m_CrystalElectronics->getCalibVector()[i];
146  }
147  if (m_CrystalEnergy.isValid()) {
148  for (int i = 0; i < ECLElementNumbers::c_NCrystals; i++)
149  m_calib[i].ascale /= m_CrystalEnergy->getCalibVector()[i] * 20000.0;
150  }
151  if (m_CrystalElectronicsTime.isValid()) {
152  for (int i = 0; i < ECLElementNumbers::c_NCrystals; i++)
153  m_calib[i].tshift += m_CrystalElectronicsTime->getCalibVector()[i] * ns_per_tick;
154  }
155  if (m_CrystalTimeOffset.isValid()) {
156  for (int i = 0; i < ECLElementNumbers::c_NCrystals; i++)
157  m_calib[i].tshift += m_CrystalTimeOffset->getCalibVector()[i] * ns_per_tick;
158  }
159  if (m_CrateTimeOffset.isValid()) {
160  for (int i = 0; i < ECLElementNumbers::c_NCrystals; i++)
161  m_calib[i].tshift += m_CrateTimeOffset->getCalibVector()[i] * ns_per_tick;
162  }
163  if (m_MCTimeOffset.isValid()) {
164  for (int i = 0; i < ECLElementNumbers::c_NCrystals; i++)
165  m_calib[i].tshift += m_MCTimeOffset->getCalibVector()[i] * ns_per_tick;
166  }
168  if (m_WaveformThresholdOverride < 0) {
169  if (m_FPGAWaveform.isValid()) {
170  for (int i = 0; i < ECLElementNumbers::c_NCrystals; i++)
171  m_Awave[i] = m_FPGAWaveform->getCalibVector()[i];
172  }
173  } else {
174  //If m_WaveformThresholdOverride > 0 override ECL_FPGA_StoreWaveform;
175  for (int i = 0; i < ECLElementNumbers::c_NCrystals; i++)
176  m_Awave[i] = m_WaveformThresholdOverride * m_calib[i].ascale; // convert GeV to ADC
177  }
178 
179  if (m_HadronPulseShape)
181 
182  // Initialize channel mapper at run start to account for possible
183  // changes in ECL mapping between runs.
184  if (!m_eclMapper.initFromDB()) {
185  B2FATAL("ECLDigitizer: Can't initialize eclChannelMapper!");
186  }
187 }
188 
189 // interface to C shape fitting function function
190 void ECLDigitizerModule::shapeFitterWrapper(const int j, const int* FitA, const int ttrig,
191  int& m_lar, int& m_ltr, int& m_lq, int& m_chi) const
192 {
193  const crystallinks_t& t = m_tbl[j]; //lookup table [0,8735]
194  const fitparams_t& r = m_fitparams[t.ifunc];
195 
196  ECLShapeFit result;
197 
198  if (!m_dspDataTest) {
199  short int* id = (short int*)m_idn[t.idn].id;
200 
201  int A0 = (int) * (id + 0) - 128;
202  int Askip = (int) * (id + 1) - 128;
203 
204  int Ahard = (int) * (id + 2);
205  int k_a = (int) * ((unsigned char*)id + 26);
206  int k_b = (int) * ((unsigned char*)id + 27);
207  int k_c = (int) * ((unsigned char*)id + 28);
208  int k_16 = (int) * ((unsigned char*)id + 29);
209  int k1_chi = (int) * ((unsigned char*)id + 24);
210  int k2_chi = (int) * ((unsigned char*)id + 25);
211 
212  int chi_thres = (int) * (id + 15);
213 
214  int trg_time = ttrig;
215 
216  result = lftda_((int*)r.f, (int*)r.f1, (int*)r.fg41, (int*)r.fg43,
217  (int*)r.fg31, (int*)r.fg32, (int*)r.fg33, (int*)FitA,
218  trg_time, A0, Ahard, Askip, k_a, k_b, k_c, k_16, k1_chi,
219  k2_chi, chi_thres);
220  } else {
221  std::vector<int> adc(31);
222  for (int i = 0; i < 31; i++) adc[i] = FitA[i];
223  // NOTE: ttrig_packed is any even number in set {8*Z, 8*Z + 2, 8*Z + 4}
224  // where Z is an integer number in 0..23 range
225  // ttrig = ttrig_packed - 2 * (ttrig_packed / 8);
226  // ttrig = 6*Z + q
227  int ttrig_packed = ttrig / 6 * 8 + ttrig % 6;
228  result = ECLDspUtilities::shapeFitter(j + 1, adc, ttrig_packed);
229  }
230 
231  m_lar = result.amp;
232  m_ltr = result.time;
233  m_lq = result.quality;
234  m_chi = result.chi2;
235 
236  //== Set precision of chi^2 to be the same as in the raw data.
237  int discarded_bits = 0;
238  if ((m_chi & 0x7800000) != 0) {
239  m_chi = 0x7800000;
240  } else if ((m_chi & 0x0600000) != 0) {
241  discarded_bits = 14;
242  } else if ((m_chi & 0x0180000) != 0) {
243  discarded_bits = 12;
244  } else if ((m_chi & 0x0060000) != 0) {
245  discarded_bits = 10;
246  } else if ((m_chi & 0x0018000) != 0) {
247  discarded_bits = 8;
248  } else if ((m_chi & 0x0006000) != 0) {
249  discarded_bits = 6;
250  } else if ((m_chi & 0x0001800) != 0) {
251  discarded_bits = 4;
252  } else if ((m_chi & 0x0000600) != 0) {
253  discarded_bits = 2;
254  }
255  if (discarded_bits > 0) {
256  m_chi >>= discarded_bits;
257  m_chi <<= discarded_bits;
258  }
259 }
260 
262 {
265 
266  const double trgtick = ec.s_clock / ec.getRF() / ec.m_ntrg; // trigger tick
267  const double tscale = 2 * trgtick, toff = ec.s_clock / (2 * ec.getRF());
268 
269  // clear the storage for the event
270  memset(m_adc.data(), 0, m_adc.size()*sizeof(adccounts_t));
271 
272  const double E2GeV = 1 / Unit::GeV; // convert Geant energy units to GeV
273  const double T2us = 1 / Unit::us; // convert Geant time units to microseconds
274 
275  // emulate response for ECL hits after ADC measurements
276  for (const auto& hit : m_eclSimHits) {
277  int cellId = hit.getCellId(); // 1 .. 8736
278  int j = cellId - 1; // 0 .. 8735
279  int id = m_eclMapper.getCrateID(cellId) - 1; // 0 .. 51
280  double timeOffset = tscale * m_ttime[id] - toff;
281  double hitE = hit.getEnergyDep() * m_calib[j].ascale * E2GeV;
282  double hitTimeAve = (hit.getFlightTime() + m_calib[j].tshift + eclp->time2sensor(j, hit.getPosition())) * T2us;
283 
284  m_adc[j].energyConversion = m_calib[j].ascale * E2GeV * 20000;
285  m_adc[j].flighttime += hit.getFlightTime() * hit.getEnergyDep(); //true time weighted by energy
286  m_adc[j].timeshift += m_calib[j].tshift * hit.getEnergyDep();
287  m_adc[j].timetosensor += eclp->time2sensor(j, hit.getPosition()) * hit.getEnergyDep();
288  m_adc[j].totalHadronDep += hit.getHadronEnergyDep(); // true deposited energy hadron component (in GeV)
289  m_adc[j].totalDep += hit.getEnergyDep(); //true deposited energy (in GeV)
290 
291  if (m_HadronPulseShape == true) {
292  double hitHadronE = hit.getHadronEnergyDep() * m_calib[j].ascale * E2GeV;
293  m_adc[j].AddHit(hitE - hitHadronE, hitTimeAve + timeOffset, m_ss_HadronShapeSimulations[0]);//Gamma Component
294  m_adc[j].AddHit(hitHadronE, hitTimeAve + timeOffset, m_ss_HadronShapeSimulations[1]); //Hadron Component
295  } else {
296  m_adc[j].AddHit(hitE, hitTimeAve + timeOffset, m_ss[m_tbl[j].iss]);
297  }
298  }
299 
300  // add only background hits
301  for (const auto& hit : m_eclHits) {
302  if (hit.getBackgroundTag() == BackgroundMetaData::bg_none) continue;
303  int cellId = hit.getCellId(); // 1 .. 8736
304  int j = cellId - 1; // 0 .. 8735
305  int id = m_eclMapper.getCrateID(cellId) - 1; // 0 .. 51
306  double timeOffset = tscale * m_ttime[id] - toff;
307  double hitE = hit.getEnergyDep() * m_calib[j].ascale * E2GeV;
308  double hitTimeAve = (hit.getTimeAve() + m_calib[j].tshift) * T2us;
309  m_adc[j].AddHit(hitE, hitTimeAve + timeOffset, m_ss[m_tbl[j].iss]);
310  }
311 
312  // internuclear counter effect -- charged particle crosses diode and produces signal
313  if (m_inter) {
314  // ionisation energy in Si is I = 3.6x10^-6 MeV for electron-hole pair
315  // 5000 pairs in the diode per 1 MeV deposited in the crystal attached to the diode
316  // conversion factor to get equvalent energy deposition in the crystal to sum up it with deposition in crystal
317  const double diodeEdep2crystalEdep = E2GeV * (1 / (5000 * 3.6e-6));
318  for (const auto& hit : m_eclDiodeHits) {
319  int cellId = hit.getCellId(); // 1 .. 8736
320  int j = cellId - 1; // 0 .. 8735
321  int id = m_eclMapper.getCrateID(cellId) - 1; // 0 .. 51
322  double timeOffset = tscale * m_ttime[id] - toff;
323  double hitE = hit.getEnergyDep() * m_calib[j].ascale * diodeEdep2crystalEdep;
324  double hitTimeAve = (hit.getTimeAve() + m_calib[j].tshift) * T2us;
325 
326  adccounts_t& a = m_adc[j];
327  // cout << "internuclearcountereffect " << j << " " << hit.getEnergyDep() << " " << hit.getTimeAve() << " " << a.total << endl;
328  // for (int i = 0; i < ec.m_nsmp; i++) cout << i << " " << a.c[i] << endl;
329  if (m_HadronPulseShape) {
330  a.AddHit(hitE, hitTimeAve + timeOffset, m_ss_HadronShapeSimulations[2]); // diode component
331  } else {
332  a.AddHit(hitE, hitTimeAve + timeOffset, m_ss[1]); // m_ss[1] is the sampled diode response
333  }
334  // for (int i = 0; i < ec.m_nsmp; i++) cout << i << " " << a.c[i] << endl;
335  }
336  }
337 
338  if (m_calibration) {
339  // This has been added by Alex Bobrov for calibration
340  // of covariance matrix artificially generate 100 MeV in time for each crystal
341  double hitE = 0.1, hitTimeAve = 0.0;
342  for (int j = 0; j < ec.m_nch; j++) {
343  int cellId = j + 1; // 1 .. 8736
344  int id = m_eclMapper.getCrateID(cellId) - 1; // 0 .. 51
345  double timeOffset = tscale * m_ttime[id] - toff;
346  m_adc[j].AddHit(hitE, hitTimeAve + timeOffset, m_ss[m_tbl[j].iss]);
347  }
348  }
349 }
350 
352 {
354  float z[ec.m_nsmp], AdcNoise[ec.m_nsmp]; // buffers with electronic noise
355  // Noise generation
356  for (int i = 0; i < ec.m_nsmp; i++) z[i] = gRandom->Gaus(0, 1);
357  m_noise[m_tbl[J].inoise].generateCorrelatedNoise(z, AdcNoise);
358  for (int i = 0; i < ec.m_nsmp; i++) FitA[i] = 20 * AdcNoise[i] + 3000;
359 }
360 
362 {
364  BitStream out(ec.m_nch * ec.m_nsmp * 18 / 32);
365  out.putNBits(m_compAlgo & 0xff, 8);
366  ECLCompress* comp = selectAlgo(m_compAlgo & 0xff);
367  if (comp == nullptr)
368  B2FATAL("Unknown compression algorithm: " << m_compAlgo);
369 
370  int FitA[ec.m_nsmp]; // buffer for the waveform fitter
371  // loop over entire calorimeter
372  for (int j = 0; j < ec.m_nch; j++) {
373  adccounts_t& a = m_adc[j];
375  for (int i = 0; i < ec.m_nsmp; i++) {
376  int A = 20000 * a.c[i] + FitA[i];
377  FitA[i] = max(0, min(A, (1 << 18) - 1));
378  }
379  comp->compress(out, FitA);
380  }
381  out.resize();
382 
383  ECLWaveforms* wf = new ECLWaveforms;
384  m_eclWaveforms.assign(wf);
385 
386  std::swap(out.getStore(), wf->getStore());
387  if (comp) delete comp;
388 }
389 
391 {
393 
394  // make relation between cellid and eclhits
395  struct ch_t {int cell, id;};
396  vector<ch_t> hitmap;
397  for (const auto& hit : m_eclHits) {
398  int j = hit.getCellId() - 1; //0~8735
399  if (hit.getBackgroundTag() == BackgroundMetaData::bg_none) hitmap.push_back({j, hit.getArrayIndex()});
400  // cout<<"C:"<<hit.getBackgroundTag()<<" "<<hit.getCellId()<<" "<<hit.getEnergyDep()<<" "<<hit.getTimeAve()<<endl;
401  }
402 
403  bool isBGOverlay = m_eclWaveforms.isValid(), isTrigTime = false;
404  BitStream out;
405  ECLCompress* comp = nullptr;
406 
407  // check background overlay
408  if (isBGOverlay) {
409  std::swap(out.getStore(), m_eclWaveforms->getStore());
410  out.setPos(0);
411  unsigned int compAlgo = out.getNBits(8);
412  comp = selectAlgo(compAlgo & 0x7f);
413  if (comp == nullptr)
414  B2FATAL("Unknown compression algorithm: " << compAlgo);
415  isTrigTime = compAlgo >> 7; // crate trigger times are stored and retrived
416  if (isTrigTime) {
417  for (int i = 0; i < ECL::ECL_CRATES; i++) {
418  unsigned char t = out.getNBits(7); // [0..72)
419  m_ttime[i] = t;
420  }
421  }
422  }
423 
424  if (!m_trigTime || !isTrigTime) { // reproduce previous logic -- one trigger time for all crates
425  int DeltaT = gRandom->Uniform(0, double(ec.m_ntrg) * 0.5); // trigger decision can come in any time [0..72)
426  for (int id = 0; id < ECL::ECL_CRATES; id++) m_ttime[id] = DeltaT;
427  }
428 
429  int triggerTag0 = m_EventMetaData->getEvent();
430  for (int id = 0; id < ECL::ECL_CRATES; id++) {
431  auto eclTrig = m_eclTrigs.appendNew();
432  int triggerPhase0 = 2 * (m_ttime[id] + m_ttime[id] / 3);
433  eclTrig->setTrigId(id);
434  eclTrig->setTimeTrig(triggerPhase0);
435  eclTrig->setTrigTag(triggerTag0);
436  }
437 
438  shapeSignals();
439 
440  // We want to produce background waveforms in simulation first than
441  // dump to a disk, read from the disk to test before real data
442  if (m_waveformMaker) { makeWaveforms(); return; }
443 
444  int FitA[ec.m_nsmp]; // buffer for the waveform fitter
445 
446  // loop over entire calorimeter
447  for (int j = 0; j < ec.m_nch; j++) {
448  adccounts_t& a = m_adc[j];
449 
450  //normalize the MC true arrival times
451  if (m_adc[j].totalDep > 0) {
452  m_adc[j].flighttime /= m_adc[j].totalDep;
453  m_adc[j].timeshift /= m_adc[j].totalDep;
454  m_adc[j].timetosensor /= m_adc[j].totalDep;
455  }
456 
457  // if background waveform is here there is no need to generate
458  // electronic noise since it is already in the waveform
459  if (isBGOverlay) {
460  comp->uncompress(out, FitA);
461  } else {
462  // Signal amplitude should be above 100 keV
463  if (a.total < 0.0001) continue;
465  }
466 
467  for (int i = 0; i < ec.m_nsmp; i++) {
468  int A = 20000 * a.c[i] + FitA[i];
469  FitA[i] = max(0, min(A, (1 << 18) - 1));
470  }
471 
472  int energyFit = 0; // fit output : Amplitude 18 bits
473  int tFit = 0; // fit output : T_ave 12 bits
474  int qualityFit = 0; // fit output : quality 2 bits
475  int chi = 0; // fit output : chi square it is not available in the experiment
476 
477  int id = m_eclMapper.getCrateID(j + 1) - 1; // 0 .. 51
478  int ttrig = 2 * m_ttime[id];
479 
480  shapeFitterWrapper(j, FitA, ttrig, energyFit, tFit, qualityFit, chi);
481 
482  if (energyFit > m_ADCThreshold) {
483  int CellId = j + 1;
484 
485  //note energyFit and m_Awave is in ADC units
486  if (energyFit > m_Awave[CellId - 1]) {
487  //only save waveforms above ADC threshold
488  const auto eclDsp = m_eclDsps.appendNew();
489  eclDsp->setCellId(CellId);
490  eclDsp->setDspA(FitA);
491  }
492 
493  // only store extra MC info if requested and above threshold
495  const auto eclDspWithExtraMCInfo = m_eclDspsWithExtraMCInfo.appendNew();
496  eclDspWithExtraMCInfo->setCellId(CellId);
497  eclDspWithExtraMCInfo->setDspA(FitA);
498  eclDspWithExtraMCInfo->setEnergyDep(a.totalDep);
499  eclDspWithExtraMCInfo->setHadronEnergyDep(a.totalHadronDep);
500  eclDspWithExtraMCInfo->setFlightTime(a.flighttime);
501  eclDspWithExtraMCInfo->setTimeShift(a.timeshift);
502  eclDspWithExtraMCInfo->setTimeToSensor(a.timetosensor);
503  eclDspWithExtraMCInfo->setEnergyConversion(a.energyConversion * 20000);
504  }
505 
506  const auto eclDigit = m_eclDigits.appendNew();
507  eclDigit->setCellId(CellId); // cellId in range from 1 to 8736
508  eclDigit->setAmp(energyFit); // E (GeV) = energyFit/20000;
509  eclDigit->setTimeFit(tFit); // t0 (us)= (1520 - m_ltr)*24.*12/508/(3072/2) ;
510  eclDigit->setQuality(qualityFit);
511  if (qualityFit == 2)
512  eclDigit->setChi(chi);
513  else eclDigit->setChi(0);
514  for (const auto& hit : hitmap)
515  if (hit.cell == j) eclDigit->addRelationTo(m_eclHits[hit.id]);
516 
517  // set relation to DspWithExtraInfo
518  for (auto& DspWithExtraMCInfo : m_eclDspsWithExtraMCInfo) {
519  if (eclDigit->getCellId() == DspWithExtraMCInfo.getCellId()) DspWithExtraMCInfo.addRelationTo(eclDigit);
520  }
521  }
522  } //store each crystal hit
523  if (comp) delete comp;
524 }
525 
527 {
528 }
529 
531 {
532 }
533 
535 {
536 
537  if (m_waveformParametersMC.hasChanged()) {
538 
539  m_ss_HadronShapeSimulations.resize(3);
540 
541  //read MC templates from database
542  float photonParsPSD[10];
543  float hadronParsPSD[10];
544  float diodeParsPSD[10];
545  for (int i = 0; i < 10; i++) {
546  photonParsPSD[i] = m_waveformParametersMC->getPhotonParameters()[i + 1];
547  hadronParsPSD[i] = m_waveformParametersMC->getHadronParameters()[i + 1];
548  diodeParsPSD[i] = m_waveformParametersMC->getDiodeParameters()[i + 1];
549  }
550 
551  //Initialize templates for hadron shape simulations
552  m_ss_HadronShapeSimulations[0].InitSample(photonParsPSD, m_waveformParametersMC->getPhotonParameters()[0]);
553  m_ss_HadronShapeSimulations[1].InitSample(hadronParsPSD, m_waveformParametersMC->getHadronParameters()[0]);
554  m_ss_HadronShapeSimulations[2].InitSample(diodeParsPSD, m_waveformParametersMC->getDiodeParameters()[0]);
555 
556  }
557 
558 }
559 
561 {
563 
564  TFile* rootfile = nullptr;
565  TTree* tree = nullptr;
566  TTree* tree2 = nullptr;
567  TTree* tree3 = nullptr;
568 
570  bool hasChanged = false;
571  hasChanged |= m_waveformParameters.hasChanged();
572  hasChanged |= m_algoParameters.hasChanged();
573  hasChanged |= m_noiseParameters.hasChanged();
574 
575  if (!hasChanged) return;
576 
577  tree = const_cast<TTree*>(&*m_waveformParameters);
578  tree2 = const_cast<TTree*>(&*m_algoParameters);
579  tree3 = const_cast<TTree*>(&*m_noiseParameters);
580  } else {
581  string dataFileName;
582  if (m_background) {
583  dataFileName = FileSystem::findFile("/data/ecl/ECL-WF-BG.root");
584  B2DEBUG(150, "ECLDigitizer: Reading configuration data with background from: " << dataFileName);
585  } else {
586  dataFileName = FileSystem::findFile("/data/ecl/ECL-WF.root");
587  B2DEBUG(150, "ECLDigitizer: Reading configuration data without background from: " << dataFileName);
588  }
589  assert(! dataFileName.empty());
590 
591  rootfile = new TFile(dataFileName.c_str(), "read");
592  tree = (TTree*) rootfile->Get("EclWF");
593  tree2 = (TTree*) rootfile->Get("EclAlgo");
594  tree3 = (TTree*) rootfile->Get("EclNoise");
595  }
596 
597  if (tree == 0 || tree2 == 0 || tree3 == 0) B2FATAL("Data not found");
598 
599  m_tbl.resize(ec.m_nch);
600 
601  const int maxncellid = 512;
602  int ncellId;
603  vector<int> cellId(maxncellid);//[ncellId] buffer for crystal identification number
604 
605  tree->SetBranchAddress("ncellId", &ncellId);
606  tree->SetBranchAddress("cellId", cellId.data());
607 
608  vector<int> eclWaveformDataTable(ec.m_nch);
609  for (Long64_t j = 0, jmax = tree->GetEntries(); j < jmax; j++) {
610  tree->GetEntry(j);
611  assert(ncellId <= maxncellid);
612  for (int i = 0; i < ncellId; ++i)
613  eclWaveformDataTable[cellId[i] - 1] = j;
614  }
615  B2DEBUG(150, "ECLDigitizer: " << tree->GetEntries() << " sets of wave form covariance matricies will be used.");
616 
617  ECLWFAlgoParams* algo = new ECLWFAlgoParams;
618  tree2->SetBranchAddress("Algopars", &algo);
619  tree2->SetBranchAddress("ncellId", &ncellId);
620  tree2->SetBranchAddress("cellId", cellId.data());
621  Long64_t jmax2 = tree2->GetEntries();
622  vector<ECLWFAlgoParams> eclWFAlgoParams;
623  eclWFAlgoParams.reserve(jmax2);
624  for (Long64_t j = 0; j < jmax2; j++) {
625  tree2->GetEntry(j);
626  assert(ncellId <= maxncellid);
627  eclWFAlgoParams.push_back(*algo);
628  for (int i = 0; i < ncellId; ++i)
629  m_tbl[cellId[i] - 1].idn = j;
630  }
631  if (algo) delete algo;
632  B2DEBUG(150, "ECLDigitizer: " << eclWFAlgoParams.size() << " parameter sets of fitting algorithm were read.");
633 
634  ECLNoiseData* noise = new ECLNoiseData;
635  tree3->SetBranchAddress("NoiseM", &noise);
636  tree3->SetBranchAddress("ncellId", &ncellId);
637  tree3->SetBranchAddress("cellId", cellId.data());
638 
639  Long64_t jmax3 = tree3->GetEntries();
640  m_noise.reserve(jmax3);
641  for (Long64_t j = 0; j < jmax3; j++) {
642  tree3->GetEntry(j);
643  assert(ncellId <= maxncellid);
644  m_noise.push_back(*noise);
645  if (ncellId == 0) {
646  for (int i = 0; i < ec.m_nch; i++)
647  m_tbl[i].inoise = 0;
648  break;
649  } else {
650  for (int i = 0; i < ncellId; ++i)
651  m_tbl[cellId[i] - 1].inoise = j;
652  }
653  }
654  if (noise) delete noise;
655  B2DEBUG(150, "ECLDigitizer: " << eclWFAlgoParams.size() << " noise matricies were loaded.");
656 
657  // repack fitting algorithm parameters
658  m_idn.resize(eclWFAlgoParams.size());
659  for (int i = 0, imax = eclWFAlgoParams.size(); i < imax; i++)
660  repack(eclWFAlgoParams[i], m_idn[i]);
661 
662  vector<uint_pair_t> pairIdx;
663  for (int i = 0; i < ec.m_nch; i++) {
664  unsigned int wfIdx = eclWaveformDataTable[i];
665  unsigned int algoIdx = m_tbl[i].idn;
666  uint_pair_t p(wfIdx, algoIdx);
667  vector<uint_pair_t>::iterator ip = find(pairIdx.begin(), pairIdx.end(), p);
668  if (ip != pairIdx.end()) { // crystal i already have the same parameters as before
669  m_tbl[i].ifunc = ip - pairIdx.begin();
670  } else { // new combination of parameters
671  m_tbl[i].ifunc = pairIdx.size();
672  pairIdx.push_back(p);
673  }
674  }
675 
676  // now we know how many distinct elements of the (Waveform x AlgoParams) matrix should be
677  m_fitparams.resize(pairIdx.size());
678 
679  ECLWaveformData* eclWFData = new ECLWaveformData;
680  tree->SetBranchAddress("CovarianceM", &eclWFData);
681  tree->SetBranchStatus("ncellId", 0); // do not read it at the moment
682  tree->SetBranchStatus("cellId", 0); // do not read it at the moment
683 
684  // fill parameters for each (Waveform x AlgoParams) parameter set
685  for (unsigned int ip = 0; ip < pairIdx.size(); ip++) {
686  const uint_pair_t& p = pairIdx[ip];
687  tree->GetEntry(p.first); // retrieve data to eclWFData pointer
688  getfitparams(*eclWFData, eclWFAlgoParams[p.second], m_fitparams[ip]);
689  }
690  B2DEBUG(150, "ECLDigitizer: " << m_fitparams.size() << " fitting crystals groups were created.");
691 
692  // at the moment there is only one sampled signal shape in the pool
693  // since all shaper parameters are the same for all crystals
694  m_ss.resize(2);
695  float MP[10]; eclWFData->getWaveformParArray(MP);
696  m_ss[0].InitSample(MP, 27.7221);
697  // parameters vector from ps.dat file, time offset 0.5 usec added to
698  // have peak position with parameters from ps.dat roughly in the
699  // same place as in current MC
700  // 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};
701  // m_ss[0].InitSample(crystal_params, 0.9999272*19.5216);
702  for (int i = 0; i < ec.m_nch; i++) m_tbl[i].iss = 0;
703  // one sampled diode response in the pool, parameters vector from
704  // pg.dat file, time offset 0.5 usec added to have peak position with
705  // parameters from ps.dat roughly in the same place as in current MC
706  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};
707  m_ss[1].InitSample(diode_params, 0.9569100 * 9.98822);
708 
709  B2DEBUG(150, "ECLDigitizer: " << m_ss.size() << " sampled signal templates were created.");
710 
711  if (eclWFData) delete eclWFData;
712 
713  if (!m_useWaveformParameters) rootfile->Close();
714 }
715 
717 {
718  // filling short int array
719  t.id[ 0] = eclWFAlgo.getlAT() + 128;
720  t.id[ 1] = eclWFAlgo.getsT() + 128;
721  t.id[ 2] = eclWFAlgo.gethT();
722  t.id[ 3] = 17952;
723  t.id[ 4] = 19529;
724  t.id[ 5] = 69;
725  t.id[ 6] = 0;
726  t.id[ 7] = 0;
727  t.id[ 8] = 257;
728  t.id[ 9] = -1;
729  t.id[10] = 0;
730  t.id[11] = 0;
731 
732  // filling unsigned char array
733  t.ic[12 * 2 + 0] = eclWFAlgo.getk1();
734  t.ic[12 * 2 + 1] = eclWFAlgo.getk2();
735  t.ic[13 * 2 + 0] = eclWFAlgo.getka();
736  t.ic[13 * 2 + 1] = eclWFAlgo.getkb();
737  t.ic[14 * 2 + 0] = eclWFAlgo.getkc();
738  t.ic[14 * 2 + 1] = eclWFAlgo.gety0s() - 16;
739 
740  // again, filling short int array
741  t.id[15] = eclWFAlgo.getcT();
742 }
743 
745 {
747 
748  double ssd[16][16];
749  eclWFData.getMatrix(ssd);
750  vector<double> MP(10);
751  eclWFData.getWaveformParArray(MP.data());
752 
753  // shape function parameters
754  int iff = 1 << 14; //Alex Bobrov for correct chi square
755 
756  int ia = 1 << eclWFAlgo.getka();
757  int ib = 1 << eclWFAlgo.getkb();
758  int ic = 1 << eclWFAlgo.getkc();
759 
760  double dbl_f [192][16];
761  double dbl_f1 [192][16];
762  double dbl_fg31[192][16];
763  double dbl_fg32[192][16];
764  double dbl_fg33[192][16];
765  double dbl_fg41[24][16];
766  double dbl_fg43[24][16];
767 
768  int_array_192x16_t& ref_f = p.f;
769  int_array_192x16_t& ref_f1 = p.f1;
770  int_array_192x16_t& ref_fg31 = p.fg31;
771  int_array_192x16_t& ref_fg32 = p.fg32;
772  int_array_192x16_t& ref_fg33 = p.fg33;
773  int_array_24x16_t& ref_fg41 = p.fg41;
774  int_array_24x16_t& ref_fg43 = p.fg43;
775 
776  // Dynamic normalization coefficient for shape function
777  double fm = 0;
778 
779  double unitscale = 1.0;
780  if (m_unitscale > 0) unitscale = m_unitscale;
781 
782  ShaperDSP_t dsp(MP, unitscale);
783  dsp.settimestride(ec.m_step);
784  dsp.setseedoffset(ec.m_step / ec.m_ndt);
785  dsp.settimeseed(-(ec.m_step - (ec.m_step / ec.m_ndt)));
786  vector<dd_t> f(16);
787  for (int k = 0; k < 2 * ec.m_ndt; k++, dsp.nextseed()) { // calculate fit parameters around 0 +- 1 ADC tick
788  dsp.fillvector(f);
789 
790  double g0g0 = 0, g0g1 = 0, g1g1 = 0, g0g2 = 0, g1g2 = 0, g2g2 = 0;
791  double sg0[16], sg1[16], sg2[16];
792  for (int j = 0; j < 16; j++) {
793  double g0 = 0, g1 = 0, g2 = 0;
794  for (int i = 0; i < 16; i++) {
795  if (fm < f[i].first) fm = f[i].first;
796  g0 += ssd[j][i] * f[i].first;
797  g1 += ssd[j][i] * f[i].second;
798  g2 += ssd[j][i];
799  }
800  g0g0 += g0 * f[j].first;
801  g0g1 += g1 * f[j].first;
802  g1g1 += g1 * f[j].second;
803  g0g2 += g0;
804  g1g2 += g1;
805  g2g2 += g2;
806  sg0[j] = g0;
807  sg1[j] = g1;
808  sg2[j] = g2;
809  }
810 
811  double a00 = g1g1 * g2g2 - g1g2 * g1g2;
812  double a11 = g0g0 * g2g2 - g0g2 * g0g2;
813  double a22 = g0g0 * g1g1 - g0g1 * g0g1;
814  double a01 = g1g2 * g0g2 - g0g1 * g2g2;
815  double a02 = g0g1 * g1g2 - g1g1 * g0g2;
816  double a12 = g0g1 * g0g2 - g0g0 * g1g2;
817 
818  double igg2 = 1 / (a11 * g1g1 + g0g1 * a01 + g1g2 * a12);
819 
820  // to fixed point representation
821  const double isd = 3. / 4., sd = 1 / isd ; // conversion factor (???)
822  for (int i = 0; i < 16; i++) {
823  double w = i ? 1.0 : 1. / 16.;
824 
825  dbl_f [k][i] = (f[i].first * iff * w);
826  dbl_f1 [k][i] = (f[i].second * iff * w * sd);
827 
828  double fg31 = (a00 * sg0[i] + a01 * sg1[i] + a02 * sg2[i]) * igg2;
829  double fg32 = (a01 * sg0[i] + a11 * sg1[i] + a12 * sg2[i]) * igg2;
830  double fg33 = (a02 * sg0[i] + a12 * sg1[i] + a22 * sg2[i]) * igg2;
831 
832  dbl_fg31[k][i] = (fg31 * ia * w);
833  dbl_fg32[k][i] = (fg32 * ib * w * isd);
834  dbl_fg33[k][i] = (fg33 * ic * w);
835  }
836 
837  //first approximation without time correction
838  int jk = 23 + ((48 - k) >> 2);
839  if (jk >= 0 && jk < 24 && (48 - k) % 4 == 0) {
840 
841  double igg1 = 1 / a11;
842  // to fixed point
843  for (int i = 0; i < 16; i++) {
844  double w = i ? 1.0 : 1. / 16.;
845 
846  double fg41 = (g2g2 * sg0[i] - g0g2 * sg2[i]) * igg1;
847  double fg43 = (g0g0 * sg2[i] - g0g2 * sg0[i]) * igg1;
848  dbl_fg41[jk][i] = (fg41 * ia * w);
849  dbl_fg43[jk][i] = (fg43 * ic * w);
850  }
851  }
852  }
853  // Ignore dynamically calculated normalization coefficient if
854  // unitscale is set to a static positive value.
855  if (m_unitscale > 0) fm = 1.0;
856  for (int k = 0; k < 2 * ec.m_ndt; k++) {
857  for (int j = 0; j < 16; j++) {
858  ref_f [k][j] = lrint(dbl_f [k][j] / fm);
859  ref_f1 [k][j] = lrint(dbl_f1 [k][j] / fm);
860  ref_fg31[k][j] = lrint(dbl_fg31[k][j] * fm);
861  ref_fg32[k][j] = lrint(dbl_fg32[k][j] * fm);
862  ref_fg33[k][j] = lrint(dbl_fg33[k][j]);
863  if (k >= 24) continue;
864  ref_fg41[k][j] = lrint(dbl_fg41[k][j] * fm);
865  ref_fg43[k][j] = lrint(dbl_fg43[k][j]);
866  }
867  }
868 }
bool hasChanged()
Check whether the object has changed since the last call to hasChanged of the accessor).
int m_ADCThreshold
ADC threshold for wavefom fits.
std::vector< calibration_t > m_calib
Storage for calibration constants.
double m_unitscale
Normalization coefficient for ECL signal shape.
StoreArray< ECLDsp > m_eclDsps
Generated waveforms.
void callbackHadronSignalShapes()
callback hadron signal shapes from database
void shapeSignals()
Emulate response of energy deposition in a crystal and attached photodiode and make waveforms.
ECL::ECLChannelMapper m_eclMapper
Channel Mapper.
DBObjPtr< TTree > m_algoParameters
Shape fitting algorithm parameters.
StoreObjPtr< ECLWaveforms > m_eclWaveforms
Compressed waveforms.
std::vector< algoparams_t > m_idn
Fit algorihtm parameters shared by group of crystals.
StoreArray< ECLHit > m_eclDiodeHits
Diode hits array.
std::pair< unsigned int, unsigned int > uint_pair_t
a pair of unsigned ints
bool m_dspDataTest
DSP data usage flag.
std::vector< fitparams_t > m_fitparams
Pairs of (waveform parameters, fit parameters)
virtual void initialize() override
Initialize variables
virtual void event() override
Actual digitization of all hits in the ECL.
double m_WaveformThresholdOverride
If gt 0, value will override ECL_FPGA_StoreWaveform and apply value (in GeV) as threshold for all cry...
DBObjPtr< ECLCrystalCalib > m_CrateTimeOffset
Crate time offset.
DBObjPtr< TTree > m_waveformParameters
CellID-specific signal shapes.
double m_DspWithExtraMCInfoThreshold
Energy threshold above which to store DSPs with extra information.
StoreArray< ECLTrig > m_eclTrigs
Trigger information.
virtual void endRun() override
Nothing so far.
DBObjPtr< ECLCrystalCalib > m_MCTimeOffset
MC time offset.
virtual void terminate() override
Free memory.
StoreArray< ECLSimHit > m_eclSimHits
SimHits array.
bool m_background
Module parameters.
StoreArray< ECLDigit > m_eclDigits
Output Arrays.
void getfitparams(const ECLWaveformData &, const ECLWFAlgoParams &, fitparams_t &)
load waveform fit parameters for the shapeFitter function
void shapeFitterWrapper(const int j, const int *FitA, const int m_ttrig, int &m_lar, int &m_ltr, int &m_lq, int &m_chi) const
function wrapper for waveform fit
std::vector< double > m_Awave
Storage for waveform saving thresholds.
DBObjPtr< ECLCrystalCalib > m_CrystalElectronicsTime
Crystal electronics time.
std::vector< signalsample_t > m_ss_HadronShapeSimulations
tabulated shape line for hadron shape simulations
fitparams_t::int_array_192x16_t int_array_192x16_t
weighting coefficients for time and amplitude calculation
DBObjPtr< ECLCrystalCalib > m_CrystalElectronics
Crystal electronics.
fitparams_t::int_array_24x16_t int_array_24x16_t
weighting coefficients amplitude calculation.
virtual void beginRun() override
Nothing so far.
unsigned int m_compAlgo
compression algorithm for background waveforms
bool m_HadronPulseShape
hadron pulse shape flag
DBObjPtr< ECLCrystalCalib > m_CrystalEnergy
Crystal energy.
bool m_trigTime
Use trigger time from beam background overlay.
std::vector< crystallinks_t > m_tbl
Lookup table for ECL channels.
std::vector< adccounts_t > m_adc
Storage for adc hits from entire calorimeter (8736 crystals)
StoreArray< ECLDspWithExtraMCInfo > m_eclDspsWithExtraMCInfo
Generated waveforms with extra MC information.
StoreArray< ECLHit > m_eclHits
Hits array.
bool m_useWaveformParameters
If true, use m_waveformParameters, m_algoParameters, m_noiseParameters.
DBObjPtr< ECLCrystalCalib > m_FPGAWaveform
FPGA waveform.
bool m_calibration
calibration flag
bool m_storeDspWithExtraMCInfo
DSP with extra info flag.
unsigned char m_ttime[ECL::ECL_CRATES]
storage for trigger time in each ECL.
DBObjPtr< TTree > m_noiseParameters
Electronics noise covariance matrix.
bool m_loadOnce
Always load waveform parameters at least once.
std::vector< ECLNoiseData > m_noise
parameters for correlated noise simulation
DBObjPtr< ECLDigitWaveformParametersForMC > m_waveformParametersMC
Hadron signal shapes.
void repack(const ECLWFAlgoParams &, algoparams_t &)
repack waveform fit parameters from ROOT format to plain array of unsigned short for the shapeFitter ...
void makeWaveforms()
Produce and compress waveforms for beam background overlay.
bool m_waveformMaker
produce only waveform digits
StoreObjPtr< EventMetaData > m_EventMetaData
Event metadata.
void readDSPDB()
read Shaper-DSP data from root file
void makeElectronicNoiseAndPedestal(int j, int *FitA)
fill the waveform array FitA by electronic noise and bias it for channel J [0-8735]
std::string m_eclWaveformsName
name of background waveforms storage
DBObjPtr< ECLCrystalCalib > m_CrystalTimeOffset
Crystal time offset.
bool m_inter
internuclear counter effect
std::vector< signalsample_t > m_ss
tabulated shape line
Container for constant matrix used to generate electronic noise.
Container for constant parameters used in waveform fits.
int getkc() const
getter for multipliers power of 2 for fg33,fg43
int getlAT() const
getter for the threshold to calculate time
int getsT() const
getter for the threshold to send data to collector
int getcT() const
getter for the chi2 threshold for quality bit
int getk1() const
getter for the multipliers power of 2 for f
int getk2() const
getter for multipliers power of 2 for chi2 calculation
int getkb() const
getter for multipliers power of 2 for fg32
int gethT() const
getter for the hardware threshold
int gety0s() const
getter for the start point for pedestal calculation
int getka() const
getter for multipliers power of 2 for fg31 fg41
ECLWaveformData - container for inverse covariant matrix and shape parameters for time and amplitude ...
void getMatrix(float M[16][16]) const
Getter method for all matrix as two dimentional array (floats)
void getWaveformParArray(float P[10]) const
Getter method for waveform shape parameters as one dimentional array of floats.
Class to store ECL waveforms for entire calorimeter.
Definition: ECLWaveforms.h:21
std::vector< unsigned int > & getStore()
get data
Definition: ECLWaveforms.h:27
Bit stream struct.
Definition: BitStream.h:19
bool initFromDB()
Initialize channel mapper from the conditions database.
int getCrateID(int iCOPPERNode, int iFINESSE, bool pcie40=false)
Get crate number by given COPPER node number and FINESSE number.
Abstract class (interface) for ECL waveform compression/decompression to/from the BitStream storage.
Definition: ECLCompress.h:41
virtual void uncompress(BitStream &in, int *adc)=0
Decompress the ECL waveform.
virtual void compress(BitStream &out, const int *adc)=0
Compress the ECL waveform.
static ECLShapeFit shapeFitter(int cid, std::vector< int > adc, int ttrig, bool adjusted_timing=true)
Emulate shape fitting algorithm from ShaperDSP using algorithm from ecl/utility/src/ECLDspEmulator....
The Class for ECL Geometry Parameters.
static ECLGeometryPar * Instance()
Static method to get a reference to the ECLGeometryPar instance.
double time2sensor(int cid, const G4ThreeVector &hit_pos)
function to calculate flight time to diode sensor
Singleton class to hold the ECL configuration.
static EclConfiguration & get()
return this instance
static constexpr double m_step
time between points in internal units t_{asrto}*m_rf/2.
static constexpr double s_clock
digitization clock in RF units
static constexpr int m_nch
total number of electronic channels (crystals) in calorimeter
void setBackground(bool val)
set the background flag
static double getRF()
See m_rf.
static constexpr int m_ntrg
number of trigger counts per ADC clock tick
static constexpr int m_ndt
number of points per ADC tick where signal fit procedure parameters are evaluated
static constexpr int m_nsmp
number of ADC measurements for signal fitting
Class include function that calculate electronic response from energy deposit
Definition: shaperdsp.h:26
void fillvector(std::vector< double > &) const
fill vector with response function values and its derivative
Definition: shaperdsp.cc:417
void settimeseed(double)
set initial time
Definition: shaperdsp.cc:371
void settimestride(double)
set grid step for function calculation
Definition: shaperdsp.cc:361
void nextseed()
substruct toffset to tzero
Definition: shaperdsp.cc:376
void setseedoffset(double)
set timeoffset
Definition: shaperdsp.cc:366
static std::string findFile(const std::string &path, bool silent=false)
Search for given file or directory in local or central release directory, and return absolute path if...
Definition: FileSystem.cc:148
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 setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
static const double us
[microsecond]
Definition: Unit.h:97
static const double GeV
Standard of [energy, momentum, mass].
Definition: Unit.h:51
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
const int c_NCrystals
Number of crystals.
Abstract base class for different kinds of events.
calibration constants per channel
ShaperDSP fit results from _lftda function.
a struct for the parameters of the algorithm
a struct for the fit parameters