Heavy ion Analysis Libriares
Loading...
Searching...
No Matches
CorrFitMath3DCF.cxx
1/*
2 * CorrFit1DCFSum.cxx
3 *
4 * Created on: 21 mar 2016
5 * Author: Daniel Wielanek
6 * E-mail: daniel.wielanek@gmail.com
7 * Warsaw University of Technology, Faculty of Physics
8 */
9#include "CorrFitMath3DCF.h"
10
11namespace Hal {
13 CorrFitMath3DCF(f1->GetParametersNo() + f2->GetParametersNo()) {
14 SetSubfunction(f1, 0);
15 SetSubfunction(f2, 1);
16 // SetParameterName(GetParametersNo()-1,"N");
17 Int_t par1 = f1->GetParametersNo();
18 for (int i = 0; i < f1->GetParametersNo(); i++) {
20 if (f1->IsParFixed(i)) {
21 FixParameter(i, f1->GetParMin(i));
22 } else {
23 SetParLimits(i, f1->GetParMin(i), f1->GetParMax(i));
24 }
25 }
26 for (int i = 0; i < f2->GetParametersNo(); i++) {
27 SetParameterName(i + par1, f2->GetParameterName(i));
28 if (f2->IsParFixed(i)) {
29 FixParameter(i + par1, f2->GetParMin(i));
30 } else {
31 SetParLimits(i + par1, f2->GetParMin(i), f2->GetParMax(i));
32 }
33 }
34 for (int i = 0; i < f1->GetParametersNo(); i++) {
35 for (int j = 0; j < f2->GetParametersNo(); j++) {
36 if (f1->GetParameterName(i) == f2->GetParameterName(j)) {
37 SetParameterName(i, f1->GetParameterName(i) + "_{1}");
38 SetParameterName(par1 + j, f2->GetParameterName(j) + "_{2}");
39 }
40 }
41 }
42 }
43
44 Double_t CorrFit3DCF_Sum::CalculateCF(const Double_t* x, const Double_t* params) const {
45 return Call(0, x, params) + Call(1, x, params + GetSubFunction(0)->GetParametersNo());
46 }
47
48 CorrFit3DCF_Sum::~CorrFit3DCF_Sum() {}
49
50 CorrFitMath3DCF::CorrFitMath3DCF(Int_t parameters_no, Int_t subfunct) : CorrFit3DCF(parameters_no), fSubFunctionsNo(subfunct) {
51 fSubFunction = new CorrFit3DCF*[fSubFunctionsNo];
52 }
53
54 Double_t CorrFitMath3DCF::Call(Int_t subfucnt, const Double_t* x, const Double_t* params) const {
55 return fSubFunction[subfucnt]->CalculateCF(x, params);
56 }
57
58 CorrFitMath3DCF::~CorrFitMath3DCF() { delete[] fSubFunction; }
59
61 for (int i = 0; i < GetSubunctionsNo(); i++) {
62 for (int j = 0; j < GetSubFunction(i)->GetParametersNo(); j++) {
63 GetSubFunction(i)->fTempParamsEval[j] = fTempParamsEval[j];
64 }
65 GetSubFunction(i)->ParametersChanged();
66 }
67 }
68} // namespace Hal
virtual Double_t CalculateCF(const Double_t *x, const Double_t *params) const
CorrFit3DCF_Sum(CorrFit3DCF *f1, CorrFit3DCF *f2)
virtual Double_t CalculateCF(const Double_t *x, const Double_t *params) const =0
void ParametersChanged() const
void ParametersChanged() const
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
Double_t * fTempParamsEval
Definition CorrFit.h:109
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
Bool_t IsParFixed(Int_t par) const
Definition CorrFit.cxx:121