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 **************************************************************************/
20 //-----------------------------------------------------------------------------
24 // Origin: Marian Ivanov, Uni. of Bratislava, ivanov@fmph.uniba.sk
26 // Declaration of class AliTPCRF1D
28 //-----------------------------------------------------------------------------
30 #include "AliTPCRF1D.h"
38 extern TStyle * gStyle;
40 static Double_t funGauss(Double_t *x, Double_t * par)
42 return TMath::Exp(-(x[0]*x[0])/(2*par[0]*par[0]));
45 static Double_t funCosh(Double_t *x, Double_t * par)
47 return 1/TMath::CosH(3.14159*x[0]/(2*par[0]));
50 static Double_t funGati(Double_t *x, Double_t * par)
52 //par[1] = is equal to k3
53 //par[0] is equal to pad wire distance
55 Float_t K3R=TMath::Sqrt(K3);
56 Float_t K2=(TMath::Pi()/2)*(1-K3R/2.);
57 Float_t K1=K2*K3R/(4*TMath::ATan(K3R));
58 Float_t l=x[0]/par[0];
59 Float_t tan2=TMath::TanH(K2*l);
61 Float_t res = K1*(1-tan2)/(1+K3*tan2);
68 ///////////////////////////////////////////////////////////////////////////
69 ///////////////////////////////////////////////////////////////////////////
70 ///////////////////////////////////////////////////////////////////////////
71 ///////////////////////////////////////////////////////////////////////////
77 AliTPCRF1D::AliTPCRF1D(Bool_t direct,Int_t np,Float_t step)
81 fcharge = new Float_t[fNRF];
93 AliTPCRF1D::~AliTPCRF1D()
95 if (fcharge!=0) delete fcharge;
96 if (fGRF !=0 ) fGRF->Delete();
99 Float_t AliTPCRF1D::GetRF(Float_t xin)
102 //return linear aproximation of RF
103 Float_t x = TMath::Abs((xin-fOffset)*fDSTEPM1)+fNRF/2;
108 res = fcharge[i1]*(Float_t(i1+1)-x)+fcharge[i1+1]*(x-Float_t(i1));
112 Float_t AliTPCRF1D::GetGRF(Float_t xin)
115 return fkNorm*fGRF->Eval(xin)/fInteg;
121 void AliTPCRF1D::SetParam( TF1 * GRF,Float_t padwidth,
122 Float_t kNorm, Float_t sigma)
124 fpadWidth = padwidth;
127 if (sigma==0) sigma= fpadWidth/TMath::Sqrt(12.);
129 fDSTEPM1 = 10/TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
130 sprintf(fType,"User");
135 void AliTPCRF1D::SetGauss(Float_t sigma, Float_t padWidth,
139 fpadWidth = padWidth;
141 if (fGRF !=0 ) fGRF->Delete();
142 fGRF = new TF1("fun",funGauss,-5,5,2);
145 fGRF->SetParameters(funParam);
146 fDSTEPM1 = 10./TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
147 //by default I set the step as one tenth of sigma
149 sprintf(fType,"Gauss");
152 void AliTPCRF1D::SetCosh(Float_t sigma, Float_t padWidth,
156 fpadWidth = padWidth;
158 if (fGRF !=0 ) fGRF->Delete();
159 fGRF = new TF1("fun", funCosh, -5.,5.,2);
161 fGRF->SetParameters(funParam);
163 fDSTEPM1 = 10./TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
164 //by default I set the step as one tenth of sigma
166 sprintf(fType,"Cosh");
169 void AliTPCRF1D::SetGati(Float_t K3, Float_t padDistance, Float_t padWidth,
173 fpadWidth = padWidth;
175 if (fGRF !=0 ) fGRF->Delete();
176 fGRF = new TF1("fun", funGati, -5.,5.,2);
177 funParam[0]=padDistance;
179 fGRF->SetParameters(funParam);
180 forigsigma=padDistance;
181 fDSTEPM1 = 10./TMath::Sqrt(padDistance*padDistance+fpadWidth*fpadWidth/12);
182 //by default I set the step as one tenth of sigma
184 sprintf(fType,"Gati");
187 void AliTPCRF1D::Draw(Float_t x1,Float_t x2,Int_t N)
190 TCanvas * c1 = new TCanvas("canRF","Pad response function",700,900);
192 TPad * pad1 = new TPad("pad1RF","",0.05,0.55,0.95,0.95,21);
194 TPad * pad2 = new TPad("pad2RF","",0.05,0.05,0.95,0.45,21);
197 sprintf(s,"RF response function for %1.2f cm pad width",
200 TH1F * hRFo = new TH1F("hRFo","Original charge distribution",N+1,x1,x2);
202 gStyle->SetOptFit(1);
203 gStyle->SetOptStat(0);
204 TH1F * hRFc = new TH1F("hRFc",s,N+1,x1,x2);
209 for (Float_t i = 0;i<N+1;i++)
211 x+=(x2-x1)/Float_t(N);
223 void AliTPCRF1D::Update()
226 for (Int_t i =0; i<fNRF;i++) fcharge[i] = 0;
227 if ( fGRF == 0 ) return;
228 fInteg = fGRF->Integral(-5*forigsigma,5*forigsigma,funParam,0.00001);
229 if ( fInteg == 0 ) fInteg = 1;
230 if (fDirect==kFALSE){
231 //integrate charge over pad for different distance of pad
232 for (Int_t i =0; i<fNRF;i++)
233 { //x in cm fpadWidth in cm
234 Float_t x = (Float_t)(i-fNRF/2)/fDSTEPM1;
235 Float_t x1=TMath::Max(x-fpadWidth/2,-5*forigsigma);
236 Float_t x2=TMath::Min(x+fpadWidth/2,5*forigsigma);
238 fkNorm*fGRF->Integral(x1,x2,funParam,0.0001)/fInteg;
242 for (Int_t i =0; i<fNRF;i++)
243 { //x in cm fpadWidth in cm
244 Float_t x = (Float_t)(i-fNRF/2)/fDSTEPM1;
245 fcharge[i] = fkNorm*fGRF->Eval(x);
251 for (Float_t x =-fNRF/fDSTEPM1; x<fNRF/fDSTEPM1;x+=1/fDSTEPM1)
252 { //x in cm fpadWidth in cm
253 Float_t weight = GetRF(x+fOffset);
260 fSigma = TMath::Sqrt(fSigma/sum-mean*mean);
265 void AliTPCRF1D::Streamer(TBuffer &R__b)
267 // Stream an object of class AliTPC.
269 if (R__b.IsReading()) {
270 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
271 TObject::Streamer(R__b);
272 //read pad parameters
274 //read charge parameters
283 R__b >> fPadDistance;
291 if (strncmp(fType,"User",3)==0){
296 if (strncmp(fType,"Gauss",3)==0)
297 fGRF = new TF1("fun",funGauss,-5.,5.,4);
298 if (strncmp(fType,"Cosh",3)==0)
299 fGRF = new TF1("fun",funCosh,-5.,5.,4);
300 if (strncmp(fType,"Gati",3)==0)
301 fGRF = new TF1("fun",funGati,-5.,5.,4);
304 R__b.ReadFastArray(fcharge,fNRF);
305 R__b.ReadFastArray(funParam,5);
306 if (fGRF!=0) fGRF->SetParameters(funParam);
309 R__b.WriteVersion(AliTPCRF1D::IsA());
310 TObject::Streamer(R__b);
313 //write charge parameters
322 R__b << fPadDistance;
325 //write interpolation parameters
326 if (strncmp(fType,"User",3)==0) R__b <<fGRF;
329 R__b.WriteFastArray(fcharge,fNRF);
330 R__b.WriteFastArray(funParam,5);