Heavy ion Analysis Libriares
Loading...
Searching...
No Matches
CorrFitMath1DCF.cxx
1/*
2 * CorrFitMath1D.cxx
3 *
4 * Created on: 20 lis 2017
5 * Author: Daniel Wielanek
6 * E-mail: daniel.wielanek@gmail.com
7 * Warsaw University of Technology, Faculty of Physics
8 */
9#include "CorrFitMath1DCF.h"
10
11namespace Hal {
12 CorrFitMath1DCF::CorrFitMath1DCF(Int_t parameters_no, Int_t subfunc) : CorrFit1DCF(parameters_no), fSubFunctionsNo(subfunc) {
13 fSubFunctions = new CorrFit1DCF*[subfunc];
14 }
15
16 void CorrFitMath1DCF::Draw(Option_t* draw_option) {
17 CorrFit1DCF::Draw(draw_option);
18 Int_t line_style[5] = {2, 3, 7, 9, 10};
19 for (int i = 0; i < fSubFunctionsNo; i++) {
20 if (i < 5) {
21 GetSubFunction(i)->SetLineStyle(line_style[i]);
22 } else {
23 GetSubFunction(i)->SetLineColor(i - 2);
24 }
25 GetSubFunction(i)->Draw("SAME");
26 }
27 }
28
29 Double_t CorrFitMath1DCF::Call(Int_t subfucnt, const Double_t* x, const Double_t* params) const {
30 return fSubFunctions[subfucnt]->CalculateCF(x, params);
31 }
32
33 CorrFitMath1DCF::~CorrFitMath1DCF() { delete[] fSubFunctions; }
34
35 Double_t CorrFit1DCF_Sum::CalculateCF(const Double_t* x, const Double_t* params) const {
36 return (Call(0, x, params) + Call(1, x, params + GetSubFunction(0)->GetParametersNo()));
37 }
38
40 CorrFitMath1DCF(f1->GetParametersNo() + f2->GetParametersNo(), 2) {
41 SetSubfunction(f1, 0);
42 SetSubfunction(f2, 1);
43 Int_t par1 = f1->GetParametersNo();
44 for (int i = 0; i < f1->GetParametersNo(); i++) {
46 if (f1->IsParFixed(i)) {
47 FixParameter(i, f1->GetParMin(i));
48 } else {
49 SetParLimits(i, f1->GetParMin(i), f1->GetParMax(i));
50 }
51 }
52 for (int i = 0; i < f2->GetParametersNo(); i++) {
53 SetParameterName(i + par1, f2->GetParameterName(i));
54 if (f2->IsParFixed(i)) {
55 FixParameter(i + par1, f2->GetParMin(i));
56 } else {
57 SetParLimits(i + par1, f2->GetParMin(i), f2->GetParMax(i));
58 }
59 }
60 for (int i = 0; i < f1->GetParametersNo(); i++) {
61 for (int j = 0; j < f2->GetParametersNo(); j++) {
62 if (f1->GetParameterName(i) == f2->GetParameterName(j)) {
63 SetParameterName(i, GetParameterName(i) + "_{1}");
64 SetParameterName(par1 + j, GetParameterName(j + par1) + "_{2}");
65 }
66 }
67 }
68 }
69
70 CorrFit1DCF_Sum::~CorrFit1DCF_Sum() {}
71
72 Double_t CorrFit1DCF_SumRsame::CalculateCF(const Double_t* x, const Double_t* params) const {
73 for (int i = 0; i < GetSubFunction(1)->GetParametersNo(); i++) {
74 fNewParams[i] = params[i + GetSubFunction(0)->GetParametersNo()];
75 }
76 fNewParams[fSecondRadius] = params[fFirstRadius];
77 return Call(0, x, params) + Call(1, x, fNewParams);
78 }
79
81 CorrFitMath1DCF(f1->GetParametersNo() + f2->GetParametersNo() - 1, 2) {
82 SetSubfunction(f1, 0);
83 SetSubfunction(f2, 1);
85 fNewParams = new Double_t[fSecondParamtersNo];
86 Int_t par1 = f1->GetParametersNo();
87 for (int i = 0; i < f1->GetParametersNo(); i++) {
89 if (f1->IsParFixed(i)) {
90 FixParameter(i, f1->GetParMin(i));
91 } else {
92 SetParLimits(i, f1->GetParMin(i), f1->GetParMax(i));
93 }
94 if (i == RadiusID()) { fFirstRadius = i; }
95 }
96 for (int i = 0; i < f2->GetParametersNo(); i++) {
97 if (i == RadiusID()) { fSecondRadius = i; }
98 SetParameterName(i + par1, f2->GetParameterName(i));
99 if (i == RadiusID()) {
100 if (f2->IsParFixed(i)) {
101 FixParameter(i + par1, f1->GetParMin(i));
102 } else {
103 SetParLimits(i + par1, f1->GetParMin(i), f1->GetParMax(i));
104 }
105 } else {
106 if (f2->IsParFixed(i)) {
107 FixParameter(i + par1, f2->GetParMin(i));
108 } else {
109 SetParLimits(i + par1, f2->GetParMin(i), f2->GetParMax(i));
110 }
111 }
112 }
113 for (int i = 0; i < f1->GetParametersNo(); i++) {
114 if (i == RadiusID()) continue;
115 for (int j = 0; j < f2->GetParametersNo(); j++) {
116 if (j == RadiusID()) continue;
117 if ((f1->GetParameterName(i) == f2->GetParameterName(j))) {
118 SetParameterName(i, GetParameterName(i) + "_{1}");
119 SetParameterName(par1 + j, GetParameterName(j + par1) + "_{2}");
120 }
121 }
122 }
123 }
124
125 CorrFit1DCF_SumRsame::~CorrFit1DCF_SumRsame() { delete[] fNewParams; }
126
127 Double_t CorrFit1DCF_Multi::CalculateCF(const Double_t* x, const Double_t* params) const {
128 return Call(0, x, params) * Call(1, x, params + fFirstFuncParams);
129 }
130
131 CorrFit1DCF_Multi::CorrFit1DCF_Multi(CorrFit1DCF* sign, CorrFit1DCF* back) :
132 CorrFitMath1DCF(sign->GetParametersNo() + back->GetParametersNo(), 2), fFirstFuncParams(sign->GetParametersNo()) {
133 SetSubfunction(sign, 0);
134 SetSubfunction(back, 1);
135 TString background_name = back->ClassName();
136 if (background_name == "CorrFit1DCF_Poly") { back->FixParameter(0, 1); }
137 for (int i = 0; i < GetSubFunction(0)->GetParametersNo(); i++) {
138 SetParameterName(i, GetSubFunction(0)->GetParameterName(i));
139 if (GetSubFunction(0)->IsParFixed(i)) {
140 FixParameter(i, GetSubFunction(0)->GetParMin(i));
141 } else {
142 SetParLimits(i, GetSubFunction(0)->GetParMin(i), GetSubFunction(0)->GetParMax(i));
143 }
144 }
145 for (int i = 0; i < GetSubFunction(1)->GetParametersNo(); i++) {
146 SetParameterName(i + fFirstFuncParams, GetSubFunction(1)->GetParameterName(i));
147 if (GetSubFunction(1)->IsParFixed(i)) {
148 FixParameter(i + fFirstFuncParams, GetSubFunction(1)->GetParMin(i));
149 } else {
150 SetParLimits(i + fFirstFuncParams, GetSubFunction(1)->GetParMin(i), GetSubFunction(1)->GetParMax(i));
151 }
152 }
153 }
154
155 CorrFit1DCF_Multi::~CorrFit1DCF_Multi() {}
156
159 for (int i = 0; i < GetSubfucntionsNo(); i++) {
160 GetSubFunction(i)->Check();
161 }
162 }
163} // namespace Hal
virtual Double_t CalculateCF(const Double_t *x, const Double_t *params) const
virtual Double_t CalculateCF(const Double_t *x, const Double_t *params) const
CorrFit1DCF_SumRsame(CorrFit1DCF *f1=NULL, CorrFit1DCF *f2=NULL)
CorrFit1DCF_Sum(CorrFit1DCF *f1=NULL, CorrFit1DCF *f2=NULL)
virtual Double_t CalculateCF(const Double_t *x, const Double_t *params) const
Int_t RadiusID() const
virtual Double_t CalculateCF(const Double_t *, const Double_t *) const
Definition CorrFit1DCF.h:93
virtual void Check()
Int_t GetSubfucntionsNo() const
void SetSubfunction(CorrFit1DCF *f, Int_t i)
Double_t Call(Int_t subfucnt, const Double_t *x, const Double_t *params) const
Int_t GetParametersNo() const
Definition CorrFit.h:269
TString GetParameterName(Int_t no) const
Definition CorrFit.h:280
Double_t GetParMax(Int_t par) const
Definition CorrFit.cxx:136
Double_t GetParMin(Int_t par) const
Definition CorrFit.cxx:131
void SetLineColor(Color_t color)
Definition CorrFit.h:196
void FixParameter(Int_t par, Double_t val)
Definition CorrFit.cxx:106
void SetParameterName(Int_t par, TString name)
Definition CorrFit.cxx:110
void SetParLimits(Int_t par, Double_t min, Double_t max)
Definition CorrFit.cxx:102
void SetLineStyle(Style_t style)
Definition CorrFit.h:201
Bool_t IsParFixed(Int_t par) const
Definition CorrFit.cxx:121