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