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.8 2001/10/21 19:07:24 hristov
19 Several pointers were set to zero in the default constructors to avoid memory management problems
21 Revision 1.7 2001/01/30 09:23:15 hristov
22 Streamers removed (R.Brun)
24 Revision 1.6 2000/09/07 11:23:27 kowal2
25 Improved algoritms, coding convensions applied.
27 Revision 1.5 2000/06/30 12:07:50 kowal2
28 Updated from the TPC-PreRelease branch
30 Revision 1.4.4.3 2000/06/26 07:39:42 kowal2
31 Changes to obey the coding rules
33 Revision 1.4.4.2 2000/06/25 08:38:41 kowal2
34 Splitted from AliTPCtracking
36 Revision 1.4.4.1 2000/06/14 16:48:24 kowal2
37 Parameter setting improved. Removed compiler warnings
39 Revision 1.4 2000/04/17 09:37:33 kowal2
40 removed obsolete AliTPCDigitsDisplay.C
42 Revision 1.3.8.2 2000/04/10 08:40:46 kowal2
44 Small changes by M. Ivanov, improvements of algorithms
46 Revision 1.3.8.1 2000/04/10 07:56:53 kowal2
47 Not used anymore - removed
49 Revision 1.3 1999/10/05 17:15:46 fca
50 Minor syntax for the Alpha OSF
52 Revision 1.2 1999/09/29 09:24:34 fca
53 Introduction of the Copyright and cvs Log
57 ///////////////////////////////////////////////////////////////////////////////
59 // Pad response function object in two dimesions //
60 // This class contains the basic functions for the //
61 // calculation of PRF according generic charge distribution //
62 // In Update function object calculate table of response function //
63 // in discrete x and y position //
64 // This table is used for interpolation od response function in any position //
65 // (function GetPRF) //
67 // Origin: Marian Ivanov, Uni. of Bratislava, ivanov@fmph.uniba.sk //
69 ///////////////////////////////////////////////////////////////////////////////
73 #include "AliTPCPRF2D.h"
75 #include <Riostream.h>
84 #include "TPaveText.h"
87 extern TStyle * gStyle;
89 const Double_t AliTPCPRF2D::fgkDegtoRad = 0.01745329251994;
90 const Double_t AliTPCPRF2D::fgkSQRT12=3.464101;
91 const Int_t AliTPCPRF2D::fgkNPRF = 100;
94 static Double_t funGauss2D(Double_t *x, Double_t * par)
96 //Gauss function -needde by the generic function object
97 return ( TMath::Exp(-(x[0]*x[0])/(2*par[0]*par[0]))*
98 TMath::Exp(-(x[1]*x[1])/(2*par[1]*par[1])));
102 static Double_t funCosh2D(Double_t *x, Double_t * par)
104 //Cosh function -needde by the generic function object
105 return ( 1/(TMath::CosH(3.14159*x[0]/(2*par[0]))*
106 TMath::CosH(3.14159*x[1]/(2*par[1]))));
109 static Double_t funGati2D(Double_t *x, Double_t * par)
111 //Gati function -needde by the generic function object
113 Float_t k3R=TMath::Sqrt(k3);
114 Float_t k2=(TMath::Pi()/2)*(1-k3R/2.);
115 Float_t k1=k2*k3R/(4*TMath::ATan(k3R));
116 Float_t l=x[0]/par[0];
117 Float_t tan2=TMath::TanH(k2*l);
119 Float_t res = k1*(1-tan2)/(1+k3*tan2);
120 //par[4] = is equal to k3Y
123 k2=(TMath::Pi()/2)*(1-k3R/2.);
124 k1=k2*k3R/(4*TMath::ATan(k3R));
126 tan2=TMath::TanH(k2*l);
128 res = res*k1*(1-tan2)/(1+k3*tan2);
132 ///////////////////////////////////////////////////////////////////////////
133 ///////////////////////////////////////////////////////////////////////////
135 ClassImp(AliTPCPRF2D)
137 AliTPCPRF2D::AliTPCPRF2D()
139 //default constructor for response function object
155 //chewron default values
157 SetChevron(0.2,0.0,1.0);
159 SetInterpolationType(2,0);
162 AliTPCPRF2D::~AliTPCPRF2D()
164 if (fChargeArray!=0) delete [] fChargeArray;
165 if (fGRF !=0 ) fGRF->Delete();
168 void AliTPCPRF2D::SetY(Float_t y1, Float_t y2, Int_t nYdiv)
171 //set virtual line position
172 //first and last line and number of lines
178 void AliTPCPRF2D::SetPad(Float_t width, Float_t height)
180 //set base chevron parameters
184 void AliTPCPRF2D::SetChevron(Float_t hstep,
188 //set shaping of chewron parameters
194 void AliTPCPRF2D::SetChParam(Float_t width, Float_t height,
195 Float_t hstep, Float_t shifty, Float_t fac)
197 SetPad(width,height);
198 SetChevron(hstep,shifty,fac);
202 Float_t AliTPCPRF2D::GetPRF(Float_t xin, Float_t yin)
204 //function which return pad response
205 //for the charge in distance xin
206 //return cubic aproximation of PRF or PRF at nearest virtual wire
207 if (fChargeArray==0) return 0;
208 //transform position to "wire position"
209 Float_t y=fDYtoWire*(yin-fY1);
210 if (fNYdiv == 1) y=fY1;
211 //normaly it find nearest line charge
213 Int_t i=Int_t(0.5+y);
214 if (y<0) i=Int_t(-0.5+y);
215 if ((i<0) || (i>=fNYdiv) ) return 0;
216 fcharge = &(fChargeArray[i*fNPRF]);
217 return GetPRFActiv(xin);
220 //make interpolation from more fore lines
223 if ((i<0) || (i>=fNYdiv) ) return 0;
229 fcharge =&(fChargeArray[(i-1)*fNPRF]);
230 z0 = GetPRFActiv(xin);
232 fcharge =&(fChargeArray[i*fNPRF]);
235 fcharge =&(fChargeArray[(i+1)*fNPRF]);
236 z2 = GetPRFActiv(xin);
239 fcharge =&(fChargeArray[(i+2)*fNPRF]);
240 z3 = GetPRFActiv(xin);
249 Float_t dy=y-Float_t(i);
251 res = a+b*dy+c*dy*dy+d*dy*dy*dy;
258 Float_t AliTPCPRF2D::GetPRFActiv(Float_t xin)
260 //GEt response function on given charege line
261 //return spline aproximaton
262 Float_t x = (xin*fDStepM1)+fNPRF/2;
265 if ( (i>1) && ((i+2)<fNPRF)) {
268 b = (fcharge[i+1]-fcharge[i-1])*0.5;
269 k = fcharge[i+1]-a-b;
270 l = (fcharge[i+2]-fcharge[i])*0.5-b;
273 Float_t dx=x-Float_t(i);
274 Float_t res = a+b*dx+c*dx*dx+d*dx*dx*dx;
281 Float_t AliTPCPRF2D::GetGRF(Float_t xin, Float_t yin)
283 //function which returnoriginal charge distribution
284 //this function is just normalised for fKnorm
286 return fKNorm*GetGRF()->Eval(xin,yin)/fInteg;
292 void AliTPCPRF2D::SetParam( TF2 * GRF, Float_t kNorm,
293 Float_t sigmaX, Float_t sigmaY)
295 //adjust parameters of the original charge distribution
296 //and pad size parameters
297 if (fGRF !=0 ) fGRF->Delete();
300 sprintf(fType,"User");
301 if (sigmaX ==0) sigmaX=(fWidth*(1+TMath::Abs(fK)))/fgkSQRT12;
302 if (sigmaY ==0) sigmaY=(fWidth*(1+TMath::Abs(fK)))/fgkSQRT12;
305 Double_t estimsigma =
306 TMath::Sqrt(sigmaX*sigmaX+(fWidth*fWidth*(1+TMath::Abs(fK))/12)+
307 TMath::Tan(fPadAngle*fgkDegtoRad)*TMath::Tan(fPadAngle*fgkDegtoRad)*fHeightFull*fHeightFull/12);
308 if (estimsigma < 5*sigmaX) {
309 fDStep = estimsigma/10.;
310 fNPRF = Int_t(estimsigma*8./fDStep);
314 Double_t width = fWidth*(1+TMath::Abs(fK))+TMath::Abs(TMath::Tan(fPadAngle*fgkDegtoRad))*fHeightFull;
315 fNPRF = Int_t((width+8.*sigmaX)/fDStep);
321 void AliTPCPRF2D::SetGauss(Float_t sigmaX, Float_t sigmaY,
325 // set parameters for Gauss generic charge distribution
330 sprintf(fType,"Gauss");
331 if (fGRF !=0 ) fGRF->Delete();
332 fGRF = new TF2("fun",funGauss2D,-5.,5.,-5.,5.,4);
337 funParam[3]=fHeightS;
339 fGRF->SetParameters(funParam);
340 Double_t estimsigma =
341 TMath::Sqrt(sigmaX*sigmaX+(fWidth*fWidth*(1+TMath::Abs(fK))/12)+
342 TMath::Tan(fPadAngle)*TMath::Tan(fPadAngle*fgkDegtoRad)*fHeightFull*fHeightFull/12);
343 if (estimsigma < 5*sigmaX) {
344 fDStep = estimsigma/10.;
345 fNPRF = Int_t(estimsigma*8./fDStep);
349 Double_t width = fWidth*(1+TMath::Abs(fK))+TMath::Abs(TMath::Tan(fPadAngle*fgkDegtoRad))*fHeightFull;
350 fNPRF = Int_t((width+8.*sigmaX)/fDStep);
355 void AliTPCPRF2D::SetCosh(Float_t sigmaX, Float_t sigmaY,
358 // set parameters for Cosh generic charge distribution
363 sprintf(fType,"Cosh");
364 if (fGRF !=0 ) fGRF->Delete();
365 fGRF = new TF2("fun", funCosh2D,-5.,5.,-5.,5.,4);
369 funParam[3]=fHeightS;
370 fGRF->SetParameters(funParam);
372 Double_t estimsigma = TMath::Sqrt(sigmaX*sigmaX+fWidth*fWidth*(1+TMath::Abs(fK))/12);
373 if (estimsigma < 5*sigmaX) {
374 fDStep = estimsigma/10.;
375 fNPRF = Int_t(estimsigma*8./fDStep);
379 fNPRF = Int_t((1.2*fWidth*(1+TMath::Abs(fK))+8.*sigmaX)/fDStep);
384 void AliTPCPRF2D::SetGati(Float_t K3X, Float_t K3Y,
388 // set parameters for Gati generic charge distribution
393 fPadDistance=padDistance;
394 sprintf(fType,"Gati");
395 if (fGRF !=0 ) fGRF->Delete();
396 fGRF = new TF2("fun", funGati2D,-5.,5.,-5.,5.,5);
398 funParam[0]=padDistance;
401 funParam[3]=fHeightS;
403 fGRF->SetParameters(funParam);
404 fOrigSigmaX=padDistance;
405 fOrigSigmaY=padDistance;
406 Float_t sigmaX = fOrigSigmaX;
407 Double_t estimsigma = TMath::Sqrt(sigmaX*sigmaX+fWidth*fWidth*(1+TMath::Abs(fK))/12);
408 if (estimsigma < 5*sigmaX) {
409 fDStep = estimsigma/10.;
410 fNPRF = Int_t(estimsigma*8./fDStep);
414 fNPRF = Int_t((1.2*fWidth*(1+TMath::Abs(fK))+8.*sigmaX)/fDStep);
420 void AliTPCPRF2D::Update()
423 //update fields with interpolated values for
426 if ( fGRF == 0 ) return;
427 //initialize interpolated values to 0
429 if (fChargeArray!=0) delete [] fChargeArray;
430 fChargeArray = new Float_t[fNPRF*fNYdiv];
431 fNChargeArray = fNPRF*fNYdiv;
432 for (i =0; i<fNPRF*fNYdiv;i++) fChargeArray[i] = 0;
433 //firstly calculate total integral of charge
435 ////////////////////////////////////////////////////////
436 //I'm waiting for normal integral
437 //in this moment only sum
438 Float_t x2= 4*fOrigSigmaX;
439 Float_t y2= 4*fOrigSigmaY;
440 Float_t dx = fOrigSigmaX/Float_t(fNdiv*6);
441 Float_t dy = fOrigSigmaY/Float_t(fNdiv*6);
442 Int_t nx = Int_t(0.5+x2/dx);
443 Int_t ny = Int_t(0.5+y2/dy);
447 for (ix=-nx;ix<=nx;ix++)
448 for ( iy=-ny;iy<=ny;iy++)
449 dInteg+=fGRF->Eval(Float_t(ix)*dx,Float_t(iy)*dy)*dx*dy;
450 /////////////////////////////////////////////////////
452 if ( fInteg == 0 ) fInteg = 1;
454 for (i=0; i<fNYdiv; i++){
455 if (fNYdiv == 1) fCurrentY = fY1;
457 fCurrentY = fY1+Double_t(i)*(fY2-fY1)/Double_t(fNYdiv-1);
458 fcharge = &(fChargeArray[i*fNPRF]);
461 //calculate conversion coefitient to convert position to virtual wire
462 fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1);
467 void AliTPCPRF2D::Update1()
470 //update fields with interpolated values for
471 //PRF calculation for given charge line
473 Double_t cos = TMath::Cos(fChargeAngle);
474 Double_t sin = TMath::Sin(fChargeAngle);
475 const Double_t kprec =0.00000001;
476 //integrate charge over pad for different distance of pad
477 for (i =0; i<fNPRF;i++){
478 //x in cm fWidth in cm
480 Double_t xch = fDStep * (Double_t)(i-fNPRF/2);
485 for (Double_t ym=-fHeightFull/2.-fShiftY; ym<fHeightFull/2.-kprec;ym+=fHeightS){
486 Double_t y2chev=TMath::Min((ym+fHeightS),Double_t(fHeightFull/2.)); // end of chevron step
487 Double_t y1chev= ym; //beginning of chevron step
488 Double_t y2 = TMath::Min(y2chev,fCurrentY+3.5*fOrigSigmaY);
489 Double_t y1 = TMath::Max((y1chev),Double_t(-fHeightFull/2.));
490 y1 = TMath::Max(y1chev,fCurrentY-3.5*fOrigSigmaY);
492 Double_t x0 = fWidth*(-1.-(Double_t(k)*fK))*0.5+ym*TMath::Tan(fPadAngle*fgkDegtoRad);
493 Double_t kx = Double_t(k)*(fK*fWidth)/fHeightS;
494 kx = TMath::Tan(TMath::ATan(kx))+TMath::Tan(fPadAngle*fgkDegtoRad);
496 Int_t ny = TMath::Max(Int_t(fNdiv*TMath::Exp(-(y1-fCurrentY)*(y1-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),4);
497 Double_t dy = TMath::Min(fOrigSigmaY/Double_t(ny),y2-y1);
500 //loop over different y strips with variable step size dy
501 if (y2>(y1+kprec)) for (Double_t y = y1; y<y2+kprec;){
504 ny = TMath::Max(Int_t(fNdiv*TMath::Exp(-(y-fCurrentY)*(y-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),5);
505 ndy = fOrigSigmaY/Double_t(ny);
508 if (ndy<kprec) ndy=2*kprec; //calculate new delta y
512 //calculation of x borders and initial step
513 Double_t deltay = (y-y1chev);
515 Double_t xp1 = x0+deltay*kx;
516 //x begining of pad at position y
517 Double_t xp2 =xp1+fWidth; //x end of pad at position y
518 Double_t xp3 =xp1+kx*dy; //...at position y+dy
519 Double_t xp4 =xp2+kx*dy; //..
521 Double_t x1 = TMath::Min(xp1,xp3);
522 x1 = TMath::Max(xp1,xch-3.5*fOrigSigmaX); //beging of integration
523 Double_t x2 = TMath::Max(xp2,xp4);
524 x2 = TMath::Min(xp2+dy*kx,xch+3.5*fOrigSigmaX); //end of integration
526 Int_t nx = TMath::Max(Int_t(fNdiv*TMath::Exp(-(x1-xch)*(x1-xch)/(2*fOrigSigmaX*fOrigSigmaX))*
527 TMath::Exp(-(y1-fCurrentY)*(y1-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),2);
528 Double_t dx = TMath::Min(fOrigSigmaX/Double_t(nx),x2-x1)/5.; //on the border more iteration
532 for (Double_t x = x1; x<x2+kprec ;){
534 nx = TMath::Max(Int_t(fNdiv*TMath::Exp(-(x-xch)*(x-xch)/(2*fOrigSigmaX*fOrigSigmaX))),3);
535 ndx = fOrigSigmaX/Double_t(nx);
539 if ( ( (x+dx+ndx)<TMath::Max(xp3,xp1)) || ( (x+dx+ndx)>TMath::Min(xp4,xp2))) {
542 if (ndx<kprec) ndx=2*kprec;
543 //INTEGRAL APROXIMATION
544 Double_t ddx,ddy,dddx,dddy;
546 ddy = fCurrentY-(y+dy/2.);
547 dddx = cos*ddx-sin*ddy;
548 dddy = sin*ddx+cos*ddy;
549 Double_t z0=fGRF->Eval(dddx,dddy); //middle point
553 dddx = cos*ddx-sin*ddy;
554 dddy = sin*ddx+cos*ddy;
555 Double_t z1=fGRF->Eval(dddx,dddy); //point down
558 ddy = fCurrentY-(y+dy);
559 dddx = cos*ddx-sin*ddy;
560 dddy = sin*ddx+cos*ddy;
561 Double_t z3=fGRF->Eval(dddx,dddy); //point up
564 ddy = fCurrentY-(y+dy/2.);
565 dddx = cos*ddx-sin*ddy;
566 dddy = sin*ddx+cos*ddy;
567 Double_t z2=fGRF->Eval(dddx,dddy); //point left
570 ddy = fCurrentY-(y+dy/2.);
571 dddx = cos*ddx-sin*ddy;
572 dddy = sin*ddx+cos*ddy;
573 Double_t z4=fGRF->Eval(dddx,dddy); //point right
576 if (z0<0) {z0=0;z1=0;z2=0;z3=0;z4=0;}
578 Double_t f2x= (z3+z1-2*z0)*4.;//second derivation in y
579 Double_t f2y= (z2+z4-2*z0)*4.;//second derivation in x
580 Double_t f1y= (z3-z1);
582 z = (z0+f2x/6.+f2y/6.);//second order aproxiation of integral
583 if (kx>kprec){ //positive derivation
584 if (x<(xp1+dy*kx)){ //calculate volume at left border
586 Double_t xx2 = TMath::Min(x+dx,xp1+dy*kx);
587 Double_t yy1 = y+(xx1-xp1)/kx;
588 Double_t yy2 = TMath::Min(y+(xx2-xp1)/kx,y+dy);
591 z-= z0*(y+dy-yy2)/dy; //constant part rectangle
592 z-= f1y*(xx2-xx1)*(y+dy-yy2)*(y+dy-yy2)/(2.*dx*dy);
594 z-=z0*(xx2-xx1)*(yy2-yy1)/(2*dx*dy); //constant part rectangle
597 if (x>xp2){ //calculate volume at right border
600 Double_t yy1 = y+(xx1-xp2)/kx;
601 Double_t yy2 = y+(xx2-xp2)/kx;
604 z-=z0*(yy1-y)/dy; //constant part
605 z-=f1y*(xx2-xx1)*(yy1-y)*(yy1-y)/(2*dx*dy);
607 z-=z0*(xx2-xx1)*(yy2-yy1)/(2*dx*dy); //constant part
610 if (kx<-kprec){ //negative derivation
611 if (x<(xp1+dy*kx)){ //calculate volume at left border
613 Double_t xx2 = TMath::Min(x+dx,xp3-dy/kx);
614 Double_t yy1 = y+(xx1-xp1)/kx;
615 Double_t yy2 = TMath::Max(y,yy1+(xx2-xx1)/kx); //yy2<yy1
617 z-= z0*(yy2-y)/dy; // constant part rectangle
618 z-= f1y*(xx2-xx1)*(yy2-y)*(yy2-y)/(2.*dx*dy);
619 z-=z0*(xx2-xx1)*(yy1-yy2)/(2*dx*dy); //constant part triangle
621 if (x>xp2){ //calculate volume at right border
622 Double_t xx1 = TMath::Max(x,xp2+dy*kx);
624 Double_t yy1 = TMath::Min(y+dy,y-(xp2-xx1)/kx);
625 Double_t yy2 = y-(xp2-xx2)/kx;
627 z-=z0*(yy2-y)/dy; //constant part rextangle
628 z-= f1y*(xx2-xx1)*(yy2-y)*(yy2-y)/(2.*dx*dy);
629 z-=z0*(xx2-xx1)*(yy1-yy2)/(2*dx*dy); //constant part triangle
633 if (z>0.) sumch+=fKNorm*z*dx*dy/fInteg;
642 }//step over different y
646 }//step over different points on line NPRF
649 void AliTPCPRF2D::UpdateSigma()
652 //calulate effective sigma X and sigma y of PRF
662 for (i=-1; i<=fNYdiv; i++){
663 if (fNYdiv == 1) y = fY1;
665 y = fY1+Float_t(i)*(fY2-fY1)/Float_t(fNYdiv-1);
666 for (x =-fNPRF*fDStep; x<fNPRF*fDStep;x+=fDStep)
668 //x in cm fWidth in cm
669 Float_t weight = GetPRF(x,y);
680 fSigmaX = TMath::Sqrt(fSigmaX/sum-fMeanX*fMeanX);
681 fSigmaY = TMath::Sqrt(fSigmaY/sum-fMeanY*fMeanY);
687 void AliTPCPRF2D::Streamer(TBuffer &R__b)
689 // Stream an object of class AliTPCPRF2D
691 if (R__b.IsReading()) {
693 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
694 AliTPCPRF2D::Class()->ReadBuffer(R__b, this, R__v, R__s, R__c);
696 if (strncmp(fType,"User",3)!=0){
698 if (strncmp(fType,"Gauss",3)==0)
699 fGRF = new TF2("fun",funGauss2D,-5.,5.,-5.,5.,4);
700 if (strncmp(fType,"Cosh",3)==0)
701 fGRF = new TF2("fun",funCosh2D,-5.,5.,-5.,5.,4);
702 if (strncmp(fType,"Gati",3)==0)
703 fGRF = new TF2("fun",funGati2D,-5.,5.,-5.,5.,5);
704 if (fGRF!=0) fGRF->SetParameters(funParam);
706 //calculate conversion coefitient to convert position to virtual wire
707 fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1);
710 AliTPCPRF2D::Class()->WriteBuffer(R__b,this);
715 TH1F * AliTPCPRF2D::GenerDrawXHisto(Float_t x1, Float_t x2,Float_t y)
717 //gener one dimensional hist of pad response function
721 sprintf(s,"Pad Response Function");
722 TH1F * hPRFc = new TH1F("hPRFc",s,kn+1,x1,x2);
726 for (Int_t i = 0;i<kn+1;i++)
728 x+=(x2-x1)/Float_t(kn);
732 hPRFc->SetXTitle("pad (cm)");
736 AliH2F * AliTPCPRF2D::GenerDrawHisto(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx, Int_t Ny)
739 //gener two dimensional histogram with PRF
742 sprintf(s,"Pad Response Function");
743 AliH2F * hPRFc = new AliH2F("hPRFc",s,Nx,x1,x2,Ny,y1,y2);
744 Float_t dx=(x2-x1)/Float_t(Nx);
745 Float_t dy=(y2-y1)/Float_t(Ny) ;
749 for ( Int_t i = 0;i<=Nx;i++,x+=dx){
751 for (Int_t j = 0;j<=Ny;j++,y+=dy){
753 hPRFc->SetCellContent(i,j,z);
756 hPRFc->SetXTitle("pad direction (cm)");
757 hPRFc->SetYTitle("pad row direction (cm)");
758 hPRFc->SetTitleOffset(1.5,"X");
759 hPRFc->SetTitleOffset(1.5,"Y");
764 AliH2F * AliTPCPRF2D::GenerDrawDistHisto(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx, Int_t Ny, Float_t thr)
766 //return histogram with distortion
767 const Float_t kminth=0.00001;
768 if (thr<kminth) thr=kminth;
770 sprintf(s,"COG distortion of PRF (threshold=%2.2f)",thr);
771 AliH2F * hPRFDist = new AliH2F("hDistortion",s,Nx,x1,x2,Ny,y1,y2);
772 Float_t dx=(x2-x1)/Float_t(Nx);
773 Float_t dy=(y2-y1)/Float_t(Ny) ;
776 for ( Int_t i = 0;i<=Nx;i++,x+=dx){
778 for(Int_t j = 0;j<=Ny;j++,y+=dy)
782 for (Int_t k=-3;k<=3;k++)
784 Float_t padx=Float_t(k)*fWidth;
785 z = GetPRF(x-padx,y);
793 ddx = (x-(sumx/sum));
796 if (TMath::Abs(ddx)<10) hPRFDist->SetCellContent(i,j,ddx);
800 hPRFDist->SetXTitle("pad direction (cm)");
801 hPRFDist->SetYTitle("pad row direction (cm)");
802 hPRFDist->SetTitleOffset(1.5,"X");
803 hPRFDist->SetTitleOffset(1.5,"Y");
811 void AliTPCPRF2D::DrawX(Float_t x1 ,Float_t x2,Float_t y1,Float_t y2, Int_t N)
814 //draw pad response function at interval <x1,x2> at given y position
817 TCanvas * c1 = new TCanvas("PRFX","Pad response function",700,900);
820 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
821 comment->SetTextAlign(12);
822 comment->SetFillColor(42);
823 DrawComment(comment);
827 TPad * pad2 = new TPad("pPRF","",0.05,0.22,0.95,0.95);
828 pad2->Divide(2,(N+1)/2);
830 gStyle->SetOptFit(1);
831 gStyle->SetOptStat(1);
832 for (Int_t i=0;i<N;i++){
836 else y = y1+i*(y2-y1)/Float_t(N-1);
838 TH1F * hPRFc =GenerDrawXHisto(x1, x2,y);
839 sprintf(ch,"PRF at wire position: %2.3f",y);
841 sprintf(ch,"PRF %d",i);
850 void AliTPCPRF2D::DrawPRF(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2, Int_t Nx, Int_t Ny)
854 TCanvas * c1 = new TCanvas("canPRF","Pad response function",700,900);
856 TPad * pad2 = new TPad("pad2PRF","",0.05,0.22,0.95,0.95);
858 gStyle->SetOptFit(1);
859 gStyle->SetOptStat(1);
860 TH2F * hPRFc = GenerDrawHisto(x1, x2, y1, y2, Nx,Ny);
864 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
865 comment->SetTextAlign(12);
866 comment->SetFillColor(42);
867 DrawComment(comment);
871 void AliTPCPRF2D::DrawDist(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2, Int_t Nx, Int_t Ny, Float_t thr)
874 //draw distortion of the COG method - for different threshold parameter
875 TCanvas * c1 = new TCanvas("padDistortion","COG distortion",700,900);
877 TPad * pad1 = new TPad("dist","",0.05,0.55,0.95,0.95,21);
879 TPad * pad2 = new TPad("dist","",0.05,0.22,0.95,0.53,21);
881 gStyle->SetOptFit(1);
882 gStyle->SetOptStat(0);
884 AliH2F * hPRFDist = GenerDrawDistHisto(x1, x2, y1, y2, Nx,Ny,thr);
887 hPRFDist->Draw("surf");
888 Float_t distmax =hPRFDist->GetMaximum();
889 Float_t distmin =hPRFDist->GetMinimum();
890 gStyle->SetOptStat(1);
892 TH1F * dist = hPRFDist->GetAmplitudes(distmin,distmax,distmin-1);
896 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
897 comment->SetTextAlign(12);
898 comment->SetFillColor(42);
899 DrawComment(comment);
903 void AliTPCPRF2D::DrawComment(TPaveText *comment)
906 //function to write comment to picture
909 //draw comments to picture
910 TText * title = comment->AddText("Pad Response Function parameters:");
911 title->SetTextSize(0.03);
912 sprintf(s,"Height of pad: %2.2f cm",fHeightFull);
914 sprintf(s,"Width pad: %2.2f cm",fWidth);
916 sprintf(s,"Pad Angle: %2.2f ",fPadAngle);
919 if (TMath::Abs(fK)>0.0001){
920 sprintf(s,"Height of one chevron unit h: %2.2f cm",2*fHeightS);
922 sprintf(s,"Overlap factor: %2.2f",fK);
926 if (strncmp(fType,"User",3)==0){
927 sprintf(s,"Charge distribution - user defined function %s ",fGRF->GetTitle());
929 sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
931 sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
934 if (strncmp(fType,"Gauss",3)==0){
935 sprintf(s,"Gauss charge distribution");
937 sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
939 sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
942 if (strncmp(fType,"Gati",3)==0){
943 sprintf(s,"Gati charge distribution");
945 sprintf(s,"K3X of Gati : %2.2f ",fK3X);
947 sprintf(s,"K3Y of Gati: %2.2f ",fK3Y);
949 sprintf(s,"Wire to Pad Distance: %2.2f ",fPadDistance);
952 if (strncmp(fType,"Cosh",3)==0){
953 sprintf(s,"Cosh charge distribution");
955 sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
957 sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
960 sprintf(s,"Normalisation: %2.2f ",fKNorm);