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 Revision 1.6 2000/09/07 11:23:27 kowal2
19 Improved algoritms, coding convensions applied.
21 Revision 1.5 2000/06/30 12:07:50 kowal2
22 Updated from the TPC-PreRelease branch
24 Revision 1.4.4.3 2000/06/26 07:39:42 kowal2
25 Changes to obey the coding rules
27 Revision 1.4.4.2 2000/06/25 08:38:41 kowal2
28 Splitted from AliTPCtracking
30 Revision 1.4.4.1 2000/06/14 16:48:24 kowal2
31 Parameter setting improved. Removed compiler warnings
33 Revision 1.4 2000/04/17 09:37:33 kowal2
34 removed obsolete AliTPCDigitsDisplay.C
36 Revision 1.3.8.2 2000/04/10 08:40:46 kowal2
38 Small changes by M. Ivanov, improvements of algorithms
40 Revision 1.3.8.1 2000/04/10 07:56:53 kowal2
41 Not used anymore - removed
43 Revision 1.3 1999/10/05 17:15:46 fca
44 Minor syntax for the Alpha OSF
46 Revision 1.2 1999/09/29 09:24:34 fca
47 Introduction of the Copyright and cvs Log
51 ///////////////////////////////////////////////////////////////////////////////
53 // Pad response function object in two dimesions //
54 // This class contains the basic functions for the //
55 // calculation of PRF according generic charge distribution //
56 // In Update function object calculate table of response function //
57 // in discrete x and y position //
58 // This table is used for interpolation od response function in any position //
59 // (function GetPRF) //
61 // Origin: Marian Ivanov, Uni. of Bratislava, ivanov@fmph.uniba.sk //
63 ///////////////////////////////////////////////////////////////////////////////
67 #include "AliTPCPRF2D.h"
78 #include "TPaveText.h"
81 extern TStyle * gStyle;
83 const Double_t AliTPCPRF2D::fgkDegtoRad = 0.01745329251994;
84 const Double_t AliTPCPRF2D::fgkSQRT12=3.464101;
85 const Int_t AliTPCPRF2D::fgkNPRF = 100;
88 static Double_t funGauss2D(Double_t *x, Double_t * par)
90 //Gauss function -needde by the generic function object
91 return ( TMath::Exp(-(x[0]*x[0])/(2*par[0]*par[0]))*
92 TMath::Exp(-(x[1]*x[1])/(2*par[1]*par[1])));
96 static Double_t funCosh2D(Double_t *x, Double_t * par)
98 //Cosh function -needde by the generic function object
99 return ( 1/(TMath::CosH(3.14159*x[0]/(2*par[0]))*
100 TMath::CosH(3.14159*x[1]/(2*par[1]))));
103 static Double_t funGati2D(Double_t *x, Double_t * par)
105 //Gati function -needde by the generic function object
107 Float_t k3R=TMath::Sqrt(k3);
108 Float_t k2=(TMath::Pi()/2)*(1-k3R/2.);
109 Float_t k1=k2*k3R/(4*TMath::ATan(k3R));
110 Float_t l=x[0]/par[0];
111 Float_t tan2=TMath::TanH(k2*l);
113 Float_t res = k1*(1-tan2)/(1+k3*tan2);
114 //par[4] = is equal to k3Y
117 k2=(TMath::Pi()/2)*(1-k3R/2.);
118 k1=k2*k3R/(4*TMath::ATan(k3R));
120 tan2=TMath::TanH(k2*l);
122 res = res*k1*(1-tan2)/(1+k3*tan2);
126 ///////////////////////////////////////////////////////////////////////////
127 ///////////////////////////////////////////////////////////////////////////
129 ClassImp(AliTPCPRF2D)
131 AliTPCPRF2D::AliTPCPRF2D()
133 //default constructor for response function object
148 //chewron default values
150 SetChevron(0.2,0.0,1.0);
152 SetInterpolationType(2,0);
155 AliTPCPRF2D::~AliTPCPRF2D()
157 if (fChargeArray!=0) delete [] fChargeArray;
158 if (fGRF !=0 ) fGRF->Delete();
161 void AliTPCPRF2D::SetY(Float_t y1, Float_t y2, Int_t nYdiv)
164 //set virtual line position
165 //first and last line and number of lines
171 void AliTPCPRF2D::SetPad(Float_t width, Float_t height)
173 //set base chevron parameters
177 void AliTPCPRF2D::SetChevron(Float_t hstep,
181 //set shaping of chewron parameters
187 void AliTPCPRF2D::SetChParam(Float_t width, Float_t height,
188 Float_t hstep, Float_t shifty, Float_t fac)
190 SetPad(width,height);
191 SetChevron(hstep,shifty,fac);
195 Float_t AliTPCPRF2D::GetPRF(Float_t xin, Float_t yin)
197 //function which return pad response
198 //for the charge in distance xin
199 //return cubic aproximation of PRF or PRF at nearest virtual wire
200 if (fChargeArray==0) return 0;
201 //transform position to "wire position"
202 Float_t y=fDYtoWire*(yin-fY1);
203 if (fNYdiv == 1) y=fY1;
204 //normaly it find nearest line charge
206 Int_t i=Int_t(0.5+y);
207 if (y<0) i=Int_t(-0.5+y);
208 if ((i<0) || (i>=fNYdiv) ) return 0;
209 fcharge = &(fChargeArray[i*fNPRF]);
210 return GetPRFActiv(xin);
213 //make interpolation from more fore lines
216 if ((i<0) || (i>=fNYdiv) ) return 0;
222 fcharge =&(fChargeArray[(i-1)*fNPRF]);
223 z0 = GetPRFActiv(xin);
225 fcharge =&(fChargeArray[i*fNPRF]);
228 fcharge =&(fChargeArray[(i+1)*fNPRF]);
229 z2 = GetPRFActiv(xin);
232 fcharge =&(fChargeArray[(i+2)*fNPRF]);
233 z3 = GetPRFActiv(xin);
242 Float_t dy=y-Float_t(i);
244 res = a+b*dy+c*dy*dy+d*dy*dy*dy;
251 Float_t AliTPCPRF2D::GetPRFActiv(Float_t xin)
253 //GEt response function on given charege line
254 //return spline aproximaton
255 Float_t x = (xin*fDStepM1)+fNPRF/2;
258 if ( (i>1) && ((i+2)<fNPRF)) {
261 b = (fcharge[i+1]-fcharge[i-1])*0.5;
262 k = fcharge[i+1]-a-b;
263 l = (fcharge[i+2]-fcharge[i])*0.5-b;
266 Float_t dx=x-Float_t(i);
267 Float_t res = a+b*dx+c*dx*dx+d*dx*dx*dx;
274 Float_t AliTPCPRF2D::GetGRF(Float_t xin, Float_t yin)
276 //function which returnoriginal charge distribution
277 //this function is just normalised for fKnorm
279 return fKNorm*GetGRF()->Eval(xin,yin)/fInteg;
285 void AliTPCPRF2D::SetParam( TF2 * GRF, Float_t kNorm,
286 Float_t sigmaX, Float_t sigmaY)
288 //adjust parameters of the original charge distribution
289 //and pad size parameters
290 if (fGRF !=0 ) fGRF->Delete();
293 sprintf(fType,"User");
294 if (sigmaX ==0) sigmaX=(fWidth*(1+TMath::Abs(fK)))/fgkSQRT12;
295 if (sigmaY ==0) sigmaY=(fWidth*(1+TMath::Abs(fK)))/fgkSQRT12;
298 Double_t estimsigma =
299 TMath::Sqrt(sigmaX*sigmaX+(fWidth*fWidth*(1+TMath::Abs(fK))/12)+
300 TMath::Tan(fPadAngle*fgkDegtoRad)*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);
314 void AliTPCPRF2D::SetGauss(Float_t sigmaX, Float_t sigmaY,
318 // set parameters for Gauss generic charge distribution
323 sprintf(fType,"Gauss");
324 if (fGRF !=0 ) fGRF->Delete();
325 fGRF = new TF2("fun",funGauss2D,-5.,5.,-5.,5.,4);
330 funParam[3]=fHeightS;
332 fGRF->SetParameters(funParam);
333 Double_t estimsigma =
334 TMath::Sqrt(sigmaX*sigmaX+(fWidth*fWidth*(1+TMath::Abs(fK))/12)+
335 TMath::Tan(fPadAngle)*TMath::Tan(fPadAngle*fgkDegtoRad)*fHeightFull*fHeightFull/12);
336 if (estimsigma < 5*sigmaX) {
337 fDStep = estimsigma/10.;
338 fNPRF = Int_t(estimsigma*8./fDStep);
342 Double_t width = fWidth*(1+TMath::Abs(fK))+TMath::Abs(TMath::Tan(fPadAngle*fgkDegtoRad))*fHeightFull;
343 fNPRF = Int_t((width+8.*sigmaX)/fDStep);
348 void AliTPCPRF2D::SetCosh(Float_t sigmaX, Float_t sigmaY,
351 // set parameters for Cosh generic charge distribution
356 sprintf(fType,"Cosh");
357 if (fGRF !=0 ) fGRF->Delete();
358 fGRF = new TF2("fun", funCosh2D,-5.,5.,-5.,5.,4);
362 funParam[3]=fHeightS;
363 fGRF->SetParameters(funParam);
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);
377 void AliTPCPRF2D::SetGati(Float_t K3X, Float_t K3Y,
381 // set parameters for Gati generic charge distribution
386 fPadDistance=padDistance;
387 sprintf(fType,"Gati");
388 if (fGRF !=0 ) fGRF->Delete();
389 fGRF = new TF2("fun", funGati2D,-5.,5.,-5.,5.,5);
391 funParam[0]=padDistance;
394 funParam[3]=fHeightS;
396 fGRF->SetParameters(funParam);
397 fOrigSigmaX=padDistance;
398 fOrigSigmaY=padDistance;
399 Float_t sigmaX = fOrigSigmaX;
400 Double_t estimsigma = TMath::Sqrt(sigmaX*sigmaX+fWidth*fWidth*(1+TMath::Abs(fK))/12);
401 if (estimsigma < 5*sigmaX) {
402 fDStep = estimsigma/10.;
403 fNPRF = Int_t(estimsigma*8./fDStep);
407 fNPRF = Int_t((1.2*fWidth*(1+TMath::Abs(fK))+8.*sigmaX)/fDStep);
413 void AliTPCPRF2D::Update()
416 //update fields with interpolated values for
419 if ( fGRF == 0 ) return;
420 //initialize interpolated values to 0
422 if (fChargeArray!=0) delete [] fChargeArray;
423 fChargeArray = new Float_t[fNPRF*fNYdiv];
424 fNChargeArray = fNPRF*fNYdiv;
425 for (i =0; i<fNPRF*fNYdiv;i++) fChargeArray[i] = 0;
426 //firstly calculate total integral of charge
428 ////////////////////////////////////////////////////////
429 //I'm waiting for normal integral
430 //in this moment only sum
431 Float_t x2= 4*fOrigSigmaX;
432 Float_t y2= 4*fOrigSigmaY;
433 Float_t dx = fOrigSigmaX/Float_t(fNdiv*6);
434 Float_t dy = fOrigSigmaY/Float_t(fNdiv*6);
435 Int_t nx = Int_t(0.5+x2/dx);
436 Int_t ny = Int_t(0.5+y2/dy);
440 for (ix=-nx;ix<=nx;ix++)
441 for ( iy=-ny;iy<=ny;iy++)
442 dInteg+=fGRF->Eval(Float_t(ix)*dx,Float_t(iy)*dy)*dx*dy;
443 /////////////////////////////////////////////////////
445 if ( fInteg == 0 ) fInteg = 1;
447 for (i=0; i<fNYdiv; i++){
448 if (fNYdiv == 1) fCurrentY = fY1;
450 fCurrentY = fY1+Double_t(i)*(fY2-fY1)/Double_t(fNYdiv-1);
451 fcharge = &(fChargeArray[i*fNPRF]);
454 //calculate conversion coefitient to convert position to virtual wire
455 fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1);
460 void AliTPCPRF2D::Update1()
463 //update fields with interpolated values for
464 //PRF calculation for given charge line
466 Double_t cos = TMath::Cos(fChargeAngle);
467 Double_t sin = TMath::Sin(fChargeAngle);
468 const Double_t kprec =0.00000001;
469 //integrate charge over pad for different distance of pad
470 for (i =0; i<fNPRF;i++){
471 //x in cm fWidth in cm
473 Double_t xch = fDStep * (Double_t)(i-fNPRF/2);
478 for (Double_t ym=-fHeightFull/2.-fShiftY; ym<fHeightFull/2.-kprec;ym+=fHeightS){
479 Double_t y2chev=TMath::Min((ym+fHeightS),Double_t(fHeightFull/2.)); // end of chevron step
480 Double_t y1chev= ym; //beginning of chevron step
481 Double_t y2 = TMath::Min(y2chev,fCurrentY+3.5*fOrigSigmaY);
482 Double_t y1 = TMath::Max((y1chev),Double_t(-fHeightFull/2.));
483 y1 = TMath::Max(y1chev,fCurrentY-3.5*fOrigSigmaY);
485 Double_t x0 = fWidth*(-1.-(Double_t(k)*fK))*0.5+ym*TMath::Tan(fPadAngle*fgkDegtoRad);
486 Double_t kx = Double_t(k)*(fK*fWidth)/fHeightS;
487 kx = TMath::Tan(TMath::ATan(kx))+TMath::Tan(fPadAngle*fgkDegtoRad);
489 Int_t ny = TMath::Max(Int_t(fNdiv*TMath::Exp(-(y1-fCurrentY)*(y1-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),4);
490 Double_t dy = TMath::Min(fOrigSigmaY/Double_t(ny),y2-y1);
493 //loop over different y strips with variable step size dy
494 if (y2>(y1+kprec)) for (Double_t y = y1; y<y2+kprec;){
497 ny = TMath::Max(Int_t(fNdiv*TMath::Exp(-(y-fCurrentY)*(y-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),5);
498 ndy = fOrigSigmaY/Double_t(ny);
501 if (ndy<kprec) ndy=2*kprec; //calculate new delta y
505 //calculation of x borders and initial step
506 Double_t deltay = (y-y1chev);
508 Double_t xp1 = x0+deltay*kx;
509 //x begining of pad at position y
510 Double_t xp2 =xp1+fWidth; //x end of pad at position y
511 Double_t xp3 =xp1+kx*dy; //...at position y+dy
512 Double_t xp4 =xp2+kx*dy; //..
514 Double_t x1 = TMath::Min(xp1,xp3);
515 x1 = TMath::Max(xp1,xch-3.5*fOrigSigmaX); //beging of integration
516 Double_t x2 = TMath::Max(xp2,xp4);
517 x2 = TMath::Min(xp2+dy*kx,xch+3.5*fOrigSigmaX); //end of integration
519 Int_t nx = TMath::Max(Int_t(fNdiv*TMath::Exp(-(x1-xch)*(x1-xch)/(2*fOrigSigmaX*fOrigSigmaX))*
520 TMath::Exp(-(y1-fCurrentY)*(y1-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),2);
521 Double_t dx = TMath::Min(fOrigSigmaX/Double_t(nx),x2-x1)/5.; //on the border more iteration
525 for (Double_t x = x1; x<x2+kprec ;){
527 nx = TMath::Max(Int_t(fNdiv*TMath::Exp(-(x-xch)*(x-xch)/(2*fOrigSigmaX*fOrigSigmaX))),3);
528 ndx = fOrigSigmaX/Double_t(nx);
532 if ( ( (x+dx+ndx)<TMath::Max(xp3,xp1)) || ( (x+dx+ndx)>TMath::Min(xp4,xp2))) {
535 if (ndx<kprec) ndx=2*kprec;
536 //INTEGRAL APROXIMATION
537 Double_t ddx,ddy,dddx,dddy;
539 ddy = fCurrentY-(y+dy/2.);
540 dddx = cos*ddx-sin*ddy;
541 dddy = sin*ddx+cos*ddy;
542 Double_t z0=fGRF->Eval(dddx,dddy); //middle point
546 dddx = cos*ddx-sin*ddy;
547 dddy = sin*ddx+cos*ddy;
548 Double_t z1=fGRF->Eval(dddx,dddy); //point down
551 ddy = fCurrentY-(y+dy);
552 dddx = cos*ddx-sin*ddy;
553 dddy = sin*ddx+cos*ddy;
554 Double_t z3=fGRF->Eval(dddx,dddy); //point up
557 ddy = fCurrentY-(y+dy/2.);
558 dddx = cos*ddx-sin*ddy;
559 dddy = sin*ddx+cos*ddy;
560 Double_t z2=fGRF->Eval(dddx,dddy); //point left
563 ddy = fCurrentY-(y+dy/2.);
564 dddx = cos*ddx-sin*ddy;
565 dddy = sin*ddx+cos*ddy;
566 Double_t z4=fGRF->Eval(dddx,dddy); //point right
569 if (z0<0) {z0=0;z1=0;z2=0;z3=0;z4=0;}
571 Double_t f2x= (z3+z1-2*z0)*4.;//second derivation in y
572 Double_t f2y= (z2+z4-2*z0)*4.;//second derivation in x
573 Double_t f1y= (z3-z1);
575 z = (z0+f2x/6.+f2y/6.);//second order aproxiation of integral
576 if (kx>kprec){ //positive derivation
577 if (x<(xp1+dy*kx)){ //calculate volume at left border
579 Double_t xx2 = TMath::Min(x+dx,xp1+dy*kx);
580 Double_t yy1 = y+(xx1-xp1)/kx;
581 Double_t yy2 = TMath::Min(y+(xx2-xp1)/kx,y+dy);
584 z-= z0*(y+dy-yy2)/dy; //constant part rectangle
585 z-= f1y*(xx2-xx1)*(y+dy-yy2)*(y+dy-yy2)/(2.*dx*dy);
587 z-=z0*(xx2-xx1)*(yy2-yy1)/(2*dx*dy); //constant part rectangle
590 if (x>xp2){ //calculate volume at right border
593 Double_t yy1 = y+(xx1-xp2)/kx;
594 Double_t yy2 = y+(xx2-xp2)/kx;
597 z-=z0*(yy1-y)/dy; //constant part
598 z-=f1y*(xx2-xx1)*(yy1-y)*(yy1-y)/(2*dx*dy);
600 z-=z0*(xx2-xx1)*(yy2-yy1)/(2*dx*dy); //constant part
603 if (kx<-kprec){ //negative derivation
604 if (x<(xp1+dy*kx)){ //calculate volume at left border
606 Double_t xx2 = TMath::Min(x+dx,xp3-dy/kx);
607 Double_t yy1 = y+(xx1-xp1)/kx;
608 Double_t yy2 = TMath::Max(y,yy1+(xx2-xx1)/kx); //yy2<yy1
610 z-= z0*(yy2-y)/dy; // constant part rectangle
611 z-= f1y*(xx2-xx1)*(yy2-y)*(yy2-y)/(2.*dx*dy);
612 z-=z0*(xx2-xx1)*(yy1-yy2)/(2*dx*dy); //constant part triangle
614 if (x>xp2){ //calculate volume at right border
615 Double_t xx1 = TMath::Max(x,xp2+dy*kx);
617 Double_t yy1 = TMath::Min(y+dy,y-(xp2-xx1)/kx);
618 Double_t yy2 = y-(xp2-xx2)/kx;
620 z-=z0*(yy2-y)/dy; //constant part rextangle
621 z-= f1y*(xx2-xx1)*(yy2-y)*(yy2-y)/(2.*dx*dy);
622 z-=z0*(xx2-xx1)*(yy1-yy2)/(2*dx*dy); //constant part triangle
626 if (z>0.) sumch+=fKNorm*z*dx*dy/fInteg;
635 }//step over different y
639 }//step over different points on line NPRF
642 void AliTPCPRF2D::UpdateSigma()
645 //calulate effective sigma X and sigma y of PRF
655 for (i=-1; i<=fNYdiv; i++){
656 if (fNYdiv == 1) y = fY1;
658 y = fY1+Float_t(i)*(fY2-fY1)/Float_t(fNYdiv-1);
659 for (x =-fNPRF*fDStep; x<fNPRF*fDStep;x+=fDStep)
661 //x in cm fWidth in cm
662 Float_t weight = GetPRF(x,y);
673 fSigmaX = TMath::Sqrt(fSigmaX/sum-fMeanX*fMeanX);
674 fSigmaY = TMath::Sqrt(fSigmaY/sum-fMeanY*fMeanY);
680 void AliTPCPRF2D::Streamer(TBuffer &R__b)
682 // Stream an object of class AliTPCPRF2D
684 if (R__b.IsReading()) {
686 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
687 AliTPCPRF2D::Class()->ReadBuffer(R__b, this, R__v, R__s, R__c);
689 if (strncmp(fType,"User",3)!=0){
691 if (strncmp(fType,"Gauss",3)==0)
692 fGRF = new TF2("fun",funGauss2D,-5.,5.,-5.,5.,4);
693 if (strncmp(fType,"Cosh",3)==0)
694 fGRF = new TF2("fun",funCosh2D,-5.,5.,-5.,5.,4);
695 if (strncmp(fType,"Gati",3)==0)
696 fGRF = new TF2("fun",funGati2D,-5.,5.,-5.,5.,5);
697 if (fGRF!=0) fGRF->SetParameters(funParam);
699 //calculate conversion coefitient to convert position to virtual wire
700 fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1);
703 AliTPCPRF2D::Class()->WriteBuffer(R__b,this);
708 TH1F * AliTPCPRF2D::GenerDrawXHisto(Float_t x1, Float_t x2,Float_t y)
710 //gener one dimensional hist of pad response function
714 sprintf(s,"Pad Response Function");
715 TH1F * hPRFc = new TH1F("hPRFc",s,kn+1,x1,x2);
719 for (Int_t i = 0;i<kn+1;i++)
721 x+=(x2-x1)/Float_t(kn);
725 hPRFc->SetXTitle("pad (cm)");
729 AliH2F * AliTPCPRF2D::GenerDrawHisto(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx, Int_t Ny)
732 //gener two dimensional histogram with PRF
735 sprintf(s,"Pad Response Function");
736 AliH2F * hPRFc = new AliH2F("hPRFc",s,Nx,x1,x2,Ny,y1,y2);
737 Float_t dx=(x2-x1)/Float_t(Nx);
738 Float_t dy=(y2-y1)/Float_t(Ny) ;
742 for ( Int_t i = 0;i<=Nx;i++,x+=dx){
744 for (Int_t j = 0;j<=Ny;j++,y+=dy){
746 hPRFc->SetCellContent(i,j,z);
749 hPRFc->SetXTitle("pad direction (cm)");
750 hPRFc->SetYTitle("pad row direction (cm)");
751 hPRFc->SetTitleOffset(1.5,"X");
752 hPRFc->SetTitleOffset(1.5,"Y");
757 AliH2F * AliTPCPRF2D::GenerDrawDistHisto(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx, Int_t Ny, Float_t thr)
759 //return histogram with distortion
760 const Float_t kminth=0.00001;
761 if (thr<kminth) thr=kminth;
763 sprintf(s,"COG distortion of PRF (threshold=%2.2f)",thr);
764 AliH2F * hPRFDist = new AliH2F("hDistortion",s,Nx,x1,x2,Ny,y1,y2);
765 Float_t dx=(x2-x1)/Float_t(Nx);
766 Float_t dy=(y2-y1)/Float_t(Ny) ;
769 for ( Int_t i = 0;i<=Nx;i++,x+=dx){
771 for(Int_t j = 0;j<=Ny;j++,y+=dy)
775 for (Int_t k=-3;k<=3;k++)
777 Float_t padx=Float_t(k)*fWidth;
778 z = GetPRF(x-padx,y);
786 ddx = (x-(sumx/sum));
789 if (TMath::Abs(ddx)<10) hPRFDist->SetCellContent(i,j,ddx);
793 hPRFDist->SetXTitle("pad direction (cm)");
794 hPRFDist->SetYTitle("pad row direction (cm)");
795 hPRFDist->SetTitleOffset(1.5,"X");
796 hPRFDist->SetTitleOffset(1.5,"Y");
804 void AliTPCPRF2D::DrawX(Float_t x1 ,Float_t x2,Float_t y1,Float_t y2, Int_t N)
807 //draw pad response function at interval <x1,x2> at given y position
810 TCanvas * c1 = new TCanvas("PRFX","Pad response function",700,900);
813 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
814 comment->SetTextAlign(12);
815 comment->SetFillColor(42);
816 DrawComment(comment);
820 TPad * pad2 = new TPad("pPRF","",0.05,0.22,0.95,0.95);
821 pad2->Divide(2,(N+1)/2);
823 gStyle->SetOptFit(1);
824 gStyle->SetOptStat(1);
825 for (Int_t i=0;i<N;i++){
829 else y = y1+i*(y2-y1)/Float_t(N-1);
831 TH1F * hPRFc =GenerDrawXHisto(x1, x2,y);
832 sprintf(ch,"PRF at wire position: %2.3f",y);
834 sprintf(ch,"PRF %d",i);
843 void AliTPCPRF2D::DrawPRF(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2, Int_t Nx, Int_t Ny)
847 TCanvas * c1 = new TCanvas("canPRF","Pad response function",700,900);
849 TPad * pad2 = new TPad("pad2PRF","",0.05,0.22,0.95,0.95);
851 gStyle->SetOptFit(1);
852 gStyle->SetOptStat(1);
853 TH2F * hPRFc = GenerDrawHisto(x1, x2, y1, y2, Nx,Ny);
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::DrawDist(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2, Int_t Nx, Int_t Ny, Float_t thr)
867 //draw distortion of the COG method - for different threshold parameter
868 TCanvas * c1 = new TCanvas("padDistortion","COG distortion",700,900);
870 TPad * pad1 = new TPad("dist","",0.05,0.55,0.95,0.95,21);
872 TPad * pad2 = new TPad("dist","",0.05,0.22,0.95,0.53,21);
874 gStyle->SetOptFit(1);
875 gStyle->SetOptStat(0);
877 AliH2F * hPRFDist = GenerDrawDistHisto(x1, x2, y1, y2, Nx,Ny,thr);
880 hPRFDist->Draw("surf");
881 Float_t distmax =hPRFDist->GetMaximum();
882 Float_t distmin =hPRFDist->GetMinimum();
883 gStyle->SetOptStat(1);
885 TH1F * dist = hPRFDist->GetAmplitudes(distmin,distmax,distmin-1);
889 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
890 comment->SetTextAlign(12);
891 comment->SetFillColor(42);
892 DrawComment(comment);
896 void AliTPCPRF2D::DrawComment(TPaveText *comment)
899 //function to write comment to picture
902 //draw comments to picture
903 TText * title = comment->AddText("Pad Response Function parameters:");
904 title->SetTextSize(0.03);
905 sprintf(s,"Height of pad: %2.2f cm",fHeightFull);
907 sprintf(s,"Width pad: %2.2f cm",fWidth);
909 sprintf(s,"Pad Angle: %2.2f ",fPadAngle);
912 if (TMath::Abs(fK)>0.0001){
913 sprintf(s,"Height of one chevron unit h: %2.2f cm",2*fHeightS);
915 sprintf(s,"Overlap factor: %2.2f",fK);
919 if (strncmp(fType,"User",3)==0){
920 sprintf(s,"Charge distribution - user defined function %s ",fGRF->GetTitle());
922 sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
924 sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
927 if (strncmp(fType,"Gauss",3)==0){
928 sprintf(s,"Gauss charge distribution");
930 sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
932 sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
935 if (strncmp(fType,"Gati",3)==0){
936 sprintf(s,"Gati charge distribution");
938 sprintf(s,"K3X of Gati : %2.2f ",fK3X);
940 sprintf(s,"K3Y of Gati: %2.2f ",fK3Y);
942 sprintf(s,"Wire to Pad Distance: %2.2f ",fPadDistance);
945 if (strncmp(fType,"Cosh",3)==0){
946 sprintf(s,"Cosh charge distribution");
948 sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
950 sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
953 sprintf(s,"Normalisation: %2.2f ",fKNorm);