]>
Commit | Line | Data |
---|---|---|
4c039060 | 1 | /************************************************************************** |
2 | * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * | |
3 | * * | |
4 | * Author: The ALICE Off-line Project. * | |
5 | * Contributors are mentioned in the code where appropriate. * | |
6 | * * | |
7 | * Permission to use, copy, modify and distribute this software and its * | |
8 | * documentation strictly for non-commercial purposes is hereby granted * | |
9 | * without fee, provided that the above copyright notice appears in all * | |
10 | * copies and that both the copyright notice and this permission notice * | |
11 | * appear in the supporting documentation. The authors make no claims * | |
12 | * about the suitability of this software for any purpose. It is * | |
13 | * provided "as is" without express or implied warranty. * | |
14 | **************************************************************************/ | |
15 | ||
16 | /* | |
17 | $Log$ | |
6e7b5431 | 18 | Revision 1.5 2000/06/30 12:07:50 kowal2 |
19 | Updated from the TPC-PreRelease branch | |
20 | ||
73042f01 | 21 | Revision 1.4.4.3 2000/06/26 07:39:42 kowal2 |
22 | Changes to obey the coding rules | |
23 | ||
24 | Revision 1.4.4.2 2000/06/25 08:38:41 kowal2 | |
25 | Splitted from AliTPCtracking | |
26 | ||
27 | Revision 1.4.4.1 2000/06/14 16:48:24 kowal2 | |
28 | Parameter setting improved. Removed compiler warnings | |
29 | ||
30 | Revision 1.4 2000/04/17 09:37:33 kowal2 | |
31 | removed obsolete AliTPCDigitsDisplay.C | |
32 | ||
cc80f89e | 33 | Revision 1.3.8.2 2000/04/10 08:40:46 kowal2 |
34 | ||
35 | Small changes by M. Ivanov, improvements of algorithms | |
36 | ||
37 | Revision 1.3.8.1 2000/04/10 07:56:53 kowal2 | |
38 | Not used anymore - removed | |
39 | ||
40 | Revision 1.3 1999/10/05 17:15:46 fca | |
41 | Minor syntax for the Alpha OSF | |
42 | ||
0b34885d | 43 | Revision 1.2 1999/09/29 09:24:34 fca |
44 | Introduction of the Copyright and cvs Log | |
45 | ||
4c039060 | 46 | */ |
47 | ||
8c555625 | 48 | /////////////////////////////////////////////////////////////////////////////// |
6e7b5431 | 49 | // AliTPCPRF2D - // |
8c555625 | 50 | // Pad response function object in two dimesions // |
51 | // This class contains the basic functions for the // | |
52 | // calculation of PRF according generic charge distribution // | |
53 | // In Update function object calculate table of response function // | |
54 | // in discrete x and y position // | |
55 | // This table is used for interpolation od response function in any position // | |
56 | // (function GetPRF) // | |
57 | // // | |
58 | // Origin: Marian Ivanov, Uni. of Bratislava, ivanov@fmph.uniba.sk // | |
59 | // // | |
60 | /////////////////////////////////////////////////////////////////////////////// | |
cc80f89e | 61 | |
62 | ||
8c555625 | 63 | #include "TMath.h" |
64 | #include "AliTPCPRF2D.h" | |
65 | #include "TF2.h" | |
66 | #include <iostream.h> | |
67 | #include <string.h> | |
68 | #include "TCanvas.h" | |
69 | #include "TPad.h" | |
70 | #include "TStyle.h" | |
6e7b5431 | 71 | #include "TH1.h" |
72 | #include "AliH2F.h" | |
73 | ||
74 | ||
8c555625 | 75 | #include "TPaveText.h" |
76 | #include "TText.h" | |
77 | ||
78 | extern TStyle * gStyle; | |
79 | ||
6e7b5431 | 80 | const Double_t AliTPCPRF2D::fgkDegtoRad = 0.01745329251994; |
81 | const Double_t AliTPCPRF2D::fgkSQRT12=3.464101; | |
82 | const Int_t AliTPCPRF2D::fgkNPRF = 100; | |
8c555625 | 83 | |
84 | ||
85 | static Double_t funGauss2D(Double_t *x, Double_t * par) | |
86 | { | |
cc80f89e | 87 | //Gauss function -needde by the generic function object |
8c555625 | 88 | return ( TMath::Exp(-(x[0]*x[0])/(2*par[0]*par[0]))* |
89 | TMath::Exp(-(x[1]*x[1])/(2*par[1]*par[1]))); | |
90 | ||
91 | } | |
92 | ||
93 | static Double_t funCosh2D(Double_t *x, Double_t * par) | |
94 | { | |
cc80f89e | 95 | //Cosh function -needde by the generic function object |
8c555625 | 96 | return ( 1/(TMath::CosH(3.14159*x[0]/(2*par[0]))* |
97 | TMath::CosH(3.14159*x[1]/(2*par[1])))); | |
98 | } | |
99 | ||
100 | static Double_t funGati2D(Double_t *x, Double_t * par) | |
101 | { | |
cc80f89e | 102 | //Gati function -needde by the generic function object |
73042f01 | 103 | Float_t k3=par[1]; |
104 | Float_t k3R=TMath::Sqrt(k3); | |
105 | Float_t k2=(TMath::Pi()/2)*(1-k3R/2.); | |
106 | Float_t k1=k2*k3R/(4*TMath::ATan(k3R)); | |
8c555625 | 107 | Float_t l=x[0]/par[0]; |
73042f01 | 108 | Float_t tan2=TMath::TanH(k2*l); |
8c555625 | 109 | tan2*=tan2; |
73042f01 | 110 | Float_t res = k1*(1-tan2)/(1+k3*tan2); |
8c555625 | 111 | //par[4] = is equal to k3Y |
73042f01 | 112 | k3=par[4]; |
113 | k3R=TMath::Sqrt(k3); | |
114 | k2=(TMath::Pi()/2)*(1-k3R/2.); | |
115 | k1=k2*k3R/(4*TMath::ATan(k3R)); | |
8c555625 | 116 | l=x[1]/par[0]; |
6e7b5431 | 117 | tan2=TMath::TanH(k2*l); |
8c555625 | 118 | tan2*=tan2; |
6e7b5431 | 119 | res = res*k1*(1-tan2)/(1+k3*tan2); |
8c555625 | 120 | return res; |
121 | } | |
122 | ||
8c555625 | 123 | /////////////////////////////////////////////////////////////////////////// |
124 | /////////////////////////////////////////////////////////////////////////// | |
125 | ||
126 | ClassImp(AliTPCPRF2D) | |
127 | ||
128 | AliTPCPRF2D::AliTPCPRF2D() | |
129 | { | |
cc80f89e | 130 | //default constructor for response function object |
6e7b5431 | 131 | fNChargeArray = 0; |
132 | fChargeArray = 0; | |
133 | fNPRF =fgkNPRF ; | |
8c555625 | 134 | fSigmaX = 0; |
cc80f89e | 135 | fSigmaY = 0; |
8c555625 | 136 | |
137 | fGRF = 0; | |
6e7b5431 | 138 | fKNorm = 1; |
cc80f89e | 139 | fOrigSigmaY=0; |
140 | fOrigSigmaX=0; | |
8c555625 | 141 | fNdiv = 5; |
cc80f89e | 142 | //set daault angels |
143 | fChargeAngle = 0; | |
6e7b5431 | 144 | fPadAngle = 0; |
8c555625 | 145 | //chewron default values |
146 | SetPad(0.8,0.8); | |
147 | SetChevron(0.2,0.0,1.0); | |
148 | SetY(-0.2,0.2,2); | |
6e7b5431 | 149 | SetInterpolationType(2,0); |
8c555625 | 150 | } |
151 | ||
152 | AliTPCPRF2D::~AliTPCPRF2D() | |
153 | { | |
6e7b5431 | 154 | if (fChargeArray!=0) delete [] fChargeArray; |
155 | if (fGRF !=0 ) fGRF->Delete(); | |
8c555625 | 156 | } |
157 | ||
158 | void AliTPCPRF2D::SetY(Float_t y1, Float_t y2, Int_t nYdiv) | |
159 | { | |
160 | // | |
161 | //set virtual line position | |
162 | //first and last line and number of lines | |
163 | fNYdiv = nYdiv; | |
8c555625 | 164 | fY1=y1; |
165 | fY2=y2; | |
166 | } | |
167 | ||
168 | void AliTPCPRF2D::SetPad(Float_t width, Float_t height) | |
169 | { | |
170 | //set base chevron parameters | |
171 | fHeightFull=height; | |
172 | fWidth=width; | |
173 | } | |
174 | void AliTPCPRF2D::SetChevron(Float_t hstep, | |
175 | Float_t shifty, | |
176 | Float_t fac) | |
177 | { | |
178 | //set shaping of chewron parameters | |
179 | fHeightS=hstep; | |
180 | fShiftY=shifty; | |
6e7b5431 | 181 | fK=fac; |
8c555625 | 182 | } |
183 | ||
184 | void AliTPCPRF2D::SetChParam(Float_t width, Float_t height, | |
185 | Float_t hstep, Float_t shifty, Float_t fac) | |
186 | { | |
187 | SetPad(width,height); | |
188 | SetChevron(hstep,shifty,fac); | |
189 | } | |
190 | ||
191 | ||
6e7b5431 | 192 | Float_t AliTPCPRF2D::GetPRF(Float_t xin, Float_t yin) |
8c555625 | 193 | { |
cc80f89e | 194 | //function which return pad response |
195 | //for the charge in distance xin | |
196 | //return cubic aproximation of PRF or PRF at nearest virtual wire | |
6e7b5431 | 197 | if (fChargeArray==0) return 0; |
8c555625 | 198 | //transform position to "wire position" |
199 | Float_t y=fDYtoWire*(yin-fY1); | |
200 | if (fNYdiv == 1) y=fY1; | |
201 | //normaly it find nearest line charge | |
6e7b5431 | 202 | if (fInterY ==0){ |
8c555625 | 203 | Int_t i=Int_t(0.5+y); |
204 | if (y<0) i=Int_t(-0.5+y); | |
205 | if ((i<0) || (i>=fNYdiv) ) return 0; | |
6e7b5431 | 206 | fcharge = &(fChargeArray[i*fNPRF]); |
8c555625 | 207 | return GetPRFActiv(xin); |
208 | } | |
209 | else{ | |
210 | //make interpolation from more fore lines | |
211 | Int_t i= Int_t(y); | |
6e7b5431 | 212 | Float_t res; |
8c555625 | 213 | if ((i<0) || (i>=fNYdiv) ) return 0; |
214 | Float_t z0=0; | |
215 | Float_t z1=0; | |
216 | Float_t z2=0; | |
217 | Float_t z3=0; | |
218 | if (i>0) { | |
6e7b5431 | 219 | fcharge =&(fChargeArray[(i-1)*fNPRF]); |
8c555625 | 220 | z0 = GetPRFActiv(xin); |
221 | } | |
6e7b5431 | 222 | fcharge =&(fChargeArray[i*fNPRF]); |
8c555625 | 223 | z1=GetPRFActiv(xin); |
224 | if ((i+1)<fNYdiv){ | |
6e7b5431 | 225 | fcharge =&(fChargeArray[(i+1)*fNPRF]); |
8c555625 | 226 | z2 = GetPRFActiv(xin); |
227 | } | |
228 | if ((i+2)<fNYdiv){ | |
6e7b5431 | 229 | fcharge =&(fChargeArray[(i+2)*fNPRF]); |
8c555625 | 230 | z3 = GetPRFActiv(xin); |
231 | } | |
73042f01 | 232 | Float_t a,b,c,d,k,l; |
8c555625 | 233 | a=z1; |
234 | b=(z2-z0)/2.; | |
73042f01 | 235 | k=z2-a-b; |
236 | l=(z3-z1)/2.-b; | |
237 | d=l-2*k; | |
238 | c=k-d; | |
8c555625 | 239 | Float_t dy=y-Float_t(i); |
6e7b5431 | 240 | |
241 | res = a+b*dy+c*dy*dy+d*dy*dy*dy; | |
8c555625 | 242 | return res; |
243 | } | |
cc80f89e | 244 | return 0.; |
8c555625 | 245 | } |
246 | ||
247 | ||
248 | Float_t AliTPCPRF2D::GetPRFActiv(Float_t xin) | |
249 | { | |
cc80f89e | 250 | //GEt response function on given charege line |
251 | //return spline aproximaton | |
8c555625 | 252 | Float_t x = (xin*fDStepM1)+fNPRF/2; |
253 | Int_t i = Int_t(x); | |
254 | ||
6e7b5431 | 255 | if ( (i>1) && ((i+2)<fNPRF)) { |
73042f01 | 256 | Float_t a,b,c,d,k,l; |
8c555625 | 257 | a = fcharge[i]; |
258 | b = (fcharge[i+1]-fcharge[i-1])*0.5; | |
73042f01 | 259 | k = fcharge[i+1]-a-b; |
260 | l = (fcharge[i+2]-fcharge[i])*0.5-b; | |
261 | d=l-2.*k; | |
262 | c=k-d; | |
8c555625 | 263 | Float_t dx=x-Float_t(i); |
264 | Float_t res = a+b*dx+c*dx*dx+d*dx*dx*dx; | |
265 | return res; | |
266 | } | |
267 | else return 0; | |
268 | } | |
269 | ||
270 | ||
271 | Float_t AliTPCPRF2D::GetGRF(Float_t xin, Float_t yin) | |
272 | { | |
cc80f89e | 273 | //function which returnoriginal charge distribution |
274 | //this function is just normalised for fKnorm | |
6e7b5431 | 275 | if (GetGRF() != 0 ) |
276 | return fKNorm*GetGRF()->Eval(xin,yin)/fInteg; | |
8c555625 | 277 | else |
278 | return 0.; | |
279 | } | |
280 | ||
281 | ||
282 | void AliTPCPRF2D::SetParam( TF2 * GRF, Float_t kNorm, | |
283 | Float_t sigmaX, Float_t sigmaY) | |
284 | { | |
cc80f89e | 285 | //adjust parameters of the original charge distribution |
286 | //and pad size parameters | |
8c555625 | 287 | if (fGRF !=0 ) fGRF->Delete(); |
288 | fGRF = GRF; | |
6e7b5431 | 289 | fKNorm = kNorm; |
290 | sprintf(fType,"User"); | |
291 | if (sigmaX ==0) sigmaX=(fWidth*(1+TMath::Abs(fK)))/fgkSQRT12; | |
292 | if (sigmaY ==0) sigmaY=(fWidth*(1+TMath::Abs(fK)))/fgkSQRT12; | |
cc80f89e | 293 | fOrigSigmaX=sigmaX; |
294 | fOrigSigmaY=sigmaY; | |
6e7b5431 | 295 | Double_t estimsigma = |
296 | TMath::Sqrt(sigmaX*sigmaX+(fWidth*fWidth*(1+TMath::Abs(fK))/12)+ | |
297 | TMath::Tan(fPadAngle*fgkDegtoRad)*TMath::Tan(fPadAngle*fgkDegtoRad)*fHeightFull*fHeightFull/12); | |
298 | if (estimsigma < 5*sigmaX) { | |
299 | fDStep = estimsigma/10.; | |
300 | fNPRF = Int_t(estimsigma*8./fDStep); | |
301 | } | |
302 | else{ | |
303 | fDStep = sigmaX; | |
304 | Double_t width = fWidth*(1+TMath::Abs(fK))+TMath::Abs(TMath::Tan(fPadAngle*fgkDegtoRad))*fHeightFull; | |
305 | fNPRF = Int_t((width+8.*sigmaX)/fDStep); | |
306 | }; | |
307 | ||
8c555625 | 308 | } |
309 | ||
310 | ||
311 | void AliTPCPRF2D::SetGauss(Float_t sigmaX, Float_t sigmaY, | |
312 | Float_t kNorm) | |
313 | { | |
cc80f89e | 314 | // |
315 | // set parameters for Gauss generic charge distribution | |
316 | // | |
6e7b5431 | 317 | fKNorm = kNorm; |
318 | fOrigSigmaX=sigmaX; | |
319 | fOrigSigmaY=sigmaY; | |
320 | sprintf(fType,"Gauss"); | |
8c555625 | 321 | if (fGRF !=0 ) fGRF->Delete(); |
322 | fGRF = new TF2("fun",funGauss2D,-5.,5.,-5.,5.,4); | |
6e7b5431 | 323 | |
8c555625 | 324 | funParam[0]=sigmaX; |
325 | funParam[1]=sigmaY; | |
326 | funParam[2]=fK; | |
327 | funParam[3]=fHeightS; | |
6e7b5431 | 328 | |
329 | fGRF->SetParameters(funParam); | |
330 | Double_t estimsigma = | |
331 | TMath::Sqrt(sigmaX*sigmaX+(fWidth*fWidth*(1+TMath::Abs(fK))/12)+ | |
332 | TMath::Tan(fPadAngle)*TMath::Tan(fPadAngle*fgkDegtoRad)*fHeightFull*fHeightFull/12); | |
333 | if (estimsigma < 5*sigmaX) { | |
334 | fDStep = estimsigma/10.; | |
335 | fNPRF = Int_t(estimsigma*8./fDStep); | |
336 | } | |
337 | else{ | |
338 | fDStep = sigmaX; | |
339 | Double_t width = fWidth*(1+TMath::Abs(fK))+TMath::Abs(TMath::Tan(fPadAngle*fgkDegtoRad))*fHeightFull; | |
340 | fNPRF = Int_t((width+8.*sigmaX)/fDStep); | |
341 | }; | |
342 | ||
343 | ||
8c555625 | 344 | } |
8c555625 | 345 | void AliTPCPRF2D::SetCosh(Float_t sigmaX, Float_t sigmaY, |
346 | Float_t kNorm) | |
cc80f89e | 347 | { |
348 | // set parameters for Cosh generic charge distribution | |
349 | // | |
6e7b5431 | 350 | fKNorm = kNorm; |
351 | fOrigSigmaX=sigmaX; | |
352 | fOrigSigmaY=sigmaY; | |
353 | sprintf(fType,"Cosh"); | |
8c555625 | 354 | if (fGRF !=0 ) fGRF->Delete(); |
355 | fGRF = new TF2("fun", funCosh2D,-5.,5.,-5.,5.,4); | |
356 | funParam[0]=sigmaX; | |
357 | funParam[1]=sigmaY; | |
358 | funParam[2]=fK; | |
359 | funParam[3]=fHeightS; | |
360 | fGRF->SetParameters(funParam); | |
6e7b5431 | 361 | |
362 | Double_t estimsigma = TMath::Sqrt(sigmaX*sigmaX+fWidth*fWidth*(1+TMath::Abs(fK))/12); | |
363 | if (estimsigma < 5*sigmaX) { | |
364 | fDStep = estimsigma/10.; | |
365 | fNPRF = Int_t(estimsigma*8./fDStep); | |
366 | } | |
367 | else{ | |
368 | fDStep = sigmaX; | |
369 | fNPRF = Int_t((1.2*fWidth*(1+TMath::Abs(fK))+8.*sigmaX)/fDStep); | |
370 | }; | |
371 | ||
8c555625 | 372 | } |
373 | ||
374 | void AliTPCPRF2D::SetGati(Float_t K3X, Float_t K3Y, | |
375 | Float_t padDistance, | |
376 | Float_t kNorm) | |
377 | { | |
cc80f89e | 378 | // set parameters for Gati generic charge distribution |
379 | // | |
6e7b5431 | 380 | fKNorm = kNorm; |
8c555625 | 381 | fK3X=K3X; |
382 | fK3Y=K3Y; | |
6e7b5431 | 383 | fPadDistance=padDistance; |
384 | sprintf(fType,"Gati"); | |
385 | if (fGRF !=0 ) fGRF->Delete(); | |
386 | fGRF = new TF2("fun", funGati2D,-5.,5.,-5.,5.,5); | |
387 | ||
8c555625 | 388 | funParam[0]=padDistance; |
389 | funParam[1]=K3X; | |
390 | funParam[2]=fK; | |
391 | funParam[3]=fHeightS; | |
392 | funParam[4]=K3Y; | |
393 | fGRF->SetParameters(funParam); | |
cc80f89e | 394 | fOrigSigmaX=padDistance; |
395 | fOrigSigmaY=padDistance; | |
6e7b5431 | 396 | Float_t sigmaX = fOrigSigmaX; |
397 | Double_t estimsigma = TMath::Sqrt(sigmaX*sigmaX+fWidth*fWidth*(1+TMath::Abs(fK))/12); | |
398 | if (estimsigma < 5*sigmaX) { | |
399 | fDStep = estimsigma/10.; | |
400 | fNPRF = Int_t(estimsigma*8./fDStep); | |
401 | } | |
402 | else{ | |
403 | fDStep = sigmaX; | |
404 | fNPRF = Int_t((1.2*fWidth*(1+TMath::Abs(fK))+8.*sigmaX)/fDStep); | |
405 | }; | |
8c555625 | 406 | } |
407 | ||
408 | ||
409 | ||
410 | void AliTPCPRF2D::Update() | |
411 | { | |
cc80f89e | 412 | // |
413 | //update fields with interpolated values for | |
414 | //PRF calculation | |
415 | ||
416 | if ( fGRF == 0 ) return; | |
417 | //initialize interpolated values to 0 | |
418 | Int_t i; | |
6e7b5431 | 419 | if (fChargeArray!=0) delete [] fChargeArray; |
420 | fChargeArray = new Float_t[fNPRF*fNYdiv]; | |
421 | fNChargeArray = fNPRF*fNYdiv; | |
422 | for (i =0; i<fNPRF*fNYdiv;i++) fChargeArray[i] = 0; | |
cc80f89e | 423 | //firstly calculate total integral of charge |
424 | ||
425 | //////////////////////////////////////////////////////// | |
426 | //I'm waiting for normal integral | |
427 | //in this moment only sum | |
428 | Float_t x2= 4*fOrigSigmaX; | |
429 | Float_t y2= 4*fOrigSigmaY; | |
430 | Float_t dx = fOrigSigmaX/Float_t(fNdiv*6); | |
431 | Float_t dy = fOrigSigmaY/Float_t(fNdiv*6); | |
432 | Int_t nx = Int_t(0.5+x2/dx); | |
433 | Int_t ny = Int_t(0.5+y2/dy); | |
434 | Int_t ix,iy; | |
435 | fInteg = 0; | |
436 | Double_t dInteg =0; | |
437 | for (ix=-nx;ix<=nx;ix++) | |
438 | for ( iy=-ny;iy<=ny;iy++) | |
439 | dInteg+=fGRF->Eval(Float_t(ix)*dx,Float_t(iy)*dy)*dx*dy; | |
440 | ///////////////////////////////////////////////////// | |
441 | fInteg =dInteg; | |
442 | if ( fInteg == 0 ) fInteg = 1; | |
443 | ||
444 | for (i=0; i<fNYdiv; i++){ | |
445 | if (fNYdiv == 1) fCurrentY = fY1; | |
8c555625 | 446 | else |
cc80f89e | 447 | fCurrentY = fY1+Double_t(i)*(fY2-fY1)/Double_t(fNYdiv-1); |
6e7b5431 | 448 | fcharge = &(fChargeArray[i*fNPRF]); |
8c555625 | 449 | Update1(); |
450 | } | |
cc80f89e | 451 | //calculate conversion coefitient to convert position to virtual wire |
452 | fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1); | |
453 | fDStepM1=1/fDStep; | |
454 | UpdateSigma(); | |
8c555625 | 455 | } |
456 | ||
8c555625 | 457 | void AliTPCPRF2D::Update1() |
458 | { | |
cc80f89e | 459 | // |
460 | //update fields with interpolated values for | |
461 | //PRF calculation for given charge line | |
8c555625 | 462 | Int_t i; |
cc80f89e | 463 | Double_t cos = TMath::Cos(fChargeAngle); |
464 | Double_t sin = TMath::Sin(fChargeAngle); | |
6e7b5431 | 465 | const Double_t kprec =0.00000001; |
466 | //integrate charge over pad for different distance of pad | |
467 | for (i =0; i<fNPRF;i++){ | |
468 | //x in cm fWidth in cm | |
469 | //calculate integral | |
470 | Double_t xch = fDStep * (Double_t)(i-fNPRF/2); | |
471 | fcharge[i]=0; | |
472 | Double_t k=1; | |
cc80f89e | 473 | |
6e7b5431 | 474 | |
475 | for (Double_t ym=-fHeightFull/2.-fShiftY; ym<fHeightFull/2.-kprec;ym+=fHeightS){ | |
476 | Double_t y2chev=TMath::Min((ym+fHeightS),Double_t(fHeightFull/2.)); // end of chevron step | |
477 | Double_t y1chev= ym; //beginning of chevron step | |
478 | Double_t y2 = TMath::Min(y2chev,fCurrentY+3.5*fOrigSigmaY); | |
479 | Double_t y1 = TMath::Max((y1chev),Double_t(-fHeightFull/2.)); | |
480 | y1 = TMath::Max(y1chev,fCurrentY-3.5*fOrigSigmaY); | |
481 | ||
482 | Double_t x0 = fWidth*(-1.-(Double_t(k)*fK))*0.5+ym*TMath::Tan(fPadAngle*fgkDegtoRad); | |
483 | Double_t kx = Double_t(k)*(fK*fWidth)/fHeightS; | |
484 | kx = TMath::Tan(TMath::ATan(kx))+TMath::Tan(fPadAngle*fgkDegtoRad); | |
485 | ||
486 | Int_t ny = TMath::Max(Int_t(fNdiv*TMath::Exp(-(y1-fCurrentY)*(y1-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),4); | |
487 | Double_t dy = TMath::Min(fOrigSigmaY/Double_t(ny),y2-y1); | |
488 | Double_t ndy = dy; | |
489 | ||
490 | //loop over different y strips with variable step size dy | |
491 | if (y2>(y1+kprec)) for (Double_t y = y1; y<y2+kprec;){ | |
492 | //new step SIZE | |
cc80f89e | 493 | |
6e7b5431 | 494 | ny = TMath::Max(Int_t(fNdiv*TMath::Exp(-(y-fCurrentY)*(y-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),5); |
495 | ndy = fOrigSigmaY/Double_t(ny); | |
496 | if (ndy>(y2-y-dy)) { | |
497 | ndy =y2-y-dy; | |
498 | if (ndy<kprec) ndy=2*kprec; //calculate new delta y | |
499 | } | |
500 | // | |
501 | Double_t sumch=0; | |
502 | //calculation of x borders and initial step | |
503 | Double_t deltay = (y-y1chev); | |
504 | ||
505 | Double_t xp1 = x0+deltay*kx; | |
506 | //x begining of pad at position y | |
507 | Double_t xp2 =xp1+fWidth; //x end of pad at position y | |
508 | Double_t xp3 =xp1+kx*dy; //...at position y+dy | |
509 | Double_t xp4 =xp2+kx*dy; //.. | |
8c555625 | 510 | |
6e7b5431 | 511 | Double_t x1 = TMath::Min(xp1,xp3); |
512 | x1 = TMath::Max(xp1,xch-3.5*fOrigSigmaX); //beging of integration | |
513 | Double_t x2 = TMath::Max(xp2,xp4); | |
514 | x2 = TMath::Min(xp2+dy*kx,xch+3.5*fOrigSigmaX); //end of integration | |
515 | ||
516 | Int_t nx = TMath::Max(Int_t(fNdiv*TMath::Exp(-(x1-xch)*(x1-xch)/(2*fOrigSigmaX*fOrigSigmaX))* | |
517 | TMath::Exp(-(y1-fCurrentY)*(y1-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),2); | |
518 | Double_t dx = TMath::Min(fOrigSigmaX/Double_t(nx),x2-x1)/5.; //on the border more iteration | |
519 | Double_t ndx=dx; | |
520 | ||
521 | if (x2>(x1+kprec)) { | |
522 | for (Double_t x = x1; x<x2+kprec ;){ | |
523 | //new step SIZE | |
524 | nx = TMath::Max(Int_t(fNdiv*TMath::Exp(-(x-xch)*(x-xch)/(2*fOrigSigmaX*fOrigSigmaX))),3); | |
525 | ndx = fOrigSigmaX/Double_t(nx); | |
526 | if (ndx>(x2-x-dx)) { | |
527 | ndx =x2-x-dx; | |
8c555625 | 528 | } |
6e7b5431 | 529 | if ( ( (x+dx+ndx)<TMath::Max(xp3,xp1)) || ( (x+dx+ndx)>TMath::Min(xp4,xp2))) { |
530 | ndx/=5.; | |
531 | } | |
532 | if (ndx<kprec) ndx=2*kprec; | |
533 | //INTEGRAL APROXIMATION | |
534 | Double_t ddx,ddy,dddx,dddy; | |
535 | ddx = xch-(x+dx/2.); | |
536 | ddy = fCurrentY-(y+dy/2.); | |
537 | dddx = cos*ddx-sin*ddy; | |
538 | dddy = sin*ddx+cos*ddy; | |
539 | Double_t z0=fGRF->Eval(dddx,dddy); //middle point | |
540 | ||
541 | ddx = xch-(x+dx/2.); | |
542 | ddy = fCurrentY-(y); | |
543 | dddx = cos*ddx-sin*ddy; | |
544 | dddy = sin*ddx+cos*ddy; | |
545 | Double_t z1=fGRF->Eval(dddx,dddy); //point down | |
546 | ||
547 | ddx = xch-(x+dx/2.); | |
548 | ddy = fCurrentY-(y+dy); | |
549 | dddx = cos*ddx-sin*ddy; | |
550 | dddy = sin*ddx+cos*ddy; | |
551 | Double_t z3=fGRF->Eval(dddx,dddy); //point up | |
552 | ||
553 | ddx = xch-(x); | |
554 | ddy = fCurrentY-(y+dy/2.); | |
555 | dddx = cos*ddx-sin*ddy; | |
556 | dddy = sin*ddx+cos*ddy; | |
557 | Double_t z2=fGRF->Eval(dddx,dddy); //point left | |
558 | ||
559 | ddx = xch-(x+dx); | |
560 | ddy = fCurrentY-(y+dy/2.); | |
561 | dddx = cos*ddx-sin*ddy; | |
562 | dddy = sin*ddx+cos*ddy; | |
563 | Double_t z4=fGRF->Eval(dddx,dddy); //point right | |
564 | ||
565 | ||
566 | if (z0<0) {z0=0;z1=0;z2=0;z3=0;z4=0;} | |
567 | ||
568 | Double_t f2x= (z3+z1-2*z0)*4.;//second derivation in y | |
569 | Double_t f2y= (z2+z4-2*z0)*4.;//second derivation in x | |
570 | Double_t f1y= (z3-z1); | |
571 | Double_t z ; | |
572 | z = (z0+f2x/6.+f2y/6.);//second order aproxiation of integral | |
573 | if (kx>kprec){ //positive derivation | |
574 | if (x<(xp1+dy*kx)){ //calculate volume at left border | |
575 | Double_t xx1 = x; | |
576 | Double_t xx2 = TMath::Min(x+dx,xp1+dy*kx); | |
577 | Double_t yy1 = y+(xx1-xp1)/kx; | |
578 | Double_t yy2 = TMath::Min(y+(xx2-xp1)/kx,y+dy); | |
579 | z=z0; | |
580 | if (yy2<y+dy) { | |
581 | z-= z0*(y+dy-yy2)/dy; //constant part rectangle | |
582 | z-= f1y*(xx2-xx1)*(y+dy-yy2)*(y+dy-yy2)/(2.*dx*dy); | |
583 | } | |
584 | z-=z0*(xx2-xx1)*(yy2-yy1)/(2*dx*dy); //constant part rectangle | |
585 | ||
586 | } | |
587 | if (x>xp2){ //calculate volume at right border | |
588 | Double_t xx1 = x; | |
589 | Double_t xx2 = x+dx; | |
590 | Double_t yy1 = y+(xx1-xp2)/kx; | |
591 | Double_t yy2 = y+(xx2-xp2)/kx; | |
592 | z=z0; | |
593 | //rectangle part | |
594 | z-=z0*(yy1-y)/dy; //constant part | |
595 | z-=f1y*(xx2-xx1)*(yy1-y)*(yy1-y)/(2*dx*dy); | |
596 | //triangle part | |
597 | z-=z0*(xx2-xx1)*(yy2-yy1)/(2*dx*dy); //constant part | |
598 | } | |
599 | } | |
600 | if (kx<-kprec){ //negative derivation | |
601 | if (x<(xp1+dy*kx)){ //calculate volume at left border | |
602 | Double_t xx1 = x; | |
603 | Double_t xx2 = TMath::Min(x+dx,xp3-dy/kx); | |
604 | Double_t yy1 = y+(xx1-xp1)/kx; | |
605 | Double_t yy2 = TMath::Max(y,yy1+(xx2-xx1)/kx); //yy2<yy1 | |
606 | z = z0; | |
607 | z-= z0*(yy2-y)/dy; // constant part rectangle | |
608 | z-= f1y*(xx2-xx1)*(yy2-y)*(yy2-y)/(2.*dx*dy); | |
609 | z-=z0*(xx2-xx1)*(yy1-yy2)/(2*dx*dy); //constant part triangle | |
610 | } | |
611 | if (x>xp2){ //calculate volume at right border | |
612 | Double_t xx1 = TMath::Max(x,xp2+dy*kx); | |
613 | Double_t xx2 = x+dx; | |
614 | Double_t yy1 = TMath::Min(y+dy,y-(xp2-xx1)/kx); | |
615 | Double_t yy2 = y-(xp2-xx2)/kx; | |
616 | z=z0; | |
617 | z-=z0*(yy2-y)/dy; //constant part rextangle | |
618 | z-= f1y*(xx2-xx1)*(yy2-y)*(yy2-y)/(2.*dx*dy); | |
619 | z-=z0*(xx2-xx1)*(yy1-yy2)/(2*dx*dy); //constant part triangle | |
620 | } | |
621 | } | |
622 | ||
623 | if (z>0.) sumch+=fKNorm*z*dx*dy/fInteg; | |
624 | ||
625 | x+=dx; | |
626 | dx = ndx; | |
627 | }; //loop over x | |
628 | fcharge[i]+=sumch; | |
629 | }//if x2>x1 | |
630 | y+=dy; | |
631 | dy =ndy; | |
632 | }//step over different y | |
633 | k*=-1.; | |
634 | }//step over chevron | |
cc80f89e | 635 | |
6e7b5431 | 636 | }//step over different points on line NPRF |
cc80f89e | 637 | } |
638 | ||
639 | void AliTPCPRF2D::UpdateSigma() | |
640 | { | |
641 | // | |
642 | //calulate effective sigma X and sigma y of PRF | |
643 | fMeanX = 0; | |
644 | fMeanY = 0; | |
645 | fSigmaX = 0; | |
646 | fSigmaY = 0; | |
647 | ||
8c555625 | 648 | Float_t sum =0; |
cc80f89e | 649 | Int_t i; |
650 | Float_t x,y; | |
651 | ||
652 | for (i=-1; i<=fNYdiv; i++){ | |
653 | if (fNYdiv == 1) y = fY1; | |
654 | else | |
655 | y = fY1+Float_t(i)*(fY2-fY1)/Float_t(fNYdiv-1); | |
656 | for (x =-fNPRF*fDStep; x<fNPRF*fDStep;x+=fDStep) | |
657 | { | |
658 | //x in cm fWidth in cm | |
659 | Float_t weight = GetPRF(x,y); | |
660 | fSigmaX+=x*x*weight; | |
661 | fSigmaY+=y*y*weight; | |
662 | fMeanX+=x*weight; | |
663 | fMeanY+=y*weight; | |
664 | sum+=weight; | |
8c555625 | 665 | }; |
cc80f89e | 666 | } |
8c555625 | 667 | if (sum>0){ |
cc80f89e | 668 | fMeanX/=sum; |
669 | fMeanY/=sum; | |
670 | fSigmaX = TMath::Sqrt(fSigmaX/sum-fMeanX*fMeanX); | |
671 | fSigmaY = TMath::Sqrt(fSigmaY/sum-fMeanY*fMeanY); | |
8c555625 | 672 | } |
673 | else fSigmaX=0; | |
8c555625 | 674 | } |
675 | ||
cc80f89e | 676 | |
8c555625 | 677 | void AliTPCPRF2D::Streamer(TBuffer &R__b) |
678 | { | |
679 | // Stream an object of class AliTPCPRF2D | |
680 | ||
681 | if (R__b.IsReading()) { | |
682 | Version_t R__v = R__b.ReadVersion(); if (R__v) { } | |
683 | TObject::Streamer(R__b); | |
684 | //read chewron parameters | |
8c555625 | 685 | R__b >> fHeightFull; |
686 | R__b >> fHeightS; | |
687 | R__b >> fShiftY; | |
688 | R__b >> fWidth; | |
689 | R__b >> fK; | |
cc80f89e | 690 | R__b >> fSigmaX; |
691 | R__b >> fSigmaY; | |
692 | R__b >> fMeanX; | |
693 | R__b >> fMeanY; | |
694 | //read charge parameters | |
695 | R__b.ReadFastArray(fType,5); | |
696 | R__b >> fOrigSigmaX; | |
697 | R__b >> fOrigSigmaY; | |
6e7b5431 | 698 | R__b >> fKNorm; |
8c555625 | 699 | R__b >> fK3X; |
700 | R__b >> fK3Y; | |
701 | R__b >> fPadDistance; | |
6e7b5431 | 702 | R__b >> fInteg; |
703 | //read angle parameters | |
704 | R__b >> fChargeAngle; | |
705 | R__b >> fPadAngle; | |
8c555625 | 706 | //read functions |
707 | if (fGRF!=0) { | |
cc80f89e | 708 | fGRF->Delete(); |
8c555625 | 709 | fGRF=0; |
710 | } | |
711 | if (strncmp(fType,"User",3)==0){ | |
712 | fGRF= new TF2; | |
713 | R__b>>fGRF; | |
714 | } | |
715 | if (strncmp(fType,"Gauss",3)==0) | |
716 | fGRF = new TF2("fun",funGauss2D,-5.,5.,-5.,5.,4); | |
717 | if (strncmp(fType,"Cosh",3)==0) | |
718 | fGRF = new TF2("fun",funCosh2D,-5.,5.,-5.,5.,4); | |
719 | if (strncmp(fType,"Gati",3)==0) | |
cc80f89e | 720 | fGRF = new TF2("fun",funGati2D,-5.,5.,-5.,5.,5); |
8c555625 | 721 | //read interpolation parameters |
722 | R__b >>fY1; | |
723 | R__b >>fY2; | |
724 | R__b >>fNYdiv; | |
725 | R__b >>fDStep; | |
726 | R__b >>fNPRF; | |
6e7b5431 | 727 | R__b >>fNChargeArray; |
728 | if (fChargeArray!=0) delete [] fChargeArray; | |
729 | if (fNChargeArray>0) { | |
730 | fChargeArray = new Float_t[fNChargeArray]; | |
731 | R__b.ReadFastArray(fChargeArray,fNChargeArray); | |
732 | } | |
733 | // | |
8c555625 | 734 | R__b.ReadFastArray(funParam,5); |
735 | if (fGRF!=0) fGRF->SetParameters(funParam); | |
736 | //calculate conversion coefitient to convert position to virtual wire | |
737 | fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1); | |
738 | fDStepM1=1/fDStep; | |
739 | } else { | |
740 | R__b.WriteVersion(AliTPCPRF2D::IsA()); | |
741 | TObject::Streamer(R__b); | |
742 | //write chewron parameters | |
8c555625 | 743 | R__b << fHeightFull; |
744 | R__b << fHeightS; | |
745 | R__b << fShiftY; | |
746 | R__b << fWidth; | |
747 | R__b << fK; | |
cc80f89e | 748 | R__b << fSigmaX; |
749 | R__b << fSigmaY; | |
750 | R__b << fMeanX; | |
751 | R__b << fMeanY; | |
8c555625 | 752 | //write charge parameters |
cc80f89e | 753 | R__b.WriteFastArray(fType,5); |
754 | R__b << fOrigSigmaX; | |
755 | R__b << fOrigSigmaY; | |
6e7b5431 | 756 | R__b << fKNorm; |
8c555625 | 757 | R__b << fK3X; |
758 | R__b << fK3Y; | |
759 | R__b << fPadDistance; | |
760 | R__b << fInteg; | |
6e7b5431 | 761 | //angle parameters |
762 | R__b << fChargeAngle; | |
763 | R__b << fPadAngle; | |
8c555625 | 764 | if (strncmp(fType,"User",3)==0) R__b <<fGRF; |
765 | //write interpolation parameters | |
766 | R__b <<fY1; | |
767 | R__b <<fY2; | |
768 | R__b <<fNYdiv; | |
769 | R__b <<fDStep; | |
770 | R__b <<fNPRF; | |
6e7b5431 | 771 | R__b <<fNChargeArray; |
772 | if (fNChargeArray>0) | |
773 | R__b.WriteFastArray(fChargeArray,fNPRF*fNYdiv); | |
8c555625 | 774 | R__b.WriteFastArray(funParam,5); |
775 | } | |
776 | } | |
777 | ||
778 | ||
6e7b5431 | 779 | TH1F * AliTPCPRF2D::GenerDrawXHisto(Float_t x1, Float_t x2,Float_t y) |
780 | { | |
781 | //gener one dimensional hist of pad response function | |
782 | // at position y | |
783 | char s[100]; | |
784 | const Int_t kn=200; | |
785 | sprintf(s,"Pad Response Function"); | |
786 | TH1F * hPRFc = new TH1F("hPRFc",s,kn+1,x1,x2); | |
8c555625 | 787 | Float_t x=x1; |
788 | Float_t y1; | |
8c555625 | 789 | |
6e7b5431 | 790 | for (Int_t i = 0;i<kn+1;i++) |
8c555625 | 791 | { |
6e7b5431 | 792 | x+=(x2-x1)/Float_t(kn); |
793 | y1 = GetPRF(x,y); | |
8c555625 | 794 | hPRFc->Fill(x,y1); |
795 | }; | |
6e7b5431 | 796 | hPRFc->SetXTitle("pad (cm)"); |
797 | return hPRFc; | |
798 | } | |
8c555625 | 799 | |
6e7b5431 | 800 | AliH2F * AliTPCPRF2D::GenerDrawHisto(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx, Int_t Ny) |
801 | { | |
802 | // | |
803 | //gener two dimensional histogram with PRF | |
804 | // | |
805 | char s[100]; | |
806 | sprintf(s,"Pad Response Function"); | |
807 | AliH2F * hPRFc = new AliH2F("hPRFc",s,Nx,x1,x2,Ny,y1,y2); | |
808 | Float_t dx=(x2-x1)/Float_t(Nx); | |
809 | Float_t dy=(y2-y1)/Float_t(Ny) ; | |
810 | Float_t x,y,z; | |
811 | x = x1; | |
812 | y = y1; | |
813 | for ( Int_t i = 0;i<=Nx;i++,x+=dx){ | |
814 | y=y1; | |
815 | for (Int_t j = 0;j<=Ny;j++,y+=dy){ | |
816 | z = GetPRF(x,y); | |
817 | hPRFc->SetCellContent(i,j,z); | |
818 | }; | |
819 | }; | |
820 | hPRFc->SetXTitle("pad direction (cm)"); | |
821 | hPRFc->SetYTitle("pad row direction (cm)"); | |
822 | hPRFc->SetTitleOffset(1.5,"X"); | |
823 | hPRFc->SetTitleOffset(1.5,"Y"); | |
824 | return hPRFc; | |
825 | } | |
826 | ||
827 | ||
828 | AliH2F * AliTPCPRF2D::GenerDrawDistHisto(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx, Int_t Ny, Float_t thr) | |
829 | { | |
830 | //return histogram with distortion | |
831 | const Float_t kminth=0.00001; | |
832 | if (thr<kminth) thr=kminth; | |
833 | char s[100]; | |
834 | sprintf(s,"COG distortion of PRF (threshold=%2.2f)",thr); | |
835 | AliH2F * hPRFDist = new AliH2F("hDistortion",s,Nx,x1,x2,Ny,y1,y2); | |
836 | Float_t dx=(x2-x1)/Float_t(Nx); | |
837 | Float_t dy=(y2-y1)/Float_t(Ny) ; | |
838 | Float_t x,y,z,ddx; | |
839 | x=x1; | |
840 | for ( Int_t i = 0;i<=Nx;i++,x+=dx){ | |
841 | y=y1; | |
842 | for(Int_t j = 0;j<=Ny;j++,y+=dy) | |
843 | { | |
844 | Float_t sumx=0; | |
845 | Float_t sum=0; | |
846 | for (Int_t k=-3;k<=3;k++) | |
847 | { | |
848 | Float_t padx=Float_t(k)*fWidth; | |
849 | z = GetPRF(x-padx,y); | |
850 | if (z>thr){ | |
851 | sum+=z; | |
852 | sumx+=z*padx; | |
853 | } | |
854 | }; | |
855 | if (sum>kminth) | |
856 | { | |
857 | ddx = (x-(sumx/sum)); | |
858 | } | |
859 | else ddx=-1; | |
860 | if (TMath::Abs(ddx)<10) hPRFDist->SetCellContent(i,j,ddx); | |
861 | } | |
862 | } | |
863 | ||
864 | hPRFDist->SetXTitle("pad direction (cm)"); | |
865 | hPRFDist->SetYTitle("pad row direction (cm)"); | |
866 | hPRFDist->SetTitleOffset(1.5,"X"); | |
867 | hPRFDist->SetTitleOffset(1.5,"Y"); | |
868 | return hPRFDist; | |
869 | } | |
870 | ||
871 | ||
872 | ||
873 | ||
874 | ||
875 | void AliTPCPRF2D::DrawX(Float_t x1 ,Float_t x2,Float_t y1,Float_t y2, Int_t N) | |
876 | { | |
877 | // | |
878 | //draw pad response function at interval <x1,x2> at given y position | |
879 | // | |
880 | if (N<0) return; | |
881 | TCanvas * c1 = new TCanvas("PRFX","Pad response function",700,900); | |
882 | c1->cd(); | |
883 | ||
8c555625 | 884 | TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC"); |
885 | comment->SetTextAlign(12); | |
886 | comment->SetFillColor(42); | |
6e7b5431 | 887 | DrawComment(comment); |
8c555625 | 888 | comment->Draw(); |
6e7b5431 | 889 | c1->cd(); |
890 | ||
891 | TPad * pad2 = new TPad("pPRF","",0.05,0.22,0.95,0.95); | |
892 | pad2->Divide(2,(N+1)/2); | |
893 | pad2->Draw(); | |
894 | gStyle->SetOptFit(1); | |
895 | gStyle->SetOptStat(1); | |
896 | for (Int_t i=0;i<N;i++){ | |
897 | char ch[200]; | |
898 | Float_t y; | |
899 | if (N==1) y=y1; | |
900 | else y = y1+i*(y2-y1)/Float_t(N-1); | |
901 | pad2->cd(i+1); | |
902 | TH1F * hPRFc =GenerDrawXHisto(x1, x2,y); | |
903 | sprintf(ch,"PRF at wire position: %2.3f",y); | |
904 | hPRFc->SetTitle(ch); | |
905 | sprintf(ch,"PRF %d",i); | |
906 | hPRFc->SetName(ch); | |
907 | hPRFc->Fit("gaus"); | |
908 | } | |
909 | ||
8c555625 | 910 | } |
911 | ||
912 | ||
913 | ||
6e7b5431 | 914 | void AliTPCPRF2D::DrawPRF(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2, Int_t Nx, Int_t Ny) |
8c555625 | 915 | { |
6e7b5431 | 916 | // |
917 | // | |
8c555625 | 918 | TCanvas * c1 = new TCanvas("canPRF","Pad response function",700,900); |
919 | c1->cd(); | |
6e7b5431 | 920 | TPad * pad2 = new TPad("pad2PRF","",0.05,0.22,0.95,0.95); |
921 | pad2->Draw(); | |
8c555625 | 922 | gStyle->SetOptFit(1); |
8c555625 | 923 | gStyle->SetOptStat(1); |
6e7b5431 | 924 | TH2F * hPRFc = GenerDrawHisto(x1, x2, y1, y2, Nx,Ny); |
8c555625 | 925 | pad2->cd(); |
6e7b5431 | 926 | hPRFc->Draw("surf"); |
8c555625 | 927 | c1->cd(); |
928 | TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC"); | |
929 | comment->SetTextAlign(12); | |
930 | comment->SetFillColor(42); | |
6e7b5431 | 931 | DrawComment(comment); |
8c555625 | 932 | comment->Draw(); |
933 | } | |
934 | ||
6e7b5431 | 935 | void AliTPCPRF2D::DrawDist(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2, Int_t Nx, Int_t Ny, Float_t thr) |
8c555625 | 936 | { |
6e7b5431 | 937 | // |
938 | //draw distortion of the COG method - for different threshold parameter | |
8c555625 | 939 | TCanvas * c1 = new TCanvas("padDistortion","COG distortion",700,900); |
940 | c1->cd(); | |
6e7b5431 | 941 | TPad * pad1 = new TPad("dist","",0.05,0.55,0.95,0.95,21); |
8c555625 | 942 | pad1->Draw(); |
6e7b5431 | 943 | TPad * pad2 = new TPad("dist","",0.05,0.22,0.95,0.53,21); |
8c555625 | 944 | pad2->Draw(); |
8c555625 | 945 | gStyle->SetOptFit(1); |
946 | gStyle->SetOptStat(0); | |
6e7b5431 | 947 | |
948 | AliH2F * hPRFDist = GenerDrawDistHisto(x1, x2, y1, y2, Nx,Ny,thr); | |
949 | ||
8c555625 | 950 | pad1->cd(); |
6e7b5431 | 951 | hPRFDist->Draw("surf"); |
952 | Float_t distmax =hPRFDist->GetMaximum(); | |
953 | Float_t distmin =hPRFDist->GetMinimum(); | |
954 | gStyle->SetOptStat(1); | |
8c555625 | 955 | |
6e7b5431 | 956 | TH1F * dist = hPRFDist->GetAmplitudes(distmin,distmax,distmin-1); |
957 | pad2->cd(); | |
958 | dist->Draw(); | |
8c555625 | 959 | c1->cd(); |
960 | TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC"); | |
961 | comment->SetTextAlign(12); | |
962 | comment->SetFillColor(42); | |
6e7b5431 | 963 | DrawComment(comment); |
964 | comment->Draw(); | |
965 | } | |
966 | ||
967 | void AliTPCPRF2D::DrawComment(TPaveText *comment) | |
968 | { | |
969 | // | |
970 | //function to write comment to picture | |
971 | ||
972 | char s[100]; | |
973 | //draw comments to picture | |
974 | TText * title = comment->AddText("Pad Response Function parameters:"); | |
8c555625 | 975 | title->SetTextSize(0.03); |
6e7b5431 | 976 | sprintf(s,"Height of pad: %2.2f cm",fHeightFull); |
8c555625 | 977 | comment->AddText(s); |
6e7b5431 | 978 | sprintf(s,"Width pad: %2.2f cm",fWidth); |
8c555625 | 979 | comment->AddText(s); |
6e7b5431 | 980 | sprintf(s,"Pad Angle: %2.2f ",fPadAngle); |
8c555625 | 981 | comment->AddText(s); |
8c555625 | 982 | |
6e7b5431 | 983 | if (TMath::Abs(fK)>0.0001){ |
984 | sprintf(s,"Height of one chevron unit h: %2.2f cm",2*fHeightS); | |
985 | comment->AddText(s); | |
986 | sprintf(s,"Overlap factor: %2.2f",fK); | |
987 | comment->AddText(s); | |
988 | } | |
989 | ||
990 | if (strncmp(fType,"User",3)==0){ | |
991 | sprintf(s,"Charge distribution - user defined function %s ",fGRF->GetTitle()); | |
992 | comment->AddText(s); | |
993 | sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX); | |
994 | comment->AddText(s); | |
995 | sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY); | |
996 | comment->AddText(s); | |
997 | } | |
998 | if (strncmp(fType,"Gauss",3)==0){ | |
999 | sprintf(s,"Gauss charge distribution"); | |
1000 | comment->AddText(s); | |
1001 | sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX); | |
1002 | comment->AddText(s); | |
1003 | sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY); | |
1004 | comment->AddText(s); | |
1005 | } | |
1006 | if (strncmp(fType,"Gati",3)==0){ | |
1007 | sprintf(s,"Gati charge distribution"); | |
1008 | comment->AddText(s); | |
1009 | sprintf(s,"K3X of Gati : %2.2f ",fK3X); | |
1010 | comment->AddText(s); | |
1011 | sprintf(s,"K3Y of Gati: %2.2f ",fK3Y); | |
1012 | comment->AddText(s); | |
1013 | sprintf(s,"Wire to Pad Distance: %2.2f ",fPadDistance); | |
1014 | comment->AddText(s); | |
1015 | } | |
1016 | if (strncmp(fType,"Cosh",3)==0){ | |
1017 | sprintf(s,"Cosh charge distribution"); | |
1018 | comment->AddText(s); | |
1019 | sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX); | |
1020 | comment->AddText(s); | |
1021 | sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY); | |
1022 | comment->AddText(s); | |
1023 | } | |
1024 | sprintf(s,"Normalisation: %2.2f ",fKNorm); | |
1025 | comment->AddText(s); | |
8c555625 | 1026 | } |
1027 |