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);
181 //make interpolation from more fore lines
184 if ((i<0) || (i>=fNYdiv) ) return 0;
190 fcharge =&(fChargeArray[(i-1)*fNPRF]);
191 z0 = GetPRFActiv(xin);
193 fcharge =&(fChargeArray[i*fNPRF]);
196 fcharge =&(fChargeArray[(i+1)*fNPRF]);
197 z2 = GetPRFActiv(xin);
200 fcharge =&(fChargeArray[(i+2)*fNPRF]);
201 z3 = GetPRFActiv(xin);
210 Float_t dy=y-Float_t(i);
212 res = a+b*dy+c*dy*dy+d*dy*dy*dy;
219 Float_t AliTPCPRF2D::GetPRFActiv(Float_t xin)
221 //GEt response function on given charege line
222 //return spline aproximaton
223 Float_t x = (xin*fDStepM1)+fNPRF/2;
226 if ( (i>1) && ((i+2)<fNPRF)) {
229 b = (fcharge[i+1]-fcharge[i-1])*0.5;
230 k = fcharge[i+1]-a-b;
231 l = (fcharge[i+2]-fcharge[i])*0.5-b;
234 Float_t dx=x-Float_t(i);
235 Float_t res = a+b*dx+c*dx*dx+d*dx*dx*dx;
242 Float_t AliTPCPRF2D::GetGRF(Float_t xin, Float_t yin)
244 //function which returnoriginal charge distribution
245 //this function is just normalised for fKnorm
247 return fKNorm*GetGRF()->Eval(xin,yin)/fInteg;
253 void AliTPCPRF2D::SetParam( TF2 * GRF, Float_t kNorm,
254 Float_t sigmaX, Float_t sigmaY)
256 //adjust parameters of the original charge distribution
257 //and pad size parameters
258 if (fGRF !=0 ) fGRF->Delete();
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;
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);
275 Double_t width = fWidth*(1+TMath::Abs(fK))+TMath::Abs(TMath::Tan(fPadAngle*fgkDegtoRad))*fHeightFull;
276 fNPRF = Int_t((width+8.*sigmaX)/fDStep);
282 void AliTPCPRF2D::SetGauss(Float_t sigmaX, Float_t sigmaY,
286 // set parameters for Gauss generic charge distribution
291 sprintf(fType,"Gauss");
292 if (fGRF !=0 ) fGRF->Delete();
293 fGRF = new TF2("fun",funGauss2D,-5.,5.,-5.,5.,4);
298 funParam[3]=fHeightS;
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);
310 Double_t width = fWidth*(1+TMath::Abs(fK))+TMath::Abs(TMath::Tan(fPadAngle*fgkDegtoRad))*fHeightFull;
311 fNPRF = Int_t((width+8.*sigmaX)/fDStep);
316 void AliTPCPRF2D::SetCosh(Float_t sigmaX, Float_t sigmaY,
319 // set parameters for Cosh generic charge distribution
324 sprintf(fType,"Cosh");
325 if (fGRF !=0 ) fGRF->Delete();
326 fGRF = new TF2("fun", funCosh2D,-5.,5.,-5.,5.,4);
330 funParam[3]=fHeightS;
331 fGRF->SetParameters(funParam);
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);
340 fNPRF = Int_t((1.2*fWidth*(1+TMath::Abs(fK))+8.*sigmaX)/fDStep);
345 void AliTPCPRF2D::SetGati(Float_t K3X, Float_t K3Y,
349 // set parameters for Gati generic charge distribution
354 fPadDistance=padDistance;
355 sprintf(fType,"Gati");
356 if (fGRF !=0 ) fGRF->Delete();
357 fGRF = new TF2("fun", funGati2D,-5.,5.,-5.,5.,5);
359 funParam[0]=padDistance;
362 funParam[3]=fHeightS;
364 fGRF->SetParameters(funParam);
365 fOrigSigmaX=padDistance;
366 fOrigSigmaY=padDistance;
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);
375 fNPRF = Int_t((1.2*fWidth*(1+TMath::Abs(fK))+8.*sigmaX)/fDStep);
381 void AliTPCPRF2D::Update()
384 //update fields with interpolated values for
387 if ( fGRF == 0 ) return;
388 //initialize interpolated values to 0
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;
394 //firstly calculate total integral of charge
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);
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 /////////////////////////////////////////////////////
413 if ( fInteg == 0 ) fInteg = 1;
415 for (i=0; i<fNYdiv; i++){
416 if (fNYdiv == 1) fCurrentY = fY1;
418 fCurrentY = fY1+Double_t(i)*(fY2-fY1)/Double_t(fNYdiv-1);
419 fcharge = &(fChargeArray[i*fNPRF]);
422 //calculate conversion coefitient to convert position to virtual wire
423 fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1);
428 void AliTPCPRF2D::Update1()
431 //update fields with interpolated values for
432 //PRF calculation for given charge line
434 Double_t cos = TMath::Cos(fChargeAngle);
435 Double_t sin = TMath::Sin(fChargeAngle);
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
441 Double_t xch = fDStep * (Double_t)(i-fNPRF/2);
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);
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);
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);
461 //loop over different y strips with variable step size dy
462 if (y2>(y1+kprec)) for (Double_t y = y1; y<y2+kprec;){
465 ny = TMath::Max(Int_t(fNdiv*TMath::Exp(-(y-fCurrentY)*(y-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),5);
466 ndy = fOrigSigmaY/Double_t(ny);
469 if (ndy<kprec) ndy=2*kprec; //calculate new delta y
473 //calculation of x borders and initial step
474 Double_t deltay = (y-y1chev);
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; //..
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
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
493 for (Double_t x = x1; x<x2+kprec ;){
495 nx = TMath::Max(Int_t(fNdiv*TMath::Exp(-(x-xch)*(x-xch)/(2*fOrigSigmaX*fOrigSigmaX))),3);
496 ndx = fOrigSigmaX/Double_t(nx);
500 if ( ( (x+dx+ndx)<TMath::Max(xp3,xp1)) || ( (x+dx+ndx)>TMath::Min(xp4,xp2))) {
503 if (ndx<kprec) ndx=2*kprec;
504 //INTEGRAL APROXIMATION
505 Double_t ddx,ddy,dddx,dddy;
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
514 dddx = cos*ddx-sin*ddy;
515 dddy = sin*ddx+cos*ddy;
516 Double_t z1=fGRF->Eval(dddx,dddy); //point down
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
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
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
537 if (z0<0) {z0=0;z1=0;z2=0;z3=0;z4=0;}
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);
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
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);
552 z-= z0*(y+dy-yy2)/dy; //constant part rectangle
553 z-= f1y*(xx2-xx1)*(y+dy-yy2)*(y+dy-yy2)/(2.*dx*dy);
555 z-=z0*(xx2-xx1)*(yy2-yy1)/(2*dx*dy); //constant part rectangle
558 if (x>xp2){ //calculate volume at right border
561 Double_t yy1 = y+(xx1-xp2)/kx;
562 Double_t yy2 = y+(xx2-xp2)/kx;
565 z-=z0*(yy1-y)/dy; //constant part
566 z-=f1y*(xx2-xx1)*(yy1-y)*(yy1-y)/(2*dx*dy);
568 z-=z0*(xx2-xx1)*(yy2-yy1)/(2*dx*dy); //constant part
571 if (kx<-kprec){ //negative derivation
572 if (x<(xp1+dy*kx)){ //calculate volume at left border
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
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
582 if (x>xp2){ //calculate volume at right border
583 Double_t xx1 = TMath::Max(x,xp2+dy*kx);
585 Double_t yy1 = TMath::Min(y+dy,y-(xp2-xx1)/kx);
586 Double_t yy2 = y-(xp2-xx2)/kx;
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
594 if (z>0.) sumch+=fKNorm*z*dx*dy/fInteg;
603 }//step over different y
607 }//step over different points on line NPRF
610 void AliTPCPRF2D::UpdateSigma()
613 //calulate effective sigma X and sigma y of PRF
623 for (i=-1; i<=fNYdiv; i++){
624 if (fNYdiv == 1) y = fY1;
626 y = fY1+Float_t(i)*(fY2-fY1)/Float_t(fNYdiv-1);
627 for (x =-fNPRF*fDStep; x<fNPRF*fDStep;x+=fDStep)
629 //x in cm fWidth in cm
630 Float_t weight = GetPRF(x,y);
641 fSigmaX = TMath::Sqrt(fSigmaX/sum-fMeanX*fMeanX);
642 fSigmaY = TMath::Sqrt(fSigmaY/sum-fMeanY*fMeanY);
648 void AliTPCPRF2D::Streamer(TBuffer &R__b)
650 // Stream an object of class AliTPCPRF2D
652 if (R__b.IsReading()) {
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);
657 if (strncmp(fType,"User",3)!=0){
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);
667 //calculate conversion coefitient to convert position to virtual wire
668 fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1);
671 AliTPCPRF2D::Class()->WriteBuffer(R__b,this);
676 TH1F * AliTPCPRF2D::GenerDrawXHisto(Float_t x1, Float_t x2,Float_t y)
678 //gener one dimensional hist of pad response function
682 sprintf(s,"Pad Response Function");
683 TH1F * hPRFc = new TH1F("hPRFc",s,kn+1,x1,x2);
687 for (Int_t i = 0;i<kn+1;i++)
689 x+=(x2-x1)/Float_t(kn);
693 hPRFc->SetXTitle("pad (cm)");
697 AliH2F * AliTPCPRF2D::GenerDrawHisto(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx, Int_t Ny)
700 //gener two dimensional histogram with PRF
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) ;
710 for ( Int_t i = 0;i<=Nx;i++,x+=dx){
712 for (Int_t j = 0;j<=Ny;j++,y+=dy){
714 hPRFc->SetCellContent(i,j,z);
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");
725 AliH2F * AliTPCPRF2D::GenerDrawDistHisto(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx, Int_t Ny, Float_t thr)
727 //return histogram with distortion
728 const Float_t kminth=0.00001;
729 if (thr<kminth) thr=kminth;
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) ;
737 for ( Int_t i = 0;i<=Nx;i++,x+=dx){
739 for(Int_t j = 0;j<=Ny;j++,y+=dy)
743 for (Int_t k=-3;k<=3;k++)
745 Float_t padx=Float_t(k)*fWidth;
746 z = GetPRF(x-padx,y);
754 ddx = (x-(sumx/sum));
757 if (TMath::Abs(ddx)<10) hPRFDist->SetCellContent(i,j,ddx);
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");
772 void AliTPCPRF2D::DrawX(Float_t x1 ,Float_t x2,Float_t y1,Float_t y2, Int_t N)
775 //draw pad response function at interval <x1,x2> at given y position
778 TCanvas * c1 = new TCanvas("PRFX","Pad response function",700,900);
781 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
782 comment->SetTextAlign(12);
783 comment->SetFillColor(42);
784 DrawComment(comment);
788 TPad * pad2 = new TPad("pPRF","",0.05,0.22,0.95,0.95);
789 pad2->Divide(2,(N+1)/2);
791 gStyle->SetOptFit(1);
792 gStyle->SetOptStat(1);
793 for (Int_t i=0;i<N;i++){
797 else y = y1+i*(y2-y1)/Float_t(N-1);
799 TH1F * hPRFc =GenerDrawXHisto(x1, x2,y);
800 sprintf(ch,"PRF at wire position: %2.3f",y);
802 sprintf(ch,"PRF %d",i);
811 void AliTPCPRF2D::DrawPRF(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2, Int_t Nx, Int_t Ny)
815 TCanvas * c1 = new TCanvas("canPRF","Pad response function",700,900);
817 TPad * pad2 = new TPad("pad2PRF","",0.05,0.22,0.95,0.95);
819 gStyle->SetOptFit(1);
820 gStyle->SetOptStat(1);
821 TH2F * hPRFc = GenerDrawHisto(x1, x2, y1, y2, Nx,Ny);
825 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
826 comment->SetTextAlign(12);
827 comment->SetFillColor(42);
828 DrawComment(comment);
832 void AliTPCPRF2D::DrawDist(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2, Int_t Nx, Int_t Ny, Float_t thr)
835 //draw distortion of the COG method - for different threshold parameter
836 TCanvas * c1 = new TCanvas("padDistortion","COG distortion",700,900);
838 TPad * pad1 = new TPad("dist","",0.05,0.55,0.95,0.95,21);
840 TPad * pad2 = new TPad("dist","",0.05,0.22,0.95,0.53,21);
842 gStyle->SetOptFit(1);
843 gStyle->SetOptStat(0);
845 AliH2F * hPRFDist = GenerDrawDistHisto(x1, x2, y1, y2, Nx,Ny,thr);
848 hPRFDist->Draw("surf");
849 Float_t distmax =hPRFDist->GetMaximum();
850 Float_t distmin =hPRFDist->GetMinimum();
851 gStyle->SetOptStat(1);
853 TH1F * dist = hPRFDist->GetAmplitudes(distmin,distmax,distmin-1);
857 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
858 comment->SetTextAlign(12);
859 comment->SetFillColor(42);
860 DrawComment(comment);
864 void AliTPCPRF2D::DrawComment(TPaveText *comment)
867 //function to write comment to picture
870 //draw comments to picture
871 TText * title = comment->AddText("Pad Response Function parameters:");
872 title->SetTextSize(0.03);
873 sprintf(s,"Height of pad: %2.2f cm",fHeightFull);
875 sprintf(s,"Width pad: %2.2f cm",fWidth);
877 sprintf(s,"Pad Angle: %2.2f ",fPadAngle);
880 if (TMath::Abs(fK)>0.0001){
881 sprintf(s,"Height of one chevron unit h: %2.2f cm",2*fHeightS);
883 sprintf(s,"Overlap factor: %2.2f",fK);
887 if (strncmp(fType,"User",3)==0){
888 sprintf(s,"Charge distribution - user defined function %s ",fGRF->GetTitle());
890 sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
892 sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
895 if (strncmp(fType,"Gauss",3)==0){
896 sprintf(s,"Gauss charge distribution");
898 sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
900 sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
903 if (strncmp(fType,"Gati",3)==0){
904 sprintf(s,"Gati charge distribution");
906 sprintf(s,"K3X of Gati : %2.2f ",fK3X);
908 sprintf(s,"K3Y of Gati: %2.2f ",fK3Y);
910 sprintf(s,"Wire to Pad Distance: %2.2f ",fPadDistance);
913 if (strncmp(fType,"Cosh",3)==0){
914 sprintf(s,"Cosh charge distribution");
916 sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
918 sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
921 sprintf(s,"Normalisation: %2.2f ",fKNorm);