BmnRoot
Loading...
Searching...
No Matches
BmnHodoDigi.h
Go to the documentation of this file.
1/* Copyright (C) 2021 Institute for Nuclear Research, Moscow
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Nikolay Karpushkin [committer] */
4
15#ifndef BmnHodoDigi_H
16#define BmnHodoDigi_H 1
17
18#include "BmnDigiContainerTemplate.h" // for BmnDigiContainerTemplate
19#include "BmnHodoAddress.h" // for BmnHodoAddress
20#include "BmnHodoDigit.h" // for BmnHodoDigit
21
23 : public BmnHodoDigit
25{
26
27 public:
31 : BmnHodoDigit()
33 , fNpeaks(0)
34 , fR2history{}
35 , fSignalIntegral{0.0}
36 , fCrosstalk{0.0}
37 {}
38
41
42 void reset() override final
43 {
46 fNpeaks = 0;
47 fR2history.clear();
48 fSignalIntegral = 0.0;
49 fCrosstalk = 0.0;
50 }
51
52 size_t fNpeaks = 0;
53 std::vector<float> fR2history;
56
60 virtual const char* GetClassName() override final { return "BmnHodoDigi"; }
61
65 uint32_t GetStripSide() const { return BmnHodoAddress::GetStripSide(fAddress); };
66
70 uint32_t GetGain() const { return BmnHodoAddress::GetGain(fAddress); };
71
72 const int DrawWfm()
73 {
74 TString hist_name =
75 (fNpeaks > 1)
76 ? Form("Strip%u. From deconv. Amplitude %.0f ZL %.0f TimeMax %.0f FitR2 %.4f; time [sample]; ampl "
77 "[a.u.]",
79 : Form("Strip%u. From raw. Amplitude %d ZL %d TimeMax %d FitR2 %.4f; time [sample]; ampl [a.u.]",
81 TCanvas* canvas = new TCanvas();
82 canvas->Divide(1, 2);
83 DrawWfmWithTitle(canvas, hist_name);
84 DrawR2history(canvas);
85
86 return 1;
87 }
88
89 void DrawR2history(TCanvas* canvas)
90 {
91 if (fR2history.empty())
92 return;
93 canvas->cd(2);
94
95 std::vector<float> points(fR2history.size());
96 std::iota(std::begin(points), std::end(points), 0); // Fill with 0, 1, ..., wfm.back().
97 TGraph* tgr_ptr = new TGraph(fR2history.size(), &points[0], &fR2history[0]);
98 tgr_ptr->SetTitle(Form("R2 history; iteration; R2"));
99 tgr_ptr->Draw();
100 }
101
103};
104
105#endif // BmnHodoDigi_H
virtual void reset()
Data class for Bmn digi container template.
float fFitTimeMax
Quality of waveform fit [] – good near 0.
float fFitZL
Amplitude from fit of waveform [adc counts].
float GetFitR2() const
Fit R2 quality.
void DrawWfmWithTitle(TCanvas *canvas, TString hist_name)
int fZL
Amplitude from waveform [adc counts].
float fFitAmpl
Time over threshold [adc samples].
int fTimeMax
Energy deposition from waveform [adc counts].
static uint32_t GetStripSide(uint32_t address)
Return Strip side from address.
static uint32_t GetGain(uint32_t address)
Return Gain from address.
Data class for BmnHodo digital signal processing.
Definition BmnHodoDigi.h:25
void reset() override final
Definition BmnHodoDigi.h:42
float fSignalIntegral
Definition BmnHodoDigi.h:54
uint32_t GetGain() const
Gain.
Definition BmnHodoDigi.h:70
const int DrawWfm()
Definition BmnHodoDigi.h:72
float fCrosstalk
Definition BmnHodoDigi.h:55
std::vector< float > fR2history
Definition BmnHodoDigi.h:53
size_t fNpeaks
Definition BmnHodoDigi.h:52
virtual const char * GetClassName() override final
Class name.
Definition BmnHodoDigi.h:60
uint32_t GetStripSide() const
Strip Side.
Definition BmnHodoDigi.h:65
ClassDefOverride(BmnHodoDigi, 4)
BmnHodoDigi()
Default constructor.
Definition BmnHodoDigi.h:30
void DrawR2history(TCanvas *canvas)
Definition BmnHodoDigi.h:89
uint32_t GetStripId() const