KaliVeda  1.13/01
Heavy-Ion Analysis Toolkit
KVHistoManipulator.h
Go to the documentation of this file.
1 /*
2 $Id: KVHistoManipulator.h,v 1.8 2009/04/07 14:54:15 ebonnet Exp $
3 $Revision: 1.8 $
4 $Date: 2009/04/07 14:54:15 $
5 */
6 
9 
10 #ifndef __KVHISTOMANIPULATOR_H
11 #define __KVHISTOMANIPULATOR_H
12 
13 #include "Riostream.h"
14 
15 #include "TH1.h"
16 #include "TH2.h"
17 #include "TCutG.h"
18 #include "KVList.h"
19 #include "TString.h"
20 #include "TGraph.h"
21 #include "TGraphErrors.h"
22 #include "TMath.h"
23 #include <TF1.h>
24 
25 class TCanvas;
26 class KVNumberList;
27 class TMultiGraph;
28 
38 
39 public:
40 
41  void init(void)
42  {
43 
44  }
45 
47  virtual ~KVHistoManipulator(void);
48 
49  void SetVisDebug(Bool_t on = kTRUE)
50  {
59  kVisDebug = on;
60  };
62  {
63  return kVisDebug;
64  };
65 
66  Int_t CutStatBin(TH1* hh, Int_t stat_min = -1, Int_t stat_max = -1);
67 
68  Int_t Apply_TCutG(TH2* hh, TCutG* cut, TString mode = "in");
69 
70  TH1* ScaleHisto(TH1* hh, TF1* fx, TF1* fy = NULL, Int_t nx = -1, Int_t ny = -1,
71  Double_t xmin = -1., Double_t xmax = -1., Double_t ymin = -1., Double_t ymax = -1., Option_t* norm = "");
72  TGraph* ScaleGraph(const TGraph* hh, TF1* fx = nullptr, TF1* fy = nullptr) const
73  {
75 
76  TString axis;
77  if (fx) axis = "X";
78  if (fy) axis.Append("Y");
79  if (axis == "XY") return ScaleGraph(hh, axis, *fx, *fy);
80  else if (axis == "X") return ScaleGraph(hh, axis, *fx);
81  return ScaleGraph(hh, axis, *fy);
82  }
83  TGraph* ScaleGraph(const TGraph* hh, const TString& axis, const TF1& f1, const TF1& f2 = TF1("f2", "x")) const;
84 
85  TH1* CentreeReduite(TH1* hh, Int_t nx = -1, Int_t ny = -1, Double_t xmin = -1., Double_t xmax = -1., Double_t ymin = -1., Double_t ymax = -1.);
86  TH2* CentreeReduiteX(TH2* hh, Int_t nx = -1, Double_t xmin = -1., Double_t xmax = -1.);
87  TH2* CentreeReduiteY(TH2* hh, Int_t ny = -1, Double_t ymin = -1., Double_t ymax = -1.);
88 
89  TH2* RenormaliseHisto(TH2* hh, Int_t bmin = -1, Int_t bmax = -1, TString axis = "X", Double_t valref = 1);
90  TH2* RenormaliseHisto(TH2* hh, Double_t valmin, Double_t valmax, TString axis = "X", Double_t valref = 1);
91 
92  TH1* CumulatedHisto(TH1* hh, TString direction = "C", Int_t bmin = -1, Int_t bmax = -1, Option_t* norm = "surf");
93  TH1* CumulatedHisto(TH1* hh, Double_t xmin, Double_t xmax, const TString& direction = "C", Option_t* norm = "surf");
94  TH1* GetDerivative(TH1* hh, Int_t order);
95 
96  TGraphErrors* GetMomentEvolution(TH2* hh, TString momentx, TString momenty, TString axis = "Y", Double_t stat_min = 0);
98  TGraph* LinkGraphs(TGraph* grx, TGraph* gry);
99 
100  KVList* Give_ProjectionList(TH2* hh, Double_t MinIntegral = -1, TString axis = "X");
101  KVNumberList* Saucisson(TH1* hh, Int_t ntranches = 10);
102  TH2* PermuteAxis(TH2* hh);
103  TGraph* PermuteAxis(TGraph* gr);
104  TGraphErrors* MakeGraphFrom(TProfile* pf, Bool_t Error = kTRUE);
105 
106  void DefinePattern(TH1* ob, TString titleX = "42 0.08 0.8", TString titleY = "42 0.07 1.2", TString labelX = "42 0.05 0.005", TString labelY = "42 0.05 0.006");
107  void DefinePattern(TGraph* ob, TString titleX = "42 0.08 0.8", TString titleY = "42 0.07 1.2", TString labelX = "42 0.05 0.005", TString labelY = "42 0.05 0.006");
108  void DefinePattern(TF1* ob, TString titleX = "42 0.08 0.8", TString titleY = "42 0.07 1.2", TString labelX = "42 0.05 0.005", TString labelY = "42 0.05 0.006");
109  void DefinePattern(TAxis* ax, TString title = "42 0.08 0.8", TString label = "42 0.05 0.005");
110 
111  void DefineLineStyle(TAttLine* ob, TString line);
112  void DefineMarkerStyle(TAttMarker* ob, TString marker);
113  void DefineStyle(TObject* ob, TString line, TString marker);
114 
115  void DefineTitle(TH1* ob, TString xtit, TString ytit);
116  void DefineTitle(TGraph* ob, TString xtit, TString ytit);
117  void DefineTitle(TF1* ob, TString xtit, TString ytit);
118 
119  Double_t GetX(TH1* ob, Double_t val, Double_t eps = 1.e-07, Int_t nmax = 50, Double_t xmin = -1.0, Double_t xmax = -1.0);
120  Double_t GetXWithLimits(TH1* ob, Double_t val, Double_t xmin = -1.0, Double_t xmax = -1.0, Double_t eps = 1.e-07, Int_t nmax = 50)
121  {
123  return GetX(ob, val, eps, nmax, xmin, xmax);
124  }
125  TF1* RescaleX(TH1* hist1, TH1* hist2, Int_t degree, Double_t* params,
126  Int_t npoints = -1, const Char_t* direction = "C",
127  Double_t xmin = -1, Double_t xmax = -1, Double_t qmin = 0.05, Double_t qmax = 0.95,
128  Double_t eps = 1.e-07);
129  void RescaleX(TH1* hist1, TH1* hist2, TF1* scale_func, Int_t npoints = 2,
130  const Char_t* direction = "C", Double_t xmin = -1, Double_t xmax = -1, Double_t qmin = 0.05, Double_t qmax = 0.95,
131  Double_t eps = 1.e-07);
132  TH1* MakeHistoRescaleX(TH1* hist1, TH1* hist2, Int_t degree, Double_t* params,
133  Option_t* opt = "", Int_t npoints = -1, const Char_t* direction = "C",
134  Double_t xmin = -1, Double_t xmax = -1, Double_t qmin = 0.05, Double_t qmax = 0.95,
135  Double_t eps = 1.e-07);
136  TH1* MakeHistoRescaleX(TH1* hist1, TH1* hist2, TF1* scale_func, Int_t npoints = 2,
137  Option_t* opt = "", const Char_t* direction = "C",
138  Double_t xmin = -1, Double_t xmax = -1, Double_t qmin = 0.05, Double_t qmax = 0.95,
139  Double_t eps = 1.e-07);
140 
141  Double_t GetChisquare(TH1* h1, TF1* f1, Bool_t norm = kTRUE, Bool_t err = kTRUE, Double_t* para = nullptr);
142  Double_t GetLikelihood(TH1* h1, TF1* f1, Bool_t norm = kTRUE, Double_t* para = nullptr);
143 
144  TGraph* DivideGraphs(TGraph* G1, TGraph* G2);
145  TGraph* ComputeNewGraphFrom(TGraph* g0, TGraph* g1, const TString& formula);
146  TGraph* ComputeNewGraphFrom(TList* lgr, TString formula);
147 
148  std::vector<Double_t> GetLimits(TGraph* G1);
149  std::vector<Double_t> GetLimits(TProfile* G1);
150  std::vector<Double_t> GetLimits(TMultiGraph* mgr);
151  std::vector<Double_t> GetLimits(TSeqCollection* mgr);
152  void ApplyCurrentLimitsToAllCanvas(Bool_t AlsoLog = kFALSE);
154 
155  ClassDef(KVHistoManipulator, 1) //Propose differentes operations sur les histo
156 };
157 
160 
161 #endif
int Int_t
#define R__EXTERN
R__EXTERN KVHistoManipulator * gHistoManipulator
................ global variable
char Char_t
bool Bool_t
double Double_t
const char Option_t
#define ClassDef(name, id)
float xmin
float xmax
Toolkit for various operations on histograms & graphs not provided by ROOT.
TF1 * RescaleX(TH1 *hist1, TH1 *hist2, Int_t degree, Double_t *params, Int_t npoints=-1, const Char_t *direction="C", Double_t xmin=-1, Double_t xmax=-1, Double_t qmin=0.05, Double_t qmax=0.95, Double_t eps=1.e-07)
KVHistoManipulator()
Default constructor.
void DefineMarkerStyle(TAttMarker *ob, TString marker)
TH1 * MakeHistoRescaleX(TH1 *hist1, TH1 *hist2, Int_t degree, Double_t *params, Option_t *opt="", Int_t npoints=-1, const Char_t *direction="C", Double_t xmin=-1, Double_t xmax=-1, Double_t qmin=0.05, Double_t qmax=0.95, Double_t eps=1.e-07)
Bool_t kVisDebug
= kTRUE for visual debugging
TH2 * CentreeReduiteX(TH2 *hh, Int_t nx=-1, Double_t xmin=-1., Double_t xmax=-1.)
std::vector< Double_t > GetLimits(TGraph *G1)
TGraph * ComputeNewGraphFrom(TGraph *g0, TGraph *g1, const TString &formula)
Double_t GetXWithLimits(TH1 *ob, Double_t val, Double_t xmin=-1.0, Double_t xmax=-1.0, Double_t eps=1.e-07, Int_t nmax=50)
TH2 * CentreeReduiteY(TH2 *hh, Int_t ny=-1, Double_t ymin=-1., Double_t ymax=-1.)
void DefineStyle(TObject *ob, TString line, TString marker)
Double_t GetX(TH1 *ob, Double_t val, Double_t eps=1.e-07, Int_t nmax=50, Double_t xmin=-1.0, Double_t xmax=-1.0)
Double_t GetChisquare(TH1 *h1, TF1 *f1, Bool_t norm=kTRUE, Bool_t err=kTRUE, Double_t *para=nullptr)
KVNumberList * Saucisson(TH1 *hh, Int_t ntranches=10)
TGraph * ScaleGraph(const TGraph *hh, TF1 *fx=nullptr, TF1 *fy=nullptr) const
TGraph * LinkGraphs(TGraph *grx, TGraph *gry)
void DefinePattern(TH1 *ob, TString titleX="42 0.08 0.8", TString titleY="42 0.07 1.2", TString labelX="42 0.05 0.005", TString labelY="42 0.05 0.006")
virtual ~KVHistoManipulator(void)
void SetVisDebug(Bool_t on=kTRUE)
TH2 * RenormaliseHisto(TH2 *hh, Int_t bmin=-1, Int_t bmax=-1, TString axis="X", Double_t valref=1)
Int_t CutStatBin(TH1 *hh, Int_t stat_min=-1, Int_t stat_max=-1)
TGraphErrors * GetMomentEvolution(TH2 *hh, TString momentx, TString momenty, TString axis="Y", Double_t stat_min=0)
void DefineTitle(TH1 *ob, TString xtit, TString ytit)
TGraph * ExtractMeanAndSigmaFromProfile(TProfile *pf, TGraph *&sigma)
void ApplyCurrentLimitsToAllCanvas(Bool_t AlsoLog=kFALSE)
TH1 * CentreeReduite(TH1 *hh, Int_t nx=-1, Int_t ny=-1, Double_t xmin=-1., Double_t xmax=-1., Double_t ymin=-1., Double_t ymax=-1.)
TCanvas * fVDCanvas
used for visual debugging
TH1 * ScaleHisto(TH1 *hh, TF1 *fx, TF1 *fy=NULL, Int_t nx=-1, Int_t ny=-1, Double_t xmin=-1., Double_t xmax=-1., Double_t ymin=-1., Double_t ymax=-1., Option_t *norm="")
TGraph * DivideGraphs(TGraph *G1, TGraph *G2)
Double_t GetLikelihood(TH1 *h1, TF1 *f1, Bool_t norm=kTRUE, Double_t *para=nullptr)
KVList * Give_ProjectionList(TH2 *hh, Double_t MinIntegral=-1, TString axis="X")
TGraphErrors * MakeGraphFrom(TProfile *pf, Bool_t Error=kTRUE)
TH1 * CumulatedHisto(TH1 *hh, TString direction="C", Int_t bmin=-1, Int_t bmax=-1, Option_t *norm="surf")
TH1 * GetDerivative(TH1 *hh, Int_t order)
void DefineLineStyle(TAttLine *ob, TString line)
Double_t GetCorrelationFactor(TH2 *hh)
Bool_t IsVisDebug() const
Int_t Apply_TCutG(TH2 *hh, TCutG *cut, TString mode="in")
Extended TList class which owns its objects by default.
Definition: KVList.h:27
Strings used to represent a set of ranges of values.
Definition: KVNumberList.h:83
TString & Append(char c, Ssiz_t rep=1)