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