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