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.7 2001/01/30 09:23:15 hristov
19 Streamers removed (R.Brun)
21 Revision 1.6 2000/09/07 11:23:27 kowal2
22 Improved algoritms, coding convensions applied.
24 Revision 1.5 2000/06/30 12:07:50 kowal2
25 Updated from the TPC-PreRelease branch
27 Revision 1.4.4.3 2000/06/26 07:39:42 kowal2
28 Changes to obey the coding rules
30 Revision 1.4.4.2 2000/06/25 08:38:41 kowal2
31 Splitted from AliTPCtracking
33 Revision 1.4.4.1 2000/06/14 16:48:24 kowal2
34 Parameter setting improved. Removed compiler warnings
36 Revision 1.4 2000/04/17 09:37:33 kowal2
37 removed obsolete AliTPCDigitsDisplay.C
39 Revision 1.3.8.2 2000/04/10 08:40:46 kowal2
41 Small changes by M. Ivanov, improvements of algorithms
43 Revision 1.3.8.1 2000/04/10 07:56:53 kowal2
44 Not used anymore - removed
46 Revision 1.3 1999/10/05 17:15:46 fca
47 Minor syntax for the Alpha OSF
49 Revision 1.2 1999/09/29 09:24:34 fca
50 Introduction of the Copyright and cvs Log
54 ///////////////////////////////////////////////////////////////////////////////
56 // Pad response function object in two dimesions //
57 // This class contains the basic functions for the //
58 // calculation of PRF according generic charge distribution //
59 // In Update function object calculate table of response function //
60 // in discrete x and y position //
61 // This table is used for interpolation od response function in any position //
62 // (function GetPRF) //
64 // Origin: Marian Ivanov, Uni. of Bratislava, ivanov@fmph.uniba.sk //
66 ///////////////////////////////////////////////////////////////////////////////
70 #include "AliTPCPRF2D.h"
81 #include "TPaveText.h"
84 extern TStyle * gStyle;
86 const Double_t AliTPCPRF2D::fgkDegtoRad = 0.01745329251994;
87 const Double_t AliTPCPRF2D::fgkSQRT12=3.464101;
88 const Int_t AliTPCPRF2D::fgkNPRF = 100;
91 static Double_t funGauss2D(Double_t *x, Double_t * par)
93 //Gauss function -needde by the generic function object
94 return ( TMath::Exp(-(x[0]*x[0])/(2*par[0]*par[0]))*
95 TMath::Exp(-(x[1]*x[1])/(2*par[1]*par[1])));
99 static Double_t funCosh2D(Double_t *x, Double_t * par)
101 //Cosh function -needde by the generic function object
102 return ( 1/(TMath::CosH(3.14159*x[0]/(2*par[0]))*
103 TMath::CosH(3.14159*x[1]/(2*par[1]))));
106 static Double_t funGati2D(Double_t *x, Double_t * par)
108 //Gati function -needde by the generic function object
110 Float_t k3R=TMath::Sqrt(k3);
111 Float_t k2=(TMath::Pi()/2)*(1-k3R/2.);
112 Float_t k1=k2*k3R/(4*TMath::ATan(k3R));
113 Float_t l=x[0]/par[0];
114 Float_t tan2=TMath::TanH(k2*l);
116 Float_t res = k1*(1-tan2)/(1+k3*tan2);
117 //par[4] = is equal to k3Y
120 k2=(TMath::Pi()/2)*(1-k3R/2.);
121 k1=k2*k3R/(4*TMath::ATan(k3R));
123 tan2=TMath::TanH(k2*l);
125 res = res*k1*(1-tan2)/(1+k3*tan2);
129 ///////////////////////////////////////////////////////////////////////////
130 ///////////////////////////////////////////////////////////////////////////
132 ClassImp(AliTPCPRF2D)
134 AliTPCPRF2D::AliTPCPRF2D()
136 //default constructor for response function object
152 //chewron default values
154 SetChevron(0.2,0.0,1.0);
156 SetInterpolationType(2,0);
159 AliTPCPRF2D::~AliTPCPRF2D()
161 if (fChargeArray!=0) delete [] fChargeArray;
162 if (fGRF !=0 ) fGRF->Delete();
165 void AliTPCPRF2D::SetY(Float_t y1, Float_t y2, Int_t nYdiv)
168 //set virtual line position
169 //first and last line and number of lines
175 void AliTPCPRF2D::SetPad(Float_t width, Float_t height)
177 //set base chevron parameters
181 void AliTPCPRF2D::SetChevron(Float_t hstep,
185 //set shaping of chewron parameters
191 void AliTPCPRF2D::SetChParam(Float_t width, Float_t height,
192 Float_t hstep, Float_t shifty, Float_t fac)
194 SetPad(width,height);
195 SetChevron(hstep,shifty,fac);
199 Float_t AliTPCPRF2D::GetPRF(Float_t xin, Float_t yin)
201 //function which return pad response
202 //for the charge in distance xin
203 //return cubic aproximation of PRF or PRF at nearest virtual wire
204 if (fChargeArray==0) return 0;
205 //transform position to "wire position"
206 Float_t y=fDYtoWire*(yin-fY1);
207 if (fNYdiv == 1) y=fY1;
208 //normaly it find nearest line charge
210 Int_t i=Int_t(0.5+y);
211 if (y<0) i=Int_t(-0.5+y);
212 if ((i<0) || (i>=fNYdiv) ) return 0;
213 fcharge = &(fChargeArray[i*fNPRF]);
214 return GetPRFActiv(xin);
217 //make interpolation from more fore lines
220 if ((i<0) || (i>=fNYdiv) ) return 0;
226 fcharge =&(fChargeArray[(i-1)*fNPRF]);
227 z0 = GetPRFActiv(xin);
229 fcharge =&(fChargeArray[i*fNPRF]);
232 fcharge =&(fChargeArray[(i+1)*fNPRF]);
233 z2 = GetPRFActiv(xin);
236 fcharge =&(fChargeArray[(i+2)*fNPRF]);
237 z3 = GetPRFActiv(xin);
246 Float_t dy=y-Float_t(i);
248 res = a+b*dy+c*dy*dy+d*dy*dy*dy;
255 Float_t AliTPCPRF2D::GetPRFActiv(Float_t xin)
257 //GEt response function on given charege line
258 //return spline aproximaton
259 Float_t x = (xin*fDStepM1)+fNPRF/2;
262 if ( (i>1) && ((i+2)<fNPRF)) {
265 b = (fcharge[i+1]-fcharge[i-1])*0.5;
266 k = fcharge[i+1]-a-b;
267 l = (fcharge[i+2]-fcharge[i])*0.5-b;
270 Float_t dx=x-Float_t(i);
271 Float_t res = a+b*dx+c*dx*dx+d*dx*dx*dx;
278 Float_t AliTPCPRF2D::GetGRF(Float_t xin, Float_t yin)
280 //function which returnoriginal charge distribution
281 //this function is just normalised for fKnorm
283 return fKNorm*GetGRF()->Eval(xin,yin)/fInteg;
289 void AliTPCPRF2D::SetParam( TF2 * GRF, Float_t kNorm,
290 Float_t sigmaX, Float_t sigmaY)
292 //adjust parameters of the original charge distribution
293 //and pad size parameters
294 if (fGRF !=0 ) fGRF->Delete();
297 sprintf(fType,"User");
298 if (sigmaX ==0) sigmaX=(fWidth*(1+TMath::Abs(fK)))/fgkSQRT12;
299 if (sigmaY ==0) sigmaY=(fWidth*(1+TMath::Abs(fK)))/fgkSQRT12;
302 Double_t estimsigma =
303 TMath::Sqrt(sigmaX*sigmaX+(fWidth*fWidth*(1+TMath::Abs(fK))/12)+
304 TMath::Tan(fPadAngle*fgkDegtoRad)*TMath::Tan(fPadAngle*fgkDegtoRad)*fHeightFull*fHeightFull/12);
305 if (estimsigma < 5*sigmaX) {
306 fDStep = estimsigma/10.;
307 fNPRF = Int_t(estimsigma*8./fDStep);
311 Double_t width = fWidth*(1+TMath::Abs(fK))+TMath::Abs(TMath::Tan(fPadAngle*fgkDegtoRad))*fHeightFull;
312 fNPRF = Int_t((width+8.*sigmaX)/fDStep);
318 void AliTPCPRF2D::SetGauss(Float_t sigmaX, Float_t sigmaY,
322 // set parameters for Gauss generic charge distribution
327 sprintf(fType,"Gauss");
328 if (fGRF !=0 ) fGRF->Delete();
329 fGRF = new TF2("fun",funGauss2D,-5.,5.,-5.,5.,4);
334 funParam[3]=fHeightS;
336 fGRF->SetParameters(funParam);
337 Double_t estimsigma =
338 TMath::Sqrt(sigmaX*sigmaX+(fWidth*fWidth*(1+TMath::Abs(fK))/12)+
339 TMath::Tan(fPadAngle)*TMath::Tan(fPadAngle*fgkDegtoRad)*fHeightFull*fHeightFull/12);
340 if (estimsigma < 5*sigmaX) {
341 fDStep = estimsigma/10.;
342 fNPRF = Int_t(estimsigma*8./fDStep);
346 Double_t width = fWidth*(1+TMath::Abs(fK))+TMath::Abs(TMath::Tan(fPadAngle*fgkDegtoRad))*fHeightFull;
347 fNPRF = Int_t((width+8.*sigmaX)/fDStep);
352 void AliTPCPRF2D::SetCosh(Float_t sigmaX, Float_t sigmaY,
355 // set parameters for Cosh generic charge distribution
360 sprintf(fType,"Cosh");
361 if (fGRF !=0 ) fGRF->Delete();
362 fGRF = new TF2("fun", funCosh2D,-5.,5.,-5.,5.,4);
366 funParam[3]=fHeightS;
367 fGRF->SetParameters(funParam);
369 Double_t estimsigma = TMath::Sqrt(sigmaX*sigmaX+fWidth*fWidth*(1+TMath::Abs(fK))/12);
370 if (estimsigma < 5*sigmaX) {
371 fDStep = estimsigma/10.;
372 fNPRF = Int_t(estimsigma*8./fDStep);
376 fNPRF = Int_t((1.2*fWidth*(1+TMath::Abs(fK))+8.*sigmaX)/fDStep);
381 void AliTPCPRF2D::SetGati(Float_t K3X, Float_t K3Y,
385 // set parameters for Gati generic charge distribution
390 fPadDistance=padDistance;
391 sprintf(fType,"Gati");
392 if (fGRF !=0 ) fGRF->Delete();
393 fGRF = new TF2("fun", funGati2D,-5.,5.,-5.,5.,5);
395 funParam[0]=padDistance;
398 funParam[3]=fHeightS;
400 fGRF->SetParameters(funParam);
401 fOrigSigmaX=padDistance;
402 fOrigSigmaY=padDistance;
403 Float_t sigmaX = fOrigSigmaX;
404 Double_t estimsigma = TMath::Sqrt(sigmaX*sigmaX+fWidth*fWidth*(1+TMath::Abs(fK))/12);
405 if (estimsigma < 5*sigmaX) {
406 fDStep = estimsigma/10.;
407 fNPRF = Int_t(estimsigma*8./fDStep);
411 fNPRF = Int_t((1.2*fWidth*(1+TMath::Abs(fK))+8.*sigmaX)/fDStep);
417 void AliTPCPRF2D::Update()
420 //update fields with interpolated values for
423 if ( fGRF == 0 ) return;
424 //initialize interpolated values to 0
426 if (fChargeArray!=0) delete [] fChargeArray;
427 fChargeArray = new Float_t[fNPRF*fNYdiv];
428 fNChargeArray = fNPRF*fNYdiv;
429 for (i =0; i<fNPRF*fNYdiv;i++) fChargeArray[i] = 0;
430 //firstly calculate total integral of charge
432 ////////////////////////////////////////////////////////
433 //I'm waiting for normal integral
434 //in this moment only sum
435 Float_t x2= 4*fOrigSigmaX;
436 Float_t y2= 4*fOrigSigmaY;
437 Float_t dx = fOrigSigmaX/Float_t(fNdiv*6);
438 Float_t dy = fOrigSigmaY/Float_t(fNdiv*6);
439 Int_t nx = Int_t(0.5+x2/dx);
440 Int_t ny = Int_t(0.5+y2/dy);
444 for (ix=-nx;ix<=nx;ix++)
445 for ( iy=-ny;iy<=ny;iy++)
446 dInteg+=fGRF->Eval(Float_t(ix)*dx,Float_t(iy)*dy)*dx*dy;
447 /////////////////////////////////////////////////////
449 if ( fInteg == 0 ) fInteg = 1;
451 for (i=0; i<fNYdiv; i++){
452 if (fNYdiv == 1) fCurrentY = fY1;
454 fCurrentY = fY1+Double_t(i)*(fY2-fY1)/Double_t(fNYdiv-1);
455 fcharge = &(fChargeArray[i*fNPRF]);
458 //calculate conversion coefitient to convert position to virtual wire
459 fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1);
464 void AliTPCPRF2D::Update1()
467 //update fields with interpolated values for
468 //PRF calculation for given charge line
470 Double_t cos = TMath::Cos(fChargeAngle);
471 Double_t sin = TMath::Sin(fChargeAngle);
472 const Double_t kprec =0.00000001;
473 //integrate charge over pad for different distance of pad
474 for (i =0; i<fNPRF;i++){
475 //x in cm fWidth in cm
477 Double_t xch = fDStep * (Double_t)(i-fNPRF/2);
482 for (Double_t ym=-fHeightFull/2.-fShiftY; ym<fHeightFull/2.-kprec;ym+=fHeightS){
483 Double_t y2chev=TMath::Min((ym+fHeightS),Double_t(fHeightFull/2.)); // end of chevron step
484 Double_t y1chev= ym; //beginning of chevron step
485 Double_t y2 = TMath::Min(y2chev,fCurrentY+3.5*fOrigSigmaY);
486 Double_t y1 = TMath::Max((y1chev),Double_t(-fHeightFull/2.));
487 y1 = TMath::Max(y1chev,fCurrentY-3.5*fOrigSigmaY);
489 Double_t x0 = fWidth*(-1.-(Double_t(k)*fK))*0.5+ym*TMath::Tan(fPadAngle*fgkDegtoRad);
490 Double_t kx = Double_t(k)*(fK*fWidth)/fHeightS;
491 kx = TMath::Tan(TMath::ATan(kx))+TMath::Tan(fPadAngle*fgkDegtoRad);
493 Int_t ny = TMath::Max(Int_t(fNdiv*TMath::Exp(-(y1-fCurrentY)*(y1-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),4);
494 Double_t dy = TMath::Min(fOrigSigmaY/Double_t(ny),y2-y1);
497 //loop over different y strips with variable step size dy
498 if (y2>(y1+kprec)) for (Double_t y = y1; y<y2+kprec;){
501 ny = TMath::Max(Int_t(fNdiv*TMath::Exp(-(y-fCurrentY)*(y-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),5);
502 ndy = fOrigSigmaY/Double_t(ny);
505 if (ndy<kprec) ndy=2*kprec; //calculate new delta y
509 //calculation of x borders and initial step
510 Double_t deltay = (y-y1chev);
512 Double_t xp1 = x0+deltay*kx;
513 //x begining of pad at position y
514 Double_t xp2 =xp1+fWidth; //x end of pad at position y
515 Double_t xp3 =xp1+kx*dy; //...at position y+dy
516 Double_t xp4 =xp2+kx*dy; //..
518 Double_t x1 = TMath::Min(xp1,xp3);
519 x1 = TMath::Max(xp1,xch-3.5*fOrigSigmaX); //beging of integration
520 Double_t x2 = TMath::Max(xp2,xp4);
521 x2 = TMath::Min(xp2+dy*kx,xch+3.5*fOrigSigmaX); //end of integration
523 Int_t nx = TMath::Max(Int_t(fNdiv*TMath::Exp(-(x1-xch)*(x1-xch)/(2*fOrigSigmaX*fOrigSigmaX))*
524 TMath::Exp(-(y1-fCurrentY)*(y1-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),2);
525 Double_t dx = TMath::Min(fOrigSigmaX/Double_t(nx),x2-x1)/5.; //on the border more iteration
529 for (Double_t x = x1; x<x2+kprec ;){
531 nx = TMath::Max(Int_t(fNdiv*TMath::Exp(-(x-xch)*(x-xch)/(2*fOrigSigmaX*fOrigSigmaX))),3);
532 ndx = fOrigSigmaX/Double_t(nx);
536 if ( ( (x+dx+ndx)<TMath::Max(xp3,xp1)) || ( (x+dx+ndx)>TMath::Min(xp4,xp2))) {
539 if (ndx<kprec) ndx=2*kprec;
540 //INTEGRAL APROXIMATION
541 Double_t ddx,ddy,dddx,dddy;
543 ddy = fCurrentY-(y+dy/2.);
544 dddx = cos*ddx-sin*ddy;
545 dddy = sin*ddx+cos*ddy;
546 Double_t z0=fGRF->Eval(dddx,dddy); //middle point
550 dddx = cos*ddx-sin*ddy;
551 dddy = sin*ddx+cos*ddy;
552 Double_t z1=fGRF->Eval(dddx,dddy); //point down
555 ddy = fCurrentY-(y+dy);
556 dddx = cos*ddx-sin*ddy;
557 dddy = sin*ddx+cos*ddy;
558 Double_t z3=fGRF->Eval(dddx,dddy); //point up
561 ddy = fCurrentY-(y+dy/2.);
562 dddx = cos*ddx-sin*ddy;
563 dddy = sin*ddx+cos*ddy;
564 Double_t z2=fGRF->Eval(dddx,dddy); //point left
567 ddy = fCurrentY-(y+dy/2.);
568 dddx = cos*ddx-sin*ddy;
569 dddy = sin*ddx+cos*ddy;
570 Double_t z4=fGRF->Eval(dddx,dddy); //point right
573 if (z0<0) {z0=0;z1=0;z2=0;z3=0;z4=0;}
575 Double_t f2x= (z3+z1-2*z0)*4.;//second derivation in y
576 Double_t f2y= (z2+z4-2*z0)*4.;//second derivation in x
577 Double_t f1y= (z3-z1);
579 z = (z0+f2x/6.+f2y/6.);//second order aproxiation of integral
580 if (kx>kprec){ //positive derivation
581 if (x<(xp1+dy*kx)){ //calculate volume at left border
583 Double_t xx2 = TMath::Min(x+dx,xp1+dy*kx);
584 Double_t yy1 = y+(xx1-xp1)/kx;
585 Double_t yy2 = TMath::Min(y+(xx2-xp1)/kx,y+dy);
588 z-= z0*(y+dy-yy2)/dy; //constant part rectangle
589 z-= f1y*(xx2-xx1)*(y+dy-yy2)*(y+dy-yy2)/(2.*dx*dy);
591 z-=z0*(xx2-xx1)*(yy2-yy1)/(2*dx*dy); //constant part rectangle
594 if (x>xp2){ //calculate volume at right border
597 Double_t yy1 = y+(xx1-xp2)/kx;
598 Double_t yy2 = y+(xx2-xp2)/kx;
601 z-=z0*(yy1-y)/dy; //constant part
602 z-=f1y*(xx2-xx1)*(yy1-y)*(yy1-y)/(2*dx*dy);
604 z-=z0*(xx2-xx1)*(yy2-yy1)/(2*dx*dy); //constant part
607 if (kx<-kprec){ //negative derivation
608 if (x<(xp1+dy*kx)){ //calculate volume at left border
610 Double_t xx2 = TMath::Min(x+dx,xp3-dy/kx);
611 Double_t yy1 = y+(xx1-xp1)/kx;
612 Double_t yy2 = TMath::Max(y,yy1+(xx2-xx1)/kx); //yy2<yy1
614 z-= z0*(yy2-y)/dy; // constant part rectangle
615 z-= f1y*(xx2-xx1)*(yy2-y)*(yy2-y)/(2.*dx*dy);
616 z-=z0*(xx2-xx1)*(yy1-yy2)/(2*dx*dy); //constant part triangle
618 if (x>xp2){ //calculate volume at right border
619 Double_t xx1 = TMath::Max(x,xp2+dy*kx);
621 Double_t yy1 = TMath::Min(y+dy,y-(xp2-xx1)/kx);
622 Double_t yy2 = y-(xp2-xx2)/kx;
624 z-=z0*(yy2-y)/dy; //constant part rextangle
625 z-= f1y*(xx2-xx1)*(yy2-y)*(yy2-y)/(2.*dx*dy);
626 z-=z0*(xx2-xx1)*(yy1-yy2)/(2*dx*dy); //constant part triangle
630 if (z>0.) sumch+=fKNorm*z*dx*dy/fInteg;
639 }//step over different y
643 }//step over different points on line NPRF
646 void AliTPCPRF2D::UpdateSigma()
649 //calulate effective sigma X and sigma y of PRF
659 for (i=-1; i<=fNYdiv; i++){
660 if (fNYdiv == 1) y = fY1;
662 y = fY1+Float_t(i)*(fY2-fY1)/Float_t(fNYdiv-1);
663 for (x =-fNPRF*fDStep; x<fNPRF*fDStep;x+=fDStep)
665 //x in cm fWidth in cm
666 Float_t weight = GetPRF(x,y);
677 fSigmaX = TMath::Sqrt(fSigmaX/sum-fMeanX*fMeanX);
678 fSigmaY = TMath::Sqrt(fSigmaY/sum-fMeanY*fMeanY);
684 void AliTPCPRF2D::Streamer(TBuffer &R__b)
686 // Stream an object of class AliTPCPRF2D
688 if (R__b.IsReading()) {
690 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
691 AliTPCPRF2D::Class()->ReadBuffer(R__b, this, R__v, R__s, R__c);
693 if (strncmp(fType,"User",3)!=0){
695 if (strncmp(fType,"Gauss",3)==0)
696 fGRF = new TF2("fun",funGauss2D,-5.,5.,-5.,5.,4);
697 if (strncmp(fType,"Cosh",3)==0)
698 fGRF = new TF2("fun",funCosh2D,-5.,5.,-5.,5.,4);
699 if (strncmp(fType,"Gati",3)==0)
700 fGRF = new TF2("fun",funGati2D,-5.,5.,-5.,5.,5);
701 if (fGRF!=0) fGRF->SetParameters(funParam);
703 //calculate conversion coefitient to convert position to virtual wire
704 fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1);
707 AliTPCPRF2D::Class()->WriteBuffer(R__b,this);
712 TH1F * AliTPCPRF2D::GenerDrawXHisto(Float_t x1, Float_t x2,Float_t y)
714 //gener one dimensional hist of pad response function
718 sprintf(s,"Pad Response Function");
719 TH1F * hPRFc = new TH1F("hPRFc",s,kn+1,x1,x2);
723 for (Int_t i = 0;i<kn+1;i++)
725 x+=(x2-x1)/Float_t(kn);
729 hPRFc->SetXTitle("pad (cm)");
733 AliH2F * AliTPCPRF2D::GenerDrawHisto(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx, Int_t Ny)
736 //gener two dimensional histogram with PRF
739 sprintf(s,"Pad Response Function");
740 AliH2F * hPRFc = new AliH2F("hPRFc",s,Nx,x1,x2,Ny,y1,y2);
741 Float_t dx=(x2-x1)/Float_t(Nx);
742 Float_t dy=(y2-y1)/Float_t(Ny) ;
746 for ( Int_t i = 0;i<=Nx;i++,x+=dx){
748 for (Int_t j = 0;j<=Ny;j++,y+=dy){
750 hPRFc->SetCellContent(i,j,z);
753 hPRFc->SetXTitle("pad direction (cm)");
754 hPRFc->SetYTitle("pad row direction (cm)");
755 hPRFc->SetTitleOffset(1.5,"X");
756 hPRFc->SetTitleOffset(1.5,"Y");
761 AliH2F * AliTPCPRF2D::GenerDrawDistHisto(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx, Int_t Ny, Float_t thr)
763 //return histogram with distortion
764 const Float_t kminth=0.00001;
765 if (thr<kminth) thr=kminth;
767 sprintf(s,"COG distortion of PRF (threshold=%2.2f)",thr);
768 AliH2F * hPRFDist = new AliH2F("hDistortion",s,Nx,x1,x2,Ny,y1,y2);
769 Float_t dx=(x2-x1)/Float_t(Nx);
770 Float_t dy=(y2-y1)/Float_t(Ny) ;
773 for ( Int_t i = 0;i<=Nx;i++,x+=dx){
775 for(Int_t j = 0;j<=Ny;j++,y+=dy)
779 for (Int_t k=-3;k<=3;k++)
781 Float_t padx=Float_t(k)*fWidth;
782 z = GetPRF(x-padx,y);
790 ddx = (x-(sumx/sum));
793 if (TMath::Abs(ddx)<10) hPRFDist->SetCellContent(i,j,ddx);
797 hPRFDist->SetXTitle("pad direction (cm)");
798 hPRFDist->SetYTitle("pad row direction (cm)");
799 hPRFDist->SetTitleOffset(1.5,"X");
800 hPRFDist->SetTitleOffset(1.5,"Y");
808 void AliTPCPRF2D::DrawX(Float_t x1 ,Float_t x2,Float_t y1,Float_t y2, Int_t N)
811 //draw pad response function at interval <x1,x2> at given y position
814 TCanvas * c1 = new TCanvas("PRFX","Pad response function",700,900);
817 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
818 comment->SetTextAlign(12);
819 comment->SetFillColor(42);
820 DrawComment(comment);
824 TPad * pad2 = new TPad("pPRF","",0.05,0.22,0.95,0.95);
825 pad2->Divide(2,(N+1)/2);
827 gStyle->SetOptFit(1);
828 gStyle->SetOptStat(1);
829 for (Int_t i=0;i<N;i++){
833 else y = y1+i*(y2-y1)/Float_t(N-1);
835 TH1F * hPRFc =GenerDrawXHisto(x1, x2,y);
836 sprintf(ch,"PRF at wire position: %2.3f",y);
838 sprintf(ch,"PRF %d",i);
847 void AliTPCPRF2D::DrawPRF(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2, Int_t Nx, Int_t Ny)
851 TCanvas * c1 = new TCanvas("canPRF","Pad response function",700,900);
853 TPad * pad2 = new TPad("pad2PRF","",0.05,0.22,0.95,0.95);
855 gStyle->SetOptFit(1);
856 gStyle->SetOptStat(1);
857 TH2F * hPRFc = GenerDrawHisto(x1, x2, y1, y2, Nx,Ny);
861 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
862 comment->SetTextAlign(12);
863 comment->SetFillColor(42);
864 DrawComment(comment);
868 void AliTPCPRF2D::DrawDist(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2, Int_t Nx, Int_t Ny, Float_t thr)
871 //draw distortion of the COG method - for different threshold parameter
872 TCanvas * c1 = new TCanvas("padDistortion","COG distortion",700,900);
874 TPad * pad1 = new TPad("dist","",0.05,0.55,0.95,0.95,21);
876 TPad * pad2 = new TPad("dist","",0.05,0.22,0.95,0.53,21);
878 gStyle->SetOptFit(1);
879 gStyle->SetOptStat(0);
881 AliH2F * hPRFDist = GenerDrawDistHisto(x1, x2, y1, y2, Nx,Ny,thr);
884 hPRFDist->Draw("surf");
885 Float_t distmax =hPRFDist->GetMaximum();
886 Float_t distmin =hPRFDist->GetMinimum();
887 gStyle->SetOptStat(1);
889 TH1F * dist = hPRFDist->GetAmplitudes(distmin,distmax,distmin-1);
893 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
894 comment->SetTextAlign(12);
895 comment->SetFillColor(42);
896 DrawComment(comment);
900 void AliTPCPRF2D::DrawComment(TPaveText *comment)
903 //function to write comment to picture
906 //draw comments to picture
907 TText * title = comment->AddText("Pad Response Function parameters:");
908 title->SetTextSize(0.03);
909 sprintf(s,"Height of pad: %2.2f cm",fHeightFull);
911 sprintf(s,"Width pad: %2.2f cm",fWidth);
913 sprintf(s,"Pad Angle: %2.2f ",fPadAngle);
916 if (TMath::Abs(fK)>0.0001){
917 sprintf(s,"Height of one chevron unit h: %2.2f cm",2*fHeightS);
919 sprintf(s,"Overlap factor: %2.2f",fK);
923 if (strncmp(fType,"User",3)==0){
924 sprintf(s,"Charge distribution - user defined function %s ",fGRF->GetTitle());
926 sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
928 sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
931 if (strncmp(fType,"Gauss",3)==0){
932 sprintf(s,"Gauss charge distribution");
934 sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
936 sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
939 if (strncmp(fType,"Gati",3)==0){
940 sprintf(s,"Gati charge distribution");
942 sprintf(s,"K3X of Gati : %2.2f ",fK3X);
944 sprintf(s,"K3Y of Gati: %2.2f ",fK3Y);
946 sprintf(s,"Wire to Pad Distance: %2.2f ",fPadDistance);
949 if (strncmp(fType,"Cosh",3)==0){
950 sprintf(s,"Cosh charge distribution");
952 sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
954 sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
957 sprintf(s,"Normalisation: %2.2f ",fKNorm);