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.5.4.3 2000/06/26 07:39:42 kowal2
19 Changes to obey the coding rules
21 Revision 1.5.4.2 2000/06/25 08:38:41 kowal2
22 Splitted from AliTPCtracking
24 Revision 1.5.4.1 2000/06/14 16:48:24 kowal2
25 Parameter setting improved. Removed compiler warnings
27 Revision 1.5 2000/04/17 09:37:33 kowal2
28 removed obsolete AliTPCDigitsDisplay.C
30 Revision 1.4.8.2 2000/04/10 08:53:09 kowal2
35 Revision 1.4 1999/09/29 09:24:34 fca
36 Introduction of the Copyright and cvs Log
40 //-----------------------------------------------------------------------------
44 // Origin: Marian Ivanov, Uni. of Bratislava, ivanov@fmph.uniba.sk
46 // Declaration of class AliTPCRF1D
48 //-----------------------------------------------------------------------------
53 #include "AliTPCRF1D.h"
61 extern TStyle * gStyle;
63 Int_t AliTPCRF1D::fgNRF=100; //default number of interpolation points
64 Float_t AliTPCRF1D::fgRFDSTEP=0.01; //default step in cm
66 static Double_t funGauss(Double_t *x, Double_t * par)
68 //Gauss function -needde by the generic function object
69 return TMath::Exp(-(x[0]*x[0])/(2*par[0]*par[0]));
72 static Double_t funCosh(Double_t *x, Double_t * par)
74 //Cosh function -needde by the generic function object
75 return 1/TMath::CosH(3.14159*x[0]/(2*par[0]));
78 static Double_t funGati(Double_t *x, Double_t * par)
80 //Gati function -needde by the generic function object
82 Float_t k3R=TMath::Sqrt(k3);
83 Float_t k2=(TMath::Pi()/2)*(1-k3R/2.);
84 Float_t k1=k2*k3R/(4*TMath::ATan(k3R));
85 Float_t l=x[0]/par[0];
86 Float_t tan2=TMath::TanH(k2*l);
88 Float_t res = k1*(1-tan2)/(1+k3*tan2);
92 ///////////////////////////////////////////////////////////////////////////
93 ///////////////////////////////////////////////////////////////////////////
98 AliTPCRF1D::AliTPCRF1D(Bool_t direct,Int_t np,Float_t step)
100 //default constructor for response function object
104 fcharge = new Float_t[fNRF];
105 if (step>0) fDSTEPM1=1./step;
106 else fDSTEPM1 = 1./fgRFDSTEP;
115 AliTPCRF1D::AliTPCRF1D(const AliTPCRF1D &prf)
118 memcpy(this, &prf, sizeof(prf));
119 fcharge = new Float_t[fNRF];
120 memcpy(fcharge,prf.fcharge, fNRF);
121 fGRF = new TF1(*(prf.fGRF));
124 AliTPCRF1D & AliTPCRF1D::operator = (const AliTPCRF1D &prf)
127 if (fcharge) delete fcharge;
128 if (fGRF) delete fGRF;
129 memcpy(this, &prf, sizeof(prf));
130 fcharge = new Float_t[fNRF];
131 memcpy(fcharge,prf.fcharge, fNRF);
132 fGRF = new TF1(*(prf.fGRF));
138 AliTPCRF1D::~AliTPCRF1D()
141 if (fcharge!=0) delete [] fcharge;
142 if (fGRF !=0 ) fGRF->Delete();
145 Float_t AliTPCRF1D::GetRF(Float_t xin)
147 //function which return response
148 //for the charge in distance xin
149 //return linear aproximation of RF
150 Float_t x = TMath::Abs((xin-fOffset)*fDSTEPM1)+fNRF/2;
155 res = fcharge[i1]*(Float_t(i1+1)-x)+fcharge[i1+1]*(x-Float_t(i1));
159 Float_t AliTPCRF1D::GetGRF(Float_t xin)
161 //function which returnoriginal charge distribution
162 //this function is just normalised for fKnorm
164 return fkNorm*fGRF->Eval(xin)/fInteg;
170 void AliTPCRF1D::SetParam( TF1 * GRF,Float_t padwidth,
171 Float_t kNorm, Float_t sigma)
173 //adjust parameters of the original charge distribution
174 //and pad size parameters
175 fpadWidth = padwidth;
178 if (sigma==0) sigma= fpadWidth/TMath::Sqrt(12.);
180 fDSTEPM1 = 10/TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
181 sprintf(fType,"User");
186 void AliTPCRF1D::SetGauss(Float_t sigma, Float_t padWidth,
190 // set parameters for Gauss generic charge distribution
192 fpadWidth = padWidth;
194 if (fGRF !=0 ) fGRF->Delete();
195 fGRF = new TF1("fun",funGauss,-5,5,2);
198 fGRF->SetParameters(funParam);
199 fDSTEPM1 = 10./TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
200 //by default I set the step as one tenth of sigma
201 sprintf(fType,"Gauss");
204 void AliTPCRF1D::SetCosh(Float_t sigma, Float_t padWidth,
208 // set parameters for Cosh generic charge distribution
210 fpadWidth = padWidth;
212 if (fGRF !=0 ) fGRF->Delete();
213 fGRF = new TF1("fun", funCosh, -5.,5.,2);
215 fGRF->SetParameters(funParam);
217 fDSTEPM1 = 10./TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
218 //by default I set the step as one tenth of sigma
219 sprintf(fType,"Cosh");
222 void AliTPCRF1D::SetGati(Float_t K3, Float_t padDistance, Float_t padWidth,
226 // set parameters for Gati generic charge distribution
228 fpadWidth = padWidth;
230 if (fGRF !=0 ) fGRF->Delete();
231 fGRF = new TF1("fun", funGati, -5.,5.,2);
232 funParam[0]=padDistance;
234 fGRF->SetParameters(funParam);
235 forigsigma=padDistance;
236 fDSTEPM1 = 10./TMath::Sqrt(padDistance*padDistance+fpadWidth*fpadWidth/12);
237 //by default I set the step as one tenth of sigma
238 sprintf(fType,"Gati");
241 void AliTPCRF1D::DrawRF(Float_t x1,Float_t x2,Int_t N)
244 //Draw prf in selected region <x1,x2> with nuber of diviision = n
247 TCanvas * c1 = new TCanvas("canRF","Pad response function",700,900);
249 TPad * pad1 = new TPad("pad1RF","",0.05,0.55,0.95,0.95,21);
251 TPad * pad2 = new TPad("pad2RF","",0.05,0.05,0.95,0.45,21);
254 sprintf(s,"RF response function for %1.2f cm pad width",
257 TH1F * hRFo = new TH1F("hRFo","Original charge distribution",N+1,x1,x2);
259 gStyle->SetOptFit(1);
260 gStyle->SetOptStat(0);
261 TH1F * hRFc = new TH1F("hRFc",s,N+1,x1,x2);
266 for (Float_t i = 0;i<N+1;i++)
268 x+=(x2-x1)/Float_t(N);
280 void AliTPCRF1D::Update()
283 //update fields with interpolated values for
286 //at the begining initialize to 0
287 for (Int_t i =0; i<fNRF;i++) fcharge[i] = 0;
288 if ( fGRF == 0 ) return;
289 fInteg = fGRF->Integral(-5*forigsigma,5*forigsigma,funParam,0.00001);
290 if ( fInteg == 0 ) fInteg = 1;
291 if (fDirect==kFALSE){
292 //integrate charge over pad for different distance of pad
293 for (Int_t i =0; i<fNRF;i++)
294 { //x in cm fpadWidth in cm
295 Float_t x = (Float_t)(i-fNRF/2)/fDSTEPM1;
296 Float_t x1=TMath::Max(x-fpadWidth/2,-5*forigsigma);
297 Float_t x2=TMath::Min(x+fpadWidth/2,5*forigsigma);
299 fkNorm*fGRF->Integral(x1,x2,funParam,0.0001)/fInteg;
303 for (Int_t i =0; i<fNRF;i++)
304 { //x in cm fpadWidth in cm
305 Float_t x = (Float_t)(i-fNRF/2)/fDSTEPM1;
306 fcharge[i] = fkNorm*fGRF->Eval(x);
312 for (Float_t x =-fNRF/fDSTEPM1; x<fNRF/fDSTEPM1;x+=1/fDSTEPM1)
313 { //x in cm fpadWidth in cm
314 Float_t weight = GetRF(x+fOffset);
321 fSigma = TMath::Sqrt(fSigma/sum-mean*mean);
326 void AliTPCRF1D::Streamer(TBuffer &R__b)
328 // Stream an object of class AliTPCRF1D.
330 if (R__b.IsReading()) {
331 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
332 TObject::Streamer(R__b);
333 //read pad parameters
335 //read charge parameters
343 R__b >> dummy; //dummuy instead fK3X in previous
344 R__b >> fPadDistance;
352 if (strncmp(fType,"User",3)==0){
357 if (strncmp(fType,"Gauss",3)==0)
358 fGRF = new TF1("fun",funGauss,-5.,5.,4);
359 if (strncmp(fType,"Cosh",3)==0)
360 fGRF = new TF1("fun",funCosh,-5.,5.,4);
361 if (strncmp(fType,"Gati",3)==0)
362 fGRF = new TF1("fun",funGati,-5.,5.,4);
365 R__b.ReadFastArray(fcharge,fNRF);
366 R__b.ReadFastArray(funParam,5);
367 if (fGRF!=0) fGRF->SetParameters(funParam);
370 R__b.WriteVersion(AliTPCRF1D::IsA());
371 TObject::Streamer(R__b);
374 //write charge parameters
382 R__b << dummy; //dummuy instead fK3X in previouis
383 R__b << fPadDistance;
386 //write interpolation parameters
387 if (strncmp(fType,"User",3)==0) R__b <<fGRF;
390 R__b.WriteFastArray(fcharge,fNRF);
391 R__b.WriteFastArray(funParam,5);