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