]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCPRF2D.cxx
58e739d5990a87f67250e699a8ef15a344032b0f
[u/mrichter/AliRoot.git] / TPC / AliTPCPRF2D.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /*
17 $Log$
18 */
19
20 ///////////////////////////////////////////////////////////////////////////////
21 //  AliTPCPRF2D -                                                                         //
22 //  Pad response function object in two dimesions                            //
23 //  This class contains the basic functions for the                          //
24 //  calculation of PRF according generic charge distribution                 //
25 //  In Update function object calculate table of response function           //
26 //  in discrete x and y position                                             //
27 // This table is used for interpolation od response function in any position //
28 // (function GetPRF)                                                          //
29 //                                                                           // 
30 //  Origin: Marian Ivanov, Uni. of Bratislava, ivanov@fmph.uniba.sk          //
31 //                                                                           //
32 ///////////////////////////////////////////////////////////////////////////////
33 #include "TMath.h"
34 #include "AliTPCPRF2D.h"
35 #include "TF2.h"
36 #include <iostream.h>
37 #include <string.h>
38 #include "TCanvas.h"
39 #include "TPad.h"
40 #include "TStyle.h"
41 #include "TH1.h"
42 #include "TH2.h"
43 #include "TPaveText.h"
44 #include "TText.h"
45
46 extern TStyle * gStyle;
47
48 static const Float_t sqrt12=3.46;
49 static const Int_t   NPRF = 100;
50
51
52 static Double_t funGauss2D(Double_t *x, Double_t * par)
53
54   return ( TMath::Exp(-(x[0]*x[0])/(2*par[0]*par[0]))*
55            TMath::Exp(-(x[1]*x[1])/(2*par[1]*par[1])));
56
57 }
58
59 static Double_t funCosh2D(Double_t *x, Double_t * par)
60 {
61   return ( 1/(TMath::CosH(3.14159*x[0]/(2*par[0]))*
62            TMath::CosH(3.14159*x[1]/(2*par[1]))));
63 }    
64
65 static Double_t funGati2D(Double_t *x, Double_t * par)
66 {
67   //par[1] = is equal to k3X
68   //par[0] is equal to pad wire distance
69   Float_t K3=par[1];
70   Float_t K3R=TMath::Sqrt(K3);
71   Float_t K2=(TMath::Pi()/2)*(1-K3R/2.);
72   Float_t K1=K2*K3R/(4*TMath::ATan(K3R));
73   Float_t l=x[0]/par[0];
74   Float_t tan2=TMath::TanH(K2*l);
75   tan2*=tan2;
76   Float_t res = K1*(1-tan2)/(1+K3*tan2);
77  //par[4] = is equal to k3Y
78   K3=par[4];
79   K3R=TMath::Sqrt(K3);
80   K2=(TMath::Pi()/2)*(1-K3R/2.);
81   K1=K2*K3R/(4*TMath::ATan(K3R));
82   l=x[1]/par[0];
83   tan2=TMath::TanH(K2*l);
84   tan2*=tan2;
85   res = res*K1*(1-tan2)/(1+K3*tan2);  
86   return res;  
87 }   
88
89
90 ///////////////////////////////////////////////////////////////////////////
91 ///////////////////////////////////////////////////////////////////////////
92 ///////////////////////////////////////////////////////////////////////////
93 ///////////////////////////////////////////////////////////////////////////
94
95 ClassImp(AliTPCPRF2D)
96
97 AliTPCPRF2D::AliTPCPRF2D()
98 {
99   ffcharge = 0;
100   fNPRF =NPRF ;
101   fSigmaX = 0;
102
103   fGRF = 0;
104   fkNorm = 1;
105   forigsigmaY=0;
106   forigsigmaX=0;
107   fNdiv = 5;
108   //chewron default values   
109   SetPad(0.8,0.8);
110   SetChevron(0.2,0.0,1.0);
111   SetY(-0.2,0.2,2);
112   // SetGauss(0.22,0.22,1);  
113 }
114
115 AliTPCPRF2D::~AliTPCPRF2D()
116 {
117   if (ffcharge!=0) delete [] ffcharge;
118   if (fGRF !=0 ) fGRF->Delete();
119 }
120
121 void AliTPCPRF2D::SetY(Float_t y1, Float_t y2, Int_t nYdiv)
122 {
123   //
124   //set virtual line position
125   //first and last line and number of lines
126   fNYdiv = nYdiv;
127   if (ffcharge!=0) delete [] ffcharge;
128   ffcharge = new Float_t[fNPRF*fNYdiv];
129   fY1=y1;
130   fY2=y2;
131 }
132
133 void AliTPCPRF2D::SetPad(Float_t width, Float_t height)
134 {
135   //set base chevron parameters
136  fHeightFull=height;
137  fWidth=width;
138 }
139 void AliTPCPRF2D::SetChevron(Float_t hstep, 
140                         Float_t shifty, 
141                         Float_t fac)
142 {
143   //set shaping of chewron parameters
144   fHeightS=hstep;
145   fShiftY=shifty;
146   fK=fWidth*fac/hstep;
147 }
148
149 void AliTPCPRF2D::SetChParam(Float_t width, Float_t height,
150                   Float_t hstep, Float_t shifty, Float_t fac)
151 {
152   SetPad(width,height);
153   SetChevron(hstep,shifty,fac);
154 }
155
156
157 Float_t AliTPCPRF2D::GetPRF(Float_t xin, Float_t yin, Bool_t inter)
158 {
159   if (ffcharge==0) return 0;
160   //  Float_t y=Float_t(fNYdiv-1)*(yin-fY1)/(fY2-fY1);
161   //transform position to "wire position"
162   Float_t y=fDYtoWire*(yin-fY1);
163   if (fNYdiv == 1) y=fY1;
164   //normaly it find nearest line charge
165   if (inter ==kFALSE){   
166     Int_t i=Int_t(0.5+y);
167     if (y<0) i=Int_t(-0.5+y);
168     if ((i<0) || (i>=fNYdiv) ) return 0;
169     fcharge   = &(ffcharge[i*fNPRF]);
170     return GetPRFActiv(xin);
171   }
172   else{
173     //make interpolation from more fore lines
174     Int_t i= Int_t(y);
175     if ((i<0) || (i>=fNYdiv) ) return 0;
176     Float_t z0=0;
177     Float_t z1=0;
178     Float_t z2=0;
179     Float_t z3=0;
180     if (i>0) {
181       fcharge =&(ffcharge[(i-1)*fNPRF]);
182       z0 = GetPRFActiv(xin);
183     }
184     fcharge =&(ffcharge[i*fNPRF]);
185     z1=GetPRFActiv(xin);
186     if ((i+1)<fNYdiv){
187       fcharge =&(ffcharge[(i+1)*fNPRF]);
188       z2 = GetPRFActiv(xin);
189     }
190     if ((i+2)<fNYdiv){
191       fcharge =&(ffcharge[(i+2)*fNPRF]);
192       z3 = GetPRFActiv(xin);
193     }
194     Float_t a,b,c,d,K,L;
195     a=z1;
196     b=(z2-z0)/2.;
197     K=z2-a-b;
198     L=(z3-z1)/2.-b;
199     d=L-2*K;
200     c=K-d;
201     Float_t dy=y-Float_t(i);
202         Float_t res = a+b*dy+c*dy*dy+d*dy*dy*dy;  
203         //Float_t res = z1*(1-dy)+z2*dy;
204     return res;            
205   }        
206   return 0.;
207
208
209
210 Float_t AliTPCPRF2D::GetPRFActiv(Float_t xin)
211 {
212   //x xin DStep unit
213   //return splaine aproximaton 
214   Float_t x = (xin*fDStepM1)+fNPRF/2;
215   Int_t i = Int_t(x);
216   
217   if  ( (i>0) && ((i+2)<fNPRF)) {
218     Float_t a,b,c,d,K,L;
219     a = fcharge[i];
220     b = (fcharge[i+1]-fcharge[i-1])*0.5; 
221     K = fcharge[i+1]-a-b;
222     L = (fcharge[i+2]-fcharge[i])*0.5-b;
223     d=L-2.*K;
224     c=K-d;
225     Float_t dx=x-Float_t(i);
226     Float_t res = a+b*dx+c*dx*dx+d*dx*dx*dx;  
227     return res;
228   }
229   else return 0;
230 }
231
232
233 Float_t  AliTPCPRF2D::GetGRF(Float_t xin, Float_t yin)
234 {  
235   if (fGRF != 0 ) 
236     return fkNorm*fGRF->Eval(xin,yin)/fInteg;
237       else
238     return 0.;
239 }
240
241    
242 void AliTPCPRF2D::SetParam( TF2 * GRF,  Float_t kNorm, 
243                        Float_t sigmaX, Float_t sigmaY)
244 {
245    if (fGRF !=0 ) fGRF->Delete();
246    fGRF = GRF;
247    fkNorm = kNorm;
248    if (sigmaX ==0) sigmaX=(fWidth+fK*fHeightS)/sqrt12;
249    if (sigmaY ==0) sigmaY=(fWidth+fK*fHeightS)/sqrt12;
250    forigsigmaX=sigmaX; 
251    forigsigmaY=sigmaY; 
252    fDStep = TMath::Sqrt(sigmaX*sigmaX+fWidth*fWidth/6.)/10.; 
253    //   Update();   
254   sprintf(fType,"User");
255 }
256   
257
258 void AliTPCPRF2D::SetGauss(Float_t sigmaX, Float_t sigmaY,
259                       Float_t kNorm)
260 {
261   fkNorm = kNorm;
262   if (fGRF !=0 ) fGRF->Delete();
263   fGRF = new TF2("fun",funGauss2D,-5.,5.,-5.,5.,4);
264   funParam[0]=sigmaX;
265   funParam[1]=sigmaY;  
266   funParam[2]=fK;
267   funParam[3]=fHeightS;    
268   forigsigmaX=sigmaX;
269   forigsigmaY=sigmaY;
270   fGRF->SetParameters(funParam);
271   fDStep = TMath::Sqrt(sigmaX*sigmaX+fWidth*fWidth/6.)/10.; 
272   //by default I set the step as one tenth of sigma
273   //Update();
274   sprintf(fType,"Gauss");
275 }
276
277 void AliTPCPRF2D::SetCosh(Float_t sigmaX, Float_t sigmaY,
278                      Float_t kNorm)
279 {
280   fkNorm = kNorm;
281   if (fGRF !=0 ) fGRF->Delete();
282   fGRF = new TF2("fun", funCosh2D,-5.,5.,-5.,5.,4);   
283   funParam[0]=sigmaX;
284   funParam[1]=sigmaY;
285   funParam[2]=fK;  
286   funParam[3]=fHeightS;
287   fGRF->SetParameters(funParam);
288   forigsigmaX=sigmaX;
289   forigsigmaY=sigmaY;
290   fDStep = TMath::Sqrt(sigmaX*sigmaX+fWidth*fWidth/6.)/10.; 
291   //by default I set the step as one tenth of sigma
292   //Update();
293   sprintf(fType,"Cosh");
294 }
295
296 void AliTPCPRF2D::SetGati(Float_t K3X, Float_t K3Y,
297                      Float_t padDistance,
298                      Float_t kNorm)
299 {
300   fkNorm = kNorm;
301   if (fGRF !=0 ) fGRF->Delete();
302   fGRF = new TF2("fun", funGati2D,-5.,5.,-5.,5.,5);  
303   fK3X=K3X;
304   fK3Y=K3Y;
305   fPadDistance=padDistance;
306   funParam[0]=padDistance;
307   funParam[1]=K3X;
308   funParam[2]=fK;  
309   funParam[3]=fHeightS;
310   funParam[4]=K3Y;
311   fGRF->SetParameters(funParam);
312   forigsigmaX=padDistance;
313   forigsigmaY=padDistance;
314   fDStep = TMath::Sqrt(padDistance*padDistance+fWidth*fWidth/6.)/10.; 
315   //by default I set the step as one tenth of sigma
316   //Update();
317   sprintf(fType,"Gati");
318 }
319
320
321
322 void AliTPCPRF2D::Update()
323 {
324   for (Int_t i=0; i<fNYdiv; i++){
325     if (fNYdiv == 1) fActualY = fY1;
326     else
327       fActualY = fY1+Float_t(i)*(fY2-fY1)/Float_t(fNYdiv-1);
328     fcharge   = &(ffcharge[i*fNPRF]);
329     Update1();
330   }
331 }
332
333
334
335 void AliTPCPRF2D::Update1()
336 {
337   //initialize to 0
338   
339
340   Int_t i;
341   Float_t x;
342   for (i =0; i<fNPRF;i++)  fcharge[i] = 0;
343   if ( fGRF == 0 ) return;
344   ////////////////////////////////////////////////////////
345   //I'm waiting for normal integral
346   //in this moment only sum
347   Float_t x2=  4*forigsigmaX;
348   Float_t y2=  4*forigsigmaY;
349   Float_t dx = forigsigmaX/Float_t(fNdiv*6);
350   Float_t dy = forigsigmaY/Float_t(fNdiv*6);  
351   fInteg  = 0;
352   for (x=0.;x<x2;x+=dx)
353     for (Float_t y=0;y<y2;y+=dy) fInteg+=fGRF->Eval(x,y)*dx*dy;
354   fInteg*=4;
355   /////////////////////////////////////////////////////
356       
357   
358   if ( fInteg == 0 ) fInteg = 1; 
359   
360     //integrate charge over pad for different distance of pad
361     for (i =0; i<fNPRF;i++)
362       {      //x in cm fWidth in cm
363         //calculate integral 
364         Float_t xch = fDStep * (Float_t)(i-fNPRF/2);
365         Float_t k=1;
366         fcharge[i]=0;
367         for (Float_t y=-fHeightFull/2.-fShiftY;
368              y<fHeightFull/2.;y+=fHeightS){
369           Float_t y2=TMath::Min((y+fHeightS),Float_t(fHeightFull/2.));
370           Float_t y1=TMath::Max((y),Float_t(-fHeightFull/2.));
371           Float_t x1;
372         
373           if (k>0) 
374             x1 = (y2-y1)*fK-(fWidth+fK*fHeightS)/2.;      
375           else
376             x1 =-(fWidth+fK*fHeightS)/2. ;        
377           Float_t x2=x1+fWidth;
378
379           if (y2>y1) {
380             
381             if ((x2-x1)*fNdiv<forigsigmaX) dx=(x2-x1);
382             else{
383               dx= forigsigmaX/Float_t(fNdiv);
384               dx = (x2-x1)/Float_t(Int_t(3+(x2-x1)/dx));          
385             }       
386             Float_t dy;
387             if ((y2-y1)*fNdiv<forigsigmaY) dy=(y2-y1);
388             else{             
389               dy= forigsigmaY/Float_t(fNdiv);
390               dy = (y2-y1)/Float_t(Int_t(3+(y2-y1)/dy));
391             }
392
393             for (x=x1;x<x2;x+=dx)
394               for (Float_t y=y1;y<y2;y+=dy){
395                 if ( (y>(fActualY-(4.0*forigsigmaY))) &&
396                      (y<(fActualY+(4.0*forigsigmaY)))){
397                   Float_t xt=x-k*fK*(y-y1); 
398                   if ((TMath::Abs(xch-xt)<4*forigsigmaX)){
399                     
400                     Float_t z0=fGRF->Eval(xch-(xt+dx/2.),fActualY-(y+dy/2.));
401                     
402                     Float_t z1=fGRF->Eval(xch-(xt+dx/2.),fActualY-y);
403                     Float_t z2=fGRF->Eval(xch-xt,fActualY-(y+dy/2.));
404                     Float_t z3=fGRF->Eval(xch-(xt-dx/2.),fActualY-y);
405                     Float_t z4=fGRF->Eval(xch-xt,fActualY-(y-dy/2.));
406                     if (z0<0) z0=0;
407                     if (z1<0) z1=0;
408                     if (z2<0) z2=0;
409                     if (z3<0) z3=0;
410                     if (z4<0) z4=0;
411                     
412                     //        Float_t a=(z1-z3)/2;
413                     //        Float_t b=(z2-z4)/2;
414                     Float_t c= (z3+z1-2*z0)/2.;
415                     Float_t d= (z2+z4-2*z0)/2.;
416                     Float_t z= (z0+c/12.+d/12.);                                
417                     
418                     //Float_t z= fGRF->Eval(xch-xt,fActualY-y);
419                     if (z>0.)         fcharge[i]+=z*dx*dy/fInteg;             
420                   }
421                 }
422               }
423           }
424           k*=-1;
425         }
426       };   
427   
428   fSigmaX = 0; 
429   Float_t sum =0;
430   Float_t mean=0;
431   for (x =-fNPRF*fDStep; x<fNPRF*fDStep;x+=fDStep)
432     {      //x in cm fWidth in cm
433       Float_t weight = GetPRFActiv(x);
434       fSigmaX+=x*x*weight; 
435       mean+=x*weight;
436       sum+=weight;
437     };  
438   if (sum>0){
439     mean/=sum;
440     fSigmaX = TMath::Sqrt(fSigmaX/sum-mean*mean);   
441   }
442   else fSigmaX=0; 
443   //calculate conversion coefitient to convert position to virtual wire
444   fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1);
445   fDStepM1=1/fDStep;
446 }
447
448 void AliTPCPRF2D::Streamer(TBuffer &R__b)
449 {
450    // Stream an object of class AliTPCPRF2D
451
452    if (R__b.IsReading()) {
453       Version_t R__v = R__b.ReadVersion(); if (R__v) { }
454       TObject::Streamer(R__b);     
455       //read chewron parameters
456       R__b >> fSigmaX;
457       R__b >> fHeightFull;
458       R__b >> fHeightS;
459       R__b >> fShiftY;
460       R__b >> fWidth;
461       R__b >> fK;
462       R__b >> fActualY;
463       //read charge parameters
464       R__b >> fType[0];
465       R__b >> fType[1];
466       R__b >> fType[2];
467       R__b >> fType[3];
468       R__b >> fType[4];
469       R__b >> forigsigmaX;
470       R__b >> forigsigmaY;
471       R__b >> fkNorm;
472       R__b >> fK3X;
473       R__b >> fK3Y;
474       R__b >> fPadDistance;
475       R__b >> fInteg;
476       
477       //read functions
478       if (fGRF!=0) { 
479         delete [] fGRF;  
480         fGRF=0;
481       }
482       if (strncmp(fType,"User",3)==0){
483         fGRF= new TF2;
484         R__b>>fGRF;   
485       }
486       if (strncmp(fType,"Gauss",3)==0) 
487         fGRF = new TF2("fun",funGauss2D,-5.,5.,-5.,5.,4);
488       if (strncmp(fType,"Cosh",3)==0) 
489         fGRF = new TF2("fun",funCosh2D,-5.,5.,-5.,5.,4);
490        if (strncmp(fType,"Gati",3)==0) 
491         fGRF = new TF2("fun",funGati2D,-5.,5.,-5.,5.,5);
492       
493       //read interpolation parameters
494       R__b >>fY1;
495       R__b >>fY2;
496       R__b >>fNYdiv;  
497       R__b >>fDStep;  
498       R__b >>fNPRF;
499       if (ffcharge!=0) delete [] ffcharge;
500       ffcharge = new Float_t[fNPRF*fNYdiv];
501       R__b.ReadFastArray(ffcharge,fNPRF*fNYdiv); 
502       R__b.ReadFastArray(funParam,5); 
503       if (fGRF!=0) fGRF->SetParameters(funParam);
504       //calculate conversion coefitient to convert position to virtual wire
505       fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1);
506       fDStepM1=1/fDStep;
507    } else {
508       R__b.WriteVersion(AliTPCPRF2D::IsA());
509       TObject::Streamer(R__b);      
510       //write chewron parameters
511       R__b << fSigmaX;
512       R__b << fHeightFull;
513       R__b << fHeightS;
514       R__b << fShiftY;
515       R__b << fWidth;
516       R__b << fK;
517       R__b << fActualY;
518       //write charge parameters
519       R__b << fType[0];
520       R__b << fType[1];
521       R__b << fType[2];
522       R__b << fType[3];
523       R__b << fType[4];
524
525       R__b << forigsigmaX;
526       R__b << forigsigmaY;
527       R__b << fkNorm;
528       R__b << fK3X;
529       R__b << fK3Y;
530       R__b << fPadDistance;  
531       R__b << fInteg;
532
533       if (strncmp(fType,"User",3)==0)   R__b <<fGRF;         
534       //write interpolation parameters
535       R__b <<fY1;
536       R__b <<fY2;
537       R__b <<fNYdiv;   
538       R__b <<fDStep;
539       R__b <<fNPRF;    
540       R__b.WriteFastArray(ffcharge,fNPRF*fNYdiv); 
541       R__b.WriteFastArray(funParam,5); 
542    }
543 }
544
545
546
547
548 void AliTPCPRF2D::DrawX(Float_t x1 ,Float_t x2,Float_t y, Bool_t inter)
549
550   if (fGRF==0) return ;
551   const Int_t N=100;
552   char s[100];
553   TCanvas  * c1 = new TCanvas("canPRF","Pad response function",700,900);
554   c1->cd();
555   TPad * pad1 = new TPad("pad1PRF","",0.05,0.61,0.95,0.97,21);
556   pad1->Draw();
557   TPad * pad2 = new TPad("pad2PRF","",0.05,0.22,0.95,0.60,21);
558   pad2->Draw();
559
560   //  pad1->cd();  
561   //pad2->cd();
562   gStyle->SetOptFit(1);
563   gStyle->SetOptStat(0); 
564   sprintf(s,"PRF response function for chevron pad");  
565   TH1F * hPRFc = new TH1F("hPRFc",s,N+1,x1,x2);
566   Float_t x=x1;
567   Float_t y1;
568   //  Float_t y2;
569
570   for (Float_t i = 0;i<N+1;i++)
571     {
572       x+=(x2-x1)/Float_t(N);
573       y1 = GetPRF(x,y,inter);
574       hPRFc->Fill(x,y1);
575     };
576
577   pad1->cd();
578   fGRF->SetRange(x1,x1,x2,x2); 
579   fGRF->SetNpx(25);
580   fGRF->SetNpy(25); 
581   fGRF->Draw("lego2");
582   // hPRFo->Fit("gaus");
583   gStyle->SetOptStat(1); 
584   pad2->cd();
585   hPRFc->Fit("gaus");
586   c1->cd(); 
587   TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
588   comment->SetTextAlign(12);
589   comment->SetFillColor(42);
590   TText *title = comment->AddText("Chevron pad parameters:");
591   title->SetTextSize(0.03);
592   sprintf(s,"Full height of pad:  %2.2f",fHeightFull);
593   comment->AddText(s);
594   sprintf(s,"Height of one chevron unit h:  %2.2f cm",2*fHeightS);
595   comment->AddText(s);
596   sprintf(s,"Width of one chevron unit  w:  %2.2f cm",fWidth);
597   comment->AddText(s);
598   sprintf(s,"Overlap factor:  %2.2f",fK*fHeightS/fWidth);
599   comment->AddText(s);
600   sprintf(s,"Y position:  %2.2f ",y);
601   comment->AddText(s);
602   sprintf(s,"Sigma x of original distribution: %2.2f ",forigsigmaX);
603   comment->AddText(s);  
604   sprintf(s,"Sigma y of original distribution: %2.2f ",forigsigmaY);
605   comment->AddText(s);    
606   sprintf(s,"Type of original distribution: %s ",fType);
607   comment->AddText(s); 
608   comment->Draw();
609 }
610
611
612
613 void AliTPCPRF2D::Draw(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2, 
614                   Bool_t inter, Int_t Nx, Int_t Ny)
615
616   char s[100];
617   if (fGRF==0) return ;
618   TCanvas  * c1 = new TCanvas("canPRF","Pad response function",700,900);
619   c1->cd();
620   TPad * pad1 = new TPad("pad1PRF","",0.05,0.61,0.95,0.97,21);
621   pad1->Draw();
622   TPad * pad2 = new TPad("pad2PRF","",0.05,0.22,0.95,0.60,21);
623   pad2->Draw();
624
625   //  pad1->cd();  
626   //pad2->cd();
627   gStyle->SetOptFit(1);
628   gStyle->SetOptStat(0); 
629   sprintf(s,"PRF response function for chevron pad");  
630   TH2F * hPRFc = new TH2F("hPRFc",s,Nx+1,x1,x2,Ny+1,y1,y2);
631   Float_t dx=(x2-x1)/Float_t(Nx);
632   Float_t dy=(y2-y1)/Float_t(Ny) ;
633   Float_t x,y,z;
634   //  Float_t y2;
635   for ( x = x1;x<=x2;x+=dx){
636     for(y = y1;y<=y2;y+=dy)
637       {
638         z = GetPRF(x,y,inter);
639         hPRFc->Fill(x,y,z);
640       };
641   }
642   pad1->cd();
643   fGRF->SetRange(x1,y1,x2,y2); 
644   fGRF->SetNpx(25);
645   fGRF->SetNpy(25); 
646   fGRF->Draw("lego2");
647   // hPRFo->Fit("gaus");
648   gStyle->SetOptStat(1); 
649   pad2->cd();
650   hPRFc->Draw("lego2");
651   c1->cd(); 
652   TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
653   comment->SetTextAlign(12);
654   comment->SetFillColor(42);
655   TText *title = comment->AddText("Chevron pad parameters:");
656   title->SetTextSize(0.03);
657   sprintf(s,"Full height of pad:  %2.2f",fHeightFull);
658   comment->AddText(s);
659   sprintf(s,"Height of one chevron unit h:  %2.2f cm",2*fHeightS);
660   comment->AddText(s);
661   sprintf(s,"Width of one chevron unit  w:  %2.2f cm",fWidth);
662   comment->AddText(s);
663   sprintf(s,"Overlap factor:  %2.2f",fK*fHeightS/fWidth);
664   comment->AddText(s); 
665   sprintf(s,"Sigma x of original distribution: %2.2f ",forigsigmaX);
666   comment->AddText(s);  
667   sprintf(s,"Sigma y of original distribution: %2.2f ",forigsigmaY);
668   comment->AddText(s);    
669   sprintf(s,"Type of original distribution: %s ",fType);
670   comment->AddText(s); 
671   comment->Draw();
672 }
673
674 void AliTPCPRF2D::DrawDist(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2, 
675                   Bool_t inter, Int_t Nx, Int_t Ny, Float_t thr)
676
677   const Float_t minth=0.00001;
678   if (thr<minth) thr=minth;
679   char s[100];
680   if (fGRF==0) return ;
681   TCanvas  * c1 = new TCanvas("padDistortion","COG distortion",700,900);
682   c1->cd();
683   TPad * pad1 = new TPad("CHARGE","",0.05,0.61,0.95,0.97,21);
684   pad1->Draw();
685   TPad * pad2 = new TPad("dist","",0.05,0.22,0.95,0.60,21);
686   pad2->Draw();
687
688   //  pad1->cd();  
689   //pad2->cd();
690   gStyle->SetOptFit(1);
691   gStyle->SetOptStat(0); 
692   sprintf(s,"COG distortion (threshold=%2.2f)",thr);  
693   TH2F * hPRFDist = new TH2F("hDistortion",s,Nx+1,x1,x2,Ny+1,y1,y2);
694   Float_t dx=(x2-x1)/Float_t(Nx);
695   Float_t dy=(y2-y1)/Float_t(Ny) ;
696   Float_t x,y,z,ddx;
697   //  Float_t y2;
698   for ( x = x1;x<(x2+dx/2.);x+=dx)
699     for(y = y1;y<=(y2+dx/2.);y+=dy)
700       {
701         Float_t sumx=0;
702         Float_t sum=0;
703         for (Float_t padx=-fWidth;padx<(fWidth*1.1);padx+=fWidth)
704           {         
705             z = GetPRF(x-padx,y,inter);
706             if (z>thr){
707               sum+=z;
708               sumx+=z*padx;
709             }   
710           };    
711         if (sum>minth)  
712           {
713             ddx = (x-(sumx/sum));
714           }
715         else ddx=-1;
716         if (TMath::Abs(ddx)<10)         hPRFDist->Fill(x,y,ddx);
717       }
718   pad1->cd();
719   fGRF->SetRange(x1,y1,x2,y2); 
720   fGRF->SetNpx(25);
721   fGRF->SetNpy(25); 
722   fGRF->Draw("lego2");
723   // hPRFo->Fit("gaus");
724   //  gStyle->SetOptStat(1); 
725   pad2->cd();
726   hPRFDist->Draw("lego2");
727   
728   c1->cd(); 
729   TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
730   comment->SetTextAlign(12);
731   comment->SetFillColor(42);
732   //  TText *title = comment->AddText("Distortion of COG method");
733   //  title->SetTextSize(0.03);
734   TText * title = comment->AddText("Chevron pad parameters:");
735   title->SetTextSize(0.03);
736   sprintf(s,"Full height of pad:  %2.2f",fHeightFull);
737   comment->AddText(s);
738   sprintf(s,"Height of one chevron unit h:  %2.2f cm",2*fHeightS);
739   comment->AddText(s);
740   sprintf(s,"Width of one chevron unit  w:  %2.2f cm",fWidth);
741   comment->AddText(s);
742   sprintf(s,"Overlap factor:  %2.2f",fK*fHeightS/fWidth);
743   comment->AddText(s); 
744   sprintf(s,"Sigma x of original distribution: %2.2f ",forigsigmaX);
745   comment->AddText(s);  
746   sprintf(s,"Sigma y of original distribution: %2.2f ",forigsigmaY);
747   comment->AddText(s);    
748   sprintf(s,"Type of original distribution: %s ",fType);
749   comment->AddText(s); 
750   comment->Draw();
751   
752 }
753