GSI Object Oriented Online Offline (Go4)  GO4-6.2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
TGo4FitData.cxx
Go to the documentation of this file.
1 // $Id: TGo4FitData.cxx 3490 2022-01-21 14:19:46Z linev $
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 "TGo4FitData.h"
15 
16 #include <iostream>
17 
18 #include "TBuffer.h"
19 #include "TMath.h"
20 #include "TH1.h"
21 #include "TH2.h"
22 #include "TH3.h"
23 #include "TGraph.h"
24 #include "TClass.h"
25 
26 #include "TGo4FitAxisTrans.h"
27 
28 
30  fiDataType(0), fbUseBinScale(kFALSE), fiTakeSigmasFrom(1), fdSigmaValue(1.), fdExcludeLessThen(0.),
31  fxAxisTrans()
32 {
34 }
35 
36 TGo4FitData::TGo4FitData(const char* iName, const char* iTitle, Int_t iDataType, Bool_t AddAmpl) :
37  TGo4FitComponent(iName,iTitle), fiDataType(iDataType),
38  fbUseBinScale(kFALSE), fiTakeSigmasFrom(1), fdSigmaValue(1.), fdExcludeLessThen(0.),
39  fxAxisTrans()
40 {
41 
43 
44  if (AddAmpl)
45  NewAmplitude("Ampl", 1.0, kTRUE);
46 
47  fxAxisTrans.SetOwner(kTRUE);
48 }
49 
51 {
53 }
54 
55 
57 {
58  Int_t oldnum = GetNumberOfTransSlots();
59  if ( (nslots<0) || (nslots == oldnum) ) return kFALSE;
60 
61  if (oldnum<nslots)
62  for(Int_t n=oldnum;n<nslots;n++) {
63  TString name("Trans");
64  name+=n;
65  fxAxisTrans.Add(new TGo4FitSlot(name.Data(),"Axis transformation", this, TGo4FitAxisTrans::Class(), kFALSE));
66  }
67  else
68  for (Int_t n=oldnum;n>nslots;n--) {
69  TObject* slot = fxAxisTrans.Last();
70  fxAxisTrans.Remove(slot);
71  fxAxisTrans.Compress();
72  delete slot;
73  }
75  return kTRUE;
76 }
77 
79 {
80  return (nslot>=0) && (nslot<=fxAxisTrans.GetLast()) ? dynamic_cast<TGo4FitSlot*> (fxAxisTrans[nslot]) : 0;
81 }
82 
84 {
85  TGo4FitSlot* slot = GetAxisTransSlot(nslot);
86  return (slot==0) ? 0 : dynamic_cast<TGo4FitAxisTrans*> (slot->GetObject());
87 }
88 
89 void TGo4FitData::SetAxisTrans(Int_t nslot, TGo4FitAxisTrans *Trans, Bool_t TransOwned)
90 {
91  if (nslot<0) return;
92  if(nslot>=GetNumberOfTransSlots())
93  SetNumberOfTransSlots(nslot+1);
94  ((TGo4FitSlot*) (fxAxisTrans[nslot]))->SetObject(Trans,TransOwned);
95 }
96 
97 void TGo4FitData::AddAxisTrans(TGo4FitAxisTrans* Trans, Bool_t TransOwned) {
98  Int_t nslot = GetNumberOfTransSlots();
99  SetNumberOfTransSlots(nslot+1);
100  ((TGo4FitSlot*) (fxAxisTrans[nslot]))->SetObject(Trans,TransOwned);
101 }
102 
103 void TGo4FitData::SetAxisTransNeeded(Int_t nslot, Bool_t iNeeded) {
104  if( iNeeded && (nslot>=GetNumberOfTransSlots()))
105  SetNumberOfTransSlots(nslot+1);
106  if ((nslot>=0) && (nslot<GetNumberOfTransSlots()))
107  ((TGo4FitSlot*) (fxAxisTrans[nslot]))->SetNeeded(iNeeded);
108 }
109 
111 {
112  if (GetUseBinScale() || (GetExcludeLessThen()>0)) return kTRUE;
113  for (Int_t n=0;n<GetNumberOfTransSlots();n++)
114  if (GetAxisTrans(n)) return kTRUE;
115  return kFALSE;
116 }
117 
118 TObject* TGo4FitData::CreateDrawObject(const char* ObjName)
119 {
120  TGo4FitDataIter* iter = MakeIter();
121  if (iter==0) return 0;
122  TObject* obj = iter->CreateDrawObject(ObjName);
123  delete iter;
124  return obj;
125 }
126 
127 Bool_t TGo4FitData::Initialize(Int_t UseBuffers)
128 {
129  TGo4FitDataIter *iter = MakeIter();
130  if (iter == 0)
131  return kFALSE;
132 
133  fiBinsSize = iter->CountPoints(kTRUE);
134 
135  fiIndexesSize = iter->IndexesSize();
136  fiScalesSize = iter->ScalesSize();
137 
138  Bool_t use = ((UseBuffers < 0) && GetUseBuffers()) || (UseBuffers > 0);
139 
140  if (use)
141  for (Int_t n = 0; n < GetNumberOfTransSlots(); n++) {
142  TGo4FitAxisTrans *trans = GetAxisTrans(n);
143  if (trans && !trans->IsAllParsFixed()) {
144  use = kFALSE;
145  break;
146  }
147  }
148 
149  if (use) {
150 
151  fxValues = new Double_t[fiBinsSize];
152  fxStandDev = new Double_t[fiBinsSize];
153  fxBinsResult = new Double_t[fiBinsSize];
154 
155  if (iter->HasIndexes())
156  fxFullIndex = new Int_t[fiBinsSize * fiIndexesSize];
157  fxFullScale = new Double_t[fiBinsSize * fiScalesSize];
158  if (iter->HasWidths())
159  fxFullWidth = new Double_t[fiBinsSize * fiScalesSize];
160 
161  Int_t nbin = 0;
162  if (iter->Reset())
163  do {
164 
165  fxValues[nbin] = iter->Value();
166  fxStandDev[nbin] = iter->StandardDeviation();
167 
168  if (fxFullIndex)
169  for (Int_t n = 0; n < fiIndexesSize; n++)
170  fxFullIndex[nbin * fiIndexesSize + n] = iter->Indexes()[n];
171 
172  if (fxFullScale)
173  for (Int_t naxis = 0; naxis < fiScalesSize; naxis++)
174  fxFullScale[nbin * fiScalesSize + naxis] =
175  iter->Scales()[naxis];
176 
177  if (fxFullWidth && iter->HasWidths())
178  for (Int_t naxis = 0; naxis < fiScalesSize; naxis++)
179  fxFullWidth[nbin * fiScalesSize + naxis] =
180  iter->Widths()[naxis];
181 
182  nbin++;
183  } while (iter->Next());
184  }
185 
186  delete iter;
187 
188  return kTRUE;
189 }
190 
192 {
194 }
195 
197 {
198  fiBinsSize = 0;
199  fiIndexesSize = 0;
200  fiScalesSize = 0;
201 
202  fxValues = 0;
203  fxStandDev = 0;
204  fxBinsResult = 0;
205 
206  fxFullScale = 0;
207  fxFullWidth = 0;
208 
209  fxFullIndex = 0;
210 }
211 
213  if (fxValues) delete[] fxValues;
214 
215  if (fxStandDev) delete[] fxStandDev;
216 
217  if (fxBinsResult) delete[] fxBinsResult;
218 
219  if (fxFullIndex) delete[] fxFullIndex;
220 
221  if (fxFullScale) delete[] fxFullScale;
222 
223  if (fxFullWidth) delete[] fxFullWidth;
224 
226 }
227 
228 Bool_t TGo4FitData::DefineScaleMinMax(Int_t naxis, Double_t& min, Double_t& max) {
229  TGo4FitDataIter* iter = MakeIter();
230  if (iter==0) return kFALSE;
231  Bool_t res = kFALSE;
232  if (iter->Reset(kFALSE) && (iter->ScalesSize()<=naxis)) {
233  min = iter->Scales()[naxis]; max = min;
234  do {
235  Double_t value = iter->Scales()[naxis];
236  if (value<min) min = value; else
237  if (value>max) max = value;
238  } while (iter->Next(kFALSE));
239  res = kTRUE;
240  }
241 
242  delete iter;
243  return res;
244 }
245 
247  TGo4FitDataIter* iter = MakeIter();
248  if (iter==0) return 0;
249  Int_t res = 0;
250  if (iter->Reset(kFALSE)) res = iter->IndexesSize();
251  delete iter;
252  return res;
253 }
254 
256 {
257  TGo4FitDataIter* iter = MakeIter();
258  if (iter==0) return 0;
259 
260  Int_t res = iter->CountPoints(kTRUE);
261  delete iter;
262 
263  return res;
264 }
265 
266 const Double_t* TGo4FitData::GetScaleValues(const Int_t nbin)
267 {
268  if(fxFullScale) return &(fxFullScale[nbin*GetScalesSize()]);
269  else return 0;
270 }
271 
272 const Double_t* TGo4FitData::GetWidthValues(const Int_t nbin)
273 {
274  if(fxFullWidth) return &(fxFullWidth[nbin*GetScalesSize()]);
275  else return 0;
276 }
277 
278 const Int_t* TGo4FitData::GetFullIndex(Int_t nbin)
279 {
280  if (fxFullIndex) return &(fxFullIndex[nbin*GetIndexesSize()]);
281  else return 0;
282 }
283 
285 {
286  if (data==0) return kFALSE;
287  TGo4FitDataIter* iter = data->MakeIter();
288  if (iter==0) return kFALSE;
289 
290  Bool_t res = kFALSE;
291  if (iter->Reset(kFALSE)) res = (iter->IndexesSize()==GetIndexesSize()) && (GetIndexesSize()>0);
292  delete iter;
293 
294  return res;
295 }
296 
298 {
299 
300  if (ModelMask==0) return;
301 
302  if (BuffersAllocated())
303  for(Int_t nbin=0;nbin<GetBinsSize();nbin++) {
304  const Double_t* values = GetScaleValues(nbin);
305 
306  Bool_t res = model->CheckRangeConditions(values, GetScalesSize());
307 
308  ModelMask[nbin] = res ? 1 : 0;
309  }
310  else {
311  TGo4FitDataIter* iter = MakeIter();
312  Int_t nbin = 0;
313  if (iter->Reset()) do {
314  Bool_t res = model->CheckRangeConditions(iter->Scales(), iter->ScalesSize());
315  ModelMask[nbin] = res ? 1 : 0;
316  nbin++;
317  } while (iter->Next());
318  }
319 }
320 
321 void TGo4FitData::FillSlotList(TSeqCollection* list)
322 {
324  for(Int_t n=0;n<=fxAxisTrans.GetLast();n++)
325  list->Add(fxAxisTrans[n]);
326 }
327 
328 void TGo4FitData::Print(Option_t* option) const
329 {
330  TGo4FitComponent::Print(option);
331  std::cout << " Data type: ";
332  switch(fiDataType) {
333  case 1: std::cout << "histogram" << std::endl; break;
334  case 2: std::cout << "graph" << std::endl; break;
335  default: std::cout << fiDataType << std::endl;
336  }
337  std::cout << " Use bin scale: " << fbUseBinScale << std::endl;
338  std::cout << " Take sigmas from: " ;
339  switch(GetSigmaSource()) {
340  case 0: std::cout << "none" << std::endl; break;
341  case 1: std::cout << "data" << std::endl; break;
342  case 2: std::cout << "const value " << GetSigmaValue() << std::endl; break;
343  }
344  std::cout << " Exclude bins less then: " << GetExcludeLessThen() << std::endl;
345  std::cout << " Axis transformation data: " << std::endl;
346  fxAxisTrans.Print(option);
347 }
348 
349 void TGo4FitData::Streamer(TBuffer& b)
350 {
351  if (b.IsReading()) {
352 
353  TGo4FitData::Class()->ReadBuffer(b, this);
354 
355  for(Int_t n=0;n<=fxAxisTrans.GetLast();n++) {
356  TGo4FitSlot* dc = (TGo4FitSlot*) fxAxisTrans[n];
357  dc->SetDefaults(this, TGo4FitAxisTrans::Class());
358  }
359 
360  } else {
361  TGo4FitData::Class()->WriteBuffer(b, this);
362  }
363 }
364 
365 // *******************************************************************************************
366 
367 
369  fxIndexes(), fxScales(), fxWidths(), fdValue(0.), fdStandardDeviation(1.), fiNumPoint(0), fbReachEnd(kTRUE) {
370 }
371 
373 }
374 
375 Bool_t TGo4FitDataIter::ReserveArrays(Int_t NumDimen, Int_t NumOwnAxis, Bool_t HasWidth) {
376  TGo4FitData* data = GetData();
377  if (data==0) return kFALSE;
378 
379  fxIndexes.Set(NumDimen); fxIndexes.Reset(0);
380 
381  Int_t size = 0;
382  if (data->GetUseBinScale()) size = NumDimen;
383  else size = NumOwnAxis;
384 
385  if (size<=0) return kFALSE;
386 
387  fxScales.Set(size); fxScales.Reset(0.);
388  if (HasWidth) { fxWidths.Set(size); fxWidths.Reset(1.); }
389  else fxWidths.Set(0);
390 
391  return kTRUE;
392 }
393 
394 void TGo4FitDataIter::TransformScales(Double_t* scales)
395 {
396  TGo4FitData* data = GetData();
397  for(Int_t nslot=0;nslot<data->GetNumberOfTransSlots();nslot++) {
398  TGo4FitAxisTrans* trans = data->GetAxisTrans(nslot);
399  if (trans) trans->Transformation(scales, ScalesSize());
400  }
401 }
402 
403 Bool_t TGo4FitDataIter::ProduceScales(const Int_t* index, const Double_t* ownscales, const Double_t* ownwidths)
404 {
405  TGo4FitData* data = GetData();
406  if (data==0) return kFALSE;
407 
408  if ( (data->GetUseBinScale()) || (ownscales==0) ) {
409  if (index==0) return kFALSE;
410  Double_t add = (data->GetDataType() == TGo4FitData::dtHistogram) ? .5 : 0.;
411  for(Int_t n=0;n<fxScales.GetSize();n++)
412  fxScales[n] = index[n] + add;
413  fxWidths.Reset(1.);
414  } else {
415  for(Int_t n=0; n<fxScales.GetSize();n++)
416  fxScales[n] = ownscales[n];
417  if (ownwidths!=0)
418  for(Int_t n=0; n<fxWidths.GetSize();n++)
419  fxWidths[n] = ownwidths[n];
420  }
421 
422  if (data->GetNumberOfTransSlots()>0) {
423  if (fxWidths.GetSize()==ScalesSize()) {
424  TArrayD arr1(ScalesSize()), arr2(ScalesSize());
425  for(Int_t n=0;n<ScalesSize();n++) {
426  arr1[n] = fxScales[n]-fxWidths[n]/2.;
427  arr2[n] = fxScales[n]+fxWidths[n]/2.;
428  }
429  TransformScales(arr1.GetArray());
430  TransformScales(arr2.GetArray());
431  for(Int_t n=0;n<ScalesSize();n++)
432  fxWidths[n] = TMath::Abs(arr2[n]-arr1[n]);
433  }
434 
435  TransformScales(fxScales.GetArray());
436  }
437 
438  return kTRUE;
439 }
440 
441 Bool_t TGo4FitDataIter::NextIndex(TArrayI& Index, TArrayI& Limits) {
442  Int_t n=0;
443  while (n<Index.GetSize()) {
444  Index[n]++;
445  if (Index[n]<Limits[n]) return kTRUE;
446  Index[n] = 0; n++;
447  }
448  return kFALSE;
449 }
450 
452  TGo4FitData* data = GetData();
453  if (data==0) return kFALSE;
454  if (data->GetSigmaSource()==2) {
456  return kTRUE;
457  } else return kFALSE;
458 }
459 
461  TGo4FitData* data = GetData();
462  if (data==0) return kFALSE;
463  if (Value()<data->GetExcludeLessThen()) return kFALSE;
464  return data->CheckRangeConditions(Scales(),ScalesSize());
465 }
466 
467 Bool_t TGo4FitDataIter::Reset(Bool_t UseRanges) {
468  fbReachEnd = kTRUE;
469 
470  if (!StartReset()) return kFALSE;
471 
472  fiNumPoint = 0;
473 
474  if (!ReadCurrentPoint()) return kFALSE;
475  if (!UseRanges) { fbReachEnd = kFALSE; return kTRUE; }
476 
477  while (!CheckPointForRange()) {
478  if (!ShiftToNextPoint()) return kFALSE;
479  if (!ReadCurrentPoint()) return kFALSE;
480  }
481 
482  fbReachEnd = kFALSE;
483  return kTRUE;
484 }
485 
486 Bool_t TGo4FitDataIter::Next(Bool_t UseRanges) {
487  fiNumPoint++;
488 
489  if (fbReachEnd || (GetData()==0)) { fbReachEnd = kTRUE; return kFALSE; }
490 
491  do {
492  if (!ShiftToNextPoint()) { fbReachEnd = kTRUE; return kFALSE; }
493 
494  if (!ReadCurrentPoint()) { fbReachEnd = kTRUE; return kFALSE; }
495 
496  if (!UseRanges) return kTRUE;
497 
498  } while (!CheckPointForRange());
499 
500  return kTRUE;
501 }
502 
503 Double_t TGo4FitDataIter::xWidths() const
504 {
505  double res = 1.;
506  if(HasWidths())
507  for(int n=0;n<fxWidths.GetSize();n++)
508  res=res*fxWidths[n];
509  return res;
510 }
511 
512 Int_t TGo4FitDataIter::CountPoints(Bool_t UseRanges) {
513  if (!Reset(UseRanges)) return 0;
514  Int_t cnt=0;
515  do {
516  cnt+=1;
517  } while (Next(UseRanges)) ;
518  return cnt;
519 }
520 
521 Bool_t TGo4FitDataIter::DefineIndexesLimits(TArrayI& Limits) {
522  if (!Reset(kFALSE)) return kFALSE;
523  if (IndexesSize()<=0) return kFALSE;
524  Limits.Set(IndexesSize()); Limits.Reset(0);
525  do {
526  for(Int_t n=0;n<IndexesSize();n++)
527  if (Indexes()[n]>Limits[n]) Limits[n] = Indexes()[n];
528  } while(Next(kFALSE));
529  return kTRUE;
530 }
531 
532 TH1* TGo4FitDataIter::CreateHistogram(const char* HistoName, Bool_t UseRanges, Bool_t SetBins) {
533  TArrayI Limits;
534  if (!DefineIndexesLimits(Limits)) return 0;
535  if (!HasIndexes() || (IndexesSize()!=ScalesSize()) || !HasWidths()) return 0;
536 
537  Int_t NumDim = IndexesSize();
538  if (NumDim>3) NumDim=3;
539 
540  Double_t* dummy = 0;
541  TH1* histo = 0;
542  switch(NumDim) {
543  case 1: histo = new TH1D(HistoName, "result", Limits[0]+1, dummy); break;
544  case 2: histo = new TH2D(HistoName, "result", Limits[0]+1, dummy, Limits[1]+1, dummy); break;
545  case 3: histo = new TH3D(HistoName, "result", Limits[0]+1, dummy, Limits[1]+1, dummy, Limits[2]+1, dummy); break;
546  default: return 0;
547  }
548 
549  histo->SetDirectory(0);
550 
551  Double_t* Axises[3];
552  for (Int_t n=0;n<NumDim;n++)
553  Axises[n] = new Double_t[Limits[n]+2];
554 
555  Double_t ampl = GetData()->GetAmplValue();
556 
557  if (Reset(UseRanges)) do {
558  if (SetBins)
559  switch (NumDim) {
560  case 1: histo->SetBinContent(Indexes()[0]+1, ampl*Value()); break;
561  case 2: histo->SetBinContent(Indexes()[0]+1, Indexes()[1]+1, ampl*Value()); break;
562  case 3: histo->SetBinContent(Indexes()[0]+1, Indexes()[1]+1, Indexes()[2]+1, ampl*Value()); break;
563  }
564  for(Int_t n=0;n<NumDim;n++) {
565  Int_t indx = Indexes()[n];
566  Axises[n][indx] = Scales()[n]-Widths()[n]/2.;
567  Axises[n][indx+1] = Scales()[n]+Widths()[n]/2.;
568  }
569  } while(Next(UseRanges));
570 
571  histo->GetXaxis()->Set(Limits[0]+1,Axises[0]);
572  if (NumDim>1) histo->GetYaxis()->Set(Limits[1]+1,Axises[1]);
573  if (NumDim>2) histo->GetZaxis()->Set(Limits[2]+1,Axises[2]);
574 
575  for (Int_t n=0;n<NumDim;n++)
576  delete[] Axises[n];
577 
578  return histo;
579 }
580 
581 TGraph* TGo4FitDataIter::CreateGraph(const char* GraphName, Bool_t UseRanges, Bool_t SetBins) {
582  Int_t NumPoints = CountPoints(UseRanges);
583  if ((NumPoints<=0) || (ScalesSize()<1)) return 0;
584 
585  TGraph* gr = new TGraph(NumPoints);
586  gr->SetName(GraphName);
587  if (Reset(UseRanges)) do {
588  (gr->GetX())[Point()] = x();
589  if (SetBins)
590  (gr->GetY())[Point()] = GetData()->GetAmplValue() * Value();
591  } while(Next(UseRanges));
592 
593  return gr;
594 }
595 
596 TObject* TGo4FitDataIter::CreateDrawObject(const char* ObjName) {
597  if (!Reset(kFALSE)) return 0;
598  if (HasIndexes() && (IndexesSize()==ScalesSize()) && HasWidths()) return CreateHistogram(ObjName, kFALSE, kTRUE);
599  else return CreateGraph(ObjName, kFALSE, kTRUE);
600 }
Bool_t DefineScaleMinMax(Int_t naxis, Double_t &min, Double_t &max)
Int_t fiBinsSize
Definition: TGo4FitData.h:346
Double_t xWidths() const
TArrayD fxWidths
Definition: TGo4FitData.h:620
TH1 * CreateHistogram(const char *HistoName, Bool_t UseRanges=kFALSE, Bool_t SetBins=kFALSE)
void SetUpdateSlotList()
Bool_t ProduceScales(const Int_t *index, const Double_t *ownscales, const Double_t *ownwidths)
TArrayI fxIndexes
Definition: TGo4FitData.h:617
Double_t fdStandardDeviation
Definition: TGo4FitData.h:622
virtual void Print(Option_t *option) const
Bool_t HasWidths() const
Definition: TGo4FitData.h:499
TGo4FitSlot * SetObject(TObject *obj, Bool_t iOwned=kFALSE)
virtual Bool_t IsAnyDataTransform()
Bool_t CheckRangeConditions(const Double_t *values, Int_t numaxis)
virtual void Print(Option_t *option) const
Double_t x() const
Definition: TGo4FitData.h:484
Double_t * fxBinsResult
Definition: TGo4FitData.h:371
Double_t * fxFullScale
Definition: TGo4FitData.h:378
Bool_t SetNumberOfTransSlots(Int_t nslots)
Definition: TGo4FitData.cxx:56
Int_t fiScalesSize
Definition: TGo4FitData.h:356
const Int_t * GetFullIndex(Int_t nbin)
virtual TGo4FitDataIter * MakeIter()
Definition: TGo4FitData.h:165
virtual void Finalize()
Int_t fiDataType
Definition: TGo4FitData.h:297
TObject * CreateDrawObject(const char *ObjName)
void ReleaseAllPointers()
Bool_t DefineIndexesLimits(TArrayI &Limits)
Int_t GetSigmaSource() const
Definition: TGo4FitData.h:83
void TransformScales(Double_t *scales)
TGraph * CreateGraph(const char *GraphName, Bool_t UseRanges=kFALSE, Bool_t SetBins=kFALSE)
Bool_t NextIndex(TArrayI &Index, TArrayI &Limits)
Int_t GetNumberOfTransSlots()
Definition: TGo4FitData.h:126
TGo4FitSlot * GetAxisTransSlot(Int_t nslot)
Definition: TGo4FitData.cxx:78
virtual Bool_t Initialize(Int_t UseBuffers=-1)
Double_t GetSigmaValue() const
Definition: TGo4FitData.h:89
TGo4FitParameter * NewAmplitude(const char *Name=0, Double_t iValue=0., Bool_t IsFixed=kFALSE, Int_t AtIndx=0)
void ResetAllPoinetrs()
virtual Bool_t Reset(Bool_t UseRanges=kTRUE)
Int_t DefineBinsSize()
const Double_t * Widths() const
Definition: TGo4FitData.h:504
Bool_t IsCompatibleData(TGo4FitData *data)
Bool_t GetUseBinScale() const
Definition: TGo4FitData.h:62
Int_t Point() const
Definition: TGo4FitData.h:524
virtual ~TGo4FitDataIter()
virtual Bool_t StartReset()=0
Bool_t BuffersAllocated() const
Definition: TGo4FitData.h:238
virtual ~TGo4FitData()
Definition: TGo4FitData.cxx:50
const Double_t * GetScaleValues(const Int_t nbin)
Int_t CountPoints(Bool_t UseRanges=kTRUE)
Double_t * fxStandDev
Definition: TGo4FitData.h:366
Int_t GetBinsSize() const
Definition: TGo4FitData.h:243
virtual Bool_t Next(Bool_t UseRanges=kTRUE)
virtual TGo4FitData * GetData() const =0
TObject * CreateDrawObject(const char *ObjName)
Int_t GetDataType() const
Definition: TGo4FitData.h:57
Bool_t ReserveArrays(Int_t NumDimen, Int_t NumOwnAxis, Bool_t HasWidth)
virtual void FillSlotList(TSeqCollection *lst)
virtual void FillSlotList(TSeqCollection *list)
TGo4FitAxisTrans * GetAxisTrans(Int_t nslot)
Definition: TGo4FitData.cxx:83
Int_t DefineDimensions()
void SetAxisTrans(Int_t nslot, TGo4FitAxisTrans *Trans, Bool_t TransOwned=kFALSE)
Definition: TGo4FitData.cxx:89
Bool_t GetDeviation()
Bool_t fbUseBinScale
Definition: TGo4FitData.h:303
Bool_t HasIndexes() const
Definition: TGo4FitData.h:454
Int_t fiIndexesSize
Definition: TGo4FitData.h:351
void SetDefaults(TNamed *iOwner, TClass *iClass)
Definition: TGo4FitSlot.cxx:68
Double_t Value() const
Definition: TGo4FitData.h:514
const Int_t * Indexes() const
Definition: TGo4FitData.h:464
const Double_t * Scales() const
Definition: TGo4FitData.h:479
void AddAxisTrans(TGo4FitAxisTrans *Trans, Bool_t TransOwned=kFALSE)
Definition: TGo4FitData.cxx:97
virtual Bool_t ShiftToNextPoint()=0
Int_t GetScalesSize() const
Definition: TGo4FitData.h:249
virtual void Transformation(Double_t *scales, Int_t naxis)=0
virtual Bool_t ReadCurrentPoint()=0
Double_t GetExcludeLessThen() const
Definition: TGo4FitData.h:112
TObject * GetObject() const
Int_t IndexesSize() const
Definition: TGo4FitData.h:459
Double_t StandardDeviation() const
Definition: TGo4FitData.h:519
Bool_t CheckPointForRange()
TObjArray fxAxisTrans
Definition: TGo4FitData.h:324
void ApplyRangesForModelMask(TGo4FitComponent *model, Char_t *ModelMask)
Int_t GetIndexesSize() const
Definition: TGo4FitData.h:291
Int_t * fxFullIndex
Definition: TGo4FitData.h:392
Double_t * fxValues
Definition: TGo4FitData.h:361
const Double_t * GetWidthValues(const Int_t nbin)
TArrayD fxScales
Definition: TGo4FitData.h:619
Double_t * fxFullWidth
Definition: TGo4FitData.h:385
void SetAxisTransNeeded(Int_t nslot, Bool_t iNeeded=kFALSE)
Int_t ScalesSize() const
Definition: TGo4FitData.h:474