]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STAT/TStatToolkit.h
Removing useless flag.
[u/mrichter/AliRoot.git] / STAT / TStatToolkit.h
1 #ifndef TSTATTOOLKIT_H
2 #define TSTATTOOLKIT_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6 //
7 // some utilities which do net exist in the standard ROOT
8 //
9  
10 #include "TObject.h"
11 #include "TVectorD.h"
12 #include "TMatrixD.h"
13 //#include "TGraph2D.h"
14 //#include "TGraph.h"
15
16 class TH1F;
17 class TH1;
18 class TH2;
19 class TH3;
20 class TString;
21 class TTree;
22 class TGraph;
23 class TGraph2D;
24 class TCanvas;
25 class TMultiGraph; 
26 class TGraphErrors; 
27
28 class TStatToolkit : public TObject
29 {
30  public:
31   TStatToolkit();
32   virtual ~TStatToolkit();
33   //
34   //
35   //
36   static void    EvaluateUni(Int_t nvectors, Double_t *data, Double_t &mean, Double_t &sigma, Int_t hh);
37   static void    EvaluateUniExternal(Int_t nvectors, Double_t *data, Double_t &mean, Double_t &sigma, Int_t hh, Float_t externalfactor=1);
38   static Int_t  Freq(Int_t n, const Int_t *inlist, Int_t *outlist, Bool_t down);    
39   //
40   // HISTOGRAMS TOOLS
41   //
42   static  void TruncatedMean(const TH1 * his, TVectorD *param, Float_t down=0, Float_t up=1.0, Bool_t verbose=kFALSE);
43   static void MedianFilter(TH1 * his1D, Int_t nmedian);
44   static Bool_t  LTMHisto(TH1 * his, TVectorD &param , Float_t fraction=1);
45   //
46   static void LTM(TH1 * his, TVectorD *param=0 , Float_t fraction=1,  Bool_t verbose=kFALSE);
47   static Double_t  FitGaus(TH1* his, TVectorD *param=0, TMatrixD *matrix=0, Float_t xmin=0, Float_t xmax=0,  Bool_t verbose=kFALSE);
48   static Double_t  FitGaus(Float_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, TVectorD *param=0, TMatrixD *matrix=0, Bool_t verbose=kFALSE);
49   static Float_t  GetCOG(const Short_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, Float_t *rms=0, Float_t *sum=0);
50
51   static TGraph2D *  MakeStat2D(TH3 * his, Int_t delta0, Int_t delta1, Int_t type);
52   static TGraphErrors *  MakeStat1D(TH2 * his, Int_t deltaBin, Double_t fraction, Int_t returnType, Int_t markerStyle, Int_t markerColor);
53   //
54   // Graph tools
55   //
56   static TGraph * MakeGraphSparse(TTree * tree, const char * expr="Entry", const char * cut="1",  Int_t mstyle=25, Int_t mcolor=1, Float_t msize=-1, Float_t offset=0.0);
57   static TGraphErrors * MakeGraphErrors(TTree * tree, const char * expr="Entry", const char * cut="1",  Int_t mstyle=25, Int_t mcolor=1, Float_t msize=-1, Float_t offset=0.0);
58
59   //
60   // Fitting function
61   //
62   static TString* FitPlane(TTree * tree, const char* drawCommand, const char* formula, const char* cuts, Double_t & chi2, Int_t &npoints,  TVectorD &fitParam, TMatrixD &covMatrix, Float_t frac=-1, Int_t start=0, Int_t stop=10000000, Bool_t fix0=kFALSE);
63   static TString* FitPlaneFixed(TTree * tree, const char* drawCommand, const char* formula, const char* cuts, Double_t & chi2, Int_t &npoints,  TVectorD &fitParam, TMatrixD &covMatrix, Float_t frac=-1, Int_t start=0, Int_t stop=10000000);
64   //
65   //Linear fitter helper function
66   //
67   static TString* FitPlaneConstrain(TTree * tree, const char* drawCommand, const char* formula, const char* cuts, Double_t & chi2, Int_t &npoints,  TVectorD &fitParam, TMatrixD &covMatrix, Float_t frac=-1, Int_t start=0, Int_t stop=10000000, Double_t constrain=-1);
68   static Int_t GetFitIndex(const TString fString, const TString subString);
69  static TString FilterFit(const TString &input, const TString filter, TVectorD &vec, TMatrixD &covar);
70  static void Update1D(Double_t delta, Double_t sigma, Int_t s1, TMatrixD &param, TMatrixD &covar);
71   static void   Constrain1D(const TString &input, const TString filter, TVectorD &param, TMatrixD & covar, Double_t mean, Double_t sigma);
72   static TString  MakeFitString(const TString &input, const TVectorD &param, const TMatrixD & covar, Bool_t verbose=kFALSE);
73   //
74   // TTree function for the trneding
75   //
76   static Int_t  MakeStatAlias(TTree * tree, const char * expr, const char * cut, const char * alias);
77   static Int_t  SetStatusAlias(TTree * tree, const char * expr, const char * cut, const char * alias);
78   static TMultiGraph*  MakeStatusMultGr(TTree * tree, const char * expr, const char * cut, const char * alias, Int_t igr);  
79   static void  AddStatusPad(TCanvas* c1, Float_t padratio, Float_t bottommargin);
80   static void  DrawStatusGraphs(TObjArray* oaMultGr);
81   //
82   // TTree function for robust draw
83   //
84   static TH1* DrawHistogram(TTree * tree, const char* drawCommand, const char* cuts = "1", const char* hname = "histo", const char* htitle = "histo", Int_t nsigma = 4, Float_t fraction = 0.75);
85   //
86   // TestFunctions:
87   //
88  static  void TestGausFit(Int_t nhistos=5000);
89
90  ClassDef(TStatToolkit,0) // Various mathematical tools for physics analysis - which are not included in ROOT TMath
91  
92 };
93 #endif