Heavy ion Analysis Libriares
Loading...
Searching...
No Matches
CorrFit1DCFCumac.cxx
1/*
2 * CorrFIt1DCFCumac.cxx
3 *
4 * Created on: 13 mar 2018
5 * Author: Daniel Wielanek
6 * E-mail: daniel.wielanek@gmail.com
7 * Warsaw University of Technology, Faculty of Physics
8 */
9#include "CorrFit1DCFCumac.h"
10
11#include "Const.h"
12
13#include <TComplex.h>
14#include <TF1.h>
15#include <TMath.h>
16
17namespace Hal {
18 CorrFit1DCFCumac::CorrFit1DCFCumac() : CorrFit1DCFCumac(2) {}
19
20 CorrFit1DCFCumac::CorrFit1DCFCumac(Int_t params) : CorrFit1DCF(params), fPis(TMath::Sqrt(TMath::Pi())), fHc(.1973), fZ(0) {}
21
22 CorrFit1DCFCumac::~CorrFit1DCFCumac() {}
23
24 Double_t CorrFit1DCFCumac::F1(Double_t z) const {
25 fZ = z;
26 return Integr();
27 }
28 Double_t CorrFit1DCFCumac::F2(Double_t z) const { return (1.0 - TMath::Exp(-z * z)) / z; }
29
30 Double_t CorrFit1DCFCumac::Ff1(Double_t x) const { return TMath::Exp(x * x - fZ * fZ) / fZ; }
31
32 Double_t CorrFit1DCFCumac::Integr() const {
33 Int_t n;
34 Double_t xx, de1, a1h, ddd, abs__, eps;
35 eps = 1E-18;
36 ddd = .005;
37 n = 100;
38 de1 = fZ / n;
39 Double_t zero = 0.0;
40 return Simps(zero, fZ, de1, ddd, eps, xx, a1h, abs__);
41 }
42
43 Double_t CorrFit1DCFCumac::Simps(Double_t a1,
44 Double_t b1,
45 Double_t h1,
46 Double_t reps1,
47 Double_t aeps1,
48 Double_t x,
49 Double_t aih,
50 Double_t aiabs) const {
51 Double_t r__1;
52 Double_t result = 0;
53
54 Double_t a, b, c__, f[7], h__;
55 Int_t k;
56 Double_t p[5], s, x0, di1, di2, di3, eps, aeps, reps, delta;
57
58 r__1 = b1 - a1;
59 h__ = Sign(h1, r__1);
60 s = Sign(1.0, h__);
61 a = a1;
62 b = b1;
63 aih = 0.f;
64 aiabs = 0.f;
65 p[1] = 4.f;
66 p[3] = 4.f;
67 p[2] = 2.f;
68 p[4] = 1.f;
69 if (b - a != 0.f) {
70 goto L1;
71 } else {
72 goto L2;
73 }
74 L1:
75 reps = TMath::Abs(reps1);
76 aeps = TMath::Abs(aeps1);
77 for (k = 1; k <= 7; ++k) {
78 /* L3: */
79 f[k - 1] = 1e17f;
80 }
81 x = a;
82 c__ = 0.f;
83 f[0] = Ff1(x) / 3;
84 L4:
85 x0 = x;
86 if ((x0 + h__ * 4 - b) * s <= 0.f) {
87 goto L5;
88 } else {
89 goto L6;
90 }
91 L6:
92 h__ = (b - x0) / 4;
93 if (h__ != 0.f) {
94 goto L7;
95 } else {
96 goto L2;
97 }
98 L7:
99 for (k = 2; k <= 7; ++k) {
100 /* L8: */
101 f[k - 1] = 1e17f;
102 }
103 c__ = 1.f;
104 L5:
105 di2 = f[0];
106 di3 = TMath::Abs(f[0]);
107 for (k = 2; k <= 5; ++k) {
108 x += h__;
109 if ((x - b) * s >= 0.f) {
110 goto L24;
111 } else {
112 goto L23;
113 }
114 L24:
115 x = b;
116 L23:
117 if (f[k - 1] - 1e17f != 0.f) {
118 goto L10;
119 } else {
120 goto L11;
121 }
122 L11:
123 f[k - 1] = Ff1(x) / 3.0;
124 L10:
125 di2 += p[k - 1] * f[k - 1];
126 /* L9: */
127 di3 += p[k - 1] * (r__1 = f[k - 1], TMath::Abs(r__1));
128 }
129 di1 = (f[0] + f[2] * 4 + f[4]) * 2 * h__;
130 di2 *= h__;
131 di3 *= h__;
132 if (reps != 0.f) {
133 goto L12;
134 } else {
135 goto L13;
136 }
137 L13:
138 if (aeps != 0.f) {
139 goto L12;
140 } else {
141 goto L14;
142 }
143 L12:
144 eps = (r__1 = (aiabs + di3) * reps, TMath::Abs(r__1));
145 if (eps - aeps >= 0.f) {
146 goto L16;
147 } else {
148 goto L15;
149 }
150 L15:
151 eps = aeps;
152 L16:
153 delta = (r__1 = di2 - di1, TMath::Abs(r__1));
154 if (delta - eps >= 0.f) {
155 goto L21;
156 } else {
157 goto L20;
158 }
159 L20:
160 if (delta - eps / 8 >= 0.f) {
161 goto L14;
162 } else {
163 goto L17;
164 }
165 L17:
166 h__ *= 2;
167 f[0] = f[4];
168 f[1] = f[5];
169 f[2] = f[6];
170 for (k = 4; k <= 7; ++k) {
171 /* L19: */
172 f[k - 1] = 1e17f;
173 }
174 goto L18;
175 L14:
176 f[0] = f[4];
177 f[2] = f[5];
178 f[4] = f[6];
179 f[1] = 1e17f;
180 f[3] = 1e17f;
181 f[5] = 1e17f;
182 f[6] = 1e17f;
183 L18:
184 di1 = di2 + (di2 - di1) / 15;
185 result += di1;
186 aih += di2;
187 aiabs += di3;
188 goto L22;
189 L21:
190 h__ /= 2;
191 f[6] = f[4];
192 f[5] = f[3];
193 f[4] = f[2];
194 f[2] = f[1];
195 f[1] = 1e17f;
196 f[3] = 1e17f;
197 x = x0;
198 c__ = 0.f;
199 goto L5;
200 L22:
201 if (c__ != 0.f) {
202 goto L2;
203 } else {
204 goto L4;
205 }
206 L2:
207 return result;
208 }
209
210 Double_t CorrFit1DCFCumac::Sign(Double_t a, Double_t b) const {
211 if (b >= 0) return TMath::Abs(a);
212 return -TMath::Abs(a);
213 }
214 //*****************************************************************************//
215
216 CorrFIt1DCFCumacLamLam::CorrFIt1DCFCumacLamLam() : CorrFit1DCFCumac(8) {}
217
218 Double_t CorrFIt1DCFCumacLamLam::CalculateCF(const Double_t* x, const Double_t* params) const {
219 Double_t d0, f0, r0, ai, ak, al, an, ex, f0m, r0m, zs, den, fim, fre, fsi, r_s__, r_t__, fsq, pol, r0res, alres;
220 Double_t ret_val, r__1;
221 switch (fKinematics) {
222 case Femto::EKinematics::kPRF: ak = x[0]; break;
223 case Femto::EKinematics::kLCMS: ak = x[0] * 0.5; break;
224 default: ak = 0.0; break;
225 }
226 an = params[NormID()];
227 al = params[LambdaID()];
228 r0 = params[RadiusID()];
229 f0 = -params[ScatteringLengthID()];
230 d0 = params[EffectiveRadiusID()];
231 pol = params[LambdaPolarizationID()];
232 alres = params[ResidualAmplitudeID()];
233 r0res = params[ResidualGaussWidhtID()];
234 f0m = 1E+6;
235 r0m = 1E+6;
236 if (f0 != 0.0) { f0m = .1973 / f0 + d0 * .5 * ak * ak / .1973; }
237 if (r0 != 0.0) { r0m = .1973 / r0; }
238 den = f0m * f0m + ak * ak;
239 fre = f0m / den;
240 fim = ak / den;
241 fsq = fre * fre + fim * fim;
242 fZ = (ak * 2.0) / r0m;
243 zs = fZ * fZ;
244 ex = TMath::Exp(-zs);
245 fsi = 0.0;
246 if (params[RadiusID()] != 0.0) {
247 ai = Integr();
248 // integr_(ai, fZ_1);
249 fsi = fsq * (r0m * r0m) * (1 - d0 / (r0 * 3.5449077020000002f)) + fre * 4 * r0m * ai / 1.772453851f;
250 fsi -= fim * 2 * r0m * (1 - ex) / fZ;
251 }
252 r_t__ = (pol * pol + 3) * .25f * (1 - al * ex);
253 r_s__ = (1 - pol * pol) * .25f * (al * (ex + fsi) + 1);
254 r__1 = ak * 2.0 * r0res / .1973f;
255 ret_val = an * (r_s__ + r_t__ + alres * exp(-(r__1 * r__1)));
256 return ret_val;
257 }
258
259 CorrFIt1DCFCumacLamLam::~CorrFIt1DCFCumacLamLam() {}
260 //******************************************************************************
261 CorrFit1DCFCumacPLam::CorrFit1DCFCumacPLam() : CorrFit1DCFCumac(8) {
262 SetParameterName(SingletScatteringLengthID(), "f_{0s}");
263 SetParameterName(SingletEffectiveRadiusID(), "d_{0s}");
264 SetParameterName(TripletScatteringLenghtID(), "f_{0t}");
265 SetParameterName(TripletEffectiveRadiusID(), "d_{0t}");
266 SetParameterName(LambdaPolarizationID(), "POL");
267 FixParameter(SingletScatteringLengthID(), 2.31);
268 FixParameter(SingletEffectiveRadiusID(), 3.04);
269 FixParameter(TripletScatteringLenghtID(), 1.78);
270 FixParameter(TripletEffectiveRadiusID(), 3.22);
271 FixParameter(LambdaPolarizationID(), 0);
272 }
273
274 CorrFit1DCFCumacPLam::~CorrFit1DCFCumacPLam() {}
275
276 Double_t CorrFit1DCFCumacPLam::CalculateCF(const Double_t* x, const Double_t* params) const {
277 Double_t ak;
278 switch (fKinematics) {
279 case Femto::EKinematics::kPRF: ak = x[0]; break;
280 case Femto::EKinematics::kLCMS: ak = x[0] * 0.5; break;
281 default: ak = 0.0; break;
282 }
283
284 Double_t al = params[LambdaID()];
285 // c r0=p(2)
286 Double_t r0 = params[RadiusID()];
287 Double_t f0_s = params[SingletScatteringLengthID()];
288 Double_t d0_s = params[SingletEffectiveRadiusID()];
289 Double_t f0_t = params[TripletScatteringLenghtID()];
290 Double_t d0_t = params[TripletEffectiveRadiusID()];
291 Double_t POL = params[LambdaPolarizationID()];
292
293 Double_t r0m = fHc / r0;
294 Double_t f0m_s = (fHc / f0_s + 0.5 * ak * ak * d0_s / fHc);
295 Double_t den_s = TMath::Power(f0m_s, 2) + TMath::Power(ak, 2);
296 Double_t fre_s = f0m_s / den_s;
297 Double_t fim_s = ak / den_s;
298 Double_t fsq_s = TMath::Power(fre_s, 2) + TMath::Power(fim_s, 2);
299
300 Double_t f0m_t = (fHc / f0_t + 0.5 * ak * ak * d0_t / fHc);
301 Double_t den_t = TMath::Power(f0m_t, 2) + TMath::Power(ak, 2);
302 Double_t fre_t = f0m_t / den_t;
303 Double_t fim_t = ak / den_t;
304 Double_t fsq_t = TMath::Power(fre_t, 2) + TMath::Power(fim_t, 2);
305
306 fZ = 2. * ak / r0m;
307 Double_t zs = fZ * fZ;
308 Double_t ex = TMath::Exp(-zs);
309 Double_t ai = Integr();
310 Double_t fsi_s = 0.5 * (fsq_s * TMath::Power(r0m, 2) * (1. - d0_s / (2. * fPis * r0)) + 4. * fre_s * r0m * ai / fPis);
311 fsi_s = fsi_s - fim_s * r0m * (1. - ex) / fZ;
312 Double_t fsi_t = 0.5 * (fsq_t * TMath::Power(r0m, 2) * (1. - d0_t / (2. * fPis * r0)) + 4. * fre_t * r0m * ai / fPis);
313 fsi_t = fsi_t - fim_t * r0m * (1. - ex) / fZ;
314 Double_t R_t = .25 * (3. + TMath::Power(POL, 2)) * (1. + al * fsi_t);
315 Double_t R_s = .25 * (1. - TMath::Power(POL, 2)) * (1. + al * fsi_s);
316 return params[NormID()] * (R_t + R_s);
317 }
318 /******************************************************************************/
319 Double_t CorrFit1DCFCumacK0K0::CalculateCF(const Double_t* x, const Double_t* params) const {
320 Double_t ak = 0;
321 switch (fKinematics) {
322 case Femto::EKinematics::kPRF: ak = x[0]; break;
323 case Femto::EKinematics::kLCMS: ak = x[0] * 0.5; break;
324 default: ak = 0.0; break;
325 }
326 const Double_t mK2 = Const::KaonZeroMass() * Const::KaonZeroMass();
327 const Double_t mpi2 = Const::PionPlusMass() * Const::PionPlusMass();
328 const Double_t meta2 = 0.547862 * 0.547862;
329 const Double_t ak2 = ak * ak;
330 Double_t s = 4. * (ak2 + mK2);
331 Double_t U = TMath::Sqrt(s);
332 Double_t k1prim = TMath::Sqrt(ak * ak + mK2 - mpi2);
333 Double_t e = TMath::Sqrt(ak2 + mK2) * 2.0;
334 // Double_t k3prim =
335 // TMath::Sqrt(mpi2 * mpi2 - 2.0 * mpi2 * meta2 - 2.0 * e * e * mpi2 + meta2 * meta2 - 2.0 * e * e * meta2 + e * e * e * e)
336 // / (2.0 * e);
337
338 Double_t k2prim = TMath::Sqrt(mpi2 * mpi2 + meta2 * meta2 + s * s - 2.0 * (mpi2 * meta2 + mpi2 * s + meta2 * s)) / (2.0 * U);
339
340 TComplex num1(params[Mf0ID()] * params[Mf0ID()] - s, -params[Gamma_f0KKID()] * ak - params[Gamma_f0pipiID()] * k1prim);
341 TComplex num2(params[Ma0ID()] * params[Ma0ID()] - s, -params[Gamma_a0KKID()] * ak - params[Gamma_a0PiEtaID()] * k2prim);
342
343 TComplex f_1(params[Gamma_f0KKID()], 0);
344 TComplex f_2(params[Gamma_a0KKID()], 0);
345 f_1 = f_1 / num1;
346 f_2 = f_2 / num2;
347 TComplex fk = 0.5 * (f_1 + f_2);
348 Double_t r = params[RadiusID()] * Femto::FmToGeV();
349 Double_t qr = ak * r * 2.;
350 Double_t fkM = (fk / r).Rho();
351
352 Double_t strong =
353 params[AssymetryID()] * (fkM * fkM + 4. * fk.Re() / (fPis * r) * F1(qr) - 2.0 * fk.Im() * F2(qr) / (r)); // brakowalo pi?
354 return params[NormID()] * (1 + params[LambdaID()] * (TMath::Exp(-qr * qr) + strong));
355 }
356
357 CorrFit1DCFCumacK0K0::CorrFit1DCFCumacK0K0() : CorrFit1DCFCumac(10) {
358 SetParameterName(Mf0ID(), "m_{f_0}");
359 SetParameterName(Ma0ID(), "m_{a_0}");
360 SetParameterName(Gamma_f0KKID(), "#gamma_{f_{0K#bar{K}}");
361 SetParameterName(Gamma_a0KKID(), "#gamma_{a_{0K#bar{K}}");
362 SetParameterName(Gamma_f0pipiID(), "#gamma_{f_{0#pi#pi}");
363 SetParameterName(Gamma_a0PiEtaID(), "#gamma_{a_{0#pi#eta}}");
364 SetParameterName(AssymetryID(), "#alpha");
365 SetDefParams(3);
366 }
367
369 Double_t params[5][6] = {{0.978, 0.792, 0.199, 0.974, 0.333, 0.222},
370 {0.973, 2.763, 0.5283, 0.985, 0.4038, 0.3711},
371 {0.996, 1.305, 0.2684, 0.992, 0.5555, 0.4401},
372 {0.996, 1.305, 0.2684, 1.003, 0.8365, 0.4580},
373 {}};
374 for (int i = 0; i < 6; i++) {
375 Double_t val = 0;
376 for (int j = 0; j < 4; j++)
377 val += params[j][i];
378 params[4][i] = val / 4.0;
379 }
380
381 if (opt > 4) opt = 4;
382 FixParameter(Mf0ID(), params[opt][0]);
383 FixParameter(Gamma_f0KKID(), params[opt][1]);
384 FixParameter(Gamma_f0pipiID(), params[opt][2]);
385 FixParameter(Ma0ID(), params[opt][3]);
386 FixParameter(Gamma_a0KKID(), params[opt][4]);
387 FixParameter(Gamma_a0PiEtaID(), params[opt][5]);
389 }
390
391 //==================================================================
392 Double_t CorrFit1DCFCumacK0Kch::CalculateCF(const Double_t* x, const Double_t* params) const {
393 Double_t ak = 0;
394 switch (fKinematics) {
395 case Femto::EKinematics::kPRF: ak = x[0]; break;
396 case Femto::EKinematics::kLCMS: ak = x[0] * 0.5; break;
397 default: ak = 0.0; break;
398 }
399 const Double_t mK2 = Const::KaonZeroMass() * Const::KaonZeroMass();
400 const Double_t mpi2 = Const::PionPlusMass() * Const::PionPlusMass();
401 const Double_t meta2 = 0.547862 * 0.547862;
402 const Double_t ak2 = ak * ak;
403 Double_t s = 4. * (ak2 + mK2);
404 Double_t U = TMath::Sqrt(s);
405 Double_t k1prim = TMath::Sqrt(ak * ak + mK2 - mpi2);
406 Double_t e = TMath::Sqrt(ak2 + mK2) * 2.0;
407 Double_t k3prim =
408 TMath::Sqrt(mpi2 * mpi2 - 2.0 * mpi2 * meta2 - 2.0 * e * e * mpi2 + meta2 * meta2 - 2.0 * e * e * meta2 + e * e * e * e)
409 / (2.0 * e);
410
411 Double_t k2prim = TMath::Sqrt(mpi2 * mpi2 + meta2 * meta2 + s * s - 2.0 * (mpi2 * meta2 + mpi2 * s + meta2 * s)) / (2.0 * U);
412 double k2prim_old = k2prim;
413 k2prim = fGammaCalc.Calculate(ak);
414 std::cout << "k- " << k2prim_old << " " << k2prim << std::endl;
415 TComplex num(params[Ma0ID()] * params[Ma0ID()] - s, -params[Gamma_a0KKID()] * ak - params[Gamma_a0PiEtaID()] * k2prim);
416 TComplex f(params[Gamma_a0KKID()], 0);
417 f = f / num;
418 TComplex fk = f;
419 Double_t r = params[RadiusID()] * Femto::FmToGeV();
420 Double_t qr = ak * r * 2.;
421 Double_t fkM = (fk / r).Rho();
422
423 Double_t strong = 0.25 * (fkM * fkM + 4.0 * fk.Re() / (fPis * r) * F1(qr) - 2.0 * fk.Im() * F2(qr) / (r));
424
425 double val = params[NormID()] * (1 + params[LambdaID()] * strong);
426 std::cout << ak << " fk " << fk.Re() << " " << fk.Im() << " " << F1(qr) << std::endl;
427
428 return val;
429 }
430
431 CorrFit1DCFCumacK0Kch::CorrFit1DCFCumacK0Kch() : CorrFit1DCFCumac(10) {
432 SetParameterName(Ma0ID(), "m_{a_0}");
433 SetParameterName(Gamma_a0KKID(), "#gamma_{a_{0K#bar{K}}");
434 SetParameterName(Gamma_a0PiEtaID(), "#gamma_{a_{0#pi#eta}}");
435 SetDefParams(3);
436 }
437
439 Double_t params[5][3] = {
440 {0.985, 0.4038, 0.3711}, {0.992, 0.5555, 0.4401}, {1.003, 0.8365, 0.4580}, {0.974, 0.3330, 0.2220}, {}};
441
442 //.9467,.9698,2.763,
443 for (int i = 0; i < 3; i++) {
444 Double_t val = 0;
445 for (int j = 0; j < 4; j++)
446 val += params[j][i];
447 params[4][i] = val / 4.0;
448 }
449
450 if (opt > 4) opt = 4;
451 fGammaCalc.ReInit(Const::KaonZeroMass(), Const::KaonZeroMass(), Const::PionPlusMass(), 0.547862);
452 FixParameter(Gamma_a0KKID(), params[opt][1]);
453 FixParameter(Gamma_a0PiEtaID(), params[opt][2]);
454 FixParameter(Ma0ID(), params[opt][0]);
455 for (int i = 0; i < GetParametersNo(); i++) {
456 fTempParamsEval[i] = fParameters[i].GetMin();
457 }
458 }
459
460 Double_t CorrFit1DCFCumacDLam::CalculateCF(const Double_t* x, const Double_t* params) const {
461 Double_t ak = 0;
462 switch (fKinematics) {
463 case Femto::EKinematics::kPRF: ak = x[0]; break;
464 case Femto::EKinematics::kLCMS: ak = x[0] * 0.5; break;
465 default: ak = 0.0; break;
466 }
467
468 Double_t f0_d = params[DoubletScatteringLengthID()];
469 Double_t d0_d = params[DoubletEffectiveRadiusID()];
470 Double_t f0_q = params[QuartetScatteringLenghtID()];
471 Double_t d0_q = params[QuartetEffectiveRadiusID()];
472 Double_t r0 = Femto::FmToGeV(params[RadiusID()]);
473 Double_t q = 2.0 * ak;
474 auto component = [&](double f0, double d0) {
475 f0 = Femto::FmToGeV(f0);
476 d0 = Femto::FmToGeV(d0);
477 // TComplex tmp = TComplex(d0 * ak * ak * 0.5, -ak) - 1.0 / f0; // https://arxiv.org/pdf/2005.05012
478 TComplex tmp = TComplex(d0 * ak * ak * 0.5, -ak) + 1.0 / f0; // Yu Hu
479 auto fk = TComplex(1.0, 0) / tmp;
480 return (fk * fk).Rho() * F(d0, r0) / (2.0 * r0 * r0) + 2.0 * fk.Re() / (fPis * r0) * F1(q * r0) - fk.Im() / r0 * F2(q * r0);
481 };
482 // std::cout << "FK " << component(f0_d, d0_d) << " " << component(f0_q, d0_q) << std::endl;
483 Double_t strong = (component(f0_d, d0_d) + 2.0 * component(f0_q, d0_q)) / 3.0;
484 return params[NormID()] * (1.0 + params[LambdaID()] * strong);
485 }
486
487 CorrFit1DCFCumacDLam::CorrFit1DCFCumacDLam() : CorrFit1DCFCumac(7) {
488 SetParameterName(DoubletScatteringLengthID(), "f_{0d}");
489 SetParameterName(DoubletEffectiveRadiusID(), "d_{0d}");
490 SetParameterName(QuartetScatteringLenghtID(), "f_{0q}");
491 SetParameterName(QuartetEffectiveRadiusID(), "d_{0q}");
492 FixParameter(DoubletScatteringLengthID(), 2.31);
493 FixParameter(DoubletEffectiveRadiusID(), 3.04);
494 FixParameter(QuartetScatteringLenghtID(), 1.78);
495 FixParameter(QuartetEffectiveRadiusID(), 3.22);
496 }
497
498 Double_t CorrFit1DCFCumacDLam::F(Double_t d, Double_t r0) const { return 1.0 - d / (2.0 * fPis * r0); }
499
500
501 Double_t CorrFit1DCFCumac::Get(Double_t q, Double_t r) {
502 Double_t val[1];
503 val[0] = q;
504 Double_t params[GetParametersNo()];
505 for (int i = 0; i < GetParametersNo(); i++)
506 params[i] = GetParameter(i);
508 return CalculateCF(val, params);
509 }
510
511} // 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
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
virtual Double_t CalculateCF(const Double_t *x, const Double_t *params) const
virtual Double_t Get(Double_t q, Double_t r)
Int_t LambdaID() const
Int_t RadiusID() const
virtual Double_t CalculateCF(const Double_t *, const Double_t *) const
Definition CorrFit1DCF.h:93
Int_t NormID() const
Femto::EKinematics fKinematics
Definition CorrFitFunc.h:69
Int_t GetParametersNo() const
Definition CorrFit.h:269
Double_t GetParameter(Int_t par) const
Definition CorrFit.cxx:141
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
std::vector< FitParam > fParameters
Definition CorrFit.h:103
void ReInit(Double_t a, Double_t b, Double_t c, Double_t d)
Double_t Calculate(Double_t kstar) const