Berkeley Nuclear Data Software
SpectrumAnalysis.h
Go to the documentation of this file.
1 #ifndef _SPECTRUM_ANALYSIS_H_
2 #define _SPECTRUM_ANALYSIS_H_
6 
7 #include "TH1.h"
8 #include "TFile.h"
9 
10 #include <vector>
11 #include <map>
12 #include <string>
13 using std::map;
14 using std::vector;
15 using std::string;
16 
18 struct Peak
19 {
20  Peak();
21  double m_amplitude = 0.;
22  double m_mean = 0.;
23  double m_stdDev = 0.;
24  double m_chi2ndf = 0.;
25  double m_binWidth = 0.;
26  double m_stdDevErr = 0.;
27  double m_ampErr = 0.;
28  void printPeak();
30  bool checkResolution();
31  double m_maxRes = 2.; // percent deltaE/E
32 };
33 
35 {
36 public:
39 
42  void loadHists(string a_fname, int a_modID, double a_countTime=0);
43  void loadSimHists(string a_fname);
44  void loadBackgroundHists(string a_fname, int a_modID, double a_countTime);
48  vector<Peak> FitPeaksMan(int a_channel,
49  double a_lower,
50  double a_higher,// int a_numPeaks,
51  vector<double> a_peakMeans,
52  bool a_fullClover=false,
53  bool a_addBack=false,
54  bool a_checkPeakRes=true);
55 
57  vector<Peak> FitPeaksManPol1(int a_channel,
58  int a_lower,
59  int a_higher,
60  // int a_numPeaks,
61  vector<double> a_peakMeans,
62  bool a_checkPeakRes=true,
63  bool a_fullClover=false,
64  bool a_addBack=false);
66  void printPeakFits();
68  vector<Peak> FitPeaksErf(int a_channel,
69  double a_lower,
70  double a_higher,
71  int a_numPeaks,
72  bool a_fullClover=false,
73  bool a_addBack=false);
74 
79  map<double, std::pair<double,double>> Eu152Efficiency(int a_channel,
80  double a_activity,
81  double a_countTime,
82  double a_simBrTot=1,
83  bool a_res=false,
84  bool a_fullClover=false,
85  bool a_addBack=false,
86  bool a_draw =false);
87 
88  void calcAddBackFactor();
91  map<double, std::pair<double,double>> U235Efficiency(int a_channel,
92  double a_activity,
93  double a_activityErr,
94  double a_countTime);
95 
97  map<int, TH1*> getHists();
98  map<int,vector<TH1*>> getCloverHists();
99  map<int, vector<Peak>> getPeakFits();
100 
102  void ChannelResolution(int a_channel);
103 /*
109  void findPeaks(int a_channel, double a_minHeight, int a_numPeaks);
110 
111  void FitPeakBasic(int a_channel, double a_peakCentroid, double a_low= 5 , double a_high=5);;
112  void FitDoublePeak(int a_channel, double a_peakCentroid, double a_peakSeparation,
113  double a_relativeHeight, double a_fitWindow);
114 
115  map<string, vector<vector<double>>> m_zeroCrossings;
117  void writeToConfig();
118 */
119 
120 
121 private:
122  double m_countTime=0;
123  double m_bkgCountTime=0;
124  vector<double> m_energy;
125 
127  map<int, TH1*> m_hists;
128  map<int, TH1*> m_bkgHists;
129  map<int, vector<TH1*>> m_cloverHists;
130  map<int, vector<double>> m_calibration;
131 
133  map<int, vector<Peak>> m_peakFits;
134 };
135 #endif
Definition: SpectrumAnalysis.h:35
void ChannelResolution(int a_channel)
creates plot of resolution vs energy
Definition: SpectrumAnalysis.cpp:687
map< int, vector< TH1 * > > getCloverHists()
Definition: SpectrumAnalysis.cpp:678
map< int, vector< Peak > > getPeakFits()
Definition: SpectrumAnalysis.cpp:682
map< double, std::pair< double, double > > Eu152Efficiency(int a_channel, double a_activity, double a_countTime, double a_simBrTot=1, bool a_res=false, bool a_fullClover=false, bool a_addBack=false, bool a_draw=false)
Definition: SpectrumAnalysis.cpp:713
void loadHists(string a_fname, int a_modID, double a_countTime=0)
Definition: SpectrumAnalysis.cpp:102
vector< Peak > FitPeaksMan(int a_channel, double a_lower, double a_higher, vector< double > a_peakMeans, bool a_fullClover=false, bool a_addBack=false, bool a_checkPeakRes=true)
Definition: SpectrumAnalysis.cpp:295
void printPeakFits()
prints information about fitted peaks (mean, area, resolution, chi2/ndf)
Definition: SpectrumAnalysis.cpp:391
vector< Peak > FitPeaksErf(int a_channel, double a_lower, double a_higher, int a_numPeaks, bool a_fullClover=false, bool a_addBack=false)
Definition: SpectrumAnalysis.cpp:504
void loadSimHists(string a_fname)
Definition: SpectrumAnalysis.cpp:234
vector< Peak > FitPeaksManPol1(int a_channel, int a_lower, int a_higher, vector< double > a_peakMeans, bool a_checkPeakRes=true, bool a_fullClover=false, bool a_addBack=false)
fit peaks with 1st degree polynomial background
Definition: SpectrumAnalysis.cpp:400
void calcAddBackFactor()
Definition: SpectrumAnalysis.cpp:828
void subtractBackgroundSpectra()
Definition: SpectrumAnalysis.cpp:267
void loadBackgroundHists(string a_fname, int a_modID, double a_countTime)
Definition: SpectrumAnalysis.cpp:174
map< double, std::pair< double, double > > U235Efficiency(int a_channel, double a_activity, double a_activityErr, double a_countTime)
Definition: SpectrumAnalysis.cpp:938
map< int, TH1 * > getHists()
returns m_hists
Definition: SpectrumAnalysis.cpp:674
SpectrumAnalysis()
default constructor
Definition: SpectrumAnalysis.cpp:30
struct to hold fitted peak data
Definition: SpectrumAnalysis.h:19
Peak()
Definition: SpectrumAnalysis.cpp:14
double m_stdDevErr
Definition: SpectrumAnalysis.h:26
double m_mean
Definition: SpectrumAnalysis.h:22
bool checkResolution()
checks if resolution of fitted peak is < 1%
Definition: SpectrumAnalysis.cpp:24
double m_chi2ndf
Definition: SpectrumAnalysis.h:24
double m_amplitude
Definition: SpectrumAnalysis.h:21
double m_maxRes
Definition: SpectrumAnalysis.h:31
double m_stdDev
Definition: SpectrumAnalysis.h:23
void printPeak()
Definition: SpectrumAnalysis.cpp:16
double m_ampErr
Definition: SpectrumAnalysis.h:27
double m_binWidth
Definition: SpectrumAnalysis.h:25