GamR  0.0.0
GammaROOT
Loading...
Searching...
No Matches
AngCorr.hh
Go to the documentation of this file.
1
9
10#ifndef ANGCORR_HH
11#define ANGCORR_HH
12
13#include <vector>
14
15#include <TMath.h>
16#include <TH1.h>
17#include <TH2.h>
18#include <TF1.h>
19#include <TGraph.h>
20#include <TGraphErrors.h>
21#include <TMultiGraph.h>
22#include <Math/SpecFunc.h> /* legendre polynomials, 3j and 6j symbols */
23
24#include <nucleus/Transition.hh>
25#include <angdist/AngDist.hh>
26
27namespace GamR {
28 namespace AngDist {
29
30 class AngCorr {
34 public:
35 double spin1;
36 double spin2;
37 double spin3;
38
47 AngCorr(double spin1, double spin2, double spin3)
48 {
49 this->spin1 = spin1;
50 this->spin2 = spin2;
51 this->spin3 = spin3;
52 std::vector<double> Qks = {1.0, 1.0, 1.0, 1.0, 1.0};
53 this->Qk = new GamR::AngDist::SolidAttenuation(Qks);
54 }
56 {
57 this->spin1 = spin1;
58 this->spin2 = spin2;
59 this->spin3 = spin3;
60 this->Qk = Q;
61 }
62 double Ak(double k, double delta1, double delta2);
63 double operator() (double *x, double *p)
64 {
65 double theta = x[0]*TMath::Pi()/180.0;
66 double a2 = Ak(2.0, p[1], p[2]);
67 double a4 = Ak(4.0, p[1], p[2]);
68 return p[0]*(1+Qk->Get(2)*a2*ROOT::Math::legendre(2, cos(theta))+Qk->Get(4)*a4*ROOT::Math::legendre(4, cos(theta)));
69 }
70 int Mixed(double s1, double s2);
71 int Mixed(int i);
72 int Mixed1() { return Mixed(spin1, spin2); }
73 int Mixed2() { return Mixed(spin2, spin3); }
74 };
75
76 class AngCorrFit : public TNamed {
80 public:
81 std::vector<Color_t> fColors;
82
83 std::vector<double> spin1s;
84 std::vector<double> spin2s;
85 std::vector<double> spin3s;
86
87 std::vector<double> fAngles;
88 std::vector<double> fIntensities;
89 std::vector<double> fErrors;
90
91 std::vector<AngCorr*> fCorrs;
92 std::vector<TF1*> fCorrFuncs;
93 std::vector<TGraph*> fChiGraphs;
94
95 std::vector<double> fDelta1;
96 std::vector<double> fDelta2;
97 std::vector<std::pair<double,double> > fDelta1Err;
98 std::vector<std::pair<double,double> > fDelta2Err;
99 std::vector<double> fChiSquared;
100
101 std::vector<double> fixDelta;//[2];
102 std::vector<int> fDel;//[2];
103
105
106 public:
107 void SetData(std::vector<double> angs, std::vector<TH2D*> hists, GamR::Nucleus::Transition peak1, GamR::Nucleus::Transition peak2, TH1D *angNorm);
108 void SetSpins(std::vector<double> s1s, std::vector<double> s2s, std::vector<double> s3s) { spin1s = s1s; spin2s = s2s; spin3s = s3s; }
109 void SetQk(GamR::AngDist::SolidAttenuation *Q) { delete Qk; Qk = Q; };
110
111 AngCorrFit(std::vector<double> angs, std::vector<TH2D*> hists, GamR::Nucleus::Transition peak1, GamR::Nucleus::Transition peak2, TH1D *angNorm, const char *name, const char *title);
112 AngCorrFit(std::vector<double> angs, std::vector<TH2D*> hists, GamR::Nucleus::Transition peak1, GamR::Nucleus::Transition peak2, TH1D *angNorm);
113
114 void FixDel1(double delta) { fixDelta[0] = delta; fDel[0] = 1; }
115 void FixDel2(double delta) { fixDelta[1] = delta; fDel[1] = 1; }
116
117 int GetNPoints() {return fAngles.size(); }
118 int GetNCombs() {return fCorrs.size(); }
119 void ErasePoint(double ang);
120 double Average();
121 TGraphErrors *GetGraph();
122 TGraphErrors *GetDispGraph(double start, double interval);
123
124 void Fit();
125 TMultiGraph* GetChiGraph();
126
127 void DrawChiGraph(TCanvas *c1, double ylow, double yhigh);
128 void DrawAngCorr(TCanvas *c1, double start, double interval);
129 void PrintDelta1(int i) { std::cout << "Delta 1 = " << fDelta1[i] << " +" << fDelta1Err[i].second << "/-" << -fDelta1Err[i].first << std::endl;}
130 void PrintDelta2(int i) { std::cout << "Delta 2 = " << fDelta2[i] << " +" << fDelta2Err[i].second << "/-" << -fDelta2Err[i].first << std::endl;}
131 void Save(std::string filename, double start, double interval);
132 void SaveChiGraph(std::string filename);
133
135 };
136 }
137}
138#endif
void DrawAngCorr(TCanvas *c1, double start, double interval)
Definition AngCorr.cc:398
std::vector< double > spin3s
Definition AngCorr.hh:85
std::vector< Color_t > fColors
Definition AngCorr.hh:81
std::vector< double > fDelta2
Definition AngCorr.hh:96
std::vector< double > fAngles
Definition AngCorr.hh:87
GamR::AngDist::SolidAttenuation * Qk
Definition AngCorr.hh:104
AngCorrFit(std::vector< double > angs, std::vector< TH2D * > hists, GamR::Nucleus::Transition peak1, GamR::Nucleus::Transition peak2, TH1D *angNorm, const char *name, const char *title)
Definition AngCorr.cc:122
std::vector< double > fIntensities
Definition AngCorr.hh:88
void SetQk(GamR::AngDist::SolidAttenuation *Q)
Definition AngCorr.hh:109
void SaveChiGraph(std::string filename)
Definition AngCorr.cc:435
void FixDel2(double delta)
Definition AngCorr.hh:115
std::vector< std::pair< double, double > > fDelta2Err
Definition AngCorr.hh:98
std::vector< double > spin2s
Definition AngCorr.hh:84
void SetData(std::vector< double > angs, std::vector< TH2D * > hists, GamR::Nucleus::Transition peak1, GamR::Nucleus::Transition peak2, TH1D *angNorm)
Definition AngCorr.cc:96
std::vector< double > fixDelta
Definition AngCorr.hh:101
void FixDel1(double delta)
Definition AngCorr.hh:114
TMultiGraph * GetChiGraph()
Definition AngCorr.cc:360
std::vector< double > fChiSquared
Definition AngCorr.hh:99
void SetSpins(std::vector< double > s1s, std::vector< double > s2s, std::vector< double > s3s)
Definition AngCorr.hh:108
std::vector< AngCorr * > fCorrs
Definition AngCorr.hh:91
void ErasePoint(double ang)
Definition AngCorr.cc:129
std::vector< TGraph * > fChiGraphs
Definition AngCorr.hh:93
std::vector< int > fDel
Definition AngCorr.hh:102
std::vector< double > spin1s
Definition AngCorr.hh:83
TGraphErrors * GetDispGraph(double start, double interval)
Definition AngCorr.cc:160
std::vector< TF1 * > fCorrFuncs
Definition AngCorr.hh:92
std::vector< double > fErrors
Definition AngCorr.hh:89
std::vector< double > fDelta1
Definition AngCorr.hh:95
void Save(std::string filename, double start, double interval)
Definition AngCorr.cc:414
std::vector< std::pair< double, double > > fDelta1Err
Definition AngCorr.hh:97
void DrawChiGraph(TCanvas *c1, double ylow, double yhigh)
Definition AngCorr.cc:370
TGraphErrors * GetGraph()
Definition AngCorr.cc:149
int Mixed(double s1, double s2)
Definition AngCorr.cc:68
double Ak(double k, double delta1, double delta2)
Definition AngCorr.cc:17
AngCorr(double spin1, double spin2, double spin3, GamR::AngDist::SolidAttenuation *Q)
Definition AngCorr.hh:55
GamR::AngDist::SolidAttenuation * Qk
Definition AngCorr.hh:39
double operator()(double *x, double *p)
Definition AngCorr.hh:63
AngCorr(double spin1, double spin2, double spin3)
Definition AngCorr.hh:47
Definition Gain.cc:19