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.8 2002/02/12 16:07:43 kowal2
19 coreccted bug in SetGauss
21 Revision 1.7 2001/01/26 19:57:22 hristov
22 Major upgrade of AliRoot code
24 Revision 1.6 2000/06/30 12:07:50 kowal2
25 Updated from the TPC-PreRelease branch
27 Revision 1.5.4.3 2000/06/26 07:39:42 kowal2
28 Changes to obey the coding rules
30 Revision 1.5.4.2 2000/06/25 08:38:41 kowal2
31 Splitted from AliTPCtracking
33 Revision 1.5.4.1 2000/06/14 16:48:24 kowal2
34 Parameter setting improved. Removed compiler warnings
36 Revision 1.5 2000/04/17 09:37:33 kowal2
37 removed obsolete AliTPCDigitsDisplay.C
39 Revision 1.4.8.2 2000/04/10 08:53:09 kowal2
44 Revision 1.4 1999/09/29 09:24:34 fca
45 Introduction of the Copyright and cvs Log
49 //-----------------------------------------------------------------------------
53 // Origin: Marian Ivanov, Uni. of Bratislava, ivanov@fmph.uniba.sk
55 // Declaration of class AliTPCRF1D
57 //-----------------------------------------------------------------------------
62 #include "AliTPCRF1D.h"
64 #include <Riostream.h>
70 extern TStyle * gStyle;
72 Int_t AliTPCRF1D::fgNRF=100; //default number of interpolation points
73 Float_t AliTPCRF1D::fgRFDSTEP=0.01; //default step in cm
75 static Double_t funGauss(Double_t *x, Double_t * par)
77 //Gauss function -needde by the generic function object
78 return TMath::Exp(-(x[0]*x[0])/(2*par[0]*par[0]));
81 static Double_t funCosh(Double_t *x, Double_t * par)
83 //Cosh function -needde by the generic function object
84 return 1/TMath::CosH(3.14159*x[0]/(2*par[0]));
87 static Double_t funGati(Double_t *x, Double_t * par)
89 //Gati function -needde by the generic function object
91 Float_t k3R=TMath::Sqrt(k3);
92 Float_t k2=(TMath::Pi()/2)*(1-k3R/2.);
93 Float_t k1=k2*k3R/(4*TMath::ATan(k3R));
94 Float_t l=x[0]/par[0];
95 Float_t tan2=TMath::TanH(k2*l);
97 Float_t res = k1*(1-tan2)/(1+k3*tan2);
101 ///////////////////////////////////////////////////////////////////////////
102 ///////////////////////////////////////////////////////////////////////////
107 AliTPCRF1D::AliTPCRF1D(Bool_t direct,Int_t np,Float_t step)
109 //default constructor for response function object
113 fcharge = new Float_t[fNRF];
114 if (step>0) fDSTEPM1=1./step;
115 else fDSTEPM1 = 1./fgRFDSTEP;
124 AliTPCRF1D::AliTPCRF1D(const AliTPCRF1D &prf)
127 memcpy(this, &prf, sizeof(prf));
128 fcharge = new Float_t[fNRF];
129 memcpy(fcharge,prf.fcharge, fNRF);
130 fGRF = new TF1(*(prf.fGRF));
133 AliTPCRF1D & AliTPCRF1D::operator = (const AliTPCRF1D &prf)
136 if (fcharge) delete fcharge;
137 if (fGRF) delete fGRF;
138 memcpy(this, &prf, sizeof(prf));
139 fcharge = new Float_t[fNRF];
140 memcpy(fcharge,prf.fcharge, fNRF);
141 fGRF = new TF1(*(prf.fGRF));
147 AliTPCRF1D::~AliTPCRF1D()
150 if (fcharge!=0) delete [] fcharge;
151 if (fGRF !=0 ) fGRF->Delete();
154 Float_t AliTPCRF1D::GetRF(Float_t xin)
156 //function which return response
157 //for the charge in distance xin
158 //return linear aproximation of RF
159 Float_t x = TMath::Abs((xin-fOffset)*fDSTEPM1)+fNRF/2;
164 res = fcharge[i1]*(Float_t(i1+1)-x)+fcharge[i1+1]*(x-Float_t(i1));
168 Float_t AliTPCRF1D::GetGRF(Float_t xin)
170 //function which returnoriginal charge distribution
171 //this function is just normalised for fKnorm
173 return fkNorm*fGRF->Eval(xin)/fInteg;
179 void AliTPCRF1D::SetParam( TF1 * GRF,Float_t padwidth,
180 Float_t kNorm, Float_t sigma)
182 //adjust parameters of the original charge distribution
183 //and pad size parameters
184 fpadWidth = padwidth;
187 if (sigma==0) sigma= fpadWidth/TMath::Sqrt(12.);
189 fDSTEPM1 = 10/TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
190 sprintf(fType,"User");
195 void AliTPCRF1D::SetGauss(Float_t sigma, Float_t padWidth,
199 // set parameters for Gauss generic charge distribution
201 fpadWidth = padWidth;
203 if (fGRF !=0 ) fGRF->Delete();
204 fGRF = new TF1("fun",funGauss,-5,5,1);
207 fGRF->SetParameters(funParam);
208 fDSTEPM1 = 10./TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
209 //by default I set the step as one tenth of sigma
210 sprintf(fType,"Gauss");
213 void AliTPCRF1D::SetCosh(Float_t sigma, Float_t padWidth,
217 // set parameters for Cosh generic charge distribution
219 fpadWidth = padWidth;
221 if (fGRF !=0 ) fGRF->Delete();
222 fGRF = new TF1("fun", funCosh, -5.,5.,2);
224 fGRF->SetParameters(funParam);
226 fDSTEPM1 = 10./TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
227 //by default I set the step as one tenth of sigma
228 sprintf(fType,"Cosh");
231 void AliTPCRF1D::SetGati(Float_t K3, Float_t padDistance, Float_t padWidth,
235 // set parameters for Gati generic charge distribution
237 fpadWidth = padWidth;
239 if (fGRF !=0 ) fGRF->Delete();
240 fGRF = new TF1("fun", funGati, -5.,5.,2);
241 funParam[0]=padDistance;
243 fGRF->SetParameters(funParam);
244 forigsigma=padDistance;
245 fDSTEPM1 = 10./TMath::Sqrt(padDistance*padDistance+fpadWidth*fpadWidth/12);
246 //by default I set the step as one tenth of sigma
247 sprintf(fType,"Gati");
250 void AliTPCRF1D::DrawRF(Float_t x1,Float_t x2,Int_t N)
253 //Draw prf in selected region <x1,x2> with nuber of diviision = n
256 TCanvas * c1 = new TCanvas("canRF","Pad response function",700,900);
258 TPad * pad1 = new TPad("pad1RF","",0.05,0.55,0.95,0.95,21);
260 TPad * pad2 = new TPad("pad2RF","",0.05,0.05,0.95,0.45,21);
263 sprintf(s,"RF response function for %1.2f cm pad width",
266 TH1F * hRFo = new TH1F("hRFo","Original charge distribution",N+1,x1,x2);
268 gStyle->SetOptFit(1);
269 gStyle->SetOptStat(0);
270 TH1F * hRFc = new TH1F("hRFc",s,N+1,x1,x2);
275 for (Float_t i = 0;i<N+1;i++)
277 x+=(x2-x1)/Float_t(N);
289 void AliTPCRF1D::Update()
292 //update fields with interpolated values for
295 //at the begining initialize to 0
296 for (Int_t i =0; i<fNRF;i++) fcharge[i] = 0;
297 if ( fGRF == 0 ) return;
298 fInteg = fGRF->Integral(-5*forigsigma,5*forigsigma,funParam,0.00001);
299 if ( fInteg == 0 ) fInteg = 1;
300 if (fDirect==kFALSE){
301 //integrate charge over pad for different distance of pad
302 for (Int_t i =0; i<fNRF;i++)
303 { //x in cm fpadWidth in cm
304 Float_t x = (Float_t)(i-fNRF/2)/fDSTEPM1;
305 Float_t x1=TMath::Max(x-fpadWidth/2,-5*forigsigma);
306 Float_t x2=TMath::Min(x+fpadWidth/2,5*forigsigma);
308 fkNorm*fGRF->Integral(x1,x2,funParam,0.0001)/fInteg;
312 for (Int_t i =0; i<fNRF;i++)
313 { //x in cm fpadWidth in cm
314 Float_t x = (Float_t)(i-fNRF/2)/fDSTEPM1;
315 fcharge[i] = fkNorm*fGRF->Eval(x);
321 for (Float_t x =-fNRF/fDSTEPM1; x<fNRF/fDSTEPM1;x+=1/fDSTEPM1)
322 { //x in cm fpadWidth in cm
323 Float_t weight = GetRF(x+fOffset);
330 fSigma = TMath::Sqrt(fSigma/sum-mean*mean);
335 void AliTPCRF1D::Streamer(TBuffer &R__b)
337 // Stream an object of class AliTPCRF1D.
338 if (R__b.IsReading()) {
339 AliTPCRF1D::Class()->ReadBuffer(R__b, this);
342 if (strncmp(fType,"Gauss",3)==0) {delete fGRF; fGRF = new TF1("fun",funGauss,-5.,5.,4);}
343 if (strncmp(fType,"Cosh",3)==0) {delete fGRF; fGRF = new TF1("fun",funCosh,-5.,5.,4);}
344 if (strncmp(fType,"Gati",3)==0) {delete fGRF; fGRF = new TF1("fun",funGati,-5.,5.,4);}
345 if (fGRF) fGRF->SetParameters(funParam);
348 AliTPCRF1D::Class()->WriteBuffer(R__b, this);