KaliVeda
Toolkit for HIC analysis
Loading...
Searching...
No Matches
KVDalitzPlot.cpp
1//Created by KVClassFactory on Tue Jan 17 15:10:46 2012
2//Author: bonnet
3
4#include "KVDalitzPlot.h"
5#include "TMath.h"
6#include "TPad.h"
7#include "TCanvas.h"
8#include "KVCanvas.h"
9
11
12
13
14
15KVDalitzPlot::KVDalitzPlot(const char* name, const char* title, Bool_t ordered, Int_t nbinsx, Double_t xlow, Double_t xup, Int_t nbinsy, Double_t ylow, Double_t yup) :
16 TH2F(name, title, nbinsx, xlow, xup, nbinsy, ylow, yup)
17{
18 fOrdered = ordered;
19 fShowBorder = 1;
20 fShowCenter = 1;
21 lb1 = lb2 = lb3 = 0;
22 lc1 = lc2 = lc3 = 0;
23}
24
25
26
27
30
32{
33 // Default constructor
34}
35
36
37
38
45
47{
48 // Copy constructor
49 // This ctor is used to make a copy of an existing object (for example
50 // when a method returns an object), and it is always a good idea to
51 // implement it.
52 // If your class allocates memory in its constructor(s) then it is ESSENTIAL :-)
53
54 obj.Copy(*this);
55}
56
57
58
61
63{
64 // Destructor
65}
66
67
68
69
77
79{
80 // This method copies the current state of 'this' object into 'obj'
81 // You should add here any member variables, for example:
82 // (supposing a member variable KVDalitzPlot::fToto)
83 // CastedObj.fToto = fToto;
84 // or
85 // CastedObj.SetToto( GetToto() );
86
87 TH2F::Copy(obj);
88 //KVDalitzPlot& CastedObj = (KVDalitzPlot&)obj;
89}
90
91
92
97
99{
100 //Place the 3 observable in an equilateral triangle
101 //after normalized the sum to 1
102 //Info : the permutation is not done
103
104 if (!fOrdered) return FillMe(a1, a2, a3);
105
106 Int_t result = 0;
107 result += FillMe(a1, a2, a3);
108 result += FillMe(a2, a3, a1);
109 result += FillMe(a3, a1, a2);
110 result += FillMe(a1, a3, a2);
111 result += FillMe(a2, a1, a3);
112 result += FillMe(a3, a2, a1);
113 return result;
114}
115
116
117
119
121{
122 Double_t sum = a1 + a2 + a3;
123 if (sum > 0) {
124 Double_t a1n = a1 / sum;
125 Double_t a2n = a2 / sum;
126
127 Double_t xI = 1. / TMath::Sqrt(3) * (a1n + 2 * a2n);
128 Double_t yI = a1n;
129
130 return TH2F::Fill(xI, yI);
131 }
132 return -1;
133}
134
135
136
138
140{
141 GetXaxis()->SetNdivisions(506);
142 GetXaxis()->SetLabelFont(42);
143 GetXaxis()->SetLabelSize(0.0);
144 GetXaxis()->SetTitleSize(0.035);
146 GetXaxis()->SetTitleFont(42);
149
150 GetYaxis()->SetNdivisions(506);
151 GetYaxis()->SetLabelFont(42);
152 GetYaxis()->SetLabelSize(0.0);
153 GetYaxis()->SetTitleSize(0.035);
155 GetYaxis()->SetTitleFont(42);
158
159 TH2F::Draw("col");
160
161 SetShowBorder(1);
162 SetShowCenter(1);
163
164 gPad->GetCanvas()->SetRightMargin(0.001);
165 gPad->GetCanvas()->SetTopMargin(0.001);
166
167 TCanvas* cc = gPad->GetCanvas();
168 if (cc->InheritsFrom("KVCanvas"))((KVCanvas*)cc)->DisableClass("TLine");
169
170 // gPad->GetCanvas()->SetLeftMargin(0.001);
171 // gPad->GetCanvas()->SetBottomMargin(0.001);
172
173 gPad->SetBorderMode(0);
174 gPad->SetBorderSize(2);
175 gPad->SetLogz();
176 gPad->SetFrameBorderMode(0);
177 gPad->SetFrameLineColor(0);
178 gPad->Modified();
179 gPad->Update();
180}
181
182
183
185
187{
189
193
194 if (fShowBorder) {
195 lb1 = new TLine(0.5758319, 1., 0., 0.);
196 lb1->SetLineWidth(2);
197 lb1->Draw("same");
198 lb2 = new TLine(0.5758319, 1., 1.15, 0.);
199 lb2->SetLineWidth(2);
200 lb2->Draw("same");
201 lb3 = new TLine(0., 0., 1.15, 0.);
202 lb3->SetLineWidth(2);
203 lb3->Draw("same");
204 }
205
206 gPad->SetFrameLineColor(0);
207 gPad->Modified();
208 gPad->Update();
209}
210
211
212
214
216{
218
222
223 if (fShowCenter) {
224 lc1 = new TLine(0.575, 1, 0.575, 0.333);
225 lc1->SetLineWidth(2);
226 lc1->SetLineStyle(9);
227 lc1->Draw();
228 lc2 = new TLine(0, 0, 0.575, 0.333);
229 lc2->SetLineWidth(2);
230 lc2->SetLineStyle(9);
231 lc2->Draw();
232 lc3 = new TLine(0.575, 0.333, 1.15, 0);
233 lc3->SetLineWidth(2);
234 lc3->SetLineStyle(9);
235 lc3->Draw();
236 }
237 gPad->SetFrameLineColor(0);
238 gPad->Modified();
239 gPad->Update();
240
241}
242
243
244
250
252{
253 //Linearization of the Dalitz representation
254 //histogram is filled according to the distance from each cell
255 //to the center
256 //
257 TH1F* h1 = new TH1F(Form("distance_%s", GetName()), Form("distance_%s", GetName()), 150, 0, 1.5);
258 for (Int_t nx = 1; nx <= GetNbinsX(); nx += 1) {
259 Double_t xx = GetXaxis()->GetBinCenter(nx);
260 for (Int_t ny = 1; ny <= GetNbinsY(); ny += 1) {
261 Double_t yy = GetYaxis()->GetBinCenter(ny);
262 Double_t dist = TMath::Sqrt(TMath::Power(xx - 0.5, 2.) + TMath::Power(yy - 1. / 3., 2.));
263 h1->Fill(dist, GetBinContent(nx, ny));
264 }
265 }
266
267 return h1;
268
269}
270
271
272
273
274
int Int_t
#define SafeDelete(p)
bool Bool_t
double Double_t
const char Option_t
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
char * Form(const char *fmt,...)
#define gPad
TCanvas with mouse-controlled dynamic zoom and pan & scan.
Definition KVCanvas.h:54
Fill 3D observables in a dalitz plot ,.
void Draw(Option_t *opt="")
Bool_t fOrdered
Int_t fShowCenter
void SetShowBorder(Int_t value=1)
Int_t FillMe(Double_t a1, Double_t a2, Double_t a3)
void SetShowCenter(Int_t value=1)
KVDalitzPlot()
Default constructor.
void Copy(TObject &) const
Int_t fShowBorder
virtual ~KVDalitzPlot()
Destructor.
Int_t FillAsDalitz(Double_t a1, Double_t a2, Double_t a3)
TH1 * GetDistanceFromCenter()
virtual void SetAxisColor(Color_t color=1, Float_t alpha=1.)
virtual void SetLabelSize(Float_t size=0.04)
virtual void SetTitleFont(Style_t font=62)
virtual void SetLabelFont(Style_t font=62)
virtual void SetTitleSize(Float_t size=0.04)
virtual void SetNdivisions(Int_t n1, Int_t n2, Int_t n3, Bool_t optim=kTRUE)
virtual void SetTickLength(Float_t length=0.03)
virtual void SetLabelColor(Color_t color=1, Float_t alpha=1.)
virtual void SetLineStyle(Style_t lstyle)
virtual void SetLineWidth(Width_t lwidth)
virtual Double_t GetBinCenter(Int_t bin) const
virtual Int_t GetNbinsY() const
TAxis * GetXaxis()
virtual Int_t GetNbinsX() const
TAxis * GetYaxis()
void Draw(Option_t *option="") override
virtual Int_t Fill(const char *name, Double_t w)
void Copy(TObject &hnew) const override
virtual Double_t GetBinContent(Int_t bin) const
const char * GetName() const override
virtual void Draw(Option_t *option="")
TH1F * h1
double dist(AxisAngle const &r1, AxisAngle const &r2)
Double_t Power(Double_t x, Double_t y)
Double_t Sqrt(Double_t x)
ClassImp(TPyArg)