Belle II Software  release-08-01-10
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. More...
 
static void writeEclDsp (const char *raw_file, ECLDspData *obj)
 Convert ECLDspData from Root object to *.dat file. More...
 
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. More...
 
static void initPedestalFit ()
 Load DSP coefficients used in the pedestal fit function. More...
 
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. More...
 

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.

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:148

◆ 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.

◆ 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.

◆ 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.

◆ 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.


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