1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
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) //
28 // Origin: Marian Ivanov, Uni. of Bratislava, ivanov@fmph.uniba.sk //
30 ///////////////////////////////////////////////////////////////////////////////
34 #include "AliTPCPRF2D.h"
36 #include <Riostream.h>
45 #include "TPaveText.h"
48 extern TStyle * gStyle;
50 const Double_t AliTPCPRF2D::fgkDegtoRad = 0.01745329251994;
51 const Double_t AliTPCPRF2D::fgkSQRT12=3.464101;
52 const Int_t AliTPCPRF2D::fgkNPRF = 100;
55 static Double_t funGauss2D(Double_t *x, Double_t * par)
57 //Gauss function -needde by the generic function object
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])));
63 static Double_t funCosh2D(Double_t *x, Double_t * par)
65 //Cosh function -needde by the generic function object
66 return ( 1/(TMath::CosH(3.14159*x[0]/(2*par[0]))*
67 TMath::CosH(3.14159*x[1]/(2*par[1]))));
70 static Double_t funGati2D(Double_t *x, Double_t * par)
72 //Gati function -needde by the generic function object
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));
77 Float_t l=x[0]/par[0];
78 Float_t tan2=TMath::TanH(k2*l);
80 Float_t res = k1*(1-tan2)/(1+k3*tan2);
81 //par[4] = is equal to k3Y
84 k2=(TMath::Pi()/2)*(1-k3R/2.);
85 k1=k2*k3R/(4*TMath::ATan(k3R));
87 tan2=TMath::TanH(k2*l);
89 res = res*k1*(1-tan2)/(1+k3*tan2);
93 ///////////////////////////////////////////////////////////////////////////
94 ///////////////////////////////////////////////////////////////////////////
98 AliTPCPRF2D::AliTPCPRF2D()
100 //default constructor for response function object
116 //chewron default values
118 SetChevron(0.2,0.0,1.0);
120 SetInterpolationType(2,0);
123 AliTPCPRF2D::~AliTPCPRF2D()
125 if (fChargeArray!=0) delete [] fChargeArray;
126 if (fGRF !=0 ) fGRF->Delete();
129 void AliTPCPRF2D::SetY(Float_t y1, Float_t y2, Int_t nYdiv)
132 //set virtual line position
133 //first and last line and number of lines
139 void AliTPCPRF2D::SetPad(Float_t width, Float_t height)
141 //set base chevron parameters
145 void AliTPCPRF2D::SetChevron(Float_t hstep,
149 //set shaping of chewron parameters
155 void AliTPCPRF2D::SetChParam(Float_t width, Float_t height,
156 Float_t hstep, Float_t shifty, Float_t fac)
158 SetPad(width,height);
159 SetChevron(hstep,shifty,fac);
163 Float_t AliTPCPRF2D::GetPRF(Float_t xin, Float_t yin)
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
168 if (fChargeArray==0) return 0;
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
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;
177 fcharge = &(fChargeArray[i*fNPRF]);
178 return GetPRFActiv(xin);
180 //make interpolation from more fore lines
183 if ((i<0) || (i>=fNYdiv) ) return 0;
189 fcharge =&(fChargeArray[(i-1)*fNPRF]);
190 z0 = GetPRFActiv(xin);
192 fcharge =&(fChargeArray[i*fNPRF]);
195 fcharge =&(fChargeArray[(i+1)*fNPRF]);
196 z2 = GetPRFActiv(xin);
199 fcharge =&(fChargeArray[(i+2)*fNPRF]);
200 z3 = GetPRFActiv(xin);
209 Float_t dy=y-Float_t(i);
211 res = a+b*dy+c*dy*dy+d*dy*dy*dy;
216 Float_t AliTPCPRF2D::GetPRFActiv(Float_t xin)
218 //GEt response function on given charege line
219 //return spline aproximaton
220 Float_t x = (xin*fDStepM1)+fNPRF/2;
223 if ( (i>1) && ((i+2)<fNPRF)) {
226 b = (fcharge[i+1]-fcharge[i-1])*0.5;
227 k = fcharge[i+1]-a-b;
228 l = (fcharge[i+2]-fcharge[i])*0.5-b;
231 Float_t dx=x-Float_t(i);
232 Float_t res = a+b*dx+c*dx*dx+d*dx*dx*dx;
239 Float_t AliTPCPRF2D::GetGRF(Float_t xin, Float_t yin)
241 //function which returnoriginal charge distribution
242 //this function is just normalised for fKnorm
244 return fKNorm*GetGRF()->Eval(xin,yin)/fInteg;
250 void AliTPCPRF2D::SetParam( TF2 * GRF, Float_t kNorm,
251 Float_t sigmaX, Float_t sigmaY)
253 //adjust parameters of the original charge distribution
254 //and pad size parameters
255 if (fGRF !=0 ) fGRF->Delete();
258 sprintf(fType,"User");
259 if (sigmaX ==0) sigmaX=(fWidth*(1+TMath::Abs(fK)))/fgkSQRT12;
260 if (sigmaY ==0) sigmaY=(fWidth*(1+TMath::Abs(fK)))/fgkSQRT12;
263 Double_t estimsigma =
264 TMath::Sqrt(sigmaX*sigmaX+(fWidth*fWidth*(1+TMath::Abs(fK))/12)+
265 TMath::Tan(fPadAngle*fgkDegtoRad)*TMath::Tan(fPadAngle*fgkDegtoRad)*fHeightFull*fHeightFull/12);
266 if (estimsigma < 5*sigmaX) {
267 fDStep = estimsigma/10.;
268 fNPRF = Int_t(estimsigma*8./fDStep);
272 Double_t width = fWidth*(1+TMath::Abs(fK))+TMath::Abs(TMath::Tan(fPadAngle*fgkDegtoRad))*fHeightFull;
273 fNPRF = Int_t((width+8.*sigmaX)/fDStep);
279 void AliTPCPRF2D::SetGauss(Float_t sigmaX, Float_t sigmaY,
283 // set parameters for Gauss generic charge distribution
288 sprintf(fType,"Gauss");
289 if (fGRF !=0 ) fGRF->Delete();
290 fGRF = new TF2("funGauss2D",funGauss2D,-5.,5.,-5.,5.,4);
295 funParam[3]=fHeightS;
297 fGRF->SetParameters(funParam);
298 Double_t estimsigma =
299 TMath::Sqrt(sigmaX*sigmaX+(fWidth*fWidth*(1+TMath::Abs(fK))/12)+
300 TMath::Tan(fPadAngle)*TMath::Tan(fPadAngle*fgkDegtoRad)*fHeightFull*fHeightFull/12);
301 if (estimsigma < 5*sigmaX) {
302 fDStep = estimsigma/10.;
303 fNPRF = Int_t(estimsigma*8./fDStep);
307 Double_t width = fWidth*(1+TMath::Abs(fK))+TMath::Abs(TMath::Tan(fPadAngle*fgkDegtoRad))*fHeightFull;
308 fNPRF = Int_t((width+8.*sigmaX)/fDStep);
313 void AliTPCPRF2D::SetCosh(Float_t sigmaX, Float_t sigmaY,
316 // set parameters for Cosh generic charge distribution
321 sprintf(fType,"Cosh");
322 if (fGRF !=0 ) fGRF->Delete();
323 fGRF = new TF2("funCosh2D", funCosh2D,-5.,5.,-5.,5.,4);
327 funParam[3]=fHeightS;
328 fGRF->SetParameters(funParam);
330 Double_t estimsigma = TMath::Sqrt(sigmaX*sigmaX+fWidth*fWidth*(1+TMath::Abs(fK))/12);
331 if (estimsigma < 5*sigmaX) {
332 fDStep = estimsigma/10.;
333 fNPRF = Int_t(estimsigma*8./fDStep);
337 fNPRF = Int_t((1.2*fWidth*(1+TMath::Abs(fK))+8.*sigmaX)/fDStep);
342 void AliTPCPRF2D::SetGati(Float_t K3X, Float_t K3Y,
346 // set parameters for Gati generic charge distribution
351 fPadDistance=padDistance;
352 sprintf(fType,"Gati");
353 if (fGRF !=0 ) fGRF->Delete();
354 fGRF = new TF2("funGati2D", funGati2D,-5.,5.,-5.,5.,5);
356 funParam[0]=padDistance;
359 funParam[3]=fHeightS;
361 fGRF->SetParameters(funParam);
362 fOrigSigmaX=padDistance;
363 fOrigSigmaY=padDistance;
364 Float_t sigmaX = fOrigSigmaX;
365 Double_t estimsigma = TMath::Sqrt(sigmaX*sigmaX+fWidth*fWidth*(1+TMath::Abs(fK))/12);
366 if (estimsigma < 5*sigmaX) {
367 fDStep = estimsigma/10.;
368 fNPRF = Int_t(estimsigma*8./fDStep);
372 fNPRF = Int_t((1.2*fWidth*(1+TMath::Abs(fK))+8.*sigmaX)/fDStep);
378 void AliTPCPRF2D::Update()
381 //update fields with interpolated values for
384 if ( fGRF == 0 ) return;
385 //initialize interpolated values to 0
387 if (fChargeArray!=0) delete [] fChargeArray;
388 fChargeArray = new Float_t[fNPRF*fNYdiv];
389 fNChargeArray = fNPRF*fNYdiv;
390 for (i =0; i<fNPRF*fNYdiv;i++) fChargeArray[i] = 0;
391 //firstly calculate total integral of charge
393 ////////////////////////////////////////////////////////
394 //I'm waiting for normal integral
395 //in this moment only sum
396 Float_t x2= 4*fOrigSigmaX;
397 Float_t y2= 4*fOrigSigmaY;
398 Float_t dx = fOrigSigmaX/Float_t(fNdiv*6);
399 Float_t dy = fOrigSigmaY/Float_t(fNdiv*6);
400 Int_t nx = Int_t(0.5+x2/dx);
401 Int_t ny = Int_t(0.5+y2/dy);
405 for (ix=-nx;ix<=nx;ix++)
406 for ( iy=-ny;iy<=ny;iy++)
407 dInteg+=fGRF->Eval(Float_t(ix)*dx,Float_t(iy)*dy)*dx*dy;
408 /////////////////////////////////////////////////////
410 if ( fInteg == 0 ) fInteg = 1;
412 for (i=0; i<fNYdiv; i++){
413 if (fNYdiv == 1) fCurrentY = fY1;
415 fCurrentY = fY1+Double_t(i)*(fY2-fY1)/Double_t(fNYdiv-1);
416 fcharge = &(fChargeArray[i*fNPRF]);
419 //calculate conversion coefitient to convert position to virtual wire
420 fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1);
425 void AliTPCPRF2D::Update1()
428 //update fields with interpolated values for
429 //PRF calculation for given charge line
431 Double_t cos = TMath::Cos(fChargeAngle);
432 Double_t sin = TMath::Sin(fChargeAngle);
433 const Double_t kprec =0.00000001;
434 //integrate charge over pad for different distance of pad
435 for (i =0; i<fNPRF;i++){
436 //x in cm fWidth in cm
438 Double_t xch = fDStep * (Double_t)(i-fNPRF/2);
443 for (Double_t ym=-fHeightFull/2.-fShiftY; ym<fHeightFull/2.-kprec;ym+=fHeightS){
444 Double_t y2chev=TMath::Min((ym+fHeightS),Double_t(fHeightFull/2.)); // end of chevron step
445 Double_t y1chev= ym; //beginning of chevron step
446 Double_t y2 = TMath::Min(y2chev,fCurrentY+3.5*fOrigSigmaY);
447 Double_t y1 = TMath::Max((y1chev),Double_t(-fHeightFull/2.));
448 y1 = TMath::Max(y1chev,fCurrentY-3.5*fOrigSigmaY);
450 Double_t x0 = fWidth*(-1.-(Double_t(k)*fK))*0.5+ym*TMath::Tan(fPadAngle*fgkDegtoRad);
451 Double_t kx = Double_t(k)*(fK*fWidth)/fHeightS;
452 kx = TMath::Tan(TMath::ATan(kx))+TMath::Tan(fPadAngle*fgkDegtoRad);
454 Int_t ny = TMath::Max(Int_t(fNdiv*TMath::Exp(-(y1-fCurrentY)*(y1-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),4);
455 Double_t dy = TMath::Min(fOrigSigmaY/Double_t(ny),y2-y1);
458 //loop over different y strips with variable step size dy
459 if (y2>(y1+kprec)) for (Double_t y = y1; y<y2+kprec;){
462 ny = TMath::Max(Int_t(fNdiv*TMath::Exp(-(y-fCurrentY)*(y-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),5);
463 ndy = fOrigSigmaY/Double_t(ny);
466 if (ndy<kprec) ndy=2*kprec; //calculate new delta y
470 //calculation of x borders and initial step
471 Double_t deltay = (y-y1chev);
473 Double_t xp1 = x0+deltay*kx;
474 //x begining of pad at position y
475 Double_t xp2 =xp1+fWidth; //x end of pad at position y
476 Double_t xp3 =xp1+kx*dy; //...at position y+dy
477 Double_t xp4 =xp2+kx*dy; //..
479 Double_t x1 = TMath::Min(xp1,xp3);
480 x1 = TMath::Max(xp1,xch-3.5*fOrigSigmaX); //beging of integration
481 Double_t x2 = TMath::Max(xp2,xp4);
482 x2 = TMath::Min(xp2+dy*kx,xch+3.5*fOrigSigmaX); //end of integration
484 Int_t nx = TMath::Max(Int_t(fNdiv*TMath::Exp(-(x1-xch)*(x1-xch)/(2*fOrigSigmaX*fOrigSigmaX))*
485 TMath::Exp(-(y1-fCurrentY)*(y1-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),2);
486 Double_t dx = TMath::Min(fOrigSigmaX/Double_t(nx),x2-x1)/5.; //on the border more iteration
490 for (Double_t x = x1; x<x2+kprec ;){
492 nx = TMath::Max(Int_t(fNdiv*TMath::Exp(-(x-xch)*(x-xch)/(2*fOrigSigmaX*fOrigSigmaX))),3);
493 ndx = fOrigSigmaX/Double_t(nx);
497 if ( ( (x+dx+ndx)<TMath::Max(xp3,xp1)) || ( (x+dx+ndx)>TMath::Min(xp4,xp2))) {
500 if (ndx<kprec) ndx=2*kprec;
501 //INTEGRAL APROXIMATION
502 Double_t ddx,ddy,dddx,dddy;
504 ddy = fCurrentY-(y+dy/2.);
505 dddx = cos*ddx-sin*ddy;
506 dddy = sin*ddx+cos*ddy;
507 Double_t z0=fGRF->Eval(dddx,dddy); //middle point
511 dddx = cos*ddx-sin*ddy;
512 dddy = sin*ddx+cos*ddy;
513 Double_t z1=fGRF->Eval(dddx,dddy); //point down
516 ddy = fCurrentY-(y+dy);
517 dddx = cos*ddx-sin*ddy;
518 dddy = sin*ddx+cos*ddy;
519 Double_t z3=fGRF->Eval(dddx,dddy); //point up
522 ddy = fCurrentY-(y+dy/2.);
523 dddx = cos*ddx-sin*ddy;
524 dddy = sin*ddx+cos*ddy;
525 Double_t z2=fGRF->Eval(dddx,dddy); //point left
528 ddy = fCurrentY-(y+dy/2.);
529 dddx = cos*ddx-sin*ddy;
530 dddy = sin*ddx+cos*ddy;
531 Double_t z4=fGRF->Eval(dddx,dddy); //point right
534 if (z0<0) {z0=0;z1=0;z2=0;z3=0;z4=0;}
536 Double_t f2x= (z3+z1-2*z0)*4.;//second derivation in y
537 Double_t f2y= (z2+z4-2*z0)*4.;//second derivation in x
538 Double_t f1y= (z3-z1);
540 z = (z0+f2x/6.+f2y/6.);//second order aproxiation of integral
541 if (kx>kprec){ //positive derivation
542 if (x<(xp1+dy*kx)){ //calculate volume at left border
544 Double_t xx2 = TMath::Min(x+dx,xp1+dy*kx);
545 Double_t yy1 = y+(xx1-xp1)/kx;
546 Double_t yy2 = TMath::Min(y+(xx2-xp1)/kx,y+dy);
549 z-= z0*(y+dy-yy2)/dy; //constant part rectangle
550 z-= f1y*(xx2-xx1)*(y+dy-yy2)*(y+dy-yy2)/(2.*dx*dy);
552 z-=z0*(xx2-xx1)*(yy2-yy1)/(2*dx*dy); //constant part rectangle
555 if (x>xp2){ //calculate volume at right border
558 Double_t yy1 = y+(xx1-xp2)/kx;
559 Double_t yy2 = y+(xx2-xp2)/kx;
562 z-=z0*(yy1-y)/dy; //constant part
563 z-=f1y*(xx2-xx1)*(yy1-y)*(yy1-y)/(2*dx*dy);
565 z-=z0*(xx2-xx1)*(yy2-yy1)/(2*dx*dy); //constant part
568 if (kx<-kprec){ //negative derivation
569 if (x<(xp1+dy*kx)){ //calculate volume at left border
571 Double_t xx2 = TMath::Min(x+dx,xp3-dy/kx);
572 Double_t yy1 = y+(xx1-xp1)/kx;
573 Double_t yy2 = TMath::Max(y,yy1+(xx2-xx1)/kx); //yy2<yy1
575 z-= z0*(yy2-y)/dy; // constant part rectangle
576 z-= f1y*(xx2-xx1)*(yy2-y)*(yy2-y)/(2.*dx*dy);
577 z-=z0*(xx2-xx1)*(yy1-yy2)/(2*dx*dy); //constant part triangle
579 if (x>xp2){ //calculate volume at right border
580 Double_t xx1 = TMath::Max(x,xp2+dy*kx);
582 Double_t yy1 = TMath::Min(y+dy,y-(xp2-xx1)/kx);
583 Double_t yy2 = y-(xp2-xx2)/kx;
585 z-=z0*(yy2-y)/dy; //constant part rextangle
586 z-= f1y*(xx2-xx1)*(yy2-y)*(yy2-y)/(2.*dx*dy);
587 z-=z0*(xx2-xx1)*(yy1-yy2)/(2*dx*dy); //constant part triangle
591 if (z>0.) sumch+=fKNorm*z*dx*dy/fInteg;
600 }//step over different y
604 }//step over different points on line NPRF
607 void AliTPCPRF2D::UpdateSigma()
610 //calulate effective sigma X and sigma y of PRF
620 for (i=-1; i<=fNYdiv; i++){
621 if (fNYdiv == 1) y = fY1;
623 y = fY1+Float_t(i)*(fY2-fY1)/Float_t(fNYdiv-1);
624 for (x =-fNPRF*fDStep; x<fNPRF*fDStep;x+=fDStep)
626 //x in cm fWidth in cm
627 Float_t weight = GetPRF(x,y);
638 fSigmaX = TMath::Sqrt(fSigmaX/sum-fMeanX*fMeanX);
639 fSigmaY = TMath::Sqrt(fSigmaY/sum-fMeanY*fMeanY);
645 void AliTPCPRF2D::Streamer(TBuffer &R__b)
647 // Stream an object of class AliTPCPRF2D
649 if (R__b.IsReading()) {
651 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
652 AliTPCPRF2D::Class()->ReadBuffer(R__b, this, R__v, R__s, R__c);
654 if (strncmp(fType,"User",3)!=0){
656 if (strncmp(fType,"Gauss",3)==0)
657 fGRF = new TF2("funGauss2D",funGauss2D,-5.,5.,-5.,5.,4);
658 if (strncmp(fType,"Cosh",3)==0)
659 fGRF = new TF2("funCosh2D",funCosh2D,-5.,5.,-5.,5.,4);
660 if (strncmp(fType,"Gati",3)==0)
661 fGRF = new TF2("funGati2D",funGati2D,-5.,5.,-5.,5.,5);
662 if (fGRF!=0) fGRF->SetParameters(funParam);
664 //calculate conversion coefitient to convert position to virtual wire
665 fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1);
668 AliTPCPRF2D::Class()->WriteBuffer(R__b,this);
673 TH1F * AliTPCPRF2D::GenerDrawXHisto(Float_t x1, Float_t x2,Float_t y)
675 //gener one dimensional hist of pad response function
679 sprintf(s,"Pad Response Function");
680 TH1F * hPRFc = new TH1F("hPRFc",s,kn+1,x1,x2);
684 for (Int_t i = 0;i<kn+1;i++)
686 x+=(x2-x1)/Float_t(kn);
690 hPRFc->SetXTitle("pad (cm)");
694 AliH2F * AliTPCPRF2D::GenerDrawHisto(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx, Int_t Ny)
697 //gener two dimensional histogram with PRF
700 sprintf(s,"Pad Response Function");
701 AliH2F * hPRFc = new AliH2F("hPRFc",s,Nx,x1,x2,Ny,y1,y2);
702 Float_t dx=(x2-x1)/Float_t(Nx);
703 Float_t dy=(y2-y1)/Float_t(Ny) ;
707 for ( Int_t i = 0;i<=Nx;i++,x+=dx){
709 for (Int_t j = 0;j<=Ny;j++,y+=dy){
711 hPRFc->SetCellContent(i,j,z);
714 hPRFc->SetXTitle("pad direction (cm)");
715 hPRFc->SetYTitle("pad row direction (cm)");
716 hPRFc->SetTitleOffset(1.5,"X");
717 hPRFc->SetTitleOffset(1.5,"Y");
722 AliH2F * AliTPCPRF2D::GenerDrawDistHisto(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx, Int_t Ny, Float_t thr)
724 //return histogram with distortion
725 const Float_t kminth=0.00001;
726 if (thr<kminth) thr=kminth;
728 sprintf(s,"COG distortion of PRF (threshold=%2.2f)",thr);
729 AliH2F * hPRFDist = new AliH2F("hDistortion",s,Nx,x1,x2,Ny,y1,y2);
730 Float_t dx=(x2-x1)/Float_t(Nx);
731 Float_t dy=(y2-y1)/Float_t(Ny) ;
734 for ( Int_t i = 0;i<=Nx;i++,x+=dx){
736 for(Int_t j = 0;j<=Ny;j++,y+=dy)
740 for (Int_t k=-3;k<=3;k++)
742 Float_t padx=Float_t(k)*fWidth;
743 z = GetPRF(x-padx,y);
751 ddx = (x-(sumx/sum));
754 if (TMath::Abs(ddx)<10) hPRFDist->SetCellContent(i,j,ddx);
758 hPRFDist->SetXTitle("pad direction (cm)");
759 hPRFDist->SetYTitle("pad row direction (cm)");
760 hPRFDist->SetTitleOffset(1.5,"X");
761 hPRFDist->SetTitleOffset(1.5,"Y");
769 void AliTPCPRF2D::DrawX(Float_t x1 ,Float_t x2,Float_t y1,Float_t y2, Int_t N)
772 //draw pad response function at interval <x1,x2> at given y position
775 TCanvas * c1 = new TCanvas("PRFX","Pad response function",700,900);
778 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
779 comment->SetTextAlign(12);
780 comment->SetFillColor(42);
781 DrawComment(comment);
785 TPad * pad2 = new TPad("pPRF","",0.05,0.22,0.95,0.95);
786 pad2->Divide(2,(N+1)/2);
788 gStyle->SetOptFit(1);
789 gStyle->SetOptStat(1);
790 for (Int_t i=0;i<N;i++){
794 else y = y1+i*(y2-y1)/Float_t(N-1);
796 TH1F * hPRFc =GenerDrawXHisto(x1, x2,y);
797 sprintf(ch,"PRF at wire position: %2.3f",y);
799 sprintf(ch,"PRF %d",i);
808 void AliTPCPRF2D::DrawPRF(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2, Int_t Nx, Int_t Ny)
812 TCanvas * c1 = new TCanvas("canPRF","Pad response function",700,900);
814 TPad * pad2 = new TPad("pad2PRF","",0.05,0.22,0.95,0.95);
816 gStyle->SetOptFit(1);
817 gStyle->SetOptStat(1);
818 TH2F * hPRFc = GenerDrawHisto(x1, x2, y1, y2, Nx,Ny);
822 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
823 comment->SetTextAlign(12);
824 comment->SetFillColor(42);
825 DrawComment(comment);
829 void AliTPCPRF2D::DrawDist(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2, Int_t Nx, Int_t Ny, Float_t thr)
832 //draw distortion of the COG method - for different threshold parameter
833 TCanvas * c1 = new TCanvas("padDistortion","COG distortion",700,900);
835 TPad * pad1 = new TPad("dist","",0.05,0.55,0.95,0.95,21);
837 TPad * pad2 = new TPad("dist","",0.05,0.22,0.95,0.53,21);
839 gStyle->SetOptFit(1);
840 gStyle->SetOptStat(0);
842 AliH2F * hPRFDist = GenerDrawDistHisto(x1, x2, y1, y2, Nx,Ny,thr);
845 hPRFDist->Draw("surf");
846 Float_t distmax =hPRFDist->GetMaximum();
847 Float_t distmin =hPRFDist->GetMinimum();
848 gStyle->SetOptStat(1);
850 TH1F * dist = hPRFDist->GetAmplitudes(distmin,distmax,distmin-1);
854 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
855 comment->SetTextAlign(12);
856 comment->SetFillColor(42);
857 DrawComment(comment);
861 void AliTPCPRF2D::DrawComment(TPaveText *comment)
864 //function to write comment to picture
867 //draw comments to picture
868 TText * title = comment->AddText("Pad Response Function parameters:");
869 title->SetTextSize(0.03);
870 sprintf(s,"Height of pad: %2.2f cm",fHeightFull);
872 sprintf(s,"Width pad: %2.2f cm",fWidth);
874 sprintf(s,"Pad Angle: %2.2f ",fPadAngle);
877 if (TMath::Abs(fK)>0.0001){
878 sprintf(s,"Height of one chevron unit h: %2.2f cm",2*fHeightS);
880 sprintf(s,"Overlap factor: %2.2f",fK);
884 if (strncmp(fType,"User",3)==0){
885 sprintf(s,"Charge distribution - user defined function %s ",fGRF->GetTitle());
887 sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
889 sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
892 if (strncmp(fType,"Gauss",3)==0){
893 sprintf(s,"Gauss charge distribution");
895 sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
897 sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
900 if (strncmp(fType,"Gati",3)==0){
901 sprintf(s,"Gati charge distribution");
903 sprintf(s,"K3X of Gati : %2.2f ",fK3X);
905 sprintf(s,"K3Y of Gati: %2.2f ",fK3Y);
907 sprintf(s,"Wire to Pad Distance: %2.2f ",fPadDistance);
910 if (strncmp(fType,"Cosh",3)==0){
911 sprintf(s,"Cosh charge distribution");
913 sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
915 sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
918 sprintf(s,"Normalisation: %2.2f ",fKNorm);