Berkeley Nuclear Data Software
ActivationAnalysis.h
Go to the documentation of this file.
1 #ifndef _ACTIVATION_ANALYSIS_H_
2 #define _ACTIVATION_ANALYSIS_H_
3 // root includes
4 #include "TH1.h"
5 #include "TFile.h"
6 // c++ includes
7 #include <utility>
8 #include <string>
9 
10 
11 #include <vector>
12 #include <map>
13 #include <string>
14 using std::map;
15 using std::vector;
16 using std::string;
18 {
19  public:
20  // default constructor
22  //Takes SPE files with calibration parameters from InterSpec. Key takes "Eu", "Foil", "Bkg"
23  void readSPEFile(std::string a_fileName, std::string a_key);
24 
25  // add the start time for the beginning of irradiation
26  void setIrradiationStartTime(double a_year,
27  double a_month,
28  double a_day,
29  double a_hour,
30  double a_minute,
31  double a_second);
32  //Calculates the integrated beam current over irradiation time. Requires csv file containing beam on time.
33  void addBeamCurrentFile(std::string a_fileName);
34  //Takes distance from beam flange center to target in centimeters
35  void setBeamFlightPath(double a_flightPath);
36  //Add sTOF spectrum for flux-weighted cross section calculation
37  void addFluxHistogram(std::string a_fileName, std::string a_histName="fluxHist");
38  // path to directory where IRDFF-II.tab.txt file is located
39  void setIRDFFTablePath(std::string a_fileName);
40  //Reads the IRDFF table and gets cross section data for given reaction. Requries X-23(a,b)Y-27m) format
41  std::map<int, std::vector<double> > readIRDFFTable(std::string a_reaction);
42  // Calculate the flux weighted cross section. Requries X-23(a,b)Y-27m) format
43  void calcFluxWeightedXSec(std::string a_reaction);
44  // Analyze Eu152 spectra and computes an efficiency curve for the given calibration histogram
45  void analyzeEu152Data(double activity);
46  //Takes parameters from SpectrumAnalysis getEffCurve and calculates the efficiency at a given energy
47  std::pair<double,double> getEff(double energy);
48  // Calculate activity at the end of irradiation. Requires ../SandiaDecay/CMakeLists.txt line 12 be SHARED.
49  //if using InterSpec counts, autoFit=false and becomes the number of counts. a_time is uncertantity of the half-life
50  void calcA0(std::string a_isotope,double a_peakMean,double a_time, bool autoFit);
51 //Calculates flux. Z is z number of foil, mass is mass of foil, allIso takes specific isotope in format "Na24" "C12" or "All" for every isotope
52  //the abundancy error and mass error vectors must ordered by isotopoe so that the nth element is the mass and abundancy error of the nth isotope
53  void calcNormalizedFlux(double mass,double a_Z,std::vector<double> a_abundancies_err, std::vector<double> a_mass_err,string allIso);
54 
55 private:
56 
57  TH1F* m_effHist;
58  TH1F* m_bkgHist;
59  TH1F* m_foilHist;
60  int m_mon=0;
61  int m_day=0;
62  int m_year=0;
63  int m_hour=0;
64  int m_tj=0;
65  int m_eu=0;
66  int m_ti=0;
67  int m_bk=0;
68  double m_a0=0.0;
69  double m_xs=0.0;
70  int m_tc=0;
71  double m_counts=0;
72  double m_lin=0.0;
73  int m_minute=0;
74  int m_second=0;
75  double m_decay=0.0;
76  double m_xs_err=0.0;
77  double m_a0_err=0.0;
78  double m_offset=0.0;
79  int m_current=0;
80  double m_decay_err=0.0;
81  int m_flightPath=0;
82  double m_efficiency=0.0;
83  int m_curent_err=0;
84  double m_decayConstant=0.0;
85  std::vector<int> m_SPE;
86  std::string m_rootFile;
87  std::string m_irdffPath;
88  std::string m_bkgrootFile;
89  std::vector<int> m_bkgSPE;
90  std::vector<double> m_effParams;
91  std::vector<double> m_neutrons;
92  std::vector<double> m_neutrons_err;
93  std::vector<double> m_binsVect;
94  std::vector<double> m_widthVect;
95  time_t m_epochIrradiationStartTime;
96 };
97 #endif
Definition: ActivationAnalysis.h:18
void addBeamCurrentFile(std::string a_fileName)
Definition: ActivationAnalysis.cpp:164
void setIrradiationStartTime(double a_year, double a_month, double a_day, double a_hour, double a_minute, double a_second)
Definition: ActivationAnalysis.cpp:144
std::map< int, std::vector< double > > readIRDFFTable(std::string a_reaction)
Definition: ActivationAnalysis.cpp:191
void analyzeEu152Data(double activity)
Definition: ActivationAnalysis.cpp:394
void calcA0(std::string a_isotope, double a_peakMean, double a_time, bool autoFit)
Definition: ActivationAnalysis.cpp:431
void calcFluxWeightedXSec(std::string a_reaction)
Definition: ActivationAnalysis.cpp:295
void readSPEFile(std::string a_fileName, std::string a_key)
Definition: ActivationAnalysis.cpp:40
void setBeamFlightPath(double a_flightPath)
Definition: ActivationAnalysis.cpp:260
void addFluxHistogram(std::string a_fileName, std::string a_histName="fluxHist")
Definition: ActivationAnalysis.cpp:267
void setIRDFFTablePath(std::string a_fileName)
Definition: ActivationAnalysis.cpp:186
std::pair< double, double > getEff(double energy)
Definition: ActivationAnalysis.cpp:412
void calcNormalizedFlux(double mass, double a_Z, std::vector< double > a_abundancies_err, std::vector< double > a_mass_err, string allIso)
Definition: ActivationAnalysis.cpp:536
ActivationAnalysis()
Definition: ActivationAnalysis.cpp:35