10#include <ecl/modules/eclDigitizer/ECLDigitizerModule.h>
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>
28#include <framework/gearbox/Unit.h>
29#include <framework/logging/Logger.h>
30#include <framework/utilities/FileSystem.h>
51 m_waveformParameters(
"ECLWFParameters"),
52 m_algoParameters(
"ECLWFAlgoParams"),
53 m_noiseParameters(
"ECLWFNoiseParams")
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)",
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);
68 addParam(
"HadronPulseShapes",
m_HadronPulseShape,
"Flag to include hadron component in pulse shape construction (default: true)",
72 "If gt 0 value is applied to all crystals for waveform saving threshold. If lt 0 dbobject is used. (GeV)", -1.0);
74 "Flag to store Dsp with extra MC information in addition to normal Dsp (default: false)",
false);
76 "Threshold above with to store Dsp with extra MC information [GeV]",
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);
83 "Use ECLWF{Parameters,AlgoParams,NoiseParams} payloads",
true);
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);
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");
106 B2DEBUG(20,
"Hadron pulse shapes for ECL simulations are disabled.");
110 B2DEBUG(20,
"Diode-crossing pulse shapes for ECL simulations are enabled.");
112 B2DEBUG(20,
"Diode-crossing pulse shapes for ECL simulations are disabled.");
132 double ns_per_tick = 1.0 / (4.0 * ec.
getRF()) * 1e3;
184 B2FATAL(
"ECLDigitizer: Can't initialize eclChannelMapper!");
190 int& m_lar,
int& m_ltr,
int& m_lq,
int& m_chi)
const
198 short int*
id = (
short int*)
m_idn[t.idn].id;
200 int A0 = (
int) * (
id + 0) - 128;
201 int Askip = (int) * (
id + 1) - 128;
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);
211 int chi_thres = (int) * (
id + 15);
213 int trg_time = ttrig;
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,
220 std::vector<int> adc(31);
221 for (
int i = 0; i < 31; i++) adc[i] = FitA[i];
226 int ttrig_packed = ttrig / 6 * 8 + ttrig % 6;
232 m_lq = result.quality;
236 int discarded_bits = 0;
237 if ((m_chi & 0x7800000) != 0) {
239 }
else if ((m_chi & 0x0600000) != 0) {
241 }
else if ((m_chi & 0x0180000) != 0) {
243 }
else if ((m_chi & 0x0060000) != 0) {
245 }
else if ((m_chi & 0x0018000) != 0) {
247 }
else if ((m_chi & 0x0006000) != 0) {
249 }
else if ((m_chi & 0x0001800) != 0) {
251 }
else if ((m_chi & 0x0000600) != 0) {
254 if (discarded_bits > 0) {
255 m_chi >>= discarded_bits;
256 m_chi <<= discarded_bits;
266 const double tscale = 2 * trgtick, toff = ec.
s_clock / (2 * ec.
getRF());
276 int cellId = hit.getCellId();
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;
283 m_adc[j].energyConversion =
m_calib[j].ascale * E2GeV * 20000;
284 m_adc[j].flighttime += hit.getFlightTime() * hit.getEnergyDep();
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();
288 m_adc[j].totalDep += hit.getEnergyDep();
291 double hitHadronE = hit.getHadronEnergyDep() *
m_calib[j].ascale * E2GeV;
295 m_adc[j].AddHit(hitE, hitTimeAve + timeOffset,
m_ss[
m_tbl[j].iss]);
302 int cellId = hit.getCellId();
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]);
316 const double diodeEdep2crystalEdep = E2GeV * (1 / (5000 * 3.6e-6));
318 int cellId = hit.getCellId();
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;
331 a.AddHit(hitE, hitTimeAve + timeOffset,
m_ss[1]);
340 double hitE = 0.1, hitTimeAve = 0.0;
341 for (
int j = 0; j < ec.
m_nch; j++) {
344 double timeOffset = tscale *
m_ttime[id] - toff;
345 m_adc[j].AddHit(hitE, hitTimeAve + timeOffset,
m_ss[
m_tbl[j].iss]);
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;
367 B2FATAL(
"Unknown compression algorithm: " <<
m_compAlgo);
371 for (
int j = 0; j < ec.
m_nch; 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));
385 std::swap(out.getStore(), wf->
getStore());
386 if (comp)
delete comp;
394 struct ch_t {
int cell, id;};
397 int j = hit.getCellId() - 1;
410 unsigned int compAlgo = out.getNBits(8);
411 comp = selectAlgo(compAlgo & 0x7f);
413 B2FATAL(
"Unknown compression algorithm: " << compAlgo);
414 isTrigTime = compAlgo >> 7;
416 for (
int i = 0; i < ECL::ECL_CRATES; i++) {
417 unsigned char t = out.getNBits(7);
424 int DeltaT = gRandom->Uniform(0,
double(ec.
m_ntrg) * 0.5);
425 for (
int id = 0;
id < ECL::ECL_CRATES;
id++)
m_ttime[
id] = DeltaT;
429 for (
int id = 0;
id < ECL::ECL_CRATES;
id++) {
432 eclTrig->setTrigId(
id);
433 eclTrig->setTimeTrig(triggerPhase0);
434 eclTrig->setTrigTag(triggerTag0);
446 for (
int j = 0; j < ec.
m_nch; j++) {
450 if (
m_adc[j].totalDep > 0) {
462 if (a.total < 0.0001)
continue;
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));
485 if (energyFit >
m_Awave[CellId - 1]) {
487 const auto eclDsp =
m_eclDsps.appendNew();
488 eclDsp->setCellId(CellId);
489 eclDsp->setDspA(FitA);
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);
506 eclDigit->setCellId(CellId);
507 eclDigit->setAmp(energyFit);
508 eclDigit->setTimeFit(tFit);
509 eclDigit->setQuality(qualityFit);
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]);
518 if (eclDigit->getCellId() == DspWithExtraMCInfo.getCellId()) DspWithExtraMCInfo.addRelationTo(eclDigit);
522 if (comp)
delete comp;
541 float photonParsPSD[10];
542 float hadronParsPSD[10];
543 float diodeParsPSD[10];
544 for (
int i = 0; i < 10; i++) {
563 TFile* rootfile =
nullptr;
564 TTree* tree =
nullptr;
565 TTree* tree2 =
nullptr;
566 TTree* tree3 =
nullptr;
569 bool hasChanged =
false;
574 if (!hasChanged)
return;
583 B2DEBUG(150,
"ECLDigitizer: Reading configuration data with background from: " << dataFileName);
586 B2DEBUG(150,
"ECLDigitizer: Reading configuration data without background from: " << dataFileName);
588 assert(! dataFileName.empty());
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");
596 if (tree == 0 || tree2 == 0 || tree3 == 0) B2FATAL(
"Data not found");
600 const int maxncellid = 512;
602 vector<int> cellId(maxncellid);
604 tree->SetBranchAddress(
"ncellId", &ncellId);
605 tree->SetBranchAddress(
"cellId", cellId.data());
607 vector<int> eclWaveformDataTable(ec.
m_nch);
608 for (Long64_t j = 0, jmax = tree->GetEntries(); j < jmax; j++) {
610 assert(ncellId <= maxncellid);
611 for (
int i = 0; i < ncellId; ++i)
612 eclWaveformDataTable[cellId[i] - 1] = j;
614 B2DEBUG(150,
"ECLDigitizer: " << tree->GetEntries() <<
" sets of wave form covariance matrices will be used.");
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++) {
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;
630 if (algo)
delete algo;
631 B2DEBUG(150,
"ECLDigitizer: " << eclWFAlgoParams.size() <<
" parameter sets of fitting algorithm were read.");
634 tree3->SetBranchAddress(
"NoiseM", &noise);
635 tree3->SetBranchAddress(
"ncellId", &ncellId);
636 tree3->SetBranchAddress(
"cellId", cellId.data());
638 Long64_t jmax3 = tree3->GetEntries();
640 for (Long64_t j = 0; j < jmax3; j++) {
642 assert(ncellId <= maxncellid);
645 for (
int i = 0; i < ec.
m_nch; i++)
649 for (
int i = 0; i < ncellId; ++i)
650 m_tbl[cellId[i] - 1].inoise = j;
653 if (noise)
delete noise;
654 B2DEBUG(150,
"ECLDigitizer: " << eclWFAlgoParams.size() <<
" noise matrices were loaded.");
657 m_idn.resize(eclWFAlgoParams.size());
658 for (
int i = 0, imax = eclWFAlgoParams.size(); i < imax; i++)
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;
666 vector<uint_pair_t>::iterator ip = find(pairIdx.begin(), pairIdx.end(), p);
667 if (ip != pairIdx.end()) {
668 m_tbl[i].ifunc = ip - pairIdx.begin();
670 m_tbl[i].ifunc = pairIdx.size();
671 pairIdx.push_back(p);
679 tree->SetBranchAddress(
"CovarianceM", &eclWFData);
680 tree->SetBranchStatus(
"ncellId", 0);
681 tree->SetBranchStatus(
"cellId", 0);
684 for (
unsigned int ip = 0; ip < pairIdx.size(); ip++) {
686 tree->GetEntry(p.first);
689 B2DEBUG(150,
"ECLDigitizer: " <<
m_fitparams.size() <<
" fitting crystals groups were created.");
695 m_ss[0].InitSample(MP, 27.7221);
701 for (
int i = 0; i < ec.
m_nch; i++)
m_tbl[i].iss = 0;
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);
708 B2DEBUG(150,
"ECLDigitizer: " <<
m_ss.size() <<
" sampled signal templates were created.");
710 if (eclWFData)
delete eclWFData;
718 t.id[ 0] = eclWFAlgo.
getlAT() + 128;
719 t.id[ 1] = eclWFAlgo.
getsT() + 128;
720 t.id[ 2] = eclWFAlgo.
gethT();
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;
740 t.id[15] = eclWFAlgo.
getcT();
749 vector<double> MP(10);
755 int ia = 1 << eclWFAlgo.
getka();
756 int ib = 1 << eclWFAlgo.
getkb();
757 int ic = 1 << eclWFAlgo.
getkc();
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];
778 double unitscale = 1.0;
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;
799 g0g0 += g0 * f[j].first;
800 g0g1 += g1 * f[j].first;
801 g1g1 += g1 * f[j].second;
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;
817 double igg2 = 1 / (a11 * g1g1 + g0g1 * a01 + g1g2 * a12);
820 const double isd = 3. / 4., sd = 1 / isd ;
821 for (
int i = 0; i < 16; i++) {
822 double w = i ? 1.0 : 1. / 16.;
824 dbl_f [k][i] = (f[i].first * iff * w);
825 dbl_f1 [k][i] = (f[i].second * iff * w * sd);
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;
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);
837 int jk = 23 + ((48 - k) >> 2);
838 if (jk >= 0 && jk < 24 && (48 - k) % 4 == 0) {
840 double igg1 = 1 / a11;
842 for (
int i = 0; i < 16; i++) {
843 double w = i ? 1.0 : 1. / 16.;
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);
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]);
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 algorithm 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)
ECLDigitizerModule()
Constructor.
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
~ECLDigitizerModule()
Destructor.
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
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.
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
void fillvector(std::vector< double > &) const
fill vector with response function values and its derivative
void settimeseed(double)
set initial time
void settimestride(double)
set grid step for function calculation
void nextseed()
substruct toffset to tzero
void setseedoffset(double)
set timeoffset
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...
void setDescription(const std::string &description)
Sets the description of the module.
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
static const double us
[microsecond]
static const double GeV
Standard of [energy, momentum, mass].
void addParam(const std::string &name, T ¶mVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
const int c_NCrystals
Number of crystals.
Abstract base class for different kinds of events.
calibration constants per channel
Indices in arrays with info on ECL channels.
ShaperDSP fit results from _lftda function.
a struct for the ADC count
a struct for the parameters of the algorithm
a struct for the fit parameters