Berkeley Nuclear Data Software
NeutronAngularDistAna.h
Go to the documentation of this file.
1 //ROOT includes
2 #include "TH1.h"
3 #include "TH2.h"
4 #include "TMatrixD.h"
5 #include "TFitResult.h"
6 //c++includes
7 #include <cstdint>
8 #include <vector>
9 #include <utility>
10 #include <deque>
11 #include <map>
12 using std::vector;
13 using std::map;
14 using std::pair;
15 // externals
16 // #include "../../BasicSupport/external/eigen/Eigen/Dense"
17 #include <Eigen/Dense>
19 
21 struct EffData
22 {
23  EffData();
24  vector<double> m_lowBinEdges;
25  vector<double> m_effs;
26  vector<double> m_effsErr;
27  double m_angle;
29  vector<double> getEfficiency(double a_lowE, double a_highE);
31  void drawEff();
32 };
36 {
37 public:
41  void readEfficiencyData(std::string a_effFile);
43  void uniBinEfficiency(int a_detID, int a_nBins, double a_lowE=1);
46  TH2F* drawEffCorrAngDist(std::string a_expFile,
47  std::string a_bkgFile);
54  void addExperimentalTOFData(std::string a_fileName);
57  // stores scint/scint in m_tofHistBKGOrgOrg
58  void addExperimentalTOFBKGData(std::string a_fileName);
64  void setBackgroundScaler(double a_bkgScaler,
65  bool a_scintScint=false);
68  std::vector<double> getEnerDistTOF(int a_detID,
69  double a_tof,
70  bool a_draw=false,
71  bool a_useEff=true,
72  double a_tofWidth=8);
75  std::vector<double> getEnerDistTOF(double a_angle,
76  double a_angleWidth,
77  double a_tof,
78  bool a_draw=false,
79  bool a_useEff=true,
80  double a_tofWidth = 8);
82  void drawRFExpBkg(int a_detID, bool a_drawSub=false);
83 
85  void plotEnerDistTOF(int a_detID);
86  void plotEnerDistTOF(double a_angle,
87  double a_angleWidth);
90  map<double, std::pair<double,double>> plotEnerDistTOF(int a_detID,
91  Eigen::VectorXf a_calcYield,
92  double a_normalization = -1);
94  void plot2DEnerDistTOF(int a_detID);
98  void getAngDistTOF(int a_tof,
99  bool a_normalize=false,
100  int a_tofWidth=8);
102  EffData getEffData(int a_detID);
106  TH2F* getScintScintHist(int a_detID,
107  bool a_normEff=false);
108  private:
109  map<int, EffData> m_detEff;
110  map<int, TH2*> m_tofHists;
111  map<int, TH2*> m_tofHistsBKG;
112  map<int, TH2*> m_tofHistsOrgOrg;
113  map<int, TH2*> m_tofHistsBKGOrgOrg;
114 
117  map<int,map<double,vector<double>>> m_maxwellFitRes;
119  double m_bkgScaler;
120  double m_bkgScalerOrgOrg;
121 };
122 
123 // [0]*2*sqrt(x/3.14)*pow(1./[1],1.5)*exp(-x/[1])+[2]
Definition: NeutronAngularDistAna.h:36
EffData getEffData(int a_detID)
returns efficiency data
Definition: NeutronAngularDistAna.cpp:937
void addExperimentalTOFBKGData(std::string a_fileName)
Definition: NeutronAngularDistAna.cpp:257
void uniBinEfficiency(int a_detID, int a_nBins, double a_lowE=1)
draw efficiency with uniform a_lowE-10 MeV binning
Definition: NeutronAngularDistAna.cpp:164
std::vector< double > getEnerDistTOF(int a_detID, double a_tof, bool a_draw=false, bool a_useEff=true, double a_tofWidth=8)
Definition: NeutronAngularDistAna.cpp:322
void addExperimentalTOFData(std::string a_fileName)
Definition: NeutronAngularDistAna.cpp:234
void plot2DEnerDistTOF(int a_detID)
uses fits to draw time-since-last-rf vs neutron energy
Definition: NeutronAngularDistAna.cpp:822
void getAngDistTOF(int a_tof, bool a_normalize=false, int a_tofWidth=8)
Definition: NeutronAngularDistAna.cpp:848
void readEfficiencyData(std::string a_effFile)
Definition: NeutronAngularDistAna.cpp:93
TH2F * getScintScintHist(int a_detID, bool a_normEff=false)
Definition: NeutronAngularDistAna.cpp:948
TH2F * drawEffCorrAngDist(std::string a_expFile, std::string a_bkgFile)
Definition: NeutronAngularDistAna.cpp:194
void setBackgroundScaler(double a_bkgScaler, bool a_scintScint=false)
Definition: NeutronAngularDistAna.cpp:280
void plotEnerDistTOF(int a_detID)
plots integral and temperature from maxwellian fit across TOF domain
Definition: NeutronAngularDistAna.cpp:517
NeutronAngularDistAna()
Definition: NeutronAngularDistAna.cpp:86
void drawRFExpBkg(int a_detID, bool a_drawSub=false)
draw x-projection of m_tofHists histograms for data and background
Definition: NeutronAngularDistAna.cpp:290
J. Gordon.
Definition: NeutronAngularDistAna.h:22
vector< double > m_effsErr
Definition: NeutronAngularDistAna.h:26
EffData()
Definition: NeutronAngularDistAna.cpp:23
double m_angle
Definition: NeutronAngularDistAna.h:27
vector< double > getEfficiency(double a_lowE, double a_highE)
calculate the efficiency for a given range of neutron energy
Definition: NeutronAngularDistAna.cpp:28
void drawEff()
draws the efficiency
Definition: NeutronAngularDistAna.cpp:64
vector< double > m_lowBinEdges
Definition: NeutronAngularDistAna.h:24
vector< double > m_effs
Definition: NeutronAngularDistAna.h:25