GSI Object Oriented Online Offline (Go4)  GO4-6.3.0
TGo4Fitter.cxx
Go to the documentation of this file.
1 // $Id$
2 //-----------------------------------------------------------------------
3 // The GSI Online Offline Object Oriented (Go4) Project
4 // Experiment Data Processing at EE department, GSI
5 //-----------------------------------------------------------------------
6 // Copyright (C) 2000- GSI Helmholtzzentrum fuer Schwerionenforschung GmbH
7 // Planckstr. 1, 64291 Darmstadt, Germany
8 // Contact: http://go4.gsi.de
9 //-----------------------------------------------------------------------
10 // This software can be used under the license agreements as stated
11 // in Go4License.txt file which is part of the distribution.
12 //-----------------------------------------------------------------------
13 
14 #include "TGo4Fitter.h"
15 
16 #include <cstdlib>
17 #include <iostream>
18 #include <iomanip>
19 
20 #include "TClass.h"
21 #include "TMath.h"
22 #include "TH1.h"
23 #include "TGraph.h"
24 #include "TCanvas.h"
25 #include "TObjArray.h"
26 #include "TObjString.h"
27 #include "TROOT.h"
28 #include "THStack.h"
29 #include "TArrayD.h"
30 #include "TBuffer.h"
31 
32 #include "TGo4FitDataHistogram.h"
33 #include "TGo4FitDataGraph.h"
34 #include "TGo4FitModelGauss1.h"
35 #include "TGo4FitModelPolynom.h"
36 #include "TGo4FitAmplEstimation.h"
37 #include "TGo4FitSlot.h"
38 
40  : TGo4FitterAbstract(), fxDatas(), fxModels(), fiFitFunctionType(0), fiMemoryUsage(0), fxUserFitFunction(nullptr),
41  fxDrawObjs(nullptr)
42 {
43 }
44 
45 TGo4Fitter::TGo4Fitter(const char *iName, const char *iTitle)
46  : TGo4FitterAbstract(iName, iTitle), fxDatas(), fxModels(), fiFitFunctionType(0), fiMemoryUsage(0),
47  fxUserFitFunction(nullptr), fxDrawObjs(nullptr)
48 {
49  fxDatas.SetOwner(kTRUE);
50  fxModels.SetOwner(kTRUE);
51 }
52 
53 TGo4Fitter::TGo4Fitter(const char *iName, Int_t iFitFunctionType, Bool_t IsAddStandardActions)
54  : TGo4FitterAbstract(iName, "TGo4Fitter object"), fxDatas(), fxModels(), fiFitFunctionType(0), fiMemoryUsage(100),
55  fxUserFitFunction(nullptr), fxDrawObjs(nullptr)
56 {
57  fxDatas.SetOwner(kTRUE);
58  fxModels.SetOwner(kTRUE);
59  SetFitFunctionType(iFitFunctionType);
60  if (IsAddStandardActions)
62 }
63 
65 {
68 }
69 
70 void TGo4Fitter::SetMemoryUsage(Int_t iMemoryUsage)
71 {
72  if (iMemoryUsage < 0)
73  fiMemoryUsage = 0;
74  else if (iMemoryUsage > 100)
75  fiMemoryUsage = 100;
76  else
77  fiMemoryUsage = iMemoryUsage;
78 }
79 
81 {
83 
84  for (Int_t ndata = 0; ndata < GetNumData(); ndata++)
85  GetData(ndata)->CollectParsTo(*this);
86 
87  for (Int_t nmodel = 0; nmodel < GetNumModel(); nmodel++) {
88  TGo4FitModel *model = GetModel(nmodel);
89  for (Int_t n = 0; n < model->NumAssigments(); n++)
90  if (FindData(model->AssignmentName(n))) {
91  model->CollectParsTo(*this);
92  break;
93  }
94  }
95 }
96 
97 void TGo4Fitter::Clear(Option_t *option)
98 {
100  DeleteAllData();
101  DeleteAllModels();
102 }
103 
105 {
106  return (n >= 0) && (n < GetNumData()) ? dynamic_cast<TGo4FitData *>(fxDatas[n]) : nullptr;
107 }
108 
109 const char *TGo4Fitter::GetDataName(Int_t n)
110 {
111  return GetData(n) ? GetData(n)->GetName() : nullptr;
112 }
113 
114 TGo4FitData *TGo4Fitter::FindData(const char *DataName)
115 {
116  return !DataName ? nullptr : dynamic_cast<TGo4FitData *>(fxDatas.FindObject(DataName));
117 }
118 
120 {
121  fxDatas.Add(data);
124  return data;
125 }
126 
127 TGo4FitDataHistogram *TGo4Fitter::AddH1(const char *DataName, TH1 *histo, Bool_t Owned, Double_t lrange, Double_t rrange)
128 {
129  TGo4FitDataHistogram *data = new TGo4FitDataHistogram(DataName, histo, Owned);
130  if ((lrange < rrange) || (rrange != 0.))
131  data->SetRange(0, lrange, rrange);
132  AddData(data);
133  return data;
134 }
135 
136 TGo4FitDataHistogram *TGo4Fitter::SetH1(const char *DataName, TH1 *histo, Bool_t Owned)
137 {
138  TGo4FitDataHistogram *data = dynamic_cast<TGo4FitDataHistogram *>(FindData(DataName));
139  if (data)
140  data->SetHistogram(histo, Owned);
141  return data;
142 }
143 
144 TGo4FitDataGraph *TGo4Fitter::AddGraph(const char *DataName, TGraph *gr, Bool_t Owned, Double_t lrange, Double_t rrange)
145 {
146  TGo4FitDataGraph *data = new TGo4FitDataGraph(DataName, gr, Owned);
147  if ((lrange < rrange) || (rrange != 0.))
148  data->SetRange(0, lrange, rrange);
149  AddData(data);
150  return data;
151 }
152 
153 TGo4FitDataGraph *TGo4Fitter::SetGraph(const char *DataName, TGraph *gr, Bool_t Owned)
154 {
155  TGo4FitDataGraph *data = dynamic_cast<TGo4FitDataGraph *>(FindData(DataName));
156  if (data)
157  data->SetGraph(gr, Owned);
158  return data;
159 }
160 
162 {
163  if (comp)
164  for (Int_t n = 0; n < comp->NumSlots(); n++) {
165  TGo4FitSlot *slot = comp->GetSlot(n);
166  ClearSlot(slot, kFALSE);
167  for (Int_t n2 = 0; n2 < NumSlots(); n2++) {
168  TGo4FitSlot *slot2 = GetSlot(n2);
169  if (slot2->GetConnectedSlot() == slot)
170  slot2->ClearConnectionToSlot();
171  }
172  }
173 }
174 
175 TGo4FitData *TGo4Fitter::RemoveData(const char *DataName, Bool_t IsDel)
176 {
177  TGo4FitData *dat = FindData(DataName);
178  if (dat) {
181  fxDatas.Remove(dat);
182  if (IsDel) {
184  delete dat;
185  dat = nullptr;
186  }
187  fxDatas.Compress();
188  }
189  return dat;
190 }
191 
193 {
194  fxDatas.Delete();
195  fxDatas.Compress();
198 }
199 
201 {
202  return (n >= 0) && (n < GetNumModel()) ? dynamic_cast<TGo4FitModel *>(fxModels[n]) : nullptr;
203 }
204 
205 TGo4FitModel *TGo4Fitter::FindModel(const char *ModelName)
206 {
207  return !ModelName ? nullptr : dynamic_cast<TGo4FitModel *>(fxModels.FindObject(ModelName));
208 }
209 
211 {
212  fxModels.Add(model);
215  return model;
216 }
217 
218 TGo4FitModel *TGo4Fitter::AddModel(const char *DataName, TGo4FitModel *model)
219 {
220  model->AssignToData(DataName);
221  fxModels.Add(model);
224  return model;
225 }
226 
228 {
229  if (n < GetNumData())
230  return GetData(n);
231  else
232  return GetModel(n - GetNumData());
233 }
234 
235 void TGo4Fitter::AddPolynomX(const char *DataName, const char *NamePrefix, Int_t MaxOrder, Int_t GroupIndex,
236  Double_t lrange, Double_t rrange)
237 {
238  if (!DataName)
239  return;
240 
241  Bool_t flag = kFALSE;
242  Int_t NumTry = 0;
243  Bool_t createmodel = kFALSE;
244 
245  do {
246  TString Prefix(NamePrefix);
247  if (NumTry > 0)
248  Prefix += NumTry;
249  flag = kFALSE;
250  Bool_t findsame = kFALSE;
251 
252  for (Int_t Order = 0; Order <= MaxOrder; Order++) {
253  TString Name(Prefix);
254  Name += "_";
255  Name += Order;
256  findsame = FindModel(Name.Data());
257  if (findsame && !createmodel)
258  break;
259 
260  if (createmodel) {
261  TGo4FitModelPolynom *comp = new TGo4FitModelPolynom(Name, Order);
262  comp->SetGroupIndex(GroupIndex);
263  if ((lrange < rrange) || (rrange != 0.))
264  comp->SetRange(0, lrange, rrange);
265  AddModel(DataName, comp);
266  }
267  }
268 
269  flag = kTRUE;
270  if (findsame) {
271  NumTry++;
272  createmodel = kFALSE;
273  } else if (createmodel)
274  flag = kFALSE;
275  else
276  createmodel = kTRUE;
277 
278  } while (flag);
279 }
280 
281 void TGo4Fitter::AddPolynomX(const char *DataName, const char *NamePrefix, TArrayD &Coef, Int_t GroupIndex)
282 {
283  if (!DataName)
284  return;
285 
286  Bool_t flag = kFALSE;
287  Int_t NumTry = 0;
288  Bool_t createmodel = kFALSE;
289 
290  do {
291  Bool_t findsame = kFALSE;
292 
293  for (Int_t n = 0; n < Coef.GetSize(); n++) {
294  TString Name(NamePrefix);
295  if (NumTry > 0)
296  Name += NumTry;
297  Name += "_";
298  Name += n;
299  findsame = FindModel(Name.Data());
300  if (findsame && !createmodel)
301  break;
302 
303  if (createmodel) {
304  TGo4FitModelPolynom *model = new TGo4FitModelPolynom(Name, n);
305  model->SetAmplValue(Coef[n]);
306  model->SetGroupIndex(GroupIndex);
307  AddModel(DataName, model);
308  }
309  }
310 
311  flag = kTRUE;
312  if (findsame) {
313  NumTry++;
314  createmodel = kFALSE;
315  } else if (createmodel)
316  flag = kFALSE;
317  else
318  createmodel = kTRUE;
319 
320  } while (flag);
321 }
322 
323 void TGo4Fitter::AddPolynoms(const char *DataName, const char *NamePrefix, Int_t MaxOrder, Int_t NumAxis,
324  Int_t GroupIndex)
325 {
326  if (!DataName)
327  return;
328  TArrayD Orders(NumAxis);
329  Orders.Reset(0.);
330 
331  Bool_t flag = kFALSE;
332  Int_t NumTry = 0;
333  Bool_t createmodel = kFALSE;
334 
335  do {
336 
337  TString Prefix(NamePrefix);
338  if (NumTry > 0)
339  Prefix += NumTry;
340  flag = kFALSE;
341  Bool_t findsame = kFALSE;
342 
343  do {
344  TString Name(Prefix);
345  for (Int_t n = 0; n < NumAxis; n++) {
346  Name += "_";
347  Name += Int_t(Orders[n]);
348  }
349  findsame = FindModel(Name.Data());
350  if (findsame && !createmodel)
351  break;
352 
353  if (createmodel) {
354  TGo4FitModelPolynom *comp = new TGo4FitModelPolynom(Name, Orders);
355  comp->SetGroupIndex(GroupIndex);
356  AddModel(DataName, comp);
357  }
358 
359  Int_t nn = 0;
360  do {
361  Orders[nn] += 1.;
362  if (Orders[nn] <= MaxOrder)
363  break;
364  Orders[nn] = 0;
365  nn++;
366  } while (nn < NumAxis);
367  flag = (nn < NumAxis);
368  } while (flag);
369 
370  flag = kTRUE;
371  if (findsame) {
372  NumTry++;
373  createmodel = kFALSE;
374  } else if (createmodel)
375  flag = kFALSE;
376  else
377  createmodel = kTRUE;
378 
379  } while (flag);
380 }
381 
382 TGo4FitModelGauss1 *TGo4Fitter::AddGauss1(const char *DataName, const char *ModelName, Double_t iPosition,
383  Double_t iWidth, Double_t iAmpl, Int_t Axis)
384 {
385  TGo4FitModelGauss1 *gauss = new TGo4FitModelGauss1(ModelName, iPosition, iWidth, Axis);
386  gauss->SetAmplValue(iAmpl);
387  AddModel(DataName, gauss);
388  return gauss;
389 }
390 
391 TGo4FitModel *TGo4Fitter::CloneModel(const char *ModelName, const char *NewName)
392 {
393  TGo4FitModel *mod = FindModel(ModelName);
394  if (!mod)
395  return nullptr;
396 
397  TString newname;
398  int cnt = 0;
399 
400  do {
401  newname = NewName ? NewName : mod->GetName();
402  if (cnt > 0) {
403  newname += "_";
404  newname += cnt;
405  }
406  cnt++;
407  } while (FindModel(newname.Data()));
408 
409  TGo4FitModel *newmod = (TGo4FitModel *)mod->Clone(newname.Data());
410 
411  return AddModel(newmod);
412 }
413 
414 TGo4FitModel *TGo4Fitter::RemoveModel(const char *ModelName, Bool_t IsDel)
415 {
416  TGo4FitModel *mod = FindModel(ModelName);
417  if (mod) {
420  fxModels.Remove(mod);
421  if (IsDel) {
423  delete mod;
424  mod = nullptr;
425  }
426  fxModels.Compress();
427  }
428  return mod;
429 }
430 
431 Int_t TGo4Fitter::NumModelsAssosiatedTo(const char *DataName)
432 {
433  Int_t res = 0;
434  for (Int_t n = 0; n < GetNumModel(); n++)
435  if (GetModel(n)->IsAssignTo(DataName))
436  res++;
437  return res;
438 }
439 
440 void TGo4Fitter::DeleteModelsAssosiatedTo(const char *DataName)
441 {
442  Int_t n = 0;
443  while (n < GetNumModel()) {
444  TGo4FitModel *model = GetModel(n++);
445  if (model->IsAssignTo(DataName)) {
446  if (model->NumAssigments() == 1) {
449  fxModels.Remove(model);
450  delete model;
451  fxModels.Compress();
452  n--;
453  } else
454  model->ClearAssignmentTo(DataName);
455  }
456  }
457 }
458 
460 {
461  fxModels.Delete();
462  fxModels.Compress();
465 }
466 
467 void TGo4Fitter::AssignModelTo(const char *ModelName, const char *DataName, Double_t RatioValue, Bool_t FixRatio)
468 {
469  TGo4FitModel *model = FindModel(ModelName);
470  if (model) {
471  if (DataName)
472  model->AssignToData(DataName, RatioValue, FixRatio);
473  else {
474  model->ClearAssignments();
475  for (int n = 0; n < GetNumData(); n++)
476  model->AssignToData(GetData(n)->GetName(), RatioValue, FixRatio);
477  }
479  }
480 }
481 
482 void TGo4Fitter::ClearModelAssignmentTo(const char *ModelName, const char *DataName)
483 {
484  TGo4FitModel *model = FindModel(ModelName);
485  if (!model)
486  return;
487  if (DataName)
488  model->ClearAssignmentTo(DataName);
489  else
490  model->ClearAssignments();
492 }
493 
494 void TGo4Fitter::ChangeDataNameInAssignments(const char *oldname, const char *newname)
495 {
496  for (Int_t n = 0; n < GetNumModel(); n++)
497  GetModel(n)->ChangeDataNameInAssignments(oldname, newname);
498 }
499 
501 {
502  Int_t dbuf = -1, mbuf = -1;
503  switch (GetMemoryUsage()) {
504  case 0:
505  dbuf = 0;
506  mbuf = 0;
507  break;
508  case 1:
509  dbuf = 1;
510  mbuf = 0;
511  break;
512  case 2:
513  dbuf = 1;
514  mbuf = 1;
515  break;
516  default: dbuf = -1; mbuf = -1;
517  }
518 
519  for (Int_t i1 = 0; i1 < GetNumData(); i1++) {
520  TGo4FitData *data = GetData(i1);
521  if (!data->Initialize(dbuf))
522  return kFALSE;
523  for (Int_t i2 = 0; i2 < GetNumModel(); i2++)
525  }
526 
527  for (Int_t i2 = 0; i2 < GetNumModel(); i2++) {
528  TGo4FitModel *fModel = GetModel(i2);
529  if (!fModel->Initialize(mbuf))
530  return kFALSE;
531  }
532 
534  std::cout << " User fit function not set. Switch to least squares " << std::endl;
536  }
537 
538  return kTRUE;
539 }
540 
542 {
543  for (Int_t i = 0; i < GetNumData(); i++)
544  GetData(i)->Finalize();
545 
546  for (Int_t i = 0; i < GetNumModel(); i++)
547  GetModel(i)->Finalize();
548 }
549 
550 Double_t TGo4Fitter::PointFitFunction(Int_t FitFunctionType, Double_t value, Double_t modelvalue, Double_t standdev)
551 {
552  switch (FitFunctionType) {
553  case ff_least_squares: {
554  Double_t zn1 = (value - modelvalue);
555  return zn1 * zn1;
556  }
557  case ff_chi_square: {
558  if (standdev <= 0.)
559  return 0.;
560  Double_t zn2 = (value - modelvalue);
561  return zn2 * zn2 / standdev;
562  }
563  case ff_chi_Pearson: {
564  if (modelvalue <= 0.)
565  return 0.;
566  Double_t zn3 = (value - modelvalue);
567  return zn3 * zn3 / modelvalue;
568  }
569  case ff_chi_Neyman: {
570  Double_t zn4 = (value - modelvalue);
571  return zn4 * zn4 / ((value < 1.) ? 1. : value);
572  }
573  case ff_chi_gamma: {
574  if (value < 0.)
575  return 0.;
576  Double_t zn5 = (value + ((value < 1.) ? 0. : 1.) - modelvalue);
577  return zn5 * zn5 / (value + 1.);
578  }
579  case ff_ML_Poisson:
580  if (modelvalue <= 0.)
581  return 0.;
582  return modelvalue - value * TMath::Log(modelvalue);
583  case ff_user:
584  if (!fxUserFitFunction)
585  return 0;
586  return fxUserFitFunction(value, modelvalue, standdev);
587  default: return (value - modelvalue) * (value - modelvalue);
588  }
589 }
590 
591 Double_t TGo4Fitter::CalculateFCN(Int_t FitFunctionType, TGo4FitData *selectdata)
592 {
593  if (GetMemoryUsage() > 0)
594  RebuildAll();
595 
596  Double_t fSum = 0.;
597 
598  for (Int_t n = 0; n < GetNumData(); n++) {
599  TGo4FitData *dat = GetData(n);
600  if (selectdata && (dat != selectdata))
601  continue;
602  Double_t DataAmpl = dat->GetAmplValue();
603 
604  if (dat->BuffersAllocated()) {
605  Int_t size = dat->GetBinsSize();
606  Double_t *values = dat->GetBinsValues();
607  Double_t *res = dat->GetBinsResult();
608  Double_t *devs = dat->GetBinsDevs();
609 
610  for (Int_t nbin = 0; nbin < size; nbin++)
611  fSum +=
612  PointFitFunction(FitFunctionType, DataAmpl * values[nbin], res[nbin], DataAmpl * DataAmpl * devs[nbin]);
613  } else {
614 
615  auto iter = dat->MakeIter();
616  if (!iter->Reset())
617  continue;
618 
619  TObjArray Models;
620 
621  TArrayD Ampls(GetNumModel());
622  for (Int_t nm = 0; nm < GetNumModel(); nm++)
623  if (GetModel(nm)->IsAssignTo(dat->GetName())) {
624  TGo4FitModel *model = GetModel(nm);
625  Models.Add(model);
626  model->BeforeEval(iter->ScalesSize());
627  Ampls[Models.GetLast()] = model->GetAmplValue() * model->GetRatioValueFor(dat->GetName());
628  }
629 
630  do {
631  Double_t value = DataAmpl * iter->Value();
632  Double_t deviat = DataAmpl * DataAmpl * iter->StandardDeviation();
633 
634  Double_t modelvalue = 0.;
635  for (Int_t nm = 0; nm <= Models.GetLast(); nm++) {
636  TGo4FitModel *model = dynamic_cast<TGo4FitModel *>(Models.At(nm));
637  modelvalue += Ampls[nm] * model->EvaluateAtPoint(iter);
638  }
639 
640  fSum += PointFitFunction(FitFunctionType, value, modelvalue, deviat);
641 
642  } while (iter->Next());
643 
644  for (Int_t nm = 0; nm <= Models.GetLast(); nm++)
645  ((TGo4FitModel *)Models.At(nm))->AfterEval();
646  }
647  }
648 
649  return fSum;
650 }
651 
652 Double_t TGo4Fitter::CalculateFitFunction(Double_t *pars, Int_t FitFunctionType, const char *DataName)
653 {
654  if (FitFunctionType < 0)
655  FitFunctionType = GetFitFunctionType();
656 
657  if (pars) {
658  ExecuteDependencies(pars);
659  SetParsValues(pars);
660  }
661 
662  return CalculateFCN(FitFunctionType, FindData(DataName));
663 }
664 
665 Int_t TGo4Fitter::CalculateNDF(const char *DataName)
666 {
667  TGo4FitData *selectdata = FindData(DataName);
668 
669  Int_t NDF = 0;
670 
671  if (selectdata) {
672  NDF = selectdata->DefineBinsSize();
673  for (Int_t nm = 0; nm < GetNumModel(); nm++) {
674  TGo4FitModel *model = GetModel(nm);
675  if (model->IsAssignTo(selectdata->GetName()))
676  NDF -= model->NumFreePars();
677  }
678  } else {
679  for (Int_t n = 0; n < GetNumData(); n++)
680  NDF += GetData(n)->DefineBinsSize();
681  NDF -= NumFreePars();
682  }
683 
684  return NDF;
685 }
686 
688 {
690 }
691 
693 {
694  return CalculateNDF(nullptr);
695 }
696 
697 void TGo4Fitter::RebuildAll(Bool_t ForceBuild)
698 {
699  for (Int_t i2 = 0; i2 < GetNumModel(); i2++)
700  GetModel(i2)->RebuildShape(ForceBuild);
701 
702  for (Int_t i1 = 0; i1 < GetNumData(); i1++) {
703  TGo4FitData *data = GetData(i1);
704  if (!data->BuffersAllocated())
705  continue;
706 
707  Int_t size = data->GetBinsSize();
708  Double_t *result = data->GetBinsResult();
709  for (Int_t nbin = 0; nbin < size; nbin++)
710  result[nbin] = 0.;
711 
712  for (Int_t i2 = 0; i2 < GetNumModel(); i2++) {
713  TGo4FitModel *model = GetModel(i2);
714  if (model->IsAssignTo(data->GetName()))
715  model->AddModelToDataResult(data);
716  }
717  }
718 }
719 
720 void TGo4Fitter::EstimateAmplitudes(Int_t NumIters)
721 {
722  TGo4FitAmplEstimation abc("this", NumIters);
723  abc.DoAction(this);
724 }
725 
726 void TGo4Fitter::AddAmplEstimation(Int_t NumIters)
727 {
728  AddAction(new TGo4FitAmplEstimation("AmplEstim", NumIters));
729 }
730 
732 {
734  AddSimpleMinuit();
735 }
736 
737 void TGo4Fitter::FillSlotList(TSeqCollection *list)
738 {
740  for (Int_t i = 0; i < GetNumComp(); i++)
741  GetComp(i)->FillSlotList(list);
742 }
743 
745 {
746  return !model ? kFALSE : model->BuffersAllocated();
747 }
748 
750 {
751  return !data ? kFALSE : data->BuffersAllocated();
752 }
753 
755 {
756  return !data ? 0 : data->GetBinsSize();
757 }
758 
760 {
761  return !data ? nullptr : data->GetBinsValues();
762 }
763 
765 {
766  return !data ? nullptr : data->GetBinsDevs();
767 }
768 
770 {
771  return !data ? nullptr : data->GetBinsResult();
772 }
773 
774 Double_t *TGo4Fitter::GetModelBinsValues(TGo4FitModel *model, const char *DataName) const
775 {
776  return !model ? nullptr : model->GetModelBins(DataName);
777 }
778 
779 void TGo4Fitter::Print(Option_t *option) const
780 {
781  TString Opt(option);
782  if (Opt == "Ampls") {
783  PrintAmpls();
784  return;
785  } else if (Opt == "Pars") {
786  PrintPars();
787  return;
788  } else if (Opt == "Results") {
789  PrintResults();
790  return;
791  } else if (Opt == "Lines") {
792  PrintLines();
793  return;
794  }
795 
797  std::cout << "Fitiing function type: ";
798  switch (fiFitFunctionType) {
799  case ff_chi_square: std::cout << "ff_chi_square" << std::endl; break;
800  case ff_chi_Pearson: std::cout << "ff_chi_Pearson" << std::endl; break;
801  case ff_chi_Neyman: std::cout << "ff_chi_Neyman" << std::endl; break;
802  case ff_chi_gamma: std::cout << "ff_chi_gamma" << std::endl; break;
803  case ff_ML_Poisson: std::cout << "ff_ML_Poisson" << std::endl; break;
804  case ff_user: std::cout << "user defined" << std::endl; break;
805  default: std::cout << "ff_least_squares" << std::endl;
806  }
807  std::cout << std::endl << " LIST OF DATA OBJECTS" << std::endl;
808  fxDatas.Print(option);
809  std::cout << std::endl << " LIST OF MODEL OBJECTS" << std::endl;
810  fxModels.Print(option);
811 }
812 
813 Bool_t TGo4Fitter::CalculatesMomentums(const char *DataName, Bool_t UseRanges, Bool_t SubstractModels, Double_t &first,
814  Double_t &second)
815 {
816  TGo4FitData *data = FindData(DataName);
817  if (!data)
818  return kFALSE;
819 
820  auto iter = data->MakeIter();
821  if (!iter)
822  return kFALSE;
823 
824  Int_t size = iter->CountPoints(UseRanges);
825  if (size == 0)
826  return kFALSE;
827 
828  TObjArray Models;
829  if (SubstractModels)
830  for (Int_t nm = 0; nm < GetNumModel(); nm++)
831  if (GetModel(nm)->IsAssignTo(DataName))
832  Models.Add(GetModel(nm));
833 
834  TArrayD Ampls(Models.GetLast() + 1);
835  for (Int_t n = 0; n <= Models.GetLast(); n++) {
836  TGo4FitModel *model = (TGo4FitModel *)Models[n];
837  model->BeforeEval(iter->ScalesSize());
838  Ampls[n] = model->GetAmplValue() * model->GetRatioValueFor(DataName);
839  }
840 
841  TArrayD bins(size), scales(size);
842 
843  Int_t ppnt = 0;
844 
845  if (iter->Reset(UseRanges))
846  do {
847  Double_t value = iter->Value();
848  for (Int_t n = 0; n <= Models.GetLast(); n++) {
849  TGo4FitModel *model = (TGo4FitModel *)Models[n];
850  value -= Ampls[n] * model->EvaluateAtPoint(iter);
851  }
852  value = TMath::Abs(value);
853  bins[ppnt] = value;
854  scales[ppnt] = iter->x();
855  ppnt++;
856  } while (iter->Next(UseRanges));
857 
858  for (Int_t n = 0; n <= Models.GetLast(); n++) {
859  TGo4FitModel *model = (TGo4FitModel *)Models[n];
860  model->AfterEval();
861  }
862 
863  Int_t niter = 0;
864 
865  do {
866  Double_t sum00 = 0., sum11 = 0., sum22 = 0.;
867  for (Int_t pnt = 0; pnt < size; pnt++)
868  if ((bins[pnt] > 0.) && ((niter == 0) || (TMath::Abs(scales[pnt] - first) < second * 2.))) {
869  sum00 += bins[pnt];
870  sum11 += bins[pnt] * scales[pnt];
871  sum22 += bins[pnt] * scales[pnt] * scales[pnt];
872  }
873 
874  if (sum00 > 0.) {
875  Double_t mid = sum11 / sum00;
876  Double_t dev = TMath::Sqrt(sum22 / sum00 - mid * mid);
877 
878  if (niter > 0)
879  if ((dev / second > 0.8) && (dev / second < 1.2))
880  niter = 10;
881 
882  first = mid;
883  second = dev;
884 
885  } else {
886  niter = 10;
887  }
888 
889  } while (niter++ < 8);
890 
891  return kTRUE;
892 }
893 
894 Double_t TGo4Fitter::CalculatesIntegral(const char *DataName, const char *ModelName, Bool_t onlycounts)
895 {
896  TGo4FitData *data = FindData(DataName);
897  if (!data)
898  return 0.;
899 
900  TGo4FitModel *model = ModelName ? FindModel(ModelName) : nullptr;
901  if (ModelName && !model)
902  return 0.;
903 
904  auto iter = data->MakeIter();
905  if (!iter || !iter->Reset(kTRUE))
906  return 0.;
907 
908  double sum = 0., ampl = 1.;
909 
910  if (model) {
911  model->BeforeEval(iter->ScalesSize());
912  ampl = model->GetAmplValue() * model->GetRatioValueFor(DataName);
913  }
914 
915  do {
916  double dx = onlycounts ? 1. : iter->xWidths();
917  double value = !model ? iter->Value() : model->EvaluateAtPoint(iter);
918  sum += ampl * value * dx;
919  } while (iter->Next(kTRUE));
920 
921  if (model)
922  model->AfterEval();
923 
924  return sum;
925 }
926 
927 Double_t TGo4Fitter::CalculatesModelIntegral(const char *ModelName, Bool_t OnlyCounts)
928 {
929  TGo4FitModel *model = FindModel(ModelName);
930  if (!model)
931  return 0.;
932  return CalculatesIntegral(model->AssignmentName(0), ModelName, OnlyCounts);
933 }
934 
935 TObject *TGo4Fitter::CreateDrawObject(const char *ResName, const char *DataName, Bool_t IsModel, const char *ModelName)
936 {
937  TGo4FitData *data = FindData(DataName);
938  if (!data)
939  return nullptr;
940 
941  TObjArray Models;
942  if (IsModel) {
943  TGo4FitModel *model = FindModel(ModelName);
944  if (model) {
945  Models.Add(model);
946  } else {
947  Int_t groupindex = -1;
948 
949  if (ModelName) {
950  TString modelname(ModelName);
951  if (modelname == "Background")
952  groupindex = 0;
953  else if (modelname.Index("Group", 0, TString::kExact) == 0) {
954  modelname.Remove(0, 5);
955  char *err = nullptr;
956  groupindex = strtol(modelname.Data(), &err, 10);
957  if (err && (*err != 0))
958  groupindex = -1;
959  }
960  }
961 
962  for (Int_t nm = 0; nm < GetNumModel(); nm++) {
963  model = GetModel(nm);
964  if (model->IsAssignTo(DataName))
965  if ((groupindex < 0) || (model->GetGroupIndex() == groupindex))
966  Models.Add(model);
967  }
968  }
969  if (Models.GetLast() < 0)
970  return nullptr;
971  }
972 
973  auto iter = data->MakeIter();
974  if (!iter || !iter->Reset(kFALSE))
975  return nullptr;
976 
977  TH1 *histo = nullptr;
978  Int_t ndim = 0;
979  TGraph *gr = nullptr;
980  Bool_t UseRanges = kTRUE;
981 
982  if (iter->HasIndexes() && (iter->IndexesSize() == iter->ScalesSize()) && iter->HasWidths()) {
983  histo = iter->CreateHistogram(ResName, kFALSE, !IsModel);
984 
985  if (histo) {
986  ndim = histo->GetDimension();
987  UseRanges = Models.GetLast() >= 0;
988  }
989  } else {
990  gr = iter->CreateGraph(ResName, kTRUE, !IsModel);
991  UseRanges = kTRUE;
992  }
993 
994  if ((histo || gr) && IsModel) {
995  TArrayD Ampls(Models.GetLast() + 1);
996 
997  for (Int_t n = 0; n <= Models.GetLast(); n++) {
998  TGo4FitModel *model = (TGo4FitModel *)Models[n];
999  model->BeforeEval(iter->ScalesSize());
1000  Ampls[n] = model->GetAmplValue() * model->GetRatioValueFor(DataName);
1001  }
1002 
1003  if (iter->Reset(UseRanges))
1004  do {
1005  Double_t zn = 0;
1006  for (Int_t n = 0; n <= Models.GetLast(); n++) {
1007  TGo4FitModel *model = (TGo4FitModel *)Models[n];
1008  zn += Ampls[n] * model->EvaluateAtPoint(iter);
1009  }
1010  if (histo)
1011  switch (ndim) {
1012  case 1: histo->SetBinContent(iter->Indexes()[0] + 1, zn); break;
1013  case 2: histo->SetBinContent(iter->Indexes()[0] + 1, iter->Indexes()[1] + 1, zn); break;
1014  case 3:
1015  histo->SetBinContent(iter->Indexes()[0] + 1, iter->Indexes()[1] + 1, iter->Indexes()[2] + 1, zn);
1016  break;
1017  }
1018  if (gr) {
1019  (gr->GetX())[iter->Point()] = iter->x();
1020  (gr->GetY())[iter->Point()] = zn;
1021  }
1022 
1023  } while (iter->Next(UseRanges));
1024 
1025  for (Int_t n = 0; n <= Models.GetLast(); n++) {
1026  TGo4FitModel *model = (TGo4FitModel *)Models[n];
1027  model->AfterEval();
1028  }
1029  }
1030 
1031  if (gr && IsModel)
1032  for (Int_t n1 = 0; n1 < gr->GetN() - 1; n1++)
1033  for (Int_t n2 = n1 + 1; n2 < gr->GetN(); n2++)
1034  if ((gr->GetX())[n1] > (gr->GetX())[n2]) {
1035  Double_t xx = (gr->GetX())[n1];
1036  (gr->GetX())[n1] = (gr->GetX())[n2];
1037  (gr->GetX())[n2] = xx;
1038  Double_t yy = (gr->GetY())[n1];
1039  (gr->GetY())[n1] = (gr->GetY())[n2];
1040  (gr->GetY())[n2] = yy;
1041  }
1042 
1043  TNamed *res = nullptr;
1044  if (histo)
1045  res = histo;
1046  else
1047  res = gr;
1048  if (res) {
1049  TString title;
1050  if (IsModel)
1051  if (ModelName)
1052  title.Form("Draw of model %s", ModelName);
1053  else
1054  title.Form("Draw of full model of %s", DataName);
1055  else
1056  title.Form("Draw of data %s", DataName);
1057  res->SetTitle(title.Data());
1058  }
1059 
1060  return res;
1061 }
1062 
1063 /* valid format for draw options
1064 
1065  "" - draw all data
1066  "*" - draw all data with models
1067  "**" - draw all data with models and all components
1068  "DataName" - draw select data with model
1069  "DataName-" - draw select data without model
1070  "DataName*" - draw select data with model and components
1071  "DataName-, Background" - draw data and sum of components, belongs to background group (index = 0)
1072  "DataName-, Group1" - draw data and sum of components, belongs to group 1
1073  "ModelName" - draw model component and correspondent data
1074 */
1075 
1076 void TGo4Fitter::Draw(Option_t *option)
1077 {
1078  TString opt(option);
1079 
1080  TCanvas *fCanvas = nullptr;
1081 
1082  if ((opt.Length() > 0) && (opt[0] == '#')) {
1083  opt.Remove(0, 1);
1084  TString CanvasName;
1085  Int_t n = 0;
1086  do {
1087  CanvasName = "Canvas";
1088  CanvasName += n++;
1089  } while (gROOT->FindObject(CanvasName.Data()));
1090  fCanvas = new TCanvas(CanvasName, TString("Draw of fitter ") + GetName() + " " + opt, 3);
1091  fCanvas->cd();
1092  }
1093 
1094  Bool_t drawdata = kFALSE, drawcomp = kFALSE, drawmodel = kFALSE;
1095  TGo4FitData *selectdata = nullptr;
1096  TObjArray selectmodels;
1097  selectmodels.SetOwner(kTRUE);
1098 
1099  if (opt == "*") {
1100  opt = "";
1101  drawdata = kTRUE;
1102  }
1103 
1104  while (opt.Length() > 0) {
1105  Int_t len = opt.Index(",", 0, TString::kExact);
1106 
1107  if (len < 0)
1108  len = opt.Length();
1109  if (len == 0)
1110  break;
1111 
1112  TString optpart(opt.Data(), len);
1113  while ((optpart.Length() > 0) && (optpart[0] == ' '))
1114  optpart.Remove(0, 1);
1115 
1116  Bool_t find = kFALSE;
1117 
1118  for (Int_t n = 0; n < GetNumData(); n++)
1119  if (optpart.Index(GetDataName(n), 0, TString::kExact) == 0) {
1120  selectdata = GetData(n);
1121  drawdata = kTRUE;
1122  drawmodel = kTRUE;
1123  find = kTRUE;
1124  optpart.Remove(0, strlen(GetDataName(n)));
1125  if (optpart == "*")
1126  drawcomp = kTRUE;
1127  else if (optpart == "-")
1128  drawmodel = kFALSE;
1129  break;
1130  }
1131 
1132  if (!find) {
1133  selectmodels.Add(new TObjString(optpart));
1134  TGo4FitModel *model = FindModel(optpart.Data());
1135  if (model && !selectdata)
1136  for (Int_t n = 0; n < GetNumData(); n++)
1137  if (model->IsAssignTo(GetDataName(n)))
1138  selectdata = GetData(n);
1139  }
1140 
1141  opt.Remove(0, len);
1142  if (opt.Length() > 0)
1143  opt.Remove(0, 1);
1144  }
1145 
1146  if (!selectdata && !drawdata)
1147  if (GetNumData() > 0)
1148  selectdata = GetData(0);
1149 
1151  fxDrawObjs = new TObjArray();
1152 
1153  for (Int_t nn = 0; nn < GetNumData(); nn++) {
1154  TGo4FitData *data = GetData(nn);
1155  if (selectdata && (data != selectdata))
1156  continue;
1157 
1158  if (drawdata) {
1159  TObject *obj = CreateDrawObject(TString(data->GetName()) + "_bins", data->GetName(), kFALSE);
1160 
1161  TAttLine *line = dynamic_cast<TAttLine *>(obj);
1162  if (line) {
1163  line->SetLineColor(1);
1164  line->SetLineWidth(1);
1165  }
1166  if (obj)
1167  fxDrawObjs->Add(obj);
1168  }
1169 
1170  if (drawmodel) {
1171  TObject *mobj = CreateDrawObject(TString(data->GetName()) + "_fullmodel", data->GetName(), kTRUE);
1172 
1173  TAttLine *line = dynamic_cast<TAttLine *>(mobj);
1174  if (line) {
1175  line->SetLineColor(4);
1176  line->SetLineWidth(2);
1177  }
1178  if (mobj)
1179  fxDrawObjs->Add(mobj);
1180  }
1181 
1182  if (drawcomp)
1183  for (Int_t nmodel = 0; nmodel < GetNumModel(); nmodel++) {
1184  TGo4FitModel *model = GetModel(nmodel);
1185  if (!model->IsAssignTo(data->GetName()))
1186  continue;
1187 
1188  TObject *cobj = CreateDrawObject(TString(data->GetName()) + "_" + model->GetName(), data->GetName(), kTRUE,
1189  model->GetName());
1190 
1191  TAttLine *line = dynamic_cast<TAttLine *>(cobj);
1192  if (line) {
1193  line->SetLineColor(6);
1194  line->SetLineWidth(1);
1195  }
1196  if (cobj)
1197  fxDrawObjs->Add(cobj);
1198  }
1199 
1200  for (Int_t n = 0; n <= selectmodels.GetLast(); n++) {
1201  TString name = ((TObjString *)(selectmodels[n]))->String();
1202 
1203  TObject *cobj = CreateDrawObject(TString(data->GetName()) + "_" + name, data->GetName(), kTRUE, name.Data());
1204 
1205  TAttLine *line = dynamic_cast<TAttLine *>(cobj);
1206  if (line) {
1207  line->SetLineColor(6);
1208  line->SetLineWidth(1);
1209  }
1210  if (cobj)
1211  fxDrawObjs->Add(cobj);
1212  }
1213  }
1214 
1215  if (fxDrawObjs->GetLast() == 0) {
1216  fxDrawObjs->At(0)->Draw();
1217  } else if (fxDrawObjs->GetLast() > 0) {
1218  Bool_t allhisto = kTRUE;
1219  for (Int_t n = 0; n <= fxDrawObjs->GetLast(); n++)
1220  if (!(fxDrawObjs->At(n)->InheritsFrom(TH1::Class())))
1221  allhisto = kFALSE;
1222  if (allhisto) {
1223  THStack *stack =
1224  new THStack(TString("Stack") + "_" + fxDrawObjs->At(0)->GetName(), fxDrawObjs->At(0)->GetName());
1225  for (Int_t n = 0; n <= fxDrawObjs->GetLast(); n++)
1226  stack->Add((TH1 *)fxDrawObjs->At(n));
1227  fxDrawObjs->Clear();
1228  fxDrawObjs->Add(stack);
1229  stack->Draw("nostack");
1230  } else
1231  for (Int_t n = 0; n <= fxDrawObjs->GetLast(); n++)
1232  if (n == 0)
1233  fxDrawObjs->At(n)->Draw("A*");
1234  else
1235  fxDrawObjs->At(n)->Draw("L");
1236  }
1237 
1238  if (fCanvas)
1239  fCanvas->Update();
1240 }
1241 
1243 {
1244  std::cout << std::endl << "*** LIST OF AMPLITUDES VALUE ***" << std::endl;
1245  for (Int_t n = 0; n < GetNumComp(); n++) {
1246  TGo4FitComponent *comp = ((TGo4Fitter *)this)->GetComp(n);
1247  if (comp->GetAmplPar())
1248  std::cout << " " << comp->GetAmplFullName() << " " << comp->GetAmplValue() << " "
1249  << comp->GetAmplError() << std::endl;
1250  }
1251 }
1252 
1254 {
1255  std::cout << std::endl << " *** LIST OF LINES PARAMETERS ***" << std::endl;
1256 
1257  int MaxAxis = 0;
1258  for (Int_t n = 0; n < GetNumModel(); n++) {
1259  TGo4FitModel *m = ((TGo4Fitter *)this)->GetModel(n);
1260  if (!m)
1261  continue;
1262  Double_t zn;
1263  for (int naxis = 0; naxis < 3; naxis++)
1264  if (m->GetPosition(naxis, zn) || m->GetWidth(naxis, zn))
1265  MaxAxis = naxis;
1266  }
1267 
1268  std::cout << std::setw(10) << "Name" << std::setw(12) << "Ampl";
1269  for (Int_t naxis = 0; naxis <= MaxAxis; naxis++)
1270  std::cout << std::setw(11) << "Pos" << naxis << std::setw(11) << "Width" << naxis;
1271  std::cout << std::endl;
1272 
1273  for (Int_t n = 0; n < GetNumModel(); n++) {
1274  TGo4FitModel *m = ((TGo4Fitter *)this)->GetModel(n);
1275  if (!m)
1276  continue;
1277  std::cout << std::setw(10) << m->GetName() << std::setw(12) << m->GetAmplValue();
1278 
1279  for (int naxis = 0; naxis <= MaxAxis; naxis++) {
1280  Double_t pos, width;
1281  std::cout << std::setw(12);
1282  if (m->GetPosition(naxis, pos))
1283  std::cout << pos;
1284  else
1285  std::cout << "---";
1286  std::cout << std::setw(12);
1287  if (m->GetWidth(naxis, width))
1288  std::cout << width;
1289  else
1290  std::cout << "---";
1291  }
1292  std::cout << std::endl;
1293  }
1294 }
1295 
1297 {
1298  if (fxDrawObjs) {
1299  lst.AddAll(fxDrawObjs);
1300  delete fxDrawObjs;
1301  fxDrawObjs = nullptr;
1302  }
1303 }
1304 
1306 {
1307  if (fxDrawObjs) {
1308  for (Int_t n = 0; n <= fxDrawObjs->GetLast(); n++)
1309  gROOT->Add(fxDrawObjs->At(n));
1310  delete fxDrawObjs;
1311  fxDrawObjs = nullptr;
1312  }
1313 }
1314 
1315 TString TGo4Fitter::FindNextName(const char *Head, Int_t start, Bool_t isModel)
1316 {
1317  TString name;
1318  Int_t n = start;
1319  do {
1320  name.Form("%s%d", Head, n++);
1321  } while (isModel ? FindModel(name.Data()) != nullptr : FindData(name.Data()) != nullptr);
1322  return name;
1323 }
1324 
1325 void TGo4Fitter::Streamer(TBuffer &b)
1326 {
1327  if (b.IsReading()) {
1328  TGo4Fitter::Class()->ReadBuffer(b, this);
1332  } else {
1334  TGo4Fitter::Class()->WriteBuffer(b, this);
1335  }
1336 }
Int_t GetNumData() const
Definition: TGo4Fitter.h:123
void SetUpdateSlotList()
void CollectParsTo(TGo4FitParsList &list) override
friend class TGo4FitAmplEstimation
Definition: TGo4Fitter.h:385
Double_t CalculatesIntegral(const char *DataName, const char *ModelName=nullptr, Bool_t OnlyCounts=kFALSE)
Definition: TGo4Fitter.cxx:894
void PrintPars() const
TGo4FitDataGraph * AddGraph(const char *DataName, TGraph *gr, Bool_t Owned=kFALSE, Double_t lrange=0., Double_t rrange=0.)
Definition: TGo4Fitter.cxx:144
void RebuildAll(Bool_t ForceBuild=kFALSE)
Definition: TGo4Fitter.cxx:697
void SetMemoryUsage(Int_t iMemoryUsage)
Definition: TGo4Fitter.cxx:70
TObjArray * fxDrawObjs
Definition: TGo4Fitter.h:445
void Print(Option_t *option="") const override
void SetHistogram(TH1 *iHistogram, Bool_t iHistogramOwned=kFALSE)
virtual void Finalize()
TGo4FitSlot * GetSlot(Int_t nslot)
void ClearAssignments()
TGo4FitComponent * GetComp(Int_t n)
Definition: TGo4Fitter.cxx:227
void MoveDrawObjectsToROOT()
void SetFitFunctionType(Int_t iFitFunctionType)
Definition: TGo4Fitter.h:80
void AddAmplEstimation(Int_t NumIters=1)
Definition: TGo4Fitter.cxx:726
const char * GetAmplFullName()
virtual void Finalize()
Double_t * GetBinsDevs() const
Definition: TGo4FitData.h:271
void EstimateAmplitudes(Int_t NumIters=1)
Definition: TGo4Fitter.cxx:720
TGo4FitModel * FindModel(const char *ModelName)
Definition: TGo4Fitter.cxx:205
void ClearModelAssignmentTo(const char *ModelName, const char *DataName=nullptr)
Definition: TGo4Fitter.cxx:482
Double_t CalculateFitFunction(Double_t *pars=nullptr, Int_t FitFunctionType=-1, const char *DataName=nullptr)
Definition: TGo4Fitter.cxx:652
void AddPolynoms(const char *DataName, const char *NamePrefix, Int_t MaxOrder=1, Int_t NumAxis=1, Int_t GroupIndex=0)
Definition: TGo4Fitter.cxx:323
void ClearConnectionToSlot()
Definition: TGo4FitSlot.cxx:83
void AssignToData(const char *DataName, Double_t RatioValue=1., Bool_t FixRatio=kFALSE)
Bool_t InitFitterData() override
Definition: TGo4Fitter.cxx:500
Int_t NumModelsAssosiatedTo(const char *DataName)
Definition: TGo4Fitter.cxx:431
void AddAction(TGo4FitterAction *Action)
Double_t * GetDataBinsDevs(TGo4FitData *data) const
Definition: TGo4Fitter.cxx:764
Double_t * GetBinsResult() const
Definition: TGo4FitData.h:277
void CheckSlotsBeforeDelete(TGo4FitComponent *comp)
Definition: TGo4Fitter.cxx:161
Double_t DoCalculation() override
Definition: TGo4Fitter.cxx:687
void AddPolynomX(const char *DataName, const char *NamePrefix, Int_t MaxOrder=1, Int_t GroupIndex=0, Double_t lrange=0., Double_t rrange=0.)
Definition: TGo4Fitter.cxx:235
void DoAction(TGo4FitterAbstract *Fitter) override
virtual Bool_t Initialize(Int_t UseBuffers=-1)
void ClearAssignmentTo(const char *DataName)
Double_t PointFitFunction(Int_t FitFunctionType, Double_t value, Double_t modelvalue, Double_t standdev)
Definition: TGo4Fitter.cxx:550
void AssignModelTo(const char *ModelName, const char *DataName, Double_t RatioValue=1., Bool_t FixRatio=kFALSE)
Definition: TGo4Fitter.cxx:467
void ClearSlot(TGo4FitSlot *slot, Bool_t NonOwned)
void SetAmplValue(Double_t iAmpl)
const char * AssignmentName(Int_t n)
Definition: TGo4FitModel.h:129
void SetOwner(TNamed *iOwner)
Definition: TGo4FitNamed.h:58
TObjArray fxModels
Definition: TGo4Fitter.h:423
TGo4FitData * RemoveData(const char *DataName, Bool_t IsDel=kFALSE)
Definition: TGo4Fitter.cxx:175
Double_t * GetDataBinsResult(TGo4FitData *data) const
Definition: TGo4Fitter.cxx:769
Int_t GetDataBinsSize(TGo4FitData *data) const
Definition: TGo4Fitter.cxx:754
void DeleteAllData()
Definition: TGo4Fitter.cxx:192
TGo4FitModel * RemoveModel(const char *ModelName, Bool_t IsDel=kFALSE)
Definition: TGo4Fitter.cxx:414
Int_t fiMemoryUsage
Definition: TGo4Fitter.h:433
TUserFitFunction fxUserFitFunction
Definition: TGo4Fitter.h:443
void CollectAllPars() override
Definition: TGo4Fitter.cxx:80
Bool_t IsAssignTo(const char *DataName) const
Definition: TGo4FitModel.h:141
Int_t DefineBinsSize()
virtual Double_t EvaluateAtPoint(TGo4FitData *data, Int_t nbin, Bool_t UseRanges=kTRUE)
Int_t fiFitFunctionType
Definition: TGo4Fitter.h:428
void DeleteModelsAssosiatedTo(const char *DataName)
Definition: TGo4Fitter.cxx:440
void PrintLines() const
void PrepareSlotsForWriting()
TGo4FitData * AddData(TGo4FitData *d)
Definition: TGo4Fitter.cxx:119
virtual Bool_t GetPosition(Int_t naxis, Double_t &pos)
void ChangeDataNameInAssignments(const char *oldname, const char *newname)
void SetGraph(TGraph *iGraph, Bool_t iGraphOwned=kFALSE)
Double_t * GetDataBinsValues(TGo4FitData *data) const
Definition: TGo4Fitter.cxx:759
void FillSlotList(TSeqCollection *) override
Definition: TGo4Fitter.cxx:737
TGo4FitSlot * GetConnectedSlot() const
Definition: TGo4FitSlot.h:123
void CheckDuplicatesOnSlot()
TGo4FitModelGauss1 * AddGauss1(const char *DataName, const char *ModelName, Double_t iPosition, Double_t iWidth, Double_t iAmpl=1., Int_t Axis=0)
Definition: TGo4Fitter.cxx:382
void ChangeDataNameInAssignments(const char *oldname, const char *newname)
Definition: TGo4Fitter.cxx:494
void SetRange(Int_t naxis, Double_t min, Double_t max)
virtual void CollectAllPars()
Double_t CalculateFCN(Int_t FitFunctionType, TGo4FitData *selectdata=nullptr)
Definition: TGo4Fitter.cxx:591
TGo4FitDataGraph * SetGraph(const char *DataName, TGraph *gr, Bool_t Owned=kFALSE)
Definition: TGo4Fitter.cxx:153
void FinalizeFitterData() override
Definition: TGo4Fitter.cxx:541
TGo4FitModel * GetModel(Int_t n)
Definition: TGo4Fitter.cxx:200
Double_t CalculatesModelIntegral(const char *ModelName, Bool_t OnlyCounts=kFALSE)
Definition: TGo4Fitter.cxx:927
TObjArray fxDatas
Definition: TGo4Fitter.h:418
Int_t GetGroupIndex() const
Definition: TGo4FitModel.h:119
Int_t DoNDFCalculation() override
Definition: TGo4Fitter.cxx:692
virtual Bool_t Initialize(Int_t UseBuffers=-1)
Int_t GetNumComp() const
Definition: TGo4Fitter.h:281
virtual Bool_t BeforeEval(Int_t ndim)
void ConnectToDataIfAssigned(TGo4FitData *data)
TString FindNextName(const char *Head, Int_t start, Bool_t isModel=kTRUE)
Double_t * GetBinsValues() const
Definition: TGo4FitData.h:265
virtual void FillSlotList(TSeqCollection *lst)
virtual void AfterEval()
Definition: TGo4FitModel.h:253
Int_t GetNumModel() const
Definition: TGo4Fitter.h:178
TGo4FitData * GetData(Int_t n)
Definition: TGo4Fitter.cxx:104
void SetGroupIndex(Int_t index=-1)
Definition: TGo4FitModel.h:109
virtual Bool_t GetWidth(Int_t naxis, Double_t &width)
virtual ~TGo4Fitter()
Definition: TGo4Fitter.cxx:64
TObject * CreateDrawObject(const char *ObjName, const char *DataName, Bool_t IsModel=kFALSE, const char *ModelName=nullptr)
Definition: TGo4Fitter.cxx:935
virtual std::unique_ptr< TGo4FitDataIter > MakeIter()
Definition: TGo4FitData.h:175
Int_t NumAssigments() const
Definition: TGo4FitModel.h:124
Int_t GetFitFunctionType() const
Definition: TGo4Fitter.h:85
Bool_t DataBuffersAllocated(TGo4FitData *data)
Definition: TGo4Fitter.cxx:749
Bool_t BuffersAllocated() const
Definition: TGo4FitData.h:248
const char * GetDataName(Int_t n)
Definition: TGo4Fitter.cxx:109
Bool_t BuffersAllocated() const
TGo4FitModel * AddModel(TGo4FitModel *m)
Definition: TGo4Fitter.cxx:210
TGo4FitModel * CloneModel(const char *ModelName, const char *NewName=nullptr)
Definition: TGo4Fitter.cxx:391
Double_t * GetModelBinsValues(TGo4FitModel *model, const char *DataName) const
Definition: TGo4Fitter.cxx:774
void RebuildShape(Bool_t ForceBuild=kFALSE)
void Clear(Option_t *option="") override
TGo4FitDataHistogram * SetH1(const char *DataName, TH1 *histo, Bool_t Owned=kFALSE)
Definition: TGo4Fitter.cxx:136
Bool_t CalculatesMomentums(const char *DataName, Bool_t UseRanges, Bool_t SubstractModels, Double_t &first, Double_t &second)
Definition: TGo4Fitter.cxx:813
Double_t * GetModelBins(const char *DataName) const
TGo4FitData * FindData(const char *DataName)
Definition: TGo4Fitter.cxx:114
Bool_t AddModelToDataResult(TGo4FitData *data)
Int_t CalculateNDF(const char *DataName=nullptr)
Definition: TGo4Fitter.cxx:665
Double_t GetRatioValueFor(const char *DataName)
void ProvideLastDrawObjects(TObjArray &lst)
void Draw(Option_t *option) override
void Clear(Option_t *option="") override
Definition: TGo4Fitter.cxx:97
Bool_t ModelBuffersAllocated(TGo4FitModel *model)
Definition: TGo4Fitter.cxx:744
Int_t GetMemoryUsage() const
Definition: TGo4Fitter.h:113
void PrintAmpls() const
Int_t GetBinsSize() const
Definition: TGo4FitData.h:253
void AddStandardActions()
Definition: TGo4Fitter.cxx:731
TGo4FitParameter * GetAmplPar()
void ExecuteDependencies(Double_t *pars)
void DeleteAllModels()
Definition: TGo4Fitter.cxx:459
void SetParsValues(Double_t *pars)
TGo4FitDataHistogram * AddH1(const char *DataName, TH1 *histo, Bool_t Owned=kFALSE, Double_t lrange=0., Double_t rrange=0.)
Definition: TGo4Fitter.cxx:127
void Print(Option_t *option="") const override
Definition: TGo4Fitter.cxx:779