PandaRoot
PndEmcPSAOptimalFilterAnalyser.h
Go to the documentation of this file.
1 //****************************************************************************
2 //* This file is part of PandaRoot. *
3 //* *
4 //* PandaRoot is distributed under the terms of the *
5 //* GNU General Public License (GPL) version 3, *
6 //* copied verbatim in the file "LICENSE". *
7 //* *
8 //* Copyright (C) 2006 - 2024 FAIR GmbH and copyright holders of PandaRoot *
9 //* The copyright holders are listed in the file "COPYRIGHTHOLDERS". *
10 //* The authors are listed in the file "AUTHORS". *
11 //****************************************************************************
12 
13 #ifndef PNDEMCOPTIMALFILTERANALYSER_HH
14 #define PNDEMCOPTIMALFILTERANALYSER_HH
15 
16 #include "PndEmcAbsPSA.h"
17 #include <string>
18 #include <vector>
19 
21 
22 //#define MAX_NUMBER_OptimalFilter_HITS 100
23 
25 
26 public:
27  PndEmcPSAOptimalFilterAnalyser(Double_t cf_frac, Double_t cf_tap, const Double_t* cf_tcorr, Int_t cf_n,
28  const Double_t* of_a, const Double_t* of_b, Int_t of_m, Int_t of_b0, PndEmcAbsPulseshape* shape, Double_t threshold, Double_t baseline = 0., Int_t verbose = 1);
31  void setBaseline(Double_t baseline) { fBaseline = baseline; }
32  void setVerbose(Int_t verbose) { fVerbose = verbose; }
33 
34  Int_t Process(const PndEmcWaveform *waveform);
35  void GetHit(Int_t i, Double_t &energy, Double_t &time);
36  void Reset();
37 
38 private:
39  Double_t fBaseline;
40 
41  Double_t fCF_fraction;
42  Int_t fCF_tap;
43  Int_t fCF_nsamples;
44  Double_t* fCF_tcorr;
45 
46  Int_t fOF_b0;
47  Int_t fOF_m;
48  Double_t* fOF_a;
49  Double_t* fOF_b;
50 
51  Double_t fThreshold;
52  std::vector<Double_t> fDigiAmplitude;
53  std::vector<Double_t> fDigiTime;
54 
55  PndEmcAbsPulseshape* fPulse;
56  Double_t fSampleRate;
57  Double_t fTimeStamp;
58 
59  Int_t fVerbose;
60 
61  void analyse(std::vector<Double_t>& signal, std::vector<Double_t>& baseline, Int_t& start_position); // return position for the next signal processing
62  void zero_crossing(const std::vector<Double_t>& wf, Int_t& start_position, Double_t& time, Int_t& index);
63  void generate_baseline(std::vector<Double_t>& baseline, Double_t amplitude, Double_t time);
64  void subtract_baseline(std::vector<Double_t>& signal, const std::vector<Double_t>& baseline);
65 
66 
67  ClassDef(PndEmcPSAOptimalFilterAnalyser, 1);
68 };
69 
70 #endif
void GetHit(Int_t i, Double_t &energy, Double_t &time)
Get energy and time of hit.
unsigned int i
Definition: P4_F32vec4.h:33
Int_t Process(const PndEmcWaveform *waveform)
Find Hits in Waveform.
represents a simulated waveform in an emc crystal
Baseclass for pulseshapeanalysis ( featureextraction )
Definition: PndEmcAbsPSA.h:33
void Reset()
reset found hits
pulseshape interface