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.5 2000/06/30 12:07:50 kowal2
19 Updated from the TPC-PreRelease branch
21 Revision 1.4.4.3 2000/06/26 07:39:42 kowal2
22 Changes to obey the coding rules
24 Revision 1.4.4.2 2000/06/25 08:38:41 kowal2
25 Splitted from AliTPCtracking
27 Revision 1.4.4.1 2000/06/14 16:48:24 kowal2
28 Parameter setting improved. Removed compiler warnings
30 Revision 1.4 2000/04/17 09:37:33 kowal2
31 removed obsolete AliTPCDigitsDisplay.C
33 Revision 1.3.8.2 2000/04/10 08:40:46 kowal2
35 Small changes by M. Ivanov, improvements of algorithms
37 Revision 1.3.8.1 2000/04/10 07:56:53 kowal2
38 Not used anymore - removed
40 Revision 1.3 1999/10/05 17:15:46 fca
41 Minor syntax for the Alpha OSF
43 Revision 1.2 1999/09/29 09:24:34 fca
44 Introduction of the Copyright and cvs Log
48 ///////////////////////////////////////////////////////////////////////////////
50 // Pad response function object in two dimesions //
51 // This class contains the basic functions for the //
52 // calculation of PRF according generic charge distribution //
53 // In Update function object calculate table of response function //
54 // in discrete x and y position //
55 // This table is used for interpolation od response function in any position //
56 // (function GetPRF) //
58 // Origin: Marian Ivanov, Uni. of Bratislava, ivanov@fmph.uniba.sk //
60 ///////////////////////////////////////////////////////////////////////////////
64 #include "AliTPCPRF2D.h"
75 #include "TPaveText.h"
78 extern TStyle * gStyle;
80 const Double_t AliTPCPRF2D::fgkDegtoRad = 0.01745329251994;
81 const Double_t AliTPCPRF2D::fgkSQRT12=3.464101;
82 const Int_t AliTPCPRF2D::fgkNPRF = 100;
85 static Double_t funGauss2D(Double_t *x, Double_t * par)
87 //Gauss function -needde by the generic function object
88 return ( TMath::Exp(-(x[0]*x[0])/(2*par[0]*par[0]))*
89 TMath::Exp(-(x[1]*x[1])/(2*par[1]*par[1])));
93 static Double_t funCosh2D(Double_t *x, Double_t * par)
95 //Cosh function -needde by the generic function object
96 return ( 1/(TMath::CosH(3.14159*x[0]/(2*par[0]))*
97 TMath::CosH(3.14159*x[1]/(2*par[1]))));
100 static Double_t funGati2D(Double_t *x, Double_t * par)
102 //Gati function -needde by the generic function object
104 Float_t k3R=TMath::Sqrt(k3);
105 Float_t k2=(TMath::Pi()/2)*(1-k3R/2.);
106 Float_t k1=k2*k3R/(4*TMath::ATan(k3R));
107 Float_t l=x[0]/par[0];
108 Float_t tan2=TMath::TanH(k2*l);
110 Float_t res = k1*(1-tan2)/(1+k3*tan2);
111 //par[4] = is equal to k3Y
114 k2=(TMath::Pi()/2)*(1-k3R/2.);
115 k1=k2*k3R/(4*TMath::ATan(k3R));
117 tan2=TMath::TanH(k2*l);
119 res = res*k1*(1-tan2)/(1+k3*tan2);
123 ///////////////////////////////////////////////////////////////////////////
124 ///////////////////////////////////////////////////////////////////////////
126 ClassImp(AliTPCPRF2D)
128 AliTPCPRF2D::AliTPCPRF2D()
130 //default constructor for response function object
145 //chewron default values
147 SetChevron(0.2,0.0,1.0);
149 SetInterpolationType(2,0);
152 AliTPCPRF2D::~AliTPCPRF2D()
154 if (fChargeArray!=0) delete [] fChargeArray;
155 if (fGRF !=0 ) fGRF->Delete();
158 void AliTPCPRF2D::SetY(Float_t y1, Float_t y2, Int_t nYdiv)
161 //set virtual line position
162 //first and last line and number of lines
168 void AliTPCPRF2D::SetPad(Float_t width, Float_t height)
170 //set base chevron parameters
174 void AliTPCPRF2D::SetChevron(Float_t hstep,
178 //set shaping of chewron parameters
184 void AliTPCPRF2D::SetChParam(Float_t width, Float_t height,
185 Float_t hstep, Float_t shifty, Float_t fac)
187 SetPad(width,height);
188 SetChevron(hstep,shifty,fac);
192 Float_t AliTPCPRF2D::GetPRF(Float_t xin, Float_t yin)
194 //function which return pad response
195 //for the charge in distance xin
196 //return cubic aproximation of PRF or PRF at nearest virtual wire
197 if (fChargeArray==0) return 0;
198 //transform position to "wire position"
199 Float_t y=fDYtoWire*(yin-fY1);
200 if (fNYdiv == 1) y=fY1;
201 //normaly it find nearest line charge
203 Int_t i=Int_t(0.5+y);
204 if (y<0) i=Int_t(-0.5+y);
205 if ((i<0) || (i>=fNYdiv) ) return 0;
206 fcharge = &(fChargeArray[i*fNPRF]);
207 return GetPRFActiv(xin);
210 //make interpolation from more fore lines
213 if ((i<0) || (i>=fNYdiv) ) return 0;
219 fcharge =&(fChargeArray[(i-1)*fNPRF]);
220 z0 = GetPRFActiv(xin);
222 fcharge =&(fChargeArray[i*fNPRF]);
225 fcharge =&(fChargeArray[(i+1)*fNPRF]);
226 z2 = GetPRFActiv(xin);
229 fcharge =&(fChargeArray[(i+2)*fNPRF]);
230 z3 = GetPRFActiv(xin);
239 Float_t dy=y-Float_t(i);
241 res = a+b*dy+c*dy*dy+d*dy*dy*dy;
248 Float_t AliTPCPRF2D::GetPRFActiv(Float_t xin)
250 //GEt response function on given charege line
251 //return spline aproximaton
252 Float_t x = (xin*fDStepM1)+fNPRF/2;
255 if ( (i>1) && ((i+2)<fNPRF)) {
258 b = (fcharge[i+1]-fcharge[i-1])*0.5;
259 k = fcharge[i+1]-a-b;
260 l = (fcharge[i+2]-fcharge[i])*0.5-b;
263 Float_t dx=x-Float_t(i);
264 Float_t res = a+b*dx+c*dx*dx+d*dx*dx*dx;
271 Float_t AliTPCPRF2D::GetGRF(Float_t xin, Float_t yin)
273 //function which returnoriginal charge distribution
274 //this function is just normalised for fKnorm
276 return fKNorm*GetGRF()->Eval(xin,yin)/fInteg;
282 void AliTPCPRF2D::SetParam( TF2 * GRF, Float_t kNorm,
283 Float_t sigmaX, Float_t sigmaY)
285 //adjust parameters of the original charge distribution
286 //and pad size parameters
287 if (fGRF !=0 ) fGRF->Delete();
290 sprintf(fType,"User");
291 if (sigmaX ==0) sigmaX=(fWidth*(1+TMath::Abs(fK)))/fgkSQRT12;
292 if (sigmaY ==0) sigmaY=(fWidth*(1+TMath::Abs(fK)))/fgkSQRT12;
295 Double_t estimsigma =
296 TMath::Sqrt(sigmaX*sigmaX+(fWidth*fWidth*(1+TMath::Abs(fK))/12)+
297 TMath::Tan(fPadAngle*fgkDegtoRad)*TMath::Tan(fPadAngle*fgkDegtoRad)*fHeightFull*fHeightFull/12);
298 if (estimsigma < 5*sigmaX) {
299 fDStep = estimsigma/10.;
300 fNPRF = Int_t(estimsigma*8./fDStep);
304 Double_t width = fWidth*(1+TMath::Abs(fK))+TMath::Abs(TMath::Tan(fPadAngle*fgkDegtoRad))*fHeightFull;
305 fNPRF = Int_t((width+8.*sigmaX)/fDStep);
311 void AliTPCPRF2D::SetGauss(Float_t sigmaX, Float_t sigmaY,
315 // set parameters for Gauss generic charge distribution
320 sprintf(fType,"Gauss");
321 if (fGRF !=0 ) fGRF->Delete();
322 fGRF = new TF2("fun",funGauss2D,-5.,5.,-5.,5.,4);
327 funParam[3]=fHeightS;
329 fGRF->SetParameters(funParam);
330 Double_t estimsigma =
331 TMath::Sqrt(sigmaX*sigmaX+(fWidth*fWidth*(1+TMath::Abs(fK))/12)+
332 TMath::Tan(fPadAngle)*TMath::Tan(fPadAngle*fgkDegtoRad)*fHeightFull*fHeightFull/12);
333 if (estimsigma < 5*sigmaX) {
334 fDStep = estimsigma/10.;
335 fNPRF = Int_t(estimsigma*8./fDStep);
339 Double_t width = fWidth*(1+TMath::Abs(fK))+TMath::Abs(TMath::Tan(fPadAngle*fgkDegtoRad))*fHeightFull;
340 fNPRF = Int_t((width+8.*sigmaX)/fDStep);
345 void AliTPCPRF2D::SetCosh(Float_t sigmaX, Float_t sigmaY,
348 // set parameters for Cosh generic charge distribution
353 sprintf(fType,"Cosh");
354 if (fGRF !=0 ) fGRF->Delete();
355 fGRF = new TF2("fun", funCosh2D,-5.,5.,-5.,5.,4);
359 funParam[3]=fHeightS;
360 fGRF->SetParameters(funParam);
362 Double_t estimsigma = TMath::Sqrt(sigmaX*sigmaX+fWidth*fWidth*(1+TMath::Abs(fK))/12);
363 if (estimsigma < 5*sigmaX) {
364 fDStep = estimsigma/10.;
365 fNPRF = Int_t(estimsigma*8./fDStep);
369 fNPRF = Int_t((1.2*fWidth*(1+TMath::Abs(fK))+8.*sigmaX)/fDStep);
374 void AliTPCPRF2D::SetGati(Float_t K3X, Float_t K3Y,
378 // set parameters for Gati generic charge distribution
383 fPadDistance=padDistance;
384 sprintf(fType,"Gati");
385 if (fGRF !=0 ) fGRF->Delete();
386 fGRF = new TF2("fun", funGati2D,-5.,5.,-5.,5.,5);
388 funParam[0]=padDistance;
391 funParam[3]=fHeightS;
393 fGRF->SetParameters(funParam);
394 fOrigSigmaX=padDistance;
395 fOrigSigmaY=padDistance;
396 Float_t sigmaX = fOrigSigmaX;
397 Double_t estimsigma = TMath::Sqrt(sigmaX*sigmaX+fWidth*fWidth*(1+TMath::Abs(fK))/12);
398 if (estimsigma < 5*sigmaX) {
399 fDStep = estimsigma/10.;
400 fNPRF = Int_t(estimsigma*8./fDStep);
404 fNPRF = Int_t((1.2*fWidth*(1+TMath::Abs(fK))+8.*sigmaX)/fDStep);
410 void AliTPCPRF2D::Update()
413 //update fields with interpolated values for
416 if ( fGRF == 0 ) return;
417 //initialize interpolated values to 0
419 if (fChargeArray!=0) delete [] fChargeArray;
420 fChargeArray = new Float_t[fNPRF*fNYdiv];
421 fNChargeArray = fNPRF*fNYdiv;
422 for (i =0; i<fNPRF*fNYdiv;i++) fChargeArray[i] = 0;
423 //firstly calculate total integral of charge
425 ////////////////////////////////////////////////////////
426 //I'm waiting for normal integral
427 //in this moment only sum
428 Float_t x2= 4*fOrigSigmaX;
429 Float_t y2= 4*fOrigSigmaY;
430 Float_t dx = fOrigSigmaX/Float_t(fNdiv*6);
431 Float_t dy = fOrigSigmaY/Float_t(fNdiv*6);
432 Int_t nx = Int_t(0.5+x2/dx);
433 Int_t ny = Int_t(0.5+y2/dy);
437 for (ix=-nx;ix<=nx;ix++)
438 for ( iy=-ny;iy<=ny;iy++)
439 dInteg+=fGRF->Eval(Float_t(ix)*dx,Float_t(iy)*dy)*dx*dy;
440 /////////////////////////////////////////////////////
442 if ( fInteg == 0 ) fInteg = 1;
444 for (i=0; i<fNYdiv; i++){
445 if (fNYdiv == 1) fCurrentY = fY1;
447 fCurrentY = fY1+Double_t(i)*(fY2-fY1)/Double_t(fNYdiv-1);
448 fcharge = &(fChargeArray[i*fNPRF]);
451 //calculate conversion coefitient to convert position to virtual wire
452 fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1);
457 void AliTPCPRF2D::Update1()
460 //update fields with interpolated values for
461 //PRF calculation for given charge line
463 Double_t cos = TMath::Cos(fChargeAngle);
464 Double_t sin = TMath::Sin(fChargeAngle);
465 const Double_t kprec =0.00000001;
466 //integrate charge over pad for different distance of pad
467 for (i =0; i<fNPRF;i++){
468 //x in cm fWidth in cm
470 Double_t xch = fDStep * (Double_t)(i-fNPRF/2);
475 for (Double_t ym=-fHeightFull/2.-fShiftY; ym<fHeightFull/2.-kprec;ym+=fHeightS){
476 Double_t y2chev=TMath::Min((ym+fHeightS),Double_t(fHeightFull/2.)); // end of chevron step
477 Double_t y1chev= ym; //beginning of chevron step
478 Double_t y2 = TMath::Min(y2chev,fCurrentY+3.5*fOrigSigmaY);
479 Double_t y1 = TMath::Max((y1chev),Double_t(-fHeightFull/2.));
480 y1 = TMath::Max(y1chev,fCurrentY-3.5*fOrigSigmaY);
482 Double_t x0 = fWidth*(-1.-(Double_t(k)*fK))*0.5+ym*TMath::Tan(fPadAngle*fgkDegtoRad);
483 Double_t kx = Double_t(k)*(fK*fWidth)/fHeightS;
484 kx = TMath::Tan(TMath::ATan(kx))+TMath::Tan(fPadAngle*fgkDegtoRad);
486 Int_t ny = TMath::Max(Int_t(fNdiv*TMath::Exp(-(y1-fCurrentY)*(y1-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),4);
487 Double_t dy = TMath::Min(fOrigSigmaY/Double_t(ny),y2-y1);
490 //loop over different y strips with variable step size dy
491 if (y2>(y1+kprec)) for (Double_t y = y1; y<y2+kprec;){
494 ny = TMath::Max(Int_t(fNdiv*TMath::Exp(-(y-fCurrentY)*(y-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),5);
495 ndy = fOrigSigmaY/Double_t(ny);
498 if (ndy<kprec) ndy=2*kprec; //calculate new delta y
502 //calculation of x borders and initial step
503 Double_t deltay = (y-y1chev);
505 Double_t xp1 = x0+deltay*kx;
506 //x begining of pad at position y
507 Double_t xp2 =xp1+fWidth; //x end of pad at position y
508 Double_t xp3 =xp1+kx*dy; //...at position y+dy
509 Double_t xp4 =xp2+kx*dy; //..
511 Double_t x1 = TMath::Min(xp1,xp3);
512 x1 = TMath::Max(xp1,xch-3.5*fOrigSigmaX); //beging of integration
513 Double_t x2 = TMath::Max(xp2,xp4);
514 x2 = TMath::Min(xp2+dy*kx,xch+3.5*fOrigSigmaX); //end of integration
516 Int_t nx = TMath::Max(Int_t(fNdiv*TMath::Exp(-(x1-xch)*(x1-xch)/(2*fOrigSigmaX*fOrigSigmaX))*
517 TMath::Exp(-(y1-fCurrentY)*(y1-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),2);
518 Double_t dx = TMath::Min(fOrigSigmaX/Double_t(nx),x2-x1)/5.; //on the border more iteration
522 for (Double_t x = x1; x<x2+kprec ;){
524 nx = TMath::Max(Int_t(fNdiv*TMath::Exp(-(x-xch)*(x-xch)/(2*fOrigSigmaX*fOrigSigmaX))),3);
525 ndx = fOrigSigmaX/Double_t(nx);
529 if ( ( (x+dx+ndx)<TMath::Max(xp3,xp1)) || ( (x+dx+ndx)>TMath::Min(xp4,xp2))) {
532 if (ndx<kprec) ndx=2*kprec;
533 //INTEGRAL APROXIMATION
534 Double_t ddx,ddy,dddx,dddy;
536 ddy = fCurrentY-(y+dy/2.);
537 dddx = cos*ddx-sin*ddy;
538 dddy = sin*ddx+cos*ddy;
539 Double_t z0=fGRF->Eval(dddx,dddy); //middle point
543 dddx = cos*ddx-sin*ddy;
544 dddy = sin*ddx+cos*ddy;
545 Double_t z1=fGRF->Eval(dddx,dddy); //point down
548 ddy = fCurrentY-(y+dy);
549 dddx = cos*ddx-sin*ddy;
550 dddy = sin*ddx+cos*ddy;
551 Double_t z3=fGRF->Eval(dddx,dddy); //point up
554 ddy = fCurrentY-(y+dy/2.);
555 dddx = cos*ddx-sin*ddy;
556 dddy = sin*ddx+cos*ddy;
557 Double_t z2=fGRF->Eval(dddx,dddy); //point left
560 ddy = fCurrentY-(y+dy/2.);
561 dddx = cos*ddx-sin*ddy;
562 dddy = sin*ddx+cos*ddy;
563 Double_t z4=fGRF->Eval(dddx,dddy); //point right
566 if (z0<0) {z0=0;z1=0;z2=0;z3=0;z4=0;}
568 Double_t f2x= (z3+z1-2*z0)*4.;//second derivation in y
569 Double_t f2y= (z2+z4-2*z0)*4.;//second derivation in x
570 Double_t f1y= (z3-z1);
572 z = (z0+f2x/6.+f2y/6.);//second order aproxiation of integral
573 if (kx>kprec){ //positive derivation
574 if (x<(xp1+dy*kx)){ //calculate volume at left border
576 Double_t xx2 = TMath::Min(x+dx,xp1+dy*kx);
577 Double_t yy1 = y+(xx1-xp1)/kx;
578 Double_t yy2 = TMath::Min(y+(xx2-xp1)/kx,y+dy);
581 z-= z0*(y+dy-yy2)/dy; //constant part rectangle
582 z-= f1y*(xx2-xx1)*(y+dy-yy2)*(y+dy-yy2)/(2.*dx*dy);
584 z-=z0*(xx2-xx1)*(yy2-yy1)/(2*dx*dy); //constant part rectangle
587 if (x>xp2){ //calculate volume at right border
590 Double_t yy1 = y+(xx1-xp2)/kx;
591 Double_t yy2 = y+(xx2-xp2)/kx;
594 z-=z0*(yy1-y)/dy; //constant part
595 z-=f1y*(xx2-xx1)*(yy1-y)*(yy1-y)/(2*dx*dy);
597 z-=z0*(xx2-xx1)*(yy2-yy1)/(2*dx*dy); //constant part
600 if (kx<-kprec){ //negative derivation
601 if (x<(xp1+dy*kx)){ //calculate volume at left border
603 Double_t xx2 = TMath::Min(x+dx,xp3-dy/kx);
604 Double_t yy1 = y+(xx1-xp1)/kx;
605 Double_t yy2 = TMath::Max(y,yy1+(xx2-xx1)/kx); //yy2<yy1
607 z-= z0*(yy2-y)/dy; // constant part rectangle
608 z-= f1y*(xx2-xx1)*(yy2-y)*(yy2-y)/(2.*dx*dy);
609 z-=z0*(xx2-xx1)*(yy1-yy2)/(2*dx*dy); //constant part triangle
611 if (x>xp2){ //calculate volume at right border
612 Double_t xx1 = TMath::Max(x,xp2+dy*kx);
614 Double_t yy1 = TMath::Min(y+dy,y-(xp2-xx1)/kx);
615 Double_t yy2 = y-(xp2-xx2)/kx;
617 z-=z0*(yy2-y)/dy; //constant part rextangle
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
623 if (z>0.) sumch+=fKNorm*z*dx*dy/fInteg;
632 }//step over different y
636 }//step over different points on line NPRF
639 void AliTPCPRF2D::UpdateSigma()
642 //calulate effective sigma X and sigma y of PRF
652 for (i=-1; i<=fNYdiv; i++){
653 if (fNYdiv == 1) y = fY1;
655 y = fY1+Float_t(i)*(fY2-fY1)/Float_t(fNYdiv-1);
656 for (x =-fNPRF*fDStep; x<fNPRF*fDStep;x+=fDStep)
658 //x in cm fWidth in cm
659 Float_t weight = GetPRF(x,y);
670 fSigmaX = TMath::Sqrt(fSigmaX/sum-fMeanX*fMeanX);
671 fSigmaY = TMath::Sqrt(fSigmaY/sum-fMeanY*fMeanY);
677 void AliTPCPRF2D::Streamer(TBuffer &R__b)
679 // Stream an object of class AliTPCPRF2D
681 if (R__b.IsReading()) {
682 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
683 TObject::Streamer(R__b);
684 //read chewron parameters
694 //read charge parameters
695 R__b.ReadFastArray(fType,5);
701 R__b >> fPadDistance;
703 //read angle parameters
704 R__b >> fChargeAngle;
711 if (strncmp(fType,"User",3)==0){
715 if (strncmp(fType,"Gauss",3)==0)
716 fGRF = new TF2("fun",funGauss2D,-5.,5.,-5.,5.,4);
717 if (strncmp(fType,"Cosh",3)==0)
718 fGRF = new TF2("fun",funCosh2D,-5.,5.,-5.,5.,4);
719 if (strncmp(fType,"Gati",3)==0)
720 fGRF = new TF2("fun",funGati2D,-5.,5.,-5.,5.,5);
721 //read interpolation parameters
727 R__b >>fNChargeArray;
728 if (fChargeArray!=0) delete [] fChargeArray;
729 if (fNChargeArray>0) {
730 fChargeArray = new Float_t[fNChargeArray];
731 R__b.ReadFastArray(fChargeArray,fNChargeArray);
734 R__b.ReadFastArray(funParam,5);
735 if (fGRF!=0) fGRF->SetParameters(funParam);
736 //calculate conversion coefitient to convert position to virtual wire
737 fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1);
740 R__b.WriteVersion(AliTPCPRF2D::IsA());
741 TObject::Streamer(R__b);
742 //write chewron parameters
752 //write charge parameters
753 R__b.WriteFastArray(fType,5);
759 R__b << fPadDistance;
762 R__b << fChargeAngle;
764 if (strncmp(fType,"User",3)==0) R__b <<fGRF;
765 //write interpolation parameters
771 R__b <<fNChargeArray;
773 R__b.WriteFastArray(fChargeArray,fNPRF*fNYdiv);
774 R__b.WriteFastArray(funParam,5);
779 TH1F * AliTPCPRF2D::GenerDrawXHisto(Float_t x1, Float_t x2,Float_t y)
781 //gener one dimensional hist of pad response function
785 sprintf(s,"Pad Response Function");
786 TH1F * hPRFc = new TH1F("hPRFc",s,kn+1,x1,x2);
790 for (Int_t i = 0;i<kn+1;i++)
792 x+=(x2-x1)/Float_t(kn);
796 hPRFc->SetXTitle("pad (cm)");
800 AliH2F * AliTPCPRF2D::GenerDrawHisto(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx, Int_t Ny)
803 //gener two dimensional histogram with PRF
806 sprintf(s,"Pad Response Function");
807 AliH2F * hPRFc = new AliH2F("hPRFc",s,Nx,x1,x2,Ny,y1,y2);
808 Float_t dx=(x2-x1)/Float_t(Nx);
809 Float_t dy=(y2-y1)/Float_t(Ny) ;
813 for ( Int_t i = 0;i<=Nx;i++,x+=dx){
815 for (Int_t j = 0;j<=Ny;j++,y+=dy){
817 hPRFc->SetCellContent(i,j,z);
820 hPRFc->SetXTitle("pad direction (cm)");
821 hPRFc->SetYTitle("pad row direction (cm)");
822 hPRFc->SetTitleOffset(1.5,"X");
823 hPRFc->SetTitleOffset(1.5,"Y");
828 AliH2F * AliTPCPRF2D::GenerDrawDistHisto(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx, Int_t Ny, Float_t thr)
830 //return histogram with distortion
831 const Float_t kminth=0.00001;
832 if (thr<kminth) thr=kminth;
834 sprintf(s,"COG distortion of PRF (threshold=%2.2f)",thr);
835 AliH2F * hPRFDist = new AliH2F("hDistortion",s,Nx,x1,x2,Ny,y1,y2);
836 Float_t dx=(x2-x1)/Float_t(Nx);
837 Float_t dy=(y2-y1)/Float_t(Ny) ;
840 for ( Int_t i = 0;i<=Nx;i++,x+=dx){
842 for(Int_t j = 0;j<=Ny;j++,y+=dy)
846 for (Int_t k=-3;k<=3;k++)
848 Float_t padx=Float_t(k)*fWidth;
849 z = GetPRF(x-padx,y);
857 ddx = (x-(sumx/sum));
860 if (TMath::Abs(ddx)<10) hPRFDist->SetCellContent(i,j,ddx);
864 hPRFDist->SetXTitle("pad direction (cm)");
865 hPRFDist->SetYTitle("pad row direction (cm)");
866 hPRFDist->SetTitleOffset(1.5,"X");
867 hPRFDist->SetTitleOffset(1.5,"Y");
875 void AliTPCPRF2D::DrawX(Float_t x1 ,Float_t x2,Float_t y1,Float_t y2, Int_t N)
878 //draw pad response function at interval <x1,x2> at given y position
881 TCanvas * c1 = new TCanvas("PRFX","Pad response function",700,900);
884 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
885 comment->SetTextAlign(12);
886 comment->SetFillColor(42);
887 DrawComment(comment);
891 TPad * pad2 = new TPad("pPRF","",0.05,0.22,0.95,0.95);
892 pad2->Divide(2,(N+1)/2);
894 gStyle->SetOptFit(1);
895 gStyle->SetOptStat(1);
896 for (Int_t i=0;i<N;i++){
900 else y = y1+i*(y2-y1)/Float_t(N-1);
902 TH1F * hPRFc =GenerDrawXHisto(x1, x2,y);
903 sprintf(ch,"PRF at wire position: %2.3f",y);
905 sprintf(ch,"PRF %d",i);
914 void AliTPCPRF2D::DrawPRF(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2, Int_t Nx, Int_t Ny)
918 TCanvas * c1 = new TCanvas("canPRF","Pad response function",700,900);
920 TPad * pad2 = new TPad("pad2PRF","",0.05,0.22,0.95,0.95);
922 gStyle->SetOptFit(1);
923 gStyle->SetOptStat(1);
924 TH2F * hPRFc = GenerDrawHisto(x1, x2, y1, y2, Nx,Ny);
928 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
929 comment->SetTextAlign(12);
930 comment->SetFillColor(42);
931 DrawComment(comment);
935 void AliTPCPRF2D::DrawDist(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2, Int_t Nx, Int_t Ny, Float_t thr)
938 //draw distortion of the COG method - for different threshold parameter
939 TCanvas * c1 = new TCanvas("padDistortion","COG distortion",700,900);
941 TPad * pad1 = new TPad("dist","",0.05,0.55,0.95,0.95,21);
943 TPad * pad2 = new TPad("dist","",0.05,0.22,0.95,0.53,21);
945 gStyle->SetOptFit(1);
946 gStyle->SetOptStat(0);
948 AliH2F * hPRFDist = GenerDrawDistHisto(x1, x2, y1, y2, Nx,Ny,thr);
951 hPRFDist->Draw("surf");
952 Float_t distmax =hPRFDist->GetMaximum();
953 Float_t distmin =hPRFDist->GetMinimum();
954 gStyle->SetOptStat(1);
956 TH1F * dist = hPRFDist->GetAmplitudes(distmin,distmax,distmin-1);
960 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
961 comment->SetTextAlign(12);
962 comment->SetFillColor(42);
963 DrawComment(comment);
967 void AliTPCPRF2D::DrawComment(TPaveText *comment)
970 //function to write comment to picture
973 //draw comments to picture
974 TText * title = comment->AddText("Pad Response Function parameters:");
975 title->SetTextSize(0.03);
976 sprintf(s,"Height of pad: %2.2f cm",fHeightFull);
978 sprintf(s,"Width pad: %2.2f cm",fWidth);
980 sprintf(s,"Pad Angle: %2.2f ",fPadAngle);
983 if (TMath::Abs(fK)>0.0001){
984 sprintf(s,"Height of one chevron unit h: %2.2f cm",2*fHeightS);
986 sprintf(s,"Overlap factor: %2.2f",fK);
990 if (strncmp(fType,"User",3)==0){
991 sprintf(s,"Charge distribution - user defined function %s ",fGRF->GetTitle());
993 sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
995 sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
998 if (strncmp(fType,"Gauss",3)==0){
999 sprintf(s,"Gauss charge distribution");
1000 comment->AddText(s);
1001 sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
1002 comment->AddText(s);
1003 sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
1004 comment->AddText(s);
1006 if (strncmp(fType,"Gati",3)==0){
1007 sprintf(s,"Gati charge distribution");
1008 comment->AddText(s);
1009 sprintf(s,"K3X of Gati : %2.2f ",fK3X);
1010 comment->AddText(s);
1011 sprintf(s,"K3Y of Gati: %2.2f ",fK3Y);
1012 comment->AddText(s);
1013 sprintf(s,"Wire to Pad Distance: %2.2f ",fPadDistance);
1014 comment->AddText(s);
1016 if (strncmp(fType,"Cosh",3)==0){
1017 sprintf(s,"Cosh charge distribution");
1018 comment->AddText(s);
1019 sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
1020 comment->AddText(s);
1021 sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
1022 comment->AddText(s);
1024 sprintf(s,"Normalisation: %2.2f ",fKNorm);
1025 comment->AddText(s);