1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 Revision 1.2 1999/09/29 09:24:34 fca
19 Introduction of the Copyright and cvs Log
23 ///////////////////////////////////////////////////////////////////////////////
25 // Pad response function object in two dimesions //
26 // This class contains the basic functions for the //
27 // calculation of PRF according generic charge distribution //
28 // In Update function object calculate table of response function //
29 // in discrete x and y position //
30 // This table is used for interpolation od response function in any position //
31 // (function GetPRF) //
33 // Origin: Marian Ivanov, Uni. of Bratislava, ivanov@fmph.uniba.sk //
35 ///////////////////////////////////////////////////////////////////////////////
37 #include "AliTPCPRF2D.h"
46 #include "TPaveText.h"
49 extern TStyle * gStyle;
51 static const Float_t sqrt12=3.46;
52 static const Int_t NPRF = 100;
55 static Double_t funGauss2D(Double_t *x, Double_t * par)
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])));
62 static Double_t funCosh2D(Double_t *x, Double_t * par)
64 return ( 1/(TMath::CosH(3.14159*x[0]/(2*par[0]))*
65 TMath::CosH(3.14159*x[1]/(2*par[1]))));
68 static Double_t funGati2D(Double_t *x, Double_t * par)
70 //par[1] = is equal to k3X
71 //par[0] is equal to pad wire distance
73 Float_t K3R=TMath::Sqrt(K3);
74 Float_t K2=(TMath::Pi()/2)*(1-K3R/2.);
75 Float_t K1=K2*K3R/(4*TMath::ATan(K3R));
76 Float_t l=x[0]/par[0];
77 Float_t tan2=TMath::TanH(K2*l);
79 Float_t res = K1*(1-tan2)/(1+K3*tan2);
80 //par[4] = is equal to k3Y
83 K2=(TMath::Pi()/2)*(1-K3R/2.);
84 K1=K2*K3R/(4*TMath::ATan(K3R));
86 tan2=TMath::TanH(K2*l);
88 res = res*K1*(1-tan2)/(1+K3*tan2);
93 ///////////////////////////////////////////////////////////////////////////
94 ///////////////////////////////////////////////////////////////////////////
95 ///////////////////////////////////////////////////////////////////////////
96 ///////////////////////////////////////////////////////////////////////////
100 AliTPCPRF2D::AliTPCPRF2D()
111 //chewron default values
113 SetChevron(0.2,0.0,1.0);
115 // SetGauss(0.22,0.22,1);
118 AliTPCPRF2D::~AliTPCPRF2D()
120 if (ffcharge!=0) delete [] ffcharge;
121 if (fGRF !=0 ) fGRF->Delete();
124 void AliTPCPRF2D::SetY(Float_t y1, Float_t y2, Int_t nYdiv)
127 //set virtual line position
128 //first and last line and number of lines
130 if (ffcharge!=0) delete [] ffcharge;
131 ffcharge = new Float_t[fNPRF*fNYdiv];
136 void AliTPCPRF2D::SetPad(Float_t width, Float_t height)
138 //set base chevron parameters
142 void AliTPCPRF2D::SetChevron(Float_t hstep,
146 //set shaping of chewron parameters
152 void AliTPCPRF2D::SetChParam(Float_t width, Float_t height,
153 Float_t hstep, Float_t shifty, Float_t fac)
155 SetPad(width,height);
156 SetChevron(hstep,shifty,fac);
160 Float_t AliTPCPRF2D::GetPRF(Float_t xin, Float_t yin, Bool_t inter)
162 if (ffcharge==0) return 0;
163 // Float_t y=Float_t(fNYdiv-1)*(yin-fY1)/(fY2-fY1);
164 //transform position to "wire position"
165 Float_t y=fDYtoWire*(yin-fY1);
166 if (fNYdiv == 1) y=fY1;
167 //normaly it find nearest line charge
169 Int_t i=Int_t(0.5+y);
170 if (y<0) i=Int_t(-0.5+y);
171 if ((i<0) || (i>=fNYdiv) ) return 0;
172 fcharge = &(ffcharge[i*fNPRF]);
173 return GetPRFActiv(xin);
176 //make interpolation from more fore lines
178 if ((i<0) || (i>=fNYdiv) ) return 0;
184 fcharge =&(ffcharge[(i-1)*fNPRF]);
185 z0 = GetPRFActiv(xin);
187 fcharge =&(ffcharge[i*fNPRF]);
190 fcharge =&(ffcharge[(i+1)*fNPRF]);
191 z2 = GetPRFActiv(xin);
194 fcharge =&(ffcharge[(i+2)*fNPRF]);
195 z3 = GetPRFActiv(xin);
204 Float_t dy=y-Float_t(i);
205 Float_t res = a+b*dy+c*dy*dy+d*dy*dy*dy;
206 //Float_t res = z1*(1-dy)+z2*dy;
212 Float_t AliTPCPRF2D::GetPRFActiv(Float_t xin)
215 //return splaine aproximaton
216 Float_t x = (xin*fDStepM1)+fNPRF/2;
219 if ( (i>0) && ((i+2)<fNPRF)) {
222 b = (fcharge[i+1]-fcharge[i-1])*0.5;
223 K = fcharge[i+1]-a-b;
224 L = (fcharge[i+2]-fcharge[i])*0.5-b;
227 Float_t dx=x-Float_t(i);
228 Float_t res = a+b*dx+c*dx*dx+d*dx*dx*dx;
235 Float_t AliTPCPRF2D::GetGRF(Float_t xin, Float_t yin)
238 return fkNorm*fGRF->Eval(xin,yin)/fInteg;
244 void AliTPCPRF2D::SetParam( TF2 * GRF, Float_t kNorm,
245 Float_t sigmaX, Float_t sigmaY)
247 if (fGRF !=0 ) fGRF->Delete();
250 if (sigmaX ==0) sigmaX=(fWidth+fK*fHeightS)/sqrt12;
251 if (sigmaY ==0) sigmaY=(fWidth+fK*fHeightS)/sqrt12;
254 fDStep = TMath::Sqrt(sigmaX*sigmaX+fWidth*fWidth/6.)/10.;
256 sprintf(fType,"User");
260 void AliTPCPRF2D::SetGauss(Float_t sigmaX, Float_t sigmaY,
264 if (fGRF !=0 ) fGRF->Delete();
265 fGRF = new TF2("fun",funGauss2D,-5.,5.,-5.,5.,4);
269 funParam[3]=fHeightS;
272 fGRF->SetParameters(funParam);
273 fDStep = TMath::Sqrt(sigmaX*sigmaX+fWidth*fWidth/6.)/10.;
274 //by default I set the step as one tenth of sigma
276 sprintf(fType,"Gauss");
279 void AliTPCPRF2D::SetCosh(Float_t sigmaX, Float_t sigmaY,
283 if (fGRF !=0 ) fGRF->Delete();
284 fGRF = new TF2("fun", funCosh2D,-5.,5.,-5.,5.,4);
288 funParam[3]=fHeightS;
289 fGRF->SetParameters(funParam);
292 fDStep = TMath::Sqrt(sigmaX*sigmaX+fWidth*fWidth/6.)/10.;
293 //by default I set the step as one tenth of sigma
295 sprintf(fType,"Cosh");
298 void AliTPCPRF2D::SetGati(Float_t K3X, Float_t K3Y,
303 if (fGRF !=0 ) fGRF->Delete();
304 fGRF = new TF2("fun", funGati2D,-5.,5.,-5.,5.,5);
307 fPadDistance=padDistance;
308 funParam[0]=padDistance;
311 funParam[3]=fHeightS;
313 fGRF->SetParameters(funParam);
314 forigsigmaX=padDistance;
315 forigsigmaY=padDistance;
316 fDStep = TMath::Sqrt(padDistance*padDistance+fWidth*fWidth/6.)/10.;
317 //by default I set the step as one tenth of sigma
319 sprintf(fType,"Gati");
324 void AliTPCPRF2D::Update()
326 for (Int_t i=0; i<fNYdiv; i++){
327 if (fNYdiv == 1) fActualY = fY1;
329 fActualY = fY1+Float_t(i)*(fY2-fY1)/Float_t(fNYdiv-1);
330 fcharge = &(ffcharge[i*fNPRF]);
337 void AliTPCPRF2D::Update1()
344 for (i =0; i<fNPRF;i++) fcharge[i] = 0;
345 if ( fGRF == 0 ) return;
346 ////////////////////////////////////////////////////////
347 //I'm waiting for normal integral
348 //in this moment only sum
349 Float_t x2= 4*forigsigmaX;
350 Float_t y2= 4*forigsigmaY;
351 Float_t dx = forigsigmaX/Float_t(fNdiv*6);
352 Float_t dy = forigsigmaY/Float_t(fNdiv*6);
354 for (x=0.;x<x2;x+=dx)
355 for (Float_t y=0;y<y2;y+=dy) fInteg+=fGRF->Eval(x,y)*dx*dy;
357 /////////////////////////////////////////////////////
360 if ( fInteg == 0 ) fInteg = 1;
362 //integrate charge over pad for different distance of pad
363 for (i =0; i<fNPRF;i++)
364 { //x in cm fWidth in cm
366 Float_t xch = fDStep * (Float_t)(i-fNPRF/2);
369 for (Float_t y=-fHeightFull/2.-fShiftY;
370 y<fHeightFull/2.;y+=fHeightS){
371 Float_t y2=TMath::Min((y+fHeightS),Float_t(fHeightFull/2.));
372 Float_t y1=TMath::Max((y),Float_t(-fHeightFull/2.));
376 x1 = (y2-y1)*fK-(fWidth+fK*fHeightS)/2.;
378 x1 =-(fWidth+fK*fHeightS)/2. ;
379 Float_t x2=x1+fWidth;
383 if ((x2-x1)*fNdiv<forigsigmaX) dx=(x2-x1);
385 dx= forigsigmaX/Float_t(fNdiv);
386 dx = (x2-x1)/Float_t(Int_t(3+(x2-x1)/dx));
389 if ((y2-y1)*fNdiv<forigsigmaY) dy=(y2-y1);
391 dy= forigsigmaY/Float_t(fNdiv);
392 dy = (y2-y1)/Float_t(Int_t(3+(y2-y1)/dy));
395 for (x=x1;x<x2;x+=dx)
396 for (Float_t y=y1;y<y2;y+=dy){
397 if ( (y>(fActualY-(4.0*forigsigmaY))) &&
398 (y<(fActualY+(4.0*forigsigmaY)))){
399 Float_t xt=x-k*fK*(y-y1);
400 if ((TMath::Abs(xch-xt)<4*forigsigmaX)){
402 Float_t z0=fGRF->Eval(xch-(xt+dx/2.),fActualY-(y+dy/2.));
404 Float_t z1=fGRF->Eval(xch-(xt+dx/2.),fActualY-y);
405 Float_t z2=fGRF->Eval(xch-xt,fActualY-(y+dy/2.));
406 Float_t z3=fGRF->Eval(xch-(xt-dx/2.),fActualY-y);
407 Float_t z4=fGRF->Eval(xch-xt,fActualY-(y-dy/2.));
414 // Float_t a=(z1-z3)/2;
415 // Float_t b=(z2-z4)/2;
416 Float_t c= (z3+z1-2*z0)/2.;
417 Float_t d= (z2+z4-2*z0)/2.;
418 Float_t z= (z0+c/12.+d/12.);
420 //Float_t z= fGRF->Eval(xch-xt,fActualY-y);
421 if (z>0.) fcharge[i]+=z*dx*dy/fInteg;
433 for (x =-fNPRF*fDStep; x<fNPRF*fDStep;x+=fDStep)
434 { //x in cm fWidth in cm
435 Float_t weight = GetPRFActiv(x);
442 fSigmaX = TMath::Sqrt(fSigmaX/sum-mean*mean);
445 //calculate conversion coefitient to convert position to virtual wire
446 fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1);
450 void AliTPCPRF2D::Streamer(TBuffer &R__b)
452 // Stream an object of class AliTPCPRF2D
454 if (R__b.IsReading()) {
455 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
456 TObject::Streamer(R__b);
457 //read chewron parameters
465 //read charge parameters
476 R__b >> fPadDistance;
484 if (strncmp(fType,"User",3)==0){
488 if (strncmp(fType,"Gauss",3)==0)
489 fGRF = new TF2("fun",funGauss2D,-5.,5.,-5.,5.,4);
490 if (strncmp(fType,"Cosh",3)==0)
491 fGRF = new TF2("fun",funCosh2D,-5.,5.,-5.,5.,4);
492 if (strncmp(fType,"Gati",3)==0)
493 fGRF = new TF2("fun",funGati2D,-5.,5.,-5.,5.,5);
495 //read interpolation parameters
501 if (ffcharge!=0) delete [] ffcharge;
502 ffcharge = new Float_t[fNPRF*fNYdiv];
503 R__b.ReadFastArray(ffcharge,fNPRF*fNYdiv);
504 R__b.ReadFastArray(funParam,5);
505 if (fGRF!=0) fGRF->SetParameters(funParam);
506 //calculate conversion coefitient to convert position to virtual wire
507 fDYtoWire=Float_t(fNYdiv-1)/(fY2-fY1);
510 R__b.WriteVersion(AliTPCPRF2D::IsA());
511 TObject::Streamer(R__b);
512 //write chewron parameters
520 //write charge parameters
532 R__b << fPadDistance;
535 if (strncmp(fType,"User",3)==0) R__b <<fGRF;
536 //write interpolation parameters
542 R__b.WriteFastArray(ffcharge,fNPRF*fNYdiv);
543 R__b.WriteFastArray(funParam,5);
550 void AliTPCPRF2D::DrawX(Float_t x1 ,Float_t x2,Float_t y, Bool_t inter)
552 if (fGRF==0) return ;
555 TCanvas * c1 = new TCanvas("canPRF","Pad response function",700,900);
557 TPad * pad1 = new TPad("pad1PRF","",0.05,0.61,0.95,0.97,21);
559 TPad * pad2 = new TPad("pad2PRF","",0.05,0.22,0.95,0.60,21);
564 gStyle->SetOptFit(1);
565 gStyle->SetOptStat(0);
566 sprintf(s,"PRF response function for chevron pad");
567 TH1F * hPRFc = new TH1F("hPRFc",s,N+1,x1,x2);
572 for (Float_t i = 0;i<N+1;i++)
574 x+=(x2-x1)/Float_t(N);
575 y1 = GetPRF(x,y,inter);
580 fGRF->SetRange(x1,x1,x2,x2);
584 // hPRFo->Fit("gaus");
585 gStyle->SetOptStat(1);
589 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
590 comment->SetTextAlign(12);
591 comment->SetFillColor(42);
592 TText *title = comment->AddText("Chevron pad parameters:");
593 title->SetTextSize(0.03);
594 sprintf(s,"Full height of pad: %2.2f",fHeightFull);
596 sprintf(s,"Height of one chevron unit h: %2.2f cm",2*fHeightS);
598 sprintf(s,"Width of one chevron unit w: %2.2f cm",fWidth);
600 sprintf(s,"Overlap factor: %2.2f",fK*fHeightS/fWidth);
602 sprintf(s,"Y position: %2.2f ",y);
604 sprintf(s,"Sigma x of original distribution: %2.2f ",forigsigmaX);
606 sprintf(s,"Sigma y of original distribution: %2.2f ",forigsigmaY);
608 sprintf(s,"Type of original distribution: %s ",fType);
615 void AliTPCPRF2D::Draw(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2,
616 Bool_t inter, Int_t Nx, Int_t Ny)
619 if (fGRF==0) return ;
620 TCanvas * c1 = new TCanvas("canPRF","Pad response function",700,900);
622 TPad * pad1 = new TPad("pad1PRF","",0.05,0.61,0.95,0.97,21);
624 TPad * pad2 = new TPad("pad2PRF","",0.05,0.22,0.95,0.60,21);
629 gStyle->SetOptFit(1);
630 gStyle->SetOptStat(0);
631 sprintf(s,"PRF response function for chevron pad");
632 TH2F * hPRFc = new TH2F("hPRFc",s,Nx+1,x1,x2,Ny+1,y1,y2);
633 Float_t dx=(x2-x1)/Float_t(Nx);
634 Float_t dy=(y2-y1)/Float_t(Ny) ;
637 for ( x = x1;x<=x2;x+=dx){
638 for(y = y1;y<=y2;y+=dy)
640 z = GetPRF(x,y,inter);
645 fGRF->SetRange(x1,y1,x2,y2);
649 // hPRFo->Fit("gaus");
650 gStyle->SetOptStat(1);
652 hPRFc->Draw("lego2");
654 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
655 comment->SetTextAlign(12);
656 comment->SetFillColor(42);
657 TText *title = comment->AddText("Chevron pad parameters:");
658 title->SetTextSize(0.03);
659 sprintf(s,"Full height of pad: %2.2f",fHeightFull);
661 sprintf(s,"Height of one chevron unit h: %2.2f cm",2*fHeightS);
663 sprintf(s,"Width of one chevron unit w: %2.2f cm",fWidth);
665 sprintf(s,"Overlap factor: %2.2f",fK*fHeightS/fWidth);
667 sprintf(s,"Sigma x of original distribution: %2.2f ",forigsigmaX);
669 sprintf(s,"Sigma y of original distribution: %2.2f ",forigsigmaY);
671 sprintf(s,"Type of original distribution: %s ",fType);
676 void AliTPCPRF2D::DrawDist(Float_t x1 ,Float_t x2,Float_t y1, Float_t y2,
677 Bool_t inter, Int_t Nx, Int_t Ny, Float_t thr)
679 const Float_t minth=0.00001;
680 if (thr<minth) thr=minth;
682 if (fGRF==0) return ;
683 TCanvas * c1 = new TCanvas("padDistortion","COG distortion",700,900);
685 TPad * pad1 = new TPad("CHARGE","",0.05,0.61,0.95,0.97,21);
687 TPad * pad2 = new TPad("dist","",0.05,0.22,0.95,0.60,21);
692 gStyle->SetOptFit(1);
693 gStyle->SetOptStat(0);
694 sprintf(s,"COG distortion (threshold=%2.2f)",thr);
695 TH2F * hPRFDist = new TH2F("hDistortion",s,Nx+1,x1,x2,Ny+1,y1,y2);
696 Float_t dx=(x2-x1)/Float_t(Nx);
697 Float_t dy=(y2-y1)/Float_t(Ny) ;
700 for ( x = x1;x<(x2+dx/2.);x+=dx)
701 for(y = y1;y<=(y2+dx/2.);y+=dy)
705 for (Float_t padx=-fWidth;padx<(fWidth*1.1);padx+=fWidth)
707 z = GetPRF(x-padx,y,inter);
715 ddx = (x-(sumx/sum));
718 if (TMath::Abs(ddx)<10) hPRFDist->Fill(x,y,ddx);
721 fGRF->SetRange(x1,y1,x2,y2);
725 // hPRFo->Fit("gaus");
726 // gStyle->SetOptStat(1);
728 hPRFDist->Draw("lego2");
731 TPaveText * comment = new TPaveText(0.05,0.02,0.95,0.20,"NDC");
732 comment->SetTextAlign(12);
733 comment->SetFillColor(42);
734 // TText *title = comment->AddText("Distortion of COG method");
735 // title->SetTextSize(0.03);
736 TText * title = comment->AddText("Chevron pad parameters:");
737 title->SetTextSize(0.03);
738 sprintf(s,"Full height of pad: %2.2f",fHeightFull);
740 sprintf(s,"Height of one chevron unit h: %2.2f cm",2*fHeightS);
742 sprintf(s,"Width of one chevron unit w: %2.2f cm",fWidth);
744 sprintf(s,"Overlap factor: %2.2f",fK*fHeightS/fWidth);
746 sprintf(s,"Sigma x of original distribution: %2.2f ",forigsigmaX);
748 sprintf(s,"Sigma y of original distribution: %2.2f ",forigsigmaY);
750 sprintf(s,"Type of original distribution: %s ",fType);