CbmRoot
Loading...
Searching...
No Matches
CbmRichRingFitterQa.cxx
Go to the documentation of this file.
1/* Copyright (C) 2009-2016 UGiessen/JINR-LIT, Giessen/Dubna
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Semen Lebedev [committer] */
4
12#include "CbmRichRingFitterQa.h"
13
14#include "CbmRichRing.h"
18#include "TCanvas.h"
19#include "TH1D.h"
20#include "TMath.h"
21#include "TMatrixD.h"
22#include "TRandom.h"
23
24#include <iostream>
25#include <vector>
26
27using std::cout;
28using std::endl;
29using std::vector;
30
32 : fhErrorA(NULL)
33 , fhErrorB(NULL)
34 , fhErrorX(NULL)
35 , fhErrorY(NULL)
36 , fhErrorPhi(NULL)
37 ,
38
39 fhA(NULL)
40 , fhB(NULL)
41 , fhX(NULL)
42 , fhY(NULL)
43 , fhPhi(NULL)
44 ,
45
46 fhRadiusErr(NULL)
47 , fhCircleXcErr(NULL)
48 , fhCircleYcErr(NULL)
49 ,
50
51 fhRadius(NULL)
52 , fhCircleXc(NULL)
53 , fhCircleYc(NULL)
54 ,
55
56 fhRadiusPool(NULL)
57 , fhCircleXcPool(NULL)
58 , fhCircleYcPool(NULL)
59{
60 fhErrorA = new TH1D("fhErrorA", "fhErrorA;dA [cm];Counter", 100, -2., 2.);
61 fhErrorB = new TH1D("fhErrorB", "fhErrorB;B [cm];Counter", 100, -2., 2.);
62 fhErrorX = new TH1D("fhErrorX", "fhErrorX;X [cm];Counter", 100, -2., 2.);
63 fhErrorY = new TH1D("fhErrorY", "fhErrorY;Y [cm];Counter", 100, -2., 2.);
64 fhErrorPhi = new TH1D("fhErrorPhi", "fhErrorPhi;d#Phi [rad];Counter", 100, -2., 2.);
65
66 fhA = new TH1D("fhA", "fhA;A [cm];Counter", 100, 5., 7.);
67 fhB = new TH1D("fhB", "fhB;B [cm];Counter", 100, 5., 7.);
68 fhX = new TH1D("fhX", "fhX;X [cm];Counter", 100, -1., 1.);
69 fhY = new TH1D("fhY", "fhY;Y [cm];Counter", 100, -1., 1.);
70 Double_t pi = TMath::Pi();
71 fhPhi = new TH1D("fhPhi", "fhPhi;#Phi [rad];Counter", 100, -pi / 2. - pi / 6., pi / 2. + pi / 6.);
72
73 // circle fitting
74 fhRadiusErr = new TH1D("fhRadiusErr", "fhRadiusErr;dR [cm];Counter", 100, -2., 2.);
75 fhCircleXcErr = new TH1D("fhCircleXcErr", "fhCircleXcErr;dXc [cm];Counter", 100, -2., 2.);
76 fhCircleYcErr = new TH1D("fhCircleYcErr", "fhCircleYcErr;dYc [cm];Counter", 100, -2., 2.);
77 fhRadius = new TH1D("fhRadius", "fhRadius;Radius [cm];Counter", 100, 4., 8.);
78 fhCircleXc = new TH1D("fhCircleXc", "fhCircleXc;Xc [cm];Counter", 100, -2., 2.);
79 fhCircleYc = new TH1D("fhCircleYc", "fhCircleYc;Yc [cm];Counter", 100, -2., 2.);
80 fhRadiusPool = new TH1D("fhRadiusPool", "fhRadiusPool;Pool R;Counter", 100, -5., 5.);
81 fhCircleXcPool = new TH1D("fhCircleXcPool", "fhCircleXcPool;Pool Xc;Counter", 100, -5., 5.);
82 fhCircleYcPool = new TH1D("fhCircleYcPool", "fhCircleYcPool;Pool Yc;Counter", 100, -5., 5.);
83}
84
86
88{
89 // Double_t maxX = 15;
90 // Double_t maxY = 15;
91 Int_t nofHits = 50;
92 Double_t A = 6.;
93 Double_t B = 6.;
94 Double_t sigmaError = 0.2;
95 CbmRichRingLight ellipse;
96 Int_t nofBadFit = 0;
97 Double_t X0 = 0.; //gRandom->Rndm()*(maxX - A);
98 Double_t Y0 = 0.; //gRandom->Rndm()* (maxY - A);
99
102
103 for (Int_t iR = 0; iR < 50000; iR++) {
104 Double_t phi = 0.; //TMath::Pi()*(6./12.); //gRandom->Rndm()*TMath::Pi() - TMath::Pi()/2.;
105 ellipse.SetXYABP(X0, Y0, A, B, phi);
106 for (Int_t iH = 0; iH < nofHits; iH++) {
107 Double_t alfa = gRandom->Rndm() * 2. * TMath::Pi();
108
109 Double_t errorX = gRandom->Gaus(0, sigmaError);
110 Double_t errorY = gRandom->Gaus(0, sigmaError);
111
112 Double_t hx = A * cos(alfa);
113 Double_t hy = B * sin(alfa);
114
115 Double_t hitXRot = hx * cos(phi) - hy * sin(phi);
116 Double_t hitYRot = hx * sin(phi) + hy * cos(phi);
117
118 CbmRichHitLight hit(hitXRot + X0 + errorX, hitYRot + Y0 + errorY);
119 ellipse.AddHit(hit);
120 }
121 // ellipse fit
122 fitEllipse->DoFit(&ellipse);
123 fhErrorA->Fill(A - ellipse.GetAaxis());
124 fhErrorB->Fill(B - ellipse.GetBaxis());
125 fhErrorX->Fill(X0 - ellipse.GetCenterX());
126 fhErrorY->Fill(Y0 - ellipse.GetCenterY());
127 fhErrorPhi->Fill(phi - ellipse.GetPhi());
128 fhA->Fill(ellipse.GetAaxis());
129 fhB->Fill(ellipse.GetBaxis());
130 fhX->Fill(ellipse.GetCenterX());
131 fhY->Fill(ellipse.GetCenterY());
132 fhPhi->Fill(ellipse.GetPhi());
133
134 // circle fit
135 fitCircle->DoFit(&ellipse);
136 TMatrixD cov(3, 3);
137 CalculateFitErrors(&ellipse, sigmaError, cov);
138 Double_t mcR = (A + B) / 2.;
139 fhRadiusErr->Fill(mcR - ellipse.GetRadius());
140 fhCircleXcErr->Fill(X0 - ellipse.GetCenterX());
141 fhCircleYcErr->Fill(Y0 - ellipse.GetCenterY());
142 fhRadius->Fill(ellipse.GetRadius());
143 fhCircleXc->Fill(ellipse.GetCenterX());
144 fhCircleYc->Fill(ellipse.GetCenterY());
145 fhRadiusPool->Fill((mcR - ellipse.GetRadius()) / sqrt(cov(2, 2)));
146 fhCircleXcPool->Fill((X0 - ellipse.GetCenterX()) / sqrt(cov(0, 0)));
147 fhCircleYcPool->Fill((Y0 - ellipse.GetCenterY()) / sqrt(cov(1, 1)));
148 } // iR
149
150 Draw();
151 cout << nofBadFit << endl;
152}
153
155{
156 TCanvas* c = new TCanvas("rich_fitter_errors", "rich_fitter_errors", 900, 600);
157 c->Divide(3, 2);
158 c->cd(1);
159 fhErrorA->Draw();
160 c->cd(2);
161 fhErrorB->Draw();
162 c->cd(3);
163 fhErrorX->Draw();
164 c->cd(4);
165 fhErrorY->Draw();
166 c->cd(5);
167 fhErrorPhi->Draw();
168 cout.precision(4);
169 cout << fhErrorA->GetMean() << " " << fhErrorA->GetRMS() << endl;
170 cout << fhErrorB->GetMean() << " " << fhErrorB->GetRMS() << endl;
171 cout << fhErrorX->GetMean() << " " << fhErrorX->GetRMS() << endl;
172 cout << fhErrorY->GetMean() << " " << fhErrorY->GetRMS() << endl;
173 cout << fhErrorPhi->GetMean() << " " << fhErrorPhi->GetRMS() << endl;
174
175 TCanvas* c2 = new TCanvas("rich_fitter_params", "rich_fitter_params", 900, 600);
176 c2->Divide(3, 2);
177 c2->cd(1);
178 fhA->Draw();
179 c2->cd(2);
180 fhB->Draw();
181 c2->cd(3);
182 fhX->Draw();
183 c2->cd(4);
184 fhY->Draw();
185 c2->cd(5);
186 fhPhi->Draw();
187
188 TCanvas* c3 = new TCanvas("rich_fitter_circle", "rich_fitter_circle", 900, 900);
189 c3->Divide(3, 3);
190 c3->cd(1);
191 fhRadiusErr->Draw();
192 c3->cd(2);
193 fhCircleXcErr->Draw();
194 c3->cd(3);
195 fhCircleYcErr->Draw();
196 c3->cd(4);
197 fhRadius->Draw();
198 c3->cd(5);
199 fhCircleXc->Draw();
200 c3->cd(6);
201 fhCircleYc->Draw();
202 c3->cd(7);
203 fhRadiusPool->Draw();
204 c3->cd(8);
205 fhCircleXcPool->Draw();
206 c3->cd(9);
207 fhCircleYcPool->Draw();
208}
209
210void CbmRichRingFitterQa::CalculateFitErrors(CbmRichRingLight* ring, Double_t sigma, TMatrixD& cov)
211{
212 //TMatrixD H3(3,3);
213 TMatrixD HY3(3, 1);
214 //TMatrixD Cov3(3,3);
215 TMatrixD HC3(3, 1);
216
217 for (Int_t i = 0; i < 3; i++) {
218 HY3(i, 0) = 0;
219 HC3(i, 0) = 0;
220 for (Int_t j = 0; j < 3; j++) {
221 //H3(i,j) = 0;
222 cov(i, j) = 0;
223 }
224 }
225
226 Double_t xc = ring->GetCenterX();
227 Double_t yc = ring->GetCenterY();
228 Double_t R = ring->GetRadius();
229 for (Int_t iHit = 0; iHit < ring->GetNofHits(); iHit++) {
230 Double_t xi = ring->GetHit(iHit).fX;
231 Double_t yi = ring->GetHit(iHit).fY;
232 Double_t ri = sqrt((xi - xc) * (xi - xc) + (yi - yc) * (yi - yc));
233 Double_t err = sigma;
234
235 Double_t f1 = (-1.0 * (xi - xc)) / (ri * err);
236 Double_t f2 = (-1.0 * (yi - yc)) / (ri * err);
237 Double_t f3 = (-1.) / err;
238 Double_t Y = (R - ri) / err;
239
240 cov(0, 0) = cov(0, 0) + f1 * f1;
241 cov(0, 1) = cov(0, 1) + f1 * f2;
242 cov(0, 2) = cov(0, 2) + f1 * f3;
243
244 cov(1, 0) = cov(0, 1);
245 cov(1, 1) = cov(1, 1) + f2 * f2;
246 cov(1, 2) = cov(1, 2) + f2 * f3;
247
248 cov(2, 0) = cov(0, 2);
249 cov(2, 1) = cov(1, 2);
250 cov(2, 2) = cov(2, 2) + f3 * f3;
251
252 HY3(0, 0) = HY3(0, 0) + Y * f1;
253 HY3(1, 0) = HY3(1, 0) + Y * f2;
254 HY3(2, 0) = HY3(2, 0) + Y * f3;
255 } // iHit
256
257 //H3.Print();
258 Double_t det = 0.0;
259 cov.Invert(&det);
260 //Cov3 = H3;
261 //H3.Print();
262 //Cov3.Print();
263
264 //HC3 = H3 * HY3;
265 //HC3.Print();
266
267 //cout << "dX0= " << HC3(0,0) << " +- " << sqrt(Cov3(0,0)) << endl;
268 // cout << "dY0= " << HC3(1,0) << " +- " << sqrt(Cov3(1,1)) << endl;
269 //cout << "dR= " << HC3(2,0) << " +- " << sqrt(Cov3(2,2)) << endl;
270
271 //cout << "dX0= " << HC3(0,0) << " +- " << sqrt(Cov3(0,0)) << endl;
272 // cout << "dY0= " << HC3(1,0) << " +- " << sqrt(Cov3(1,1)) << endl;
273 // cout << "dR= " << HC3(2,0) << " +- " << sqrt(Cov3(2,2)) << endl;
274}
275
Here the ring is fitted with the COP algorithm from A. Ayriyan/G. Ososkov.
This is the implementation of ellipse fitting using MINUIT.
Here the ring is fitted with Taubin algorithm from A. Ayriyan, G. Ososkov, N. Chernov.
ClassImp(CbmRichRingFitterQa)
Test ellipse and circle fitting on toy model.
friend fvec sqrt(const fvec &a)
friend fvec cos(const fvec &a)
friend fvec sin(const fvec &a)
const double pi
Here the ring is fitted with the COP algorithm from A. Ayriyan/G. Ososkov.
virtual void DoFit(CbmRichRingLight *ring)
Inherited from CbmRichRingFitterBase.
Here the ring is fitted with Taubin algorithm from A. Ayriyan, G. Ososkov, N. Chernov.
virtual void DoFit(CbmRichRingLight *ring)
Inherited from CbmRichRingFitterBase.
Test ellipse and circle fitting on toy model.
CbmRichRingFitterQa()
Standard constructor.
void GenerateEllipse()
Generate ellipse.
void Draw(Option_t *="")
Draw generated and fitted circle/ellipse.
virtual ~CbmRichRingFitterQa()
Destructor.
void CalculateFitErrors(CbmRichRingLight *ring, Double_t sigma, TMatrixD &cov)
void SetXYABP(float x, float y, float a, float b, float p)
Set all 5 ellipse parameters.
float GetCenterX() const
float GetAaxis() const
float GetRadius() const
float GetCenterY() const
float GetPhi() const
int GetNofHits() const
Return number of hits in ring.
float GetBaxis() const
CbmRichHitLight GetHit(int ind)
Return hit by the index.
void AddHit(CbmRichHitLight hit)
Add new hit to the ring.