AliRoot Core  ee782a0 (ee782a0)
AliTPCPRF2D.cxx
Go to the documentation of this file.
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 
28 
29 #include <Riostream.h>
30 #include <TCanvas.h>
31 #include <TClass.h>
32 #include <TF2.h>
33 #include <TH1.h>
34 #include <TMath.h>
35 #include <TPad.h>
36 #include <TPaveText.h>
37 #include <TStyle.h>
38 #include <TText.h>
39 #include <string.h>
40 #include "TBuffer.h"
41 
42 #include "AliH2F.h"
43 #include "AliTPCPRF2D.h"
44 
45 
46 extern TStyle * gStyle;
47 
48 const Double_t AliTPCPRF2D::fgkDegtoRad = 0.01745329251994;
49 const Double_t AliTPCPRF2D::fgkSQRT12=3.464101;
50 const Int_t AliTPCPRF2D::fgkNPRF = 100;
51 
52 
53 static Double_t FunGauss2D(const Double_t *const x, const Double_t *const par)
54 {
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(const Double_t *const x, const Double_t *const par)
63 {
65 
66  return ( 1/(TMath::CosH(3.14159*x[0]/(2*par[0]))*
67  TMath::CosH(3.14159*x[1]/(2*par[1]))));
68 }
69 
70 static Double_t FunGati2D(const Double_t *const x, const Double_t *const par)
71 {
73 
74  Float_t k3=par[1];
75  Float_t k3R=TMath::Sqrt(k3);
76  Float_t k2=(TMath::Pi()/2)*(1-k3R/2.);
77  Float_t k1=k2*k3R/(4*TMath::ATan(k3R));
78  Float_t l=x[0]/par[0];
79  Float_t tan2=TMath::TanH(k2*l);
80  tan2*=tan2;
81  Float_t res = k1*(1-tan2)/(1+k3*tan2);
82  //par[4] = is equal to k3Y
83  k3=par[4];
84  k3R=TMath::Sqrt(k3);
85  k2=(TMath::Pi()/2)*(1-k3R/2.);
86  k1=k2*k3R/(4*TMath::ATan(k3R));
87  l=x[1]/par[0];
88  tan2=TMath::TanH(k2*l);
89  tan2*=tan2;
90  res = res*k1*(1-tan2)/(1+k3*tan2);
91  return res;
92 }
93 
96 
98 ClassImp(AliTPCPRF2D)
100 
102  :TObject(),
103  fcharge(0),
104  fY1(0.),
105  fY2(0.),
106  fNYdiv(0),
107  fNChargeArray(0),
108  fChargeArray(0),
109  fHeightFull(0.),
110  fHeightS(0.),
111  fShiftY(0.),
112  fWidth(0.),
113  fK(0.),
114  fNPRF(0),
115  fNdiv(5),
116  fDStep(0.),
117  fKNorm(1.),
118  fInteg(0.),
119  fGRF(0),
120  fK3X(0.),
121  fK3Y(0.),
122  fPadDistance(0.),
123  fOrigSigmaX(0.),
124  fOrigSigmaY(0.),
125  fChargeAngle(0.),
126  fPadAngle(0.),
127  fSigmaX(0.),
128  fSigmaY(0.),
129  fMeanX(0.),
130  fMeanY(0.),
131  fInterX(0),
132  fInterY(0),
133  fCurrentY(0.),
134  fDYtoWire(0.),
135  fDStepM1(0.)
136 {
137  //default constructor for response function object
138 
139  fNPRF =fgkNPRF ;
140  for(Int_t i=0;i<5;i++){
141  funParam[i]=0.;
142  fType[i]=0;
143  }
144 
145 
146  //chewron default values
147  SetPad(0.8,0.8);
148  SetChevron(0.2,0.0,1.0);
149  SetY(-0.2,0.2,2);
150  SetInterpolationType(2,0);
151 }
152 
154 {
155  if (fChargeArray!=0) delete [] fChargeArray;
156  if (fGRF !=0 ) fGRF->Delete();
157 }
158 
159 void AliTPCPRF2D::SetY(Float_t y1, Float_t y2, Int_t nYdiv)
160 {
163 
164  fNYdiv = nYdiv;
165  fY1=y1;
166  fY2=y2;
167 }
168 
169 void AliTPCPRF2D::SetPad(Float_t width, Float_t height)
170 {
172 
173  fHeightFull=height;
174  fWidth=width;
175 }
176 void AliTPCPRF2D::SetChevron(Float_t hstep,
177  Float_t shifty,
178  Float_t fac)
179 {
181 
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 {
200 
201  if (fChargeArray==0) return 0;
202  //transform position to "wire position"
203  Float_t y=fDYtoWire*(yin-fY1);
204  if (fNYdiv == 1) y=fY1;
205  //normaly it find nearest line charge
206  if (fInterY ==0){
207  Int_t i=Int_t(0.5+y);
208  if (y<0) i=Int_t(-0.5+y);
209  if ((i<0) || (i>=fNYdiv) ) return 0;
210  fcharge = &(fChargeArray[i*fNPRF]);
211  return GetPRFActiv(xin);
212  }
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 
248 
249 Float_t AliTPCPRF2D::GetPRFActiv(Float_t xin)
250 {
253 
254  Float_t x = (xin*fDStepM1)+fNPRF/2;
255  Int_t i = Int_t(x);
256 
257  if ( (i>1) && ((i+2)<fNPRF)) {
258  Float_t a,b,c,d,k,l;
259  a = fcharge[i];
260  b = (fcharge[i+1]-fcharge[i-1])*0.5;
261  k = fcharge[i+1]-a-b;
262  l = (fcharge[i+2]-fcharge[i])*0.5-b;
263  d=l-2.*k;
264  c=k-d;
265  Float_t dx=x-Float_t(i);
266  Float_t res = a+b*dx+c*dx*dx+d*dx*dx*dx;
267  return res;
268  }
269  else return 0;
270 }
271 
272 
273 Float_t AliTPCPRF2D::GetGRF(Float_t xin, Float_t yin)
274 {
277 
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 *const GRF, Float_t kNorm,
286  Float_t sigmaX, Float_t sigmaY)
287 {
290 
291  if (fGRF !=0 ) fGRF->Delete();
292  fGRF = GRF;
293  fKNorm = kNorm;
294  //sprintf(fType,"User");
295  snprintf(fType,5,"User");
296  if (sigmaX ==0) sigmaX=(fWidth*(1+TMath::Abs(fK)))/fgkSQRT12;
297  if (sigmaY ==0) sigmaY=(fWidth*(1+TMath::Abs(fK)))/fgkSQRT12;
298  fOrigSigmaX=sigmaX;
299  fOrigSigmaY=sigmaY;
300  Double_t estimsigma =
301  TMath::Sqrt(sigmaX*sigmaX+(fWidth*fWidth*(1+TMath::Abs(fK))/12)+
302  TMath::Tan(fPadAngle*fgkDegtoRad)*TMath::Tan(fPadAngle*fgkDegtoRad)*fHeightFull*fHeightFull/12);
303  if (estimsigma < 5*sigmaX) {
304  fDStep = estimsigma/10.;
305  fNPRF = Int_t(estimsigma*8./fDStep);
306  }
307  else{
308  fDStep = sigmaX;
309  Double_t width = fWidth*(1+TMath::Abs(fK))+TMath::Abs(TMath::Tan(fPadAngle*fgkDegtoRad))*fHeightFull;
310  fNPRF = Int_t((width+8.*sigmaX)/fDStep);
311  };
312 
313 }
314 
315 
316 void AliTPCPRF2D::SetGauss(Float_t sigmaX, Float_t sigmaY,
317  Float_t kNorm)
318 {
320 
321  fKNorm = kNorm;
322  fOrigSigmaX=sigmaX;
323  fOrigSigmaY=sigmaY;
324  //sprintf(fType,"Gauss");
325  snprintf(fType,5,"Gauss");
326  if (fGRF !=0 ) fGRF->Delete();
327  fGRF = new TF2("FunGauss2D",FunGauss2D,-5.,5.,-5.,5.,4);
328 
329  funParam[0]=sigmaX;
330  funParam[1]=sigmaY;
331  funParam[2]=fK;
332  funParam[3]=fHeightS;
333 
334  fGRF->SetParameters(funParam);
335  Double_t estimsigma =
336  TMath::Sqrt(sigmaX*sigmaX+(fWidth*fWidth*(1+TMath::Abs(fK))/12)+
337  TMath::Tan(fPadAngle)*TMath::Tan(fPadAngle*fgkDegtoRad)*fHeightFull*fHeightFull/12);
338  if (estimsigma < 5*sigmaX) {
339  fDStep = estimsigma/10.;
340  fNPRF = Int_t(estimsigma*8./fDStep);
341  }
342  else{
343  fDStep = sigmaX;
344  Double_t width = fWidth*(1+TMath::Abs(fK))+TMath::Abs(TMath::Tan(fPadAngle*fgkDegtoRad))*fHeightFull;
345  fNPRF = Int_t((width+8.*sigmaX)/fDStep);
346  };
347 
348 
349 }
350 void AliTPCPRF2D::SetCosh(Float_t sigmaX, Float_t sigmaY,
351  Float_t kNorm)
352 {
354 
355  fKNorm = kNorm;
356  fOrigSigmaX=sigmaX;
357  fOrigSigmaY=sigmaY;
358  // sprintf(fType,"Cosh");
359  snprintf(fType,5,"Cosh");
360  if (fGRF !=0 ) fGRF->Delete();
361  fGRF = new TF2("FunCosh2D", FunCosh2D,-5.,5.,-5.,5.,4);
362  funParam[0]=sigmaX;
363  funParam[1]=sigmaY;
364  funParam[2]=fK;
365  funParam[3]=fHeightS;
366  fGRF->SetParameters(funParam);
367 
368  Double_t estimsigma = TMath::Sqrt(sigmaX*sigmaX+fWidth*fWidth*(1+TMath::Abs(fK))/12);
369  if (estimsigma < 5*sigmaX) {
370  fDStep = estimsigma/10.;
371  fNPRF = Int_t(estimsigma*8./fDStep);
372  }
373  else{
374  fDStep = sigmaX;
375  fNPRF = Int_t((1.2*fWidth*(1+TMath::Abs(fK))+8.*sigmaX)/fDStep);
376  };
377 
378 }
379 
380 void AliTPCPRF2D::SetGati(Float_t K3X, Float_t K3Y,
381  Float_t padDistance,
382  Float_t kNorm)
383 {
385 
386  fKNorm = kNorm;
387  fK3X=K3X;
388  fK3Y=K3Y;
389  fPadDistance=padDistance;
390  //sprintf(fType,"Gati");
391  snprintf(fType,5,"Gati");
392  if (fGRF !=0 ) fGRF->Delete();
393  fGRF = new TF2("FunGati2D", FunGati2D,-5.,5.,-5.,5.,5);
394 
395  funParam[0]=padDistance;
396  funParam[1]=K3X;
397  funParam[2]=fK;
398  funParam[3]=fHeightS;
399  funParam[4]=K3Y;
400  fGRF->SetParameters(funParam);
401  fOrigSigmaX=padDistance;
402  fOrigSigmaY=padDistance;
403  Float_t sigmaX = fOrigSigmaX;
404  Double_t estimsigma = TMath::Sqrt(sigmaX*sigmaX+fWidth*fWidth*(1+TMath::Abs(fK))/12);
405  if (estimsigma < 5*sigmaX) {
406  fDStep = estimsigma/10.;
407  fNPRF = Int_t(estimsigma*8./fDStep);
408  }
409  else{
410  fDStep = sigmaX;
411  fNPRF = Int_t((1.2*fWidth*(1+TMath::Abs(fK))+8.*sigmaX)/fDStep);
412  };
413 }
414 
415 
416 
418 {
421 
422  if ( fGRF == 0 ) return;
423  //initialize interpolated values to 0
424  Int_t i;
425  if (fChargeArray!=0) delete [] fChargeArray;
426  fChargeArray = new Float_t[fNPRF*fNYdiv];
428  for (i =0; i<fNPRF*fNYdiv;i++) fChargeArray[i] = 0;
429  //firstly calculate total integral of charge
430 
432  //I'm waiting for normal integral
433  //in this moment only sum
434  Float_t x2= 4*fOrigSigmaX;
435  Float_t y2= 4*fOrigSigmaY;
436  Float_t dx = fOrigSigmaX/Float_t(fNdiv*6);
437  Float_t dy = fOrigSigmaY/Float_t(fNdiv*6);
438  Int_t nx = Int_t(0.5+x2/dx);
439  Int_t ny = Int_t(0.5+y2/dy);
440  Int_t ix,iy;
441  fInteg = 0;
442  Double_t dInteg =0;
443  for (ix=-nx;ix<=nx;ix++)
444  for ( iy=-ny;iy<=ny;iy++)
445  dInteg+=fGRF->Eval(Float_t(ix)*dx,Float_t(iy)*dy)*dx*dy;
447  fInteg =dInteg;
448  if ( fInteg == 0 ) fInteg = 1;
449 
450  for (i=0; i<fNYdiv; i++){
451  if (fNYdiv == 1) fCurrentY = fY1;
452  else
453  fCurrentY = fY1+Double_t(i)*(fY2-fY1)/Double_t(fNYdiv-1);
454  fcharge = &(fChargeArray[i*fNPRF]);
455  Update1();
456  }
457  //calculate conversion coefitient to convert position to virtual wire
458  fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1);
459  fDStepM1=1/fDStep;
460  UpdateSigma();
461 }
462 
464 {
467 
468  Int_t i;
469  Double_t cos = TMath::Cos(fChargeAngle);
470  Double_t sin = TMath::Sin(fChargeAngle);
471  const Double_t kprec =0.00000001;
472  //integrate charge over pad for different distance of pad
473  for (i =0; i<fNPRF;i++){
474  //x in cm fWidth in cm
475  //calculate integral
476  Double_t xch = fDStep * (Double_t)(i-fNPRF/2);
477  fcharge[i]=0;
478  Double_t k=1;
479 
480 
481  for (Double_t ym=-fHeightFull/2.-fShiftY; ym<fHeightFull/2.-kprec;ym+=fHeightS){
482  Double_t y2chev=TMath::Min((ym+fHeightS),Double_t(fHeightFull/2.)); // end of chevron step
483  Double_t y1chev= ym; //beginning of chevron step
484  Double_t y2 = TMath::Min(y2chev,fCurrentY+3.5*fOrigSigmaY);
485  Double_t y1 = TMath::Max((y1chev),Double_t(-fHeightFull/2.));
486  y1 = TMath::Max(y1chev,fCurrentY-3.5*fOrigSigmaY);
487 
488  Double_t x0 = fWidth*(-1.-(Double_t(k)*fK))*0.5+ym*TMath::Tan(fPadAngle*fgkDegtoRad);
489  Double_t kx = Double_t(k)*(fK*fWidth)/fHeightS;
490  kx = TMath::Tan(TMath::ATan(kx))+TMath::Tan(fPadAngle*fgkDegtoRad);
491 
492  Int_t ny = TMath::Max(Int_t(fNdiv*TMath::Exp(-(y1-fCurrentY)*(y1-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),4);
493  Double_t dy = TMath::Min(fOrigSigmaY/Double_t(ny),y2-y1);
494  Double_t ndy = dy;
495 
496  //loop over different y strips with variable step size dy
497  if (y2>(y1+kprec)) for (Double_t y = y1; y<y2+kprec;){
498  //new step SIZE
499 
500  ny = TMath::Max(Int_t(fNdiv*TMath::Exp(-(y-fCurrentY)*(y-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),5);
501  ndy = fOrigSigmaY/Double_t(ny);
502  if (ndy>(y2-y-dy)) {
503  ndy =y2-y-dy;
504  if (ndy<kprec) ndy=2*kprec; //calculate new delta y
505  }
506  //
507  Double_t sumch=0;
508  //calculation of x borders and initial step
509  Double_t deltay = (y-y1chev);
510 
511  Double_t xp1 = x0+deltay*kx;
512  //x begining of pad at position y
513  Double_t xp2 =xp1+fWidth; //x end of pad at position y
514  Double_t xp3 =xp1+kx*dy; //...at position y+dy
515  Double_t xp4 =xp2+kx*dy; //..
516 
517  Double_t x1 = TMath::Min(xp1,xp3);
518  x1 = TMath::Max(xp1,xch-3.5*fOrigSigmaX); //beging of integration
519  Double_t x2 = TMath::Max(xp2,xp4);
520  x2 = TMath::Min(xp2+dy*kx,xch+3.5*fOrigSigmaX); //end of integration
521 
522  Int_t nx = TMath::Max(Int_t(fNdiv*TMath::Exp(-(x1-xch)*(x1-xch)/(2*fOrigSigmaX*fOrigSigmaX))*
523  TMath::Exp(-(y1-fCurrentY)*(y1-fCurrentY)/(2*fOrigSigmaY*fOrigSigmaY))),2);
524  Double_t dx = TMath::Min(fOrigSigmaX/Double_t(nx),x2-x1)/5.; //on the border more iteration
525  Double_t ndx=dx;
526 
527  if (x2>(x1+kprec)) {
528  for (Double_t x = x1; x<x2+kprec ;){
529  //new step SIZE
530  nx = TMath::Max(Int_t(fNdiv*TMath::Exp(-(x-xch)*(x-xch)/(2*fOrigSigmaX*fOrigSigmaX))),3);
531  ndx = fOrigSigmaX/Double_t(nx);
532  if (ndx>(x2-x-dx)) {
533  ndx =x2-x-dx;
534  }
535  if ( ( (x+dx+ndx)<TMath::Max(xp3,xp1)) || ( (x+dx+ndx)>TMath::Min(xp4,xp2))) {
536  ndx/=5.;
537  }
538  if (ndx<kprec) ndx=2*kprec;
539  //INTEGRAL APROXIMATION
540  Double_t ddx,ddy,dddx,dddy;
541  ddx = xch-(x+dx/2.);
542  ddy = fCurrentY-(y+dy/2.);
543  dddx = cos*ddx-sin*ddy;
544  dddy = sin*ddx+cos*ddy;
545  Double_t z0=fGRF->Eval(dddx,dddy); //middle point
546 
547  ddx = xch-(x+dx/2.);
548  ddy = fCurrentY-(y);
549  dddx = cos*ddx-sin*ddy;
550  dddy = sin*ddx+cos*ddy;
551  Double_t z1=fGRF->Eval(dddx,dddy); //point down
552 
553  ddx = xch-(x+dx/2.);
554  ddy = fCurrentY-(y+dy);
555  dddx = cos*ddx-sin*ddy;
556  dddy = sin*ddx+cos*ddy;
557  Double_t z3=fGRF->Eval(dddx,dddy); //point up
558 
559  ddx = xch-(x);
560  ddy = fCurrentY-(y+dy/2.);
561  dddx = cos*ddx-sin*ddy;
562  dddy = sin*ddx+cos*ddy;
563  Double_t z2=fGRF->Eval(dddx,dddy); //point left
564 
565  ddx = xch-(x+dx);
566  ddy = fCurrentY-(y+dy/2.);
567  dddx = cos*ddx-sin*ddy;
568  dddy = sin*ddx+cos*ddy;
569  Double_t z4=fGRF->Eval(dddx,dddy); //point right
570 
571 
572  if (z0<0) {z0=0;z1=0;z2=0;z3=0;z4=0;}
573 
574  Double_t f2x= (z3+z1-2*z0)*4.;//second derivation in y
575  Double_t f2y= (z2+z4-2*z0)*4.;//second derivation in x
576  Double_t f1y= (z3-z1);
577  Double_t z ;
578  z = (z0+f2x/6.+f2y/6.);//second order aproxiation of integral
579  if (kx>kprec){ //positive derivation
580  if (x<(xp1+dy*kx)){ //calculate volume at left border
581  Double_t xx1 = x;
582  Double_t xx2 = TMath::Min(x+dx,xp1+dy*kx);
583  Double_t yy1 = y+(xx1-xp1)/kx;
584  Double_t yy2 = TMath::Min(y+(xx2-xp1)/kx,y+dy);
585  z=z0;
586  if (yy2<y+dy) {
587  z-= z0*(y+dy-yy2)/dy; //constant part rectangle
588  z-= f1y*(xx2-xx1)*(y+dy-yy2)*(y+dy-yy2)/(2.*dx*dy);
589  }
590  z-=z0*(xx2-xx1)*(yy2-yy1)/(2*dx*dy); //constant part rectangle
591 
592  }
593  if (x>xp2){ //calculate volume at right border
594  Double_t xx1 = x;
595  Double_t xx2 = x+dx;
596  Double_t yy1 = y+(xx1-xp2)/kx;
597  Double_t yy2 = y+(xx2-xp2)/kx;
598  z=z0;
599  //rectangle part
600  z-=z0*(yy1-y)/dy; //constant part
601  z-=f1y*(xx2-xx1)*(yy1-y)*(yy1-y)/(2*dx*dy);
602  //triangle part
603  z-=z0*(xx2-xx1)*(yy2-yy1)/(2*dx*dy); //constant part
604  }
605  }
606  if (kx<-kprec){ //negative derivation
607  if (x<(xp1+dy*kx)){ //calculate volume at left border
608  Double_t xx1 = x;
609  Double_t xx2 = TMath::Min(x+dx,xp3-dy/kx);
610  Double_t yy1 = y+(xx1-xp1)/kx;
611  Double_t yy2 = TMath::Max(y,yy1+(xx2-xx1)/kx); //yy2<yy1
612  z = z0;
613  z-= z0*(yy2-y)/dy; // constant part rectangle
614  z-= f1y*(xx2-xx1)*(yy2-y)*(yy2-y)/(2.*dx*dy);
615  z-=z0*(xx2-xx1)*(yy1-yy2)/(2*dx*dy); //constant part triangle
616  }
617  if (x>xp2){ //calculate volume at right border
618  Double_t xx1 = TMath::Max(x,xp2+dy*kx);
619  Double_t xx2 = x+dx;
620  Double_t yy1 = TMath::Min(y+dy,y-(xp2-xx1)/kx);
621  Double_t yy2 = y-(xp2-xx2)/kx;
622  z=z0;
623  z-=z0*(yy2-y)/dy; //constant part rextangle
624  z-= f1y*(xx2-xx1)*(yy2-y)*(yy2-y)/(2.*dx*dy);
625  z-=z0*(xx2-xx1)*(yy1-yy2)/(2*dx*dy); //constant part triangle
626  }
627  }
628 
629  if (z>0.) sumch+=fKNorm*z*dx*dy/fInteg;
630 
631  x+=dx;
632  dx = ndx;
633  }; //loop over x
634  fcharge[i]+=sumch;
635  }//if x2>x1
636  y+=dy;
637  dy =ndy;
638  }//step over different y
639  k*=-1.;
640  }//step over chevron
641 
642  }//step over different points on line NPRF
643 }
644 
646 {
648 
649  fMeanX = 0;
650  fMeanY = 0;
651  fSigmaX = 0;
652  fSigmaY = 0;
653 
654  Float_t sum =0;
655  Int_t i;
656  Float_t x,y;
657 
658  for (i=-1; i<=fNYdiv; i++){
659  if (fNYdiv == 1) y = fY1;
660  else
661  y = fY1+Float_t(i)*(fY2-fY1)/Float_t(fNYdiv-1);
662  for (x =-fNPRF*fDStep; x<fNPRF*fDStep;x+=fDStep)
663  {
664  //x in cm fWidth in cm
665  Float_t weight = GetPRF(x,y);
666  fSigmaX+=x*x*weight;
667  fSigmaY+=y*y*weight;
668  fMeanX+=x*weight;
669  fMeanY+=y*weight;
670  sum+=weight;
671  };
672  }
673  if (sum>0){
674  fMeanX/=sum;
675  fMeanY/=sum;
676  fSigmaX = TMath::Sqrt(fSigmaX/sum-fMeanX*fMeanX);
677  fSigmaY = TMath::Sqrt(fSigmaY/sum-fMeanY*fMeanY);
678  }
679  else fSigmaX=0;
680 }
681 
682 
683 void AliTPCPRF2D::Streamer(TBuffer &xRuub)
684 {
686 
687  if (xRuub.IsReading()) {
688  UInt_t xRuus, xRuuc;
689  Version_t xRuuv = xRuub.ReadVersion(&xRuus, &xRuuc);
690  AliTPCPRF2D::Class()->ReadBuffer(xRuub, this, xRuuv, xRuus, xRuuc);
691  //read functions
692  if (strncmp(fType,"User",3)!=0){
693  delete fGRF;
694  if (strncmp(fType,"Gauss",3)==0)
695  fGRF = new TF2("FunGauss2D",FunGauss2D,-5.,5.,-5.,5.,4);
696  if (strncmp(fType,"Cosh",3)==0)
697  fGRF = new TF2("FunCosh2D",FunCosh2D,-5.,5.,-5.,5.,4);
698  if (strncmp(fType,"Gati",3)==0)
699  fGRF = new TF2("FunGati2D",FunGati2D,-5.,5.,-5.,5.,5);
700  if (fGRF!=0) fGRF->SetParameters(funParam);
701  }
702  //calculate conversion coefitient to convert position to virtual wire
703  fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1);
704  fDStepM1=1/fDStep;
705  } else {
706  AliTPCPRF2D::Class()->WriteBuffer(xRuub,this);
707  }
708 }
709 
710 
711 TH1F * AliTPCPRF2D::GenerDrawXHisto(Float_t x1, Float_t x2,Float_t y)
712 {
715 
716  char s[100];
717  const Int_t kn=200;
718  //sprintf(s,"Pad Response Function");
719  snprintf(s,100,"Pad Response Function");
720  TH1F * hPRFc = new TH1F("hPRFc",s,kn+1,x1,x2);
721  Float_t x=x1;
722  Float_t y1;
723 
724  for (Int_t i = 0;i<kn+1;i++)
725  {
726  x+=(x2-x1)/Float_t(kn);
727  y1 = GetPRF(x,y);
728  hPRFc->Fill(x,y1);
729  };
730  hPRFc->SetXTitle("pad (cm)");
731  return hPRFc;
732 }
733 
734 AliH2F * AliTPCPRF2D::GenerDrawHisto(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx, Int_t Ny)
735 {
737 
738  char s[100];
739  //sprintf(s,"Pad Response Function");
740  snprintf(s,100,"Pad Response Function");
741  AliH2F * hPRFc = new AliH2F("hPRFc",s,Nx,x1,x2,Ny,y1,y2);
742  Float_t dx=(x2-x1)/Float_t(Nx);
743  Float_t dy=(y2-y1)/Float_t(Ny) ;
744  Float_t x,y,z;
745  x = x1;
746  y = y1;
747  for ( Int_t i = 0;i<=Nx;i++,x+=dx){
748  y=y1;
749  for (Int_t j = 0;j<=Ny;j++,y+=dy){
750  z = GetPRF(x,y);
751  hPRFc->SetBinContent(hPRFc->GetBin(i,j),z);
752  };
753  };
754  hPRFc->SetXTitle("pad direction (cm)");
755  hPRFc->SetYTitle("pad row direction (cm)");
756  hPRFc->SetTitleOffset(1.5,"X");
757  hPRFc->SetTitleOffset(1.5,"Y");
758  return hPRFc;
759 }
760 
761 
762 AliH2F * AliTPCPRF2D::GenerDrawDistHisto(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx, Int_t Ny, Float_t thr)
763 {
765 
766  const Float_t kminth=0.00001;
767  if (thr<kminth) thr=kminth;
768  char s[100];
769  //sprintf(s,"COG distortion of PRF (threshold=%2.2f)",thr);
770  snprintf(s,100,"COG distortion of PRF (threshold=%2.2f)",thr);
771  AliH2F * hPRFDist = new AliH2F("hDistortion",s,Nx,x1,x2,Ny,y1,y2);
772  Float_t dx=(x2-x1)/Float_t(Nx);
773  Float_t dy=(y2-y1)/Float_t(Ny) ;
774  Float_t x,y,z,ddx;
775  x=x1;
776  for ( Int_t i = 0;i<=Nx;i++,x+=dx){
777  y=y1;
778  for(Int_t j = 0;j<=Ny;j++,y+=dy)
779  {
780  Float_t sumx=0;
781  Float_t sum=0;
782  for (Int_t k=-3;k<=3;k++)
783  {
784  Float_t padx=Float_t(k)*fWidth;
785  z = GetPRF(x-padx,y);
786  if (z>thr){
787  sum+=z;
788  sumx+=z*padx;
789  }
790  };
791  if (sum>kminth)
792  {
793  ddx = (x-(sumx/sum));
794  }
795  else ddx=-1;
796  if (TMath::Abs(ddx)<10) hPRFDist->SetBinContent(hPRFDist->GetBin(i,j),ddx);
797  }
798  }
799 
800  hPRFDist->SetXTitle("pad direction (cm)");
801  hPRFDist->SetYTitle("pad row direction (cm)");
802  hPRFDist->SetTitleOffset(1.5,"X");
803  hPRFDist->SetTitleOffset(1.5,"Y");
804  return hPRFDist;
805 }
806 
807 
808 
809 
810 
811 void AliTPCPRF2D::DrawX(Float_t x1 ,Float_t x2,Float_t y1,Float_t y2, Int_t N)
812 {
814 
815  if (N<0) return;
816  TCanvas * c1 = new TCanvas("PRFX","Pad response function",700,900);
817  c1->cd();
818 
819  TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
820  comment->SetTextAlign(12);
821  comment->SetFillColor(42);
822  DrawComment(comment);
823  comment->Draw();
824  c1->cd();
825 
826  TPad * pad2 = new TPad("pPRF","",0.05,0.22,0.95,0.95);
827  pad2->Divide(2,(N+1)/2);
828  pad2->Draw();
829  gStyle->SetOptFit(1);
830  gStyle->SetOptStat(1);
831  for (Int_t i=0;i<N;i++){
832  char ch[200];
833  Float_t y;
834  if (N==1) y=y1;
835  else y = y1+i*(y2-y1)/Float_t(N-1);
836  pad2->cd(i+1);
837  TH1F * hPRFc =GenerDrawXHisto(x1, x2,y);
838  //sprintf(ch,"PRF at wire position: %2.3f",y);
839  snprintf(ch,40,"PRF at wire position: %2.3f",y);
840  hPRFc->SetTitle(ch);
841  //sprintf(ch,"PRF %d",i);
842  snprintf(ch,15,"PRF %d",i);
843  hPRFc->SetName(ch);
844  hPRFc->Fit("gaus");
845  }
846 
847 }
848 
849 
850 
851 void AliTPCPRF2D::DrawPRF(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2, Int_t Nx, Int_t Ny)
852 {
854 
855  TCanvas * c1 = new TCanvas("canPRF","Pad response function",700,900);
856  c1->cd();
857  TPad * pad2 = new TPad("pad2PRF","",0.05,0.22,0.95,0.95);
858  pad2->Draw();
859  gStyle->SetOptFit(1);
860  gStyle->SetOptStat(1);
861  TH2F * hPRFc = GenerDrawHisto(x1, x2, y1, y2, Nx,Ny);
862  pad2->cd();
863  hPRFc->Draw("surf");
864  c1->cd();
865  TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
866  comment->SetTextAlign(12);
867  comment->SetFillColor(42);
868  DrawComment(comment);
869  comment->Draw();
870 }
871 
872 void AliTPCPRF2D::DrawDist(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2, Int_t Nx, Int_t Ny, Float_t thr)
873 {
875 
876  TCanvas * c1 = new TCanvas("padDistortion","COG distortion",700,900);
877  c1->cd();
878  TPad * pad1 = new TPad("dist","",0.05,0.55,0.95,0.95,21);
879  pad1->Draw();
880  TPad * pad2 = new TPad("dist","",0.05,0.22,0.95,0.53,21);
881  pad2->Draw();
882  gStyle->SetOptFit(1);
883  gStyle->SetOptStat(0);
884 
885  AliH2F * hPRFDist = GenerDrawDistHisto(x1, x2, y1, y2, Nx,Ny,thr);
886 
887  pad1->cd();
888  hPRFDist->Draw("surf");
889  Float_t distmax =hPRFDist->GetMaximum();
890  Float_t distmin =hPRFDist->GetMinimum();
891  gStyle->SetOptStat(1);
892 
893  TH1F * dist = hPRFDist->GetAmplitudes(distmin,distmax,distmin-1);
894  pad2->cd();
895  dist->Draw();
896  c1->cd();
897  TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
898  comment->SetTextAlign(12);
899  comment->SetFillColor(42);
900  DrawComment(comment);
901  comment->Draw();
902 }
903 
905 {
907 
908  char s[100];
909  //draw comments to picture
910  TText * title = comment->AddText("Pad Response Function parameters:");
911  title->SetTextSize(0.03);
912  //sprintf(s,"Height of pad: %2.2f cm",fHeightFull);
913  snprintf(s,100,"Height of pad: %2.2f cm",fHeightFull);
914  comment->AddText(s);
915  //sprintf(s,"Width pad: %2.2f cm",fWidth);
916  snprintf(s,100,"Width pad: %2.2f cm",fWidth);
917  comment->AddText(s);
918  //sprintf(s,"Pad Angle: %2.2f ",fPadAngle);
919  snprintf(s,100,"Pad Angle: %2.2f ",fPadAngle);
920  comment->AddText(s);
921 
922  if (TMath::Abs(fK)>0.0001){
923  //sprintf(s,"Height of one chevron unit h: %2.2f cm",2*fHeightS);
924  snprintf(s,100,"Height of one chevron unit h: %2.2f cm",2*fHeightS);
925  comment->AddText(s);
926  //sprintf(s,"Overlap factor: %2.2f",fK);
927  snprintf(s,100,"Overlap factor: %2.2f",fK);
928  comment->AddText(s);
929  }
930 
931  if (strncmp(fType,"User",3)==0){
932  //sprintf(s,"Charge distribution - user defined function %s ",fGRF->GetTitle());
933  snprintf(s,100,"Charge distribution - user defined function %s ",fGRF->GetTitle());
934  comment->AddText(s);
935  //sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
936  snprintf(s,100,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
937  comment->AddText(s);
938  //sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
939  snprintf(s,100,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
940  comment->AddText(s);
941  }
942  if (strncmp(fType,"Gauss",3)==0){
943  //sprintf(s,"Gauss charge distribution");
944  snprintf(s,100,"Gauss charge distribution");
945  comment->AddText(s);
946  //sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
947  snprintf(s,100,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
948  comment->AddText(s);
949  //sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
950  snprintf(s,100,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
951  comment->AddText(s);
952  }
953  if (strncmp(fType,"Gati",3)==0){
954  //sprintf(s,"Gati charge distribution");
955  snprintf(s,100,"Gati charge distribution");
956  comment->AddText(s);
957  //sprintf(s,"K3X of Gati : %2.2f ",fK3X);
958  snprintf(s,100,"K3X of Gati : %2.2f ",fK3X);
959  comment->AddText(s);
960  //sprintf(s,"K3Y of Gati: %2.2f ",fK3Y);
961  snprintf(s,100,"K3Y of Gati: %2.2f ",fK3Y);
962  comment->AddText(s);
963  //sprintf(s,"Wire to Pad Distance: %2.2f ",fPadDistance);
964  snprintf(s,100,"Wire to Pad Distance: %2.2f ",fPadDistance);
965  comment->AddText(s);
966  }
967  if (strncmp(fType,"Cosh",3)==0){
968  //sprintf(s,"Cosh charge distribution");
969  snprintf(s,100,"Cosh charge distribution");
970  comment->AddText(s);
971  //sprintf(s,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
972  snprintf(s,100,"Sigma x of charge distribution: %2.2f ",fOrigSigmaX);
973  comment->AddText(s);
974  //sprintf(s,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
975  snprintf(s,100,"Sigma y of charge distribution: %2.2f ",fOrigSigmaY);
976  comment->AddText(s);
977  }
978  //sprintf(s,"Normalisation: %2.2f ",fKNorm);
979  snprintf(s,100,"Normalisation: %2.2f ",fKNorm);
980  comment->AddText(s);
981 }
982 
TBrowser b
Definition: RunAnaESD.C:12
static Double_t FunGauss2D(const Double_t *const x, const Double_t *const par)
Definition: AliTPCPRF2D.cxx:53
Float_t fChargeAngle
&#39;angle&#39; of charge distribution refernce system to pad reference system
Definition: AliTPCPRF2D.h:106
Float_t fSigmaY
sigma Y of PAD response function
Definition: AliTPCPRF2D.h:110
Float_t fShiftY
shift of the step
Definition: AliTPCPRF2D.h:87
Float_t fY2
position of last virtual vire
Definition: AliTPCPRF2D.h:77
static Double_t FunGati2D(const Double_t *const x, const Double_t *const par)
Definition: AliTPCPRF2D.cxx:70
Float_t fK3X
KX parameter (only for Gati parametrization)
Definition: AliTPCPRF2D.h:99
Int_t fNYdiv
number of wires
Definition: AliTPCPRF2D.h:78
Float_t fDStepM1
! used in GetPRFActiv to make calculation faster
Definition: AliTPCPRF2D.h:120
TStyle * gStyle
virtual void DrawPRF(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx=20, Int_t Ny=20)
Float_t fSigmaX
sigma X of PAD response function
Definition: AliTPCPRF2D.h:109
static TString comment
Definition: ConfigCosmic.C:131
void DrawComment(TPaveText *comment)
static Double_t FunCosh2D(const Double_t *const x, const Double_t *const par)
Definition: AliTPCPRF2D.cxx:62
virtual void UpdateSigma()
char fType[5]
charge type
Definition: AliTPCPRF2D.h:117
static const Double_t fgkDegtoRad
numeric constant
Definition: AliTPCPRF2D.h:122
Float_t fOrigSigmaY
sigma of original distribution;
Definition: AliTPCPRF2D.h:104
Int_t fNdiv
number of division to calculate integral
Definition: AliTPCPRF2D.h:93
virtual void DrawX(Float_t x1, Float_t x2, Float_t y1, Float_t y2=0, Int_t N=1)
Float_t fDYtoWire
! used to make PRF calculation faster in GetPRF
Definition: AliTPCPRF2D.h:119
Float_t fHeightS
height of the one step
Definition: AliTPCPRF2D.h:86
virtual void SetGati(Float_t K3X, Float_t K3Y, Float_t padDistance, Float_t kNorm=1)
TH1F * GetAmplitudes(Float_t zmin, Float_t zmax, Float_t th=1., Float_t xmin=0, Float_t xmax=0, Float_t ymin=0, Float_t ymax=0)
Definition: AliH2F.cxx:231
Float_t fPadAngle
&#39;angle&#39; of the pad assymetry
Definition: AliTPCPRF2D.h:107
virtual ~AliTPCPRF2D()
void sum()
Int_t fNChargeArray
Definition: AliTPCPRF2D.h:79
virtual void DrawDist(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx=20, Int_t Ny=20, Float_t thr=0)
virtual void SetChParam(Float_t width, Float_t height, Float_t hstep, Float_t shifty, Float_t fac)
Float_t fOrigSigmaX
sigma of original distribution;
Definition: AliTPCPRF2D.h:103
Float_t fCurrentY
in reality we calculate PRF only for one fixed y
Definition: AliTPCPRF2D.h:118
Float_t fK
k factor of the chewron
Definition: AliTPCPRF2D.h:89
virtual void SetCosh(Float_t sigmaX, Float_t sigmaY, Float_t kNorm=1)
Float_t fMeanY
mean Y value
Definition: AliTPCPRF2D.h:112
Pad response function object in two dimesions.
Definition: AliTPCPRF2D.h:19
Float_t fPadDistance
pad anode distnce (only for Gati parametrisation)
Definition: AliTPCPRF2D.h:101
Int_t fInterY
interpolation in Y
Definition: AliTPCPRF2D.h:114
Float_t GetPRFActiv(Float_t xin)
virtual void SetPad(Float_t width, Float_t height)
TF2 * fGRF
charge distribution function
Definition: AliTPCPRF2D.h:97
virtual Float_t GetPRF(Float_t xin, Float_t yin)
virtual TF2 * GetGRF()
Definition: AliTPCPRF2D.h:26
Float_t fY1
position of first "virtual" vire
Definition: AliTPCPRF2D.h:76
AliH2F * GenerDrawDistHisto(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx=20, Int_t Ny=20, Float_t thr=0)
TH1F * GenerDrawXHisto(Float_t x1, Float_t x2, Float_t y)
Float_t fMeanX
mean X value
Definition: AliTPCPRF2D.h:111
virtual void Update()
Float_t fKNorm
normalisation factor of the charge integral
Definition: AliTPCPRF2D.h:95
Definition: AliH2F.h:16
Float_t fInteg
integral of GRF on +- infinity
Definition: AliTPCPRF2D.h:96
Float_t fK3Y
KY parameter (only for Gati parametrisation)
Definition: AliTPCPRF2D.h:100
AliH2F * GenerDrawHisto(Float_t x1, Float_t x2, Float_t y1, Float_t y2, Int_t Nx=20, Int_t Ny=20)
void Update1()
void SetParam(TF2 *const GRF, Float_t kNorm, Float_t sigmaX=0, Float_t sigmaY=0)
static const Double_t fgkSQRT12
numeric constant
Definition: AliTPCPRF2D.h:123
Float_t fHeightFull
height of the full pad
Definition: AliTPCPRF2D.h:85
virtual void SetGauss(Float_t sigmaX, Float_t sigmaY, Float_t kNorm=1)
void res(Char_t i)
Definition: Resolution.C:2
virtual void SetY(Float_t y1, Float_t y2, Int_t nYdiv)
virtual void SetChevron(Float_t hstep, Float_t shifty, Float_t fac)
static const Int_t fgkNPRF
default number of division
Definition: AliTPCPRF2D.h:124
Float_t * fcharge
! field with PRF
Definition: AliTPCPRF2D.h:75
Double_t funParam[5]
parameters of used charge function
Definition: AliTPCPRF2D.h:91
Float_t fDStep
element step for point
Definition: AliTPCPRF2D.h:94
Float_t fWidth
width of the pad
Definition: AliTPCPRF2D.h:88
Float_t * fChargeArray
pointer to array of arrays
Definition: AliTPCPRF2D.h:81
Int_t fNPRF
number of interpolations point
Definition: AliTPCPRF2D.h:92