GSI Object Oriented Online Offline (Go4) GO4-6.4.0
Loading...
Searching...
No Matches
TGo4Fitter.h
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#ifndef TGO4FITTER_H
15#define TGO4FITTER_H
16
17#include "TGo4FitterAbstract.h"
18
19class TH1;
20class TGraph;
21
22#include "TArrayD.h"
23#include "TObjArray.h"
24
26class TGo4FitData;
29class TGo4FitModel;
31
32typedef Double_t (*TUserFitFunction)(Double_t, Double_t, Double_t);
33
39 public:
40 enum { ff_least_squares = 0,
46 ff_user = 6 };
47
51 TGo4Fitter();
52
56 TGo4Fitter(const char *iName, const char *iTitle);
57
62 TGo4Fitter(const char *iName, Int_t iFitFunctionType, Bool_t IsAddStandardActions);
63
67 virtual ~TGo4Fitter();
68
81 void SetFitFunctionType(Int_t iFitFunctionType) { fiFitFunctionType = iFitFunctionType; }
82
86 Int_t GetFitFunctionType() const { return fiFitFunctionType; }
87
99
108 void SetMemoryUsage(Int_t iMemoryUsage);
109
114 Int_t GetMemoryUsage() const { return fiMemoryUsage; }
115
119 void Clear(Option_t *option = "") override;
120
124 Int_t GetNumData() const { return fxDatas.GetLast()+1; }
125
129 TGo4FitData *GetData(Int_t n);
130
134 const char *GetDataName(Int_t n);
135
139 TGo4FitData *FindData(const char *DataName);
140
145
149 TGo4FitDataHistogram* AddH1(const char *DataName, TH1 *histo, Bool_t Owned = kFALSE, Double_t lrange = 0., Double_t rrange = 0.);
150
154 TGo4FitDataHistogram* SetH1(const char *DataName, TH1 *histo, Bool_t Owned = kFALSE);
155
159 TGo4FitDataGraph* AddGraph(const char *DataName, TGraph *gr, Bool_t Owned = kFALSE, Double_t lrange = 0., Double_t rrange = 0.);
160
164 TGo4FitDataGraph* SetGraph(const char *DataName, TGraph *gr, Bool_t Owned = kFALSE);
165
169 TGo4FitData *RemoveData(const char *DataName, Bool_t IsDel = kFALSE);
170
174 void DeleteAllData();
175
179 Int_t GetNumModel() const { return fxModels.GetLast()+1; }
180
184 TGo4FitModel *GetModel(Int_t n);
185
189 TGo4FitModel *FindModel(const char *ModelName);
190
195
199 TGo4FitModel *AddModel(const char *DataName, TGo4FitModel *m);
200
209 void AddPolynomX(const char *DataName, const char *NamePrefix, Int_t MaxOrder = 1, Int_t GroupIndex = 0, Double_t lrange = 0., Double_t rrange = 0.);
210
216 void AddPolynomX(const char *DataName, const char *NamePrefix, TArrayD& Coef, Int_t GroupIndex = 0);
217
225 void AddPolynoms(const char *DataName, const char *NamePrefix, Int_t MaxOrder = 1, Int_t NumAxis = 1, Int_t GroupIndex = 0);
226
230 TGo4FitModelGauss1 *AddGauss1(const char *DataName, const char *ModelName, Double_t iPosition, Double_t iWidth, Double_t iAmpl = 1., Int_t Axis = 0);
231
235 Int_t NumModelsAssosiatedTo(const char *DataName);
236
240 TGo4FitModel *CloneModel(const char *ModelName, const char *NewName = nullptr);
241
247 TGo4FitModel *RemoveModel(const char *ModelName, Bool_t IsDel = kFALSE);
248
252 void DeleteModelsAssosiatedTo(const char *DataName);
253
257 void DeleteAllModels();
258
264 void AssignModelTo(const char *ModelName, const char *DataName, Double_t RatioValue = 1., Bool_t FixRatio = kFALSE);
265
270 void ClearModelAssignmentTo(const char *ModelName, const char *DataName = nullptr);
271
276 void ChangeDataNameInAssignments(const char *oldname, const char *newname);
277
282 Int_t GetNumComp() const { return GetNumData() + GetNumModel(); }
283
289 TGo4FitComponent *GetComp(Int_t n);
290
296 void EstimateAmplitudes(Int_t NumIters = 1);
297
301 void AddAmplEstimation(Int_t NumIters = 1);
302
307 void AddStandardActions();
308
312 Double_t CalculateFitFunction(Double_t *pars = nullptr, Int_t FitFunctionType = -1, const char *DataName = nullptr);
313
319 Int_t CalculateNDF(const char *DataName = nullptr);
320
324 void FillSlotList(TSeqCollection *) override;
325
330 void Print(Option_t *option = "") const override;
331
336 Bool_t CalculatesMomentums(const char *DataName, Bool_t UseRanges, Bool_t SubstractModels, Double_t& first, Double_t& second);
337
342 Double_t CalculatesIntegral(const char *DataName, const char *ModelName = nullptr, Bool_t OnlyCounts = kFALSE);
343
348 Double_t CalculatesModelIntegral(const char *ModelName, Bool_t OnlyCounts = kFALSE);
349
354 TObject *CreateDrawObject(const char *ObjName, const char *DataName, Bool_t IsModel = kFALSE, const char *ModelName = nullptr);
355
364 void Draw(Option_t *option) override;
365
370 void ProvideLastDrawObjects(TObjArray& lst);
371
375 void PrintAmpls() const;
376
380 void PrintLines() const;
381
382 TString FindNextName(const char *Head, Int_t start, Bool_t isModel = kTRUE);
383
384 protected:
385
387
388 // a set of function to access data buffers
389
390 Bool_t ModelBuffersAllocated(TGo4FitModel *model);
391 Bool_t DataBuffersAllocated(TGo4FitData *data);
392
393 Int_t GetDataBinsSize(TGo4FitData *data) const;
394 Double_t *GetDataBinsValues(TGo4FitData *data) const;
395 Double_t *GetDataBinsDevs(TGo4FitData *data) const;
396 Double_t *GetDataBinsResult(TGo4FitData *data) const;
397 Double_t *GetModelBinsValues(TGo4FitModel *model, const char *DataName) const;
398
404 void RebuildAll(Bool_t ForceBuild = kFALSE);
405
406 void CollectAllPars() override;
407
408 Double_t DoCalculation() override;
409 Int_t DoNDFCalculation() override;
410 Double_t CalculateFCN(Int_t FitFunctionType, TGo4FitData *selectdata = nullptr);
411 Double_t PointFitFunction(Int_t FitFunctionType, Double_t value, Double_t modelvalue, Double_t standdev);
412
413 Bool_t InitFitterData() override;
414 void FinalizeFitterData() override;
415
419 TObjArray fxDatas;
420
424 TObjArray fxModels;
425
430
435
436 private:
437
440
445
446 TObjArray *fxDrawObjs{nullptr};
447
448 ClassDefOverride(TGo4Fitter,1)
449};
450
451#endif // TGO4FITTER_H
Double_t(* TUserFitFunction)(Double_t, Double_t, Double_t)
Definition TGo4Fitter.h:32
Basic abstract class, combining common properties of data and model.
Data object, which provides access to TGraph and TGraphErrors ROOT objects.
Data objects, which provides access to generic TH1 ROOT histogram.
Basic abstract class for representing data, which should be fitted.
Definition TGo4FitData.h:39
One dimensional gaussian peak.
Basic abstract class for representing model components of fitted data.
TGo4FitterAbstract()
Default constructor.
Central class of Go4Fit package.
Definition TGo4Fitter.h:38
TObjArray fxModels
Container for model components.
Definition TGo4Fitter.h:424
void FinalizeFitterData() override
Finalize fitter data.
void MoveDrawObjectsToROOT()
@ ff_least_squares
Definition TGo4Fitter.h:40
TGo4FitDataHistogram * AddH1(const char *DataName, TH1 *histo, Bool_t Owned=kFALSE, Double_t lrange=0., Double_t rrange=0.)
Create TGo4FitDataHistogram object and adds its to fitter.
TGo4FitDataHistogram * SetH1(const char *DataName, TH1 *histo, Bool_t Owned=kFALSE)
Set histogram to existing TGo4FitDataHistogram object.
void PrintLines() const
Print amplitude, position and width for each model components.
Bool_t ModelBuffersAllocated(TGo4FitModel *model)
Bool_t DataBuffersAllocated(TGo4FitData *data)
Double_t CalculatesModelIntegral(const char *ModelName, Bool_t OnlyCounts=kFALSE)
Calculates integral of model (if ModelName is specified) if OnlyCounts specified, only sum of values ...
Int_t fiMemoryUsage
Defines use of memory during actions executions.
Definition TGo4Fitter.h:434
Int_t GetNumComp() const
Return total number of TGo4FitComponent (data and model) objects in fitter.
Definition TGo4Fitter.h:282
Double_t * GetDataBinsResult(TGo4FitData *data) const
void SetUserFitFunction(TUserFitFunction iFunc)
Set user-defined fitted function.
Definition TGo4Fitter.h:98
TGo4FitModel * GetModel(Int_t n)
Return model component with given index.
Int_t NumModelsAssosiatedTo(const char *DataName)
Counts models associated with specific data.
void Print(Option_t *option="") const override
Print containment of fitter.
TGo4FitModel * RemoveModel(const char *ModelName, Bool_t IsDel=kFALSE)
Remove model component from fitter.
void DeleteModelsAssosiatedTo(const char *DataName)
Remove models associated with specific data.
const char * GetDataName(Int_t n)
Return name of data object with given index.
TGo4FitModelGauss1 * AddGauss1(const char *DataName, const char *ModelName, Double_t iPosition, Double_t iWidth, Double_t iAmpl=1., Int_t Axis=0)
Add 1-dim gaussian model to fitter.
TGo4FitData * AddData(TGo4FitData *d)
Add data object to fitter.
TGo4FitDataGraph * SetGraph(const char *DataName, TGraph *gr, Bool_t Owned=kFALSE)
Set graph to existing TGo4FitDataGraph object.
TGo4FitData * RemoveData(const char *DataName, Bool_t IsDel=kFALSE)
Remove data object from fitter.
Bool_t InitFitterData() override
Initialize fitter data.
void EstimateAmplitudes(Int_t NumIters=1)
Estimate amplitude of all model components.
void PrintAmpls() const
Print value of all amplitude parameters.
TUserFitFunction fxUserFitFunction
Pointer on user fit function.
Definition TGo4Fitter.h:444
Double_t * GetModelBinsValues(TGo4FitModel *model, const char *DataName) const
void AddAmplEstimation(Int_t NumIters=1)
Add amplitude estimation to actions list.
friend class TGo4FitAmplEstimation
Definition TGo4Fitter.h:386
TObjArray * fxDrawObjs
Definition TGo4Fitter.h:446
Double_t CalculateFCN(Int_t FitFunctionType, TGo4FitData *selectdata=nullptr)
void SetMemoryUsage(Int_t iMemoryUsage)
Set value of memory usage.
void CheckSlotsBeforeDelete(TGo4FitComponent *comp)
Double_t PointFitFunction(Int_t FitFunctionType, Double_t value, Double_t modelvalue, Double_t standdev)
void DeleteAllData()
Delete all data objects from fitter.
void ClearModelAssignmentTo(const char *ModelName, const char *DataName=nullptr)
Remove assignment to given data (if exists).
Int_t GetNumModel() const
Return number of model component in fitter.
Definition TGo4Fitter.h:179
Int_t GetDataBinsSize(TGo4FitData *data) const
TGo4FitModel * FindModel(const char *ModelName)
Return model component with given name.
TGo4FitModel * CloneModel(const char *ModelName, const char *NewName=nullptr)
Clones specified model.
void ProvideLastDrawObjects(TObjArray &lst)
Copy pointer on drawn object after last command to specified TObjArray.
void RebuildAll(Bool_t ForceBuild=kFALSE)
Update all data objects and model components according to current parameters values.
void AddPolynomX(const char *DataName, const char *NamePrefix, Int_t MaxOrder=1, Int_t GroupIndex=0, Double_t lrange=0., Double_t rrange=0.)
Construct 1-dim polynom for specified data object for x scale.
TString FindNextName(const char *Head, Int_t start, Bool_t isModel=kTRUE)
void Draw(Option_t *option) override
Draw fitter on current canvas.
void AddStandardActions()
Add list of standard actions to fitter.
TGo4FitModel * AddModel(TGo4FitModel *m)
Add model component to fitter.
Int_t GetMemoryUsage() const
Return value of memory usage parameter.
Definition TGo4Fitter.h:114
TGo4FitData * GetData(Int_t n)
Return data object with given index.
TGo4FitData * FindData(const char *DataName)
Return data object with given name.
Double_t CalculateFitFunction(Double_t *pars=nullptr, Int_t FitFunctionType=-1, const char *DataName=nullptr)
Calculate value of fit function for given set of parameters and specified type of fit function (it ca...
Int_t GetNumData() const
Return number of data objects in fitter.
Definition TGo4Fitter.h:124
void Clear(Option_t *option="") override
Remove all data, all models and all actions.
Int_t fiFitFunctionType
Defines type of fitted function.
Definition TGo4Fitter.h:429
void AddPolynoms(const char *DataName, const char *NamePrefix, Int_t MaxOrder=1, Int_t NumAxis=1, Int_t GroupIndex=0)
Construct full polynom for specified data object.
TObjArray fxDatas
Container for data objects.
Definition TGo4Fitter.h:419
Int_t CalculateNDF(const char *DataName=nullptr)
Calculates number of degree of freedom (NDF).
Double_t DoCalculation() override
Calculates value of fit function according current values of parameters.
Bool_t CalculatesMomentums(const char *DataName, Bool_t UseRanges, Bool_t SubstractModels, Double_t &first, Double_t &second)
Calculates first and second momentum for specified data Usage of ranges and subtraction of model can ...
Double_t CalculatesIntegral(const char *DataName, const char *ModelName=nullptr, Bool_t OnlyCounts=kFALSE)
Calculates integral for data or model (if ModelName is specified) if OnlyCounts specified,...
Int_t DoNDFCalculation() override
Calculates number of dimensions of freedom; Should be implemented in inherited classes.
void CollectAllPars() override
Should collect parameters from all associated to fitter objects.
void AssignModelTo(const char *ModelName, const char *DataName, Double_t RatioValue=1., Bool_t FixRatio=kFALSE)
Assign model to specified data object.
Double_t * GetDataBinsDevs(TGo4FitData *data) const
Double_t * GetDataBinsValues(TGo4FitData *data) const
void ChangeDataNameInAssignments(const char *oldname, const char *newname)
Change data name in model component assignments.
void SetFitFunctionType(Int_t iFitFunctionType)
Set fitted function type for minimization.
Definition TGo4Fitter.h:81
void FillSlotList(TSeqCollection *) override
Collect all TGo4FitSlot objects, situated in fitter, data and models to givven TObjArray.
TGo4FitDataGraph * AddGraph(const char *DataName, TGraph *gr, Bool_t Owned=kFALSE, Double_t lrange=0., Double_t rrange=0.)
Create TGo4FitDataGraph object and adds its to fitter.
Int_t GetFitFunctionType() const
Return type of fitted function.
Definition TGo4Fitter.h:86
virtual ~TGo4Fitter()
Destructor.
TObject * CreateDrawObject(const char *ObjName, const char *DataName, Bool_t IsModel=kFALSE, const char *ModelName=nullptr)
Create object (TH1 or TGraph), which can be drawn.
TGo4Fitter()
Default constructor.
void DeleteAllModels()
Delete all model objects from fitter.
TGo4FitComponent * GetComp(Int_t n)
Return TGo4FitComponent object with given index.