10#include <ecl/utility/ECLDspUtilities.h>
13#include <ecl/dbobjects/ECLCrystalCalib.h>
14#include <ecl/dbobjects/ECLDspData.h>
15#include <ecl/mapper/ECLChannelMapper.h>
18#include <framework/database/DBArray.h>
19#include <framework/database/DBObjPtr.h>
20#include <framework/logging/Logger.h>
21#include <framework/utilities/FileSystem.h>
31const short ECLDSP_FORMAT_VERSION = 1;
41void readEclDspCoefs(FILE* fl,
void* ptr,
int word_size,
int word_count)
43 int read_items = fread(ptr, word_size, word_count, fl);
44 if (read_items != word_count) {
45 B2ERROR(
"Error reading DSP coefficients");
52 fl = fopen(filename,
"rb");
54 B2ERROR(
"Can't open file " << filename);
60 int crateNum = (boardNumber - 1) / 12 + 1;
61 int shaperNum = (boardNumber - 1) % 12;
68 if (shaperNum >= 10) {
69 if (crateNum >= 37 && crateNum <= 44) {
81 int size = fread(
id, nsiz, nsiz1, fl);
83 B2ERROR(
"Error reading header of DSP file");
86 std::vector<short int> extraData;
87 extraData.push_back(ECLDSP_FORMAT_VERSION);
88 extraData.push_back(data->getPackerVersion());
89 data->setExtraData(extraData);
91 data->setverMaj(
id[8] & 0xFF);
92 data->setverMin(
id[8] >> 8);
93 data->setkb(
id[13] >> 8);
94 data->setka(
id[13] - 256 * data->getkb());
95 data->sety0Startr(
id[14] >> 8);
96 data->setkc(
id[14] - 256 * data->gety0Startr());
97 data->setchiThresh(
id[15]);
98 data->setk2(
id[16] >> 8);
99 data->setk1(
id[16] - 256 * data->getk2());
101 data->setlAT(
id[128]);
102 data->setsT(
id[192]);
103 data->setaAT(
id[208]);
105 std::vector<short int> f(49152), f1(49152), f31(49152),
106 f32(49152), f33(49152), f41(6144), f43(6144);
108 for (
int i = 0; i < 16; ++i) {
110 readEclDspCoefs(fl, &(*(f41.begin() + i * nsiz1)), nsiz, nsiz1);
112 readEclDspCoefs(fl, &(*(f31.begin() + i * nsiz1)), nsiz, nsiz1);
113 readEclDspCoefs(fl, &(*(f32.begin() + i * nsiz1)), nsiz, nsiz1);
114 readEclDspCoefs(fl, &(*(f33.begin() + i * nsiz1)), nsiz, nsiz1);
116 readEclDspCoefs(fl, &(*(f43.begin() + i * nsiz1)), nsiz, nsiz1);
118 readEclDspCoefs(fl, &(*(f.begin() + i * nsiz1)), nsiz, nsiz1);
119 readEclDspCoefs(fl, &(*(f1.begin() + i * nsiz1)), nsiz, nsiz1);
138 const unsigned short DEFAULT_HEADER[256] {
140 0x4543, 0x4c44, 0x5350, 0x2046, 0x494c, 0x4500, 0x0000, 0x0000,
141 0xABCD, 0xffff, 0x0000, 0x0000, 0x0000, 0xABCD, 0xABCD, 0xABCD,
142 0xABCD, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
143 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
144 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
145 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
146 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
147 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
148 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD,
149 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD,
150 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
151 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
152 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
153 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
154 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
155 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
156 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD,
157 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD,
158 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
159 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
160 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
161 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
162 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
163 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
164 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD,
165 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD,
166 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD,
167 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD, 0xABCD,
168 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
169 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
170 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
171 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000, 0x0000,
175 fl = fopen(filename,
"wb");
181 unsigned short header[256];
183 int format_version = data->getExtraData()[0];
184 if (format_version > ECLDSP_FORMAT_VERSION) {
185 B2WARNING(
"Format version " << format_version <<
" is not fully supported, some data might be discarded.");
188 for (
int i = 0; i < 256; i++) {
190 header[i] = (data->getverMin() << 8) | data->getverMaj();
191 }
else if (i == 13) {
192 header[i] = (data->getkb() << 8) | data->getka();
193 }
else if (i == 14) {
194 header[i] = (data->gety0Startr() << 8) | data->getkc();
195 }
else if (i == 15) {
196 header[i] = data->getchiThresh();
197 }
else if (i == 16) {
198 header[i] = (data->getk2() << 8) | data->getk1();
199 }
else if (i >= 64 && i < 80) {
200 header[i] = data->gethT();
201 }
else if (i >= 128 && i < 144) {
202 header[i] = data->getlAT();
203 }
else if (i >= 192 && i < 208) {
204 header[i] = data->getsT();
205 }
else if (i >= 208 && i < 224) {
206 header[i] = data->getaAT();
209 int high = (DEFAULT_HEADER[i] & 0xFF00) >> 8;
210 int low = (DEFAULT_HEADER[i] & 0x00FF);
211 header[i] = (low << 8) + high;
215 int size = fwrite(header, nsiz, nsiz1, fl);
217 B2FATAL(
"Error writing header of DSP file " << filename);
220 std::vector<short int> f(49152), f1(49152), f31(49152),
221 f32(49152), f33(49152), f41(6144), f43(6144);
230 for (
int i = 0; i < 16; ++i) {
232 fwrite(&(*(f41.begin() + i * nsiz1)), nsiz, nsiz1, fl);
234 fwrite(&(*(f31.begin() + i * nsiz1)), nsiz, nsiz1, fl);
235 fwrite(&(*(f32.begin() + i * nsiz1)), nsiz, nsiz1, fl);
236 fwrite(&(*(f33.begin() + i * nsiz1)), nsiz, nsiz1, fl);
238 fwrite(&(*(f43.begin() + i * nsiz1)), nsiz, nsiz1, fl);
240 fwrite(&(*(f.begin() + i * nsiz1)), nsiz, nsiz1, fl);
241 fwrite(&(*(f1.begin() + i * nsiz1)), nsiz, nsiz1, fl);
249short int* vectorsplit(std::vector<short int>& vectorFrom,
int channel)
251 size_t size = vectorFrom.size();
252 if (size % 16) B2ERROR(
"Split is impossible!" <<
"Vector size" << size);
253 return (vectorFrom.data() + (size / 16) * (channel - 1));
257 bool adjusted_timing)
270 int dsp_group = (crate - 1) / 18;
271 int dsp_id = (crate - 1) % 18;
272 std::string payload_name =
"ECLDSPPars" + std::to_string(dsp_group);
274 const ECLDspData* data = dsp_data[dsp_id * 12 + shaper - 1];
276 std::vector<short> vec_f; data->getF(vec_f);
277 std::vector<short> vec_f1; data->getF1(vec_f1);
278 std::vector<short> vec_fg31; data->getF31(vec_fg31);
279 std::vector<short> vec_fg32; data->getF32(vec_fg32);
280 std::vector<short> vec_fg33; data->getF33(vec_fg33);
281 std::vector<short> vec_fg41; data->getF41(vec_fg41);
282 std::vector<short> vec_fg43; data->getF43(vec_fg43);
284 short* f = vectorsplit(vec_f, channel);
285 short* f1 = vectorsplit(vec_f1, channel);
286 short* fg31 = vectorsplit(vec_fg31, channel);
287 short* fg32 = vectorsplit(vec_fg32, channel);
288 short* fg33 = vectorsplit(vec_fg33, channel);
289 short* fg41 = vectorsplit(vec_fg41, channel);
290 short* fg43 = vectorsplit(vec_fg43, channel);
292 int k_a = data->getka();
293 int k_b = data->getkb();
294 int k_c = data->getkc();
295 int k_1 = data->getk1();
296 int k_2 = data->getk2();
297 int k_16 = data->gety0Startr();
298 int chi_thres = data->getchiThresh();
304 int ttrig2 = ttrig - 2 * (ttrig / 8);
311 int A0 = thr_LowAmp->getCalibVector()[cid - 1];
312 int Ahard = thr_HitThresh->getCalibVector()[cid - 1];
313 int Askip = thr_StoreDigit->getCalibVector()[cid - 1];
316 auto result = lftda_(f, f1, fg41, fg43, fg31, fg32, fg33, y, ttrig2, A0,
317 Ahard, Askip, k_a, k_b, k_c, k_16, k_1, k_2,
318 chi_thres, adjusted_timing);
326 TFile* file =
new TFile(path.c_str(),
"read");
327 if (!file->IsOpen()) {
328 B2FATAL(
"Unable to load coefficients for ECL pedestal fit");
330 TTree* tree = (TTree*)file->Get(
"dsp_coefs");
331 int nentries = tree->GetEntries();
332 float fg31_i, fg32_i;
333 tree->SetBranchAddress(
"fg31", &fg31_i);
334 tree->SetBranchAddress(
"fg32", &fg32_i);
336 for (
int i = 0; i < nentries; i++) {
359 for (
int i = 0; i < 16; i++) {
360 if (adc[i] > ped_max) {
367 if (ped_max_pos < 2 || ped_max_pos > 13) {
374 int time_index = ped_max_pos * 4 - 8;
377 for (
int iter = 0; iter < 2; iter++) {
380 for (
int j = 0; j < 16; j++) {
385 time_index -= tim * 4;
386 if (time_index > 47) time_index = 47;
387 if (time_index < 0) time_index = 0;
Class for accessing arrays of objects in the database.
Class for accessing objects in the database.
This object contains ECL DSP coefs – electromagnetic calorimeter digital signal processing coefficien...
This class provides access to ECL channel map that is either a) Loaded from the database (see ecl/dbo...
bool initFromDB()
Initialize channel mapper from the conditions database.
int getShaperChannel(int cellID)
Get number of DSP channel in the shaper by given number of CellId.
int getShaperPosition(int cellID)
Get position of the shaper in the crate by given CellId.
int getCrateID(int iCOPPERNode, int iFINESSE, bool pcie40=false)
Get crate number by given COPPER node number and FINESSE number.
static ECLDspData * readEclDsp(const char *raw_file, int boardNumber)
Convert ECLDspData from *.dat file to Root object.
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....
static int pedestal_fit_initialized
Flag indicating whether arrays fg31,fg32 are filled.
static ECLPedestalFit pedestalFit(std::vector< int > adc)
Fit pedestal part of the signal waveform (first 16 samples) This method will fit the first 16 samples...
static void writeEclDsp(const char *raw_file, ECLDspData *obj)
Convert ECLDspData from Root object to *.dat file.
static float pedfit_fg32[768]
DSP coefficients used to determine time in pedestalFit.
static float pedfit_fg31[768]
DSP coefficients used to determine amplitude in pedestalFit.
static void initPedestalFit()
Load DSP coefficients used in the pedestal fit function.
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...
Abstract base class for different kinds of events.
This struct is returned by the pedestalFit method that fits the first 16 samples of the waveform (ped...
ShaperDSP fit results from _lftda function.