StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StFcsPulseAna.h
1 //Author:David Kapukchyan
2 /*
3 @[December 17, 2022]
4 > SetFitPars can now set a function parameters for a range of in npeaks as opposed to all peaks
5 
6 @[June 25, 2022]
7  > Added Refinements to drawing and copying. Also made some changes to how the tunnel probability is used to merge the peaks.
8 
9 @[March 3, 2022]
10  > Copied StFcsPulseFit so that I can inherit from and use the PeakAna class methods from MyTools
11 */
12 
17 #ifndef STFCSPULSEANA_H
18 #define STFCSPULSEANA_H
19 
20 //ROOT Headers
21 #include "TMath.h"
22 #include "TLegend.h"
23 #include "TLegendEntry.h"
24 #include "TF1.h"
25 #include "TGraphAsymmErrors.h"
26 #include "TMultiGraph.h"
27 #include "TH1.h"
28 
29 #include "TLine.h"
30 #include "TPaveText.h"
31 #include "TPaveStats.h"
32 
33 #include "StFcsDbMaker/StFcsDbPulse.h"
34 #include "PeakAna.h"
35 
36 class StFcsPulseAna : public PeakAna
37 {
38  public:
39  StFcsPulseAna();
40  StFcsPulseAna( std::string name );
41  explicit StFcsPulseAna( TGraph* Sig, std::string name = "StFcsPulseAna");
42  StFcsPulseAna(const StFcsPulseAna& old, const char* post_name="_copy", TGraph* graph=0);
44 
45  virtual ~StFcsPulseAna();
46 
47  virtual void Copy(TObject& obj) const;
48  virtual TObject* Clone(const char* newname) const;
49 
50  virtual Int_t AnalyzeForPeak();
51  virtual Int_t AnalyzeForPeak(Double_t peak, Double_t width){ return PeakAna::AnalyzeForPeak(peak,width); }//<! same as #AnalyzeForPeak() but also sets search parameters for peak
52  virtual Int_t AnalyzeForNoisyPeak(){ return PeakAna::AnalyzeForNoisyPeak(); }//<! First calls #ConvertToAna() function and then #AnalyzeForPeak() on the converted data
53 
54  static double MaxwellBoltzmannDist(double* x, double* p);
55  static void GetMBPars(const double& xpeak, const double& xrise, const double& yh, const double& ped, double& height, double& scale );
56 
57  virtual TGraphAsymmErrors* GetData()const{return (TGraphAsymmErrors*)mG_Data;}
58 
59  const char* GetName() const {return mName.c_str();}
60  std::string& Name() {return mName;}
61  const std::string& Name() const {return mName;}
62 
64 
65  void ResetFinder();
66  void ResetBaseline();
67  void ResetSum();
68 
69  Int_t Sum(Int_t Start, Int_t End);
70  Int_t SumWindow();
71  Double_t GausFit(Int_t Start=0,Int_t End=0);
72  void SignalMBPars(double& height, double& scale);
73  double SumMB();
74  Double_t MBFit(Int_t Start=0,Int_t End=0);
75  Double_t PulseFit(Int_t Start=0, Int_t End=0);
76 
77  void SetFitPars(TF1*& func, int start=-1, int end=-1);
78 
79  static void FillAdc(TGraphAsymmErrors* g, unsigned short& counter, int Start, unsigned short* adcdata);
80  static int SumDep0(TGraphAsymmErrors* gdata, int Start, int ped=0);
81  static int SumDep0Mod(TGraphAsymmErrors* gdata, int Start, int ped=0);
82 
83  void SetFitSignal(TF1* func){mF1_SignalFit=func;}
84  void SetBaselineFit(TF1* func){mF1_BaselineFit=func;}
85 
86  void SetZS(){SetBaseline(0,0.39);}
87 
88  void AnalyzeForPedestal();
89 
90  TH1F* BaselineHist()const{return mH1_Baseline;}
91  TF1* SignalFit(const char option ='n');//Replace with a function that actually calls "fit" and then returns the function?Arguments can be values and range?A little more complicated since hard to tell when user would want to call a fit or just get the function itself?Maybe no arguments means give me function, arguments means do a refit?
92  TF1* SignalFit() const {return mF1_SignalFit;}
93  TF1* BaselineFit()const{return mF1_BaselineFit;}
94 
95  virtual StFcsPulseAna* DrawCopy(Option_t* opt="",const char* name_postfix = "_copy", TGraph* graph=0) const;
96 
97  virtual void Print(Option_t* opt="") const;
98 
99  virtual void MergeByProbability(std::vector<PeakWindow>& merged) const;
100 
101  protected:
102  void Init();
103 
105 
106  TH1F* mH1_Baseline;
107  //Main fit functions
110 
111  bool FindBaseline();
112 
113  private:
114 
115  std::string mName;
116 
117  bool mWindowSum;
118  bool mFitSum;
119 
120  //Variables to hold sum values. This is to avoid tedious recomputation
121  Int_t mSumAdc;
122  Double_t mSumAdcFit;
123 
124  ClassDef(StFcsPulseAna,1)
125 
126 };
127 
128 #endif
virtual TGraphAsymmErrors * GetData() const
Overwrite from PeakAna to type convert for StFcsWaveformFitMaker.
Definition: StFcsPulseAna.h:57
void SignalMBPars(double &height, double &scale)
Figure out the height and scale of a Maxwell-Boltzmann distribution that approximates the signal...
void SetZS()
Call this for ZS data which uses thresholds that are relevant for that. (like 0 baseline and 0...
Definition: StFcsPulseAna.h:86
StFcsDbPulse * mDbPulse
Pointer to StFcsDbPulse.
virtual TObject * Clone(const char *newname) const
Clones internal graph as opposed to just copying the pointer.
static int SumDep0Mod(TGraphAsymmErrors *gdata, int Start, int ped=0)
Test my modified sum method to DEP board.
const std::string & Name() const
Definition: StFcsPulseAna.h:61
Double_t MBFit(Int_t Start=0, Int_t End=0)
Fit a Maxwell-Boltzmann distribution to mFoundPeak and return the integral minus the baseline...
virtual void MergeByProbability(std::vector< PeakWindow > &merged) const
Overwritten from PeakAna::MergeByProbability() to change merge criteria.
void ResetSum()
Only resets variables related to finding the sum.
std::string & Name()
Definition: StFcsPulseAna.h:60
TGraph * mG_Data
TGraph that stores the x,y data.
Definition: PeakAna.h:485
TF1 * mF1_BaselineFit
Gaussian function to fit to mH1_Baseline to determine baseline.
Double_t PulseFit(Int_t Start=0, Int_t End=0)
Fit the pulse shape defined in StFcsDbPulse::multiPulseShape() to all peaks and return the integral o...
Double_t GausFit(Int_t Start=0, Int_t End=0)
Do a Gaussian fit on mFoundPeak and return the integral subtracted by the baseline.
void AnalyzeForPedestal()
Analyze graph data to determine baseline internally.
static void GetMBPars(const double &xpeak, const double &xrise, const double &yh, const double &ped, double &height, double &scale)
Get parameters for a Maxwell Boltzmann distribution from above based on the 4 const parameters...
void SetBaseline(Double_t base, Double_t sigma)
Definition: PeakAna.cxx:330
void SetFitSignal(TF1 *func)
Definition: StFcsPulseAna.h:83
static void FillAdc(TGraphAsymmErrors *g, unsigned short &counter, int Start, unsigned short *adcdata)
Needed for SumDep0, et al. It basically fills in tb vs. adc sequentially so all timebins from a given...
virtual void Copy(TObject &obj) const
Internal method use Clone instead.
virtual Int_t AnalyzeForPeak()
Main analysis method for finding peaks.
Definition: PeakAna.cxx:153
virtual StFcsPulseAna * DrawCopy(Option_t *opt="", const char *name_postfix="_copy", TGraph *graph=0) const
Clone and draw, see PeakAna::Draw() for options.
void ResetFinder()
Resets all histograms and values.
TF1 * SignalFit() const
Definition: StFcsPulseAna.h:92
Int_t SumWindow()
Sum the ADC in the peak defined by mFoundPeak and subtract the baseline.
void Init()
Initialize everything to zero except signal and background histograms.
void SetFitPars(TF1 *&func, int start=-1, int end=-1)
Set the parameters of an external TF1 function that has the form of StFcsDbPulse::multiPulseShape(), optionally only set fit paramaters for peaks from index start up to and including end.
virtual ~StFcsPulseAna()
Destructor.
bool FindBaseline()
Does Gaussian fitting to mH1_Baseline to determine baseline.
TH1F * BaselineHist() const
Definition: StFcsPulseAna.h:90
const char * GetName() const
Definition: StFcsPulseAna.h:59
virtual void Print(Option_t *opt="") const
Print class variables.
static int SumDep0(TGraphAsymmErrors *gdata, int Start, int ped=0)
Test of sum method in DEP board.
StFcsPulseAna()
Default constructor.
TF1 * mF1_SignalFit
Function to fit to the data.
void setDbPulse(StFcsDbPulse *p)
Definition: StFcsPulseAna.h:63
virtual Int_t AnalyzeForPeak()
Overwritten from PeakAna to process peak tunneling after finding all peaks.
void SetBaselineFit(TF1 *func)
Definition: StFcsPulseAna.h:84
void ResetBaseline()
Resets baseline values.
StFcsPulseAna & operator=(const StFcsPulseAna &rhs)
Assignment operator.
TF1 * BaselineFit() const
Definition: StFcsPulseAna.h:93
double SumMB()
Integral of an approximated Maxwell-Boltzmann distribution minus the baseline.
TH1F * mH1_Baseline
Histogram that holds projection of mG_Data for baseline determination.
Int_t Sum(Int_t Start, Int_t End)
Add raw ADC within a given range and subtract the baseline.
static double MaxwellBoltzmannDist(double *x, double *p)
Maxwell Boltzmann Distribution function.