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