Belle II Software development
ECLDspUtilities Class Reference

This class contains static methods to make them accessible from pyROOT. More...

#include <ECLDspUtilities.h>

Static Public Member Functions

static ECLDspDatareadEclDsp (const char *raw_file, int boardNumber)
 Convert ECLDspData from *.dat file to Root object.
 
static void writeEclDsp (const char *raw_file, ECLDspData *obj)
 Convert ECLDspData from Root object to *.dat file.
 
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.cc See ecl/examples/eclShapeFitter.py for usage example.
 
static void initPedestalFit ()
 Load DSP coefficients used in the pedestal fit function.
 
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 of the waveform and return the amplitude of the peak found in that region.
 

Private Member Functions

 ECLDspUtilities ()
 Private constructor since class only contains static methods, no need to create an instance.
 

Static Private Attributes

static int pedestal_fit_initialized = 0
 Flag indicating whether arrays fg31,fg32 are filled.
 
static float pedfit_fg31 [768] = {}
 DSP coefficients used to determine amplitude in pedestalFit.
 
static float pedfit_fg32 [768] = {}
 DSP coefficients used to determine time in pedestalFit.
 

Detailed Description

This class contains static methods to make them accessible from pyROOT.

Definition at line 41 of file ECLDspUtilities.h.

Constructor & Destructor Documentation

◆ ECLDspUtilities()

ECLDspUtilities ( )
inlineprivate

Private constructor since class only contains static methods, no need to create an instance.

Definition at line 111 of file ECLDspUtilities.h.

111{}

Member Function Documentation

◆ initPedestalFit()

void initPedestalFit ( )
static

Load DSP coefficients used in the pedestal fit function.

If it is not done explicitly, pedestalFit will do it internally when it is called the first time.

However, it is preferable to call it explicitly, in a thread-safe context.

Definition at line 323 of file ECLDspUtilities.cc.

324{
325 std::string path = FileSystem::findFile("ecl/data/ecl_pedestal_peak_fit.root");
326 TFile* file = new TFile(path.c_str(), "read");
327 if (!file->IsOpen()) {
328 B2FATAL("Unable to load coefficients for ECL pedestal fit");
329 }
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);
335 //== Load DSP coefficients used in pedestal fitting
336 for (int i = 0; i < nentries; i++) {
337 tree->GetEntry(i);
338 pedfit_fg31[i] = fg31_i;
339 pedfit_fg32[i] = fg32_i;
340 }
341 file->Close();
342
344}
static int pedestal_fit_initialized
Flag indicating whether arrays fg31,fg32 are filled.
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 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

◆ pedestalFit()

ECLPedestalFit pedestalFit ( std::vector< int >  adc)
static

Fit pedestal part of the signal waveform (first 16 samples) This method will fit the first 16 samples of the waveform and return the amplitude of the peak found in that region.

Parameters
[in]adcvector of waveform samples (size >= 16)
Returns
struct with fit results

Definition at line 346 of file ECLDspUtilities.cc.

347{
350 }
351
352 ECLPedestalFit result;
353 float amp;
354
355 //== Find maximum in the pedestal
356
357 int ped_max = 0;
358 int ped_max_pos = 0;
359 for (int i = 0; i < 16; i++) {
360 if (adc[i] > ped_max) {
361 ped_max = adc[i];
362 ped_max_pos = i;
363 }
364 }
365
366 // Skip fit if the peak position is on the edge of the fit region
367 if (ped_max_pos < 2 || ped_max_pos > 13) {
368 result.amp = -128;
369 return result;
370 }
371
372 //== Get first position estimate from maximum location
373
374 int time_index = ped_max_pos * 4 - 8;
375
376 //== Run two iterations of chi2 minimization
377 for (int iter = 0; iter < 2; iter++) {
378 amp = 0;
379 float tim = 0;
380 for (int j = 0; j < 16; j++) {
381 amp += pedfit_fg31[j + time_index * 16] * adc[j];
382 tim += pedfit_fg32[j + time_index * 16] * adc[j];
383 }
384 tim = tim / amp;
385 time_index -= tim * 4;
386 if (time_index > 47) time_index = 47;
387 if (time_index < 0) time_index = 0;
388 }
389
390 result.amp = amp;
391
392 return result;
393}
static void initPedestalFit()
Load DSP coefficients used in the pedestal fit function.
This struct is returned by the pedestalFit method that fits the first 16 samples of the waveform (ped...

◆ readEclDsp()

ECLDspData * readEclDsp ( const char *  raw_file,
int  boardNumber 
)
static

Convert ECLDspData from *.dat file to Root object.

Parameters
[in]raw_filePath to dsp??.dat file.
[in]boardNumberNumber of shaperDSP board, from 1 to 52*12
Returns
ECLDspData object

Definition at line 49 of file ECLDspUtilities.cc.

50{
51 FILE* fl;
52 fl = fopen(filename, "rb");
53 if (fl == nullptr) {
54 B2ERROR("Can't open file " << filename);
55 }
56
57 ECLDspData* data = new ECLDspData(boardNumber);
58
59 //== Do not fill data for non-existent shapers
60 int crateNum = (boardNumber - 1) / 12 + 1;
61 int shaperNum = (boardNumber - 1) % 12;
62
63 if (shaperNum >= 8) {
64 if (crateNum >= 45) {
65 fclose(fl);
66 return data;
67 }
68 if (shaperNum >= 10) {
69 if (crateNum >= 37 && crateNum <= 44) {
70 fclose(fl);
71 return data;
72 }
73 }
74 }
75
76 // Word size
77 const int nsiz = 2;
78 // Size of fragment to read (in words)
79 int nsiz1 = 256;
80 short int id[256];
81 int size = fread(id, nsiz, nsiz1, fl);
82 if (size != nsiz1) {
83 B2ERROR("Error reading header of DSP file");
84 }
85
86 std::vector<short int> extraData;
87 extraData.push_back(ECLDSP_FORMAT_VERSION);
88 extraData.push_back(data->getPackerVersion());
89 data->setExtraData(extraData);
90
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());
100 data->sethT(id[64]);
101 data->setlAT(id[128]);
102 data->setsT(id[192]);
103 data->setaAT(id[208]);
104
105 std::vector<short int> f(49152), f1(49152), f31(49152),
106 f32(49152), f33(49152), f41(6144), f43(6144);
107
108 for (int i = 0; i < 16; ++i) {
109 nsiz1 = 384;
110 readEclDspCoefs(fl, &(*(f41.begin() + i * nsiz1)), nsiz, nsiz1);
111 nsiz1 = 3072;
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);
115 nsiz1 = 384;
116 readEclDspCoefs(fl, &(*(f43.begin() + i * nsiz1)), nsiz, nsiz1);
117 nsiz1 = 3072;
118 readEclDspCoefs(fl, &(*(f.begin() + i * nsiz1)), nsiz, nsiz1);
119 readEclDspCoefs(fl, &(*(f1.begin() + i * nsiz1)), nsiz, nsiz1);
120 }
121 fclose(fl);
122
123 data->setF41(f41);
124 data->setF31(f31);
125 data->setF32(f32);
126 data->setF33(f33);
127 data->setF43(f43);
128 data->setF(f);
129 data->setF1(f1);
130
131 return data;
132}
This object contains ECL DSP coefs – electromagnetic calorimeter digital signal processing coefficien...
Definition: ECLDspData.h:33

◆ shapeFitter()

ECLShapeFit shapeFitter ( int  cid,
std::vector< int >  adc,
int  ttrig,
bool  adjusted_timing = true 
)
static

Emulate shape fitting algorithm from ShaperDSP using algorithm from ecl/utility/src/ECLDspEmulator.cc See ecl/examples/eclShapeFitter.py for usage example.

Parameters
[in]cidCellID, 1..8736
[in]adcWaveform data from ECLDsp dataobject of length 31
[in]ttrigTrigger time from ECLTrig dataobject
[in]adjusted_timingOptional. Use adjusted formula to determine fit time. Implemented in ShaperDSP firmware since exp 14. If true, algorithm will determine time near 0 with higher precision, time of low-energy hits will be one of {-4,0,4} If false, time will be one of {-32, -16, 0}

Definition at line 256 of file ECLDspUtilities.cc.

258{
259 ECLChannelMapper mapper;
260 mapper.initFromDB();
261
262 //== Get mapping data
263
264 int crate = mapper.getCrateID(cid);
265 int shaper = mapper.getShaperPosition(cid);
266 int channel = mapper.getShaperChannel(cid);
267
268 //== Get DSP data
269
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);
273 DBArray<ECLDspData> dsp_data(payload_name);
274 const ECLDspData* data = dsp_data[dsp_id * 12 + shaper - 1];
275
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);
283
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);
291
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();
299
300 //== Get ADC data
301 int* y = adc.data();
302
303 //== Get trigger time
304 int ttrig2 = ttrig - 2 * (ttrig / 8);
305
306 //== Get thresholds
307 DBObjPtr<ECLCrystalCalib> thr_LowAmp("ECL_FPGA_LowAmp");
308 DBObjPtr<ECLCrystalCalib> thr_HitThresh("ECL_FPGA_HitThresh");
309 DBObjPtr<ECLCrystalCalib> thr_StoreDigit("ECL_FPGA_StoreDigit");
310
311 int A0 = thr_LowAmp->getCalibVector()[cid - 1];
312 int Ahard = thr_HitThresh->getCalibVector()[cid - 1];
313 int Askip = thr_StoreDigit->getCalibVector()[cid - 1];
314
315 //== Perform fit
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);
319
320 return result;
321}
Class for accessing arrays of objects in the database.
Definition: DBArray.h:26
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
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.

◆ writeEclDsp()

void writeEclDsp ( const char *  raw_file,
ECLDspData obj 
)
static

Convert ECLDspData from Root object to *.dat file.

Parameters
[in]raw_filePath to dsp??.dat file to be created.
[in]objObject to be written

Definition at line 134 of file ECLDspUtilities.cc.

135{
136 // Default header for DSP file.
137 // Words '0xABCD' are overwritten with current parameters.
138 const unsigned short DEFAULT_HEADER[256] {
139 // ECLDSP FILE.....
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,
172 };
173
174 FILE* fl;
175 fl = fopen(filename, "wb");
176 // Word size
177 const int nsiz = 2;
178 // Size of fragment to write (in words)
179 int nsiz1 = 256;
180 // modification of DEFAULT_HEADER
181 unsigned short header[256];
182
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.");
186 }
187
188 for (int i = 0; i < 256; i++) {
189 if (i == 8) {
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();
207 } else {
208 // Reverse bytes for header.
209 int high = (DEFAULT_HEADER[i] & 0xFF00) >> 8;
210 int low = (DEFAULT_HEADER[i] & 0x00FF);
211 header[i] = (low << 8) + high;
212 }
213 }
214
215 int size = fwrite(header, nsiz, nsiz1, fl);
216 if (size != nsiz1) {
217 B2FATAL("Error writing header of DSP file " << filename);
218 }
219
220 std::vector<short int> f(49152), f1(49152), f31(49152),
221 f32(49152), f33(49152), f41(6144), f43(6144);
222 data->getF41(f41);
223 data->getF31(f31);
224 data->getF32(f32);
225 data->getF33(f33);
226 data->getF43(f43);
227 data->getF(f);
228 data->getF1(f1);
229
230 for (int i = 0; i < 16; ++i) {
231 nsiz1 = 384;
232 fwrite(&(*(f41.begin() + i * nsiz1)), nsiz, nsiz1, fl);
233 nsiz1 = 3072;
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);
237 nsiz1 = 384;
238 fwrite(&(*(f43.begin() + i * nsiz1)), nsiz, nsiz1, fl);
239 nsiz1 = 3072;
240 fwrite(&(*(f.begin() + i * nsiz1)), nsiz, nsiz1, fl);
241 fwrite(&(*(f1.begin() + i * nsiz1)), nsiz, nsiz1, fl);
242 }
243 fclose(fl);
244}

Member Data Documentation

◆ pedestal_fit_initialized

int pedestal_fit_initialized = 0
staticprivate

Flag indicating whether arrays fg31,fg32 are filled.

Definition at line 101 of file ECLDspUtilities.h.

◆ pedfit_fg31

float pedfit_fg31 = {}
staticprivate

DSP coefficients used to determine amplitude in pedestalFit.

Definition at line 103 of file ECLDspUtilities.h.

◆ pedfit_fg32

float pedfit_fg32 = {}
staticprivate

DSP coefficients used to determine time in pedestalFit.

Definition at line 105 of file ECLDspUtilities.h.


The documentation for this class was generated from the following files: