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 **************************************************************************/
19 //-----------------------------------------------------------------------------
23 // Origin: Marian Ivanov, Uni. of Bratislava, ivanov@fmph.uniba.sk
25 // Declaration of class AliTPCRF1D
27 //-----------------------------------------------------------------------------
31 #include <Riostream.h>
41 #include "AliTPCRF1D.h"
43 extern TStyle * gStyle;
45 Int_t AliTPCRF1D::fgNRF=100; //default number of interpolation points
46 Float_t AliTPCRF1D::fgRFDSTEP=0.01; //default step in cm
48 static Double_t funGauss(Double_t *x, Double_t * par)
50 //Gauss function -needde by the generic function object
51 return TMath::Exp(-(x[0]*x[0])/(2*par[0]*par[0]));
54 static Double_t funCosh(Double_t *x, Double_t * par)
56 //Cosh function -needde by the generic function object
57 return 1/TMath::CosH(3.14159*x[0]/(2*par[0]));
60 static Double_t funGati(Double_t *x, Double_t * par)
62 //Gati function -needde by the generic function object
64 Float_t k3R=TMath::Sqrt(k3);
65 Float_t k2=(TMath::Pi()/2)*(1-k3R/2.);
66 Float_t k1=k2*k3R/(4*TMath::ATan(k3R));
67 Float_t l=x[0]/par[0];
68 Float_t tan2=TMath::TanH(k2*l);
70 Float_t res = k1*(1-tan2)/(1+k3*tan2);
74 ///////////////////////////////////////////////////////////////////////////
75 ///////////////////////////////////////////////////////////////////////////
80 AliTPCRF1D::AliTPCRF1D(Bool_t direct,Int_t np,Float_t step)
95 //default constructor for response function object
99 fcharge = new Float_t[fNRF];
100 if (step>0) fDSTEPM1=1./step;
101 else fDSTEPM1 = 1./fgRFDSTEP;
104 AliTPCRF1D::AliTPCRF1D(const AliTPCRF1D &prf)
121 memcpy(this, &prf, sizeof(prf));
122 fcharge = new Float_t[fNRF];
123 memcpy(fcharge,prf.fcharge, fNRF);
124 fGRF = new TF1(*(prf.fGRF));
125 //PH Change the name (add 0 to the end)
126 TString s(fGRF->GetName());
128 fGRF->SetName(s.Data());
131 AliTPCRF1D & AliTPCRF1D::operator = (const AliTPCRF1D &prf)
134 if (fcharge) delete fcharge;
135 if (fGRF) delete fGRF;
136 memcpy(this, &prf, sizeof(prf));
137 fcharge = new Float_t[fNRF];
138 memcpy(fcharge,prf.fcharge, fNRF);
139 fGRF = new TF1(*(prf.fGRF));
140 //PH Change the name (add 0 to the end)
141 TString s(fGRF->GetName());
143 fGRF->SetName(s.Data());
149 AliTPCRF1D::~AliTPCRF1D()
152 if (fcharge!=0) delete [] fcharge;
153 if (fGRF !=0 ) fGRF->Delete();
156 Float_t AliTPCRF1D::GetRF(Float_t xin)
158 //function which return response
159 //for the charge in distance xin
160 //return linear aproximation of RF
161 Float_t x = (xin-fOffset)*fDSTEPM1+fNRF/2;
165 if (i1+1<fNRF &&i1>0)
166 res = fcharge[i1]*(Float_t(i1+1)-x)+fcharge[i1+1]*(x-Float_t(i1));
170 Float_t AliTPCRF1D::GetGRF(Float_t xin)
172 //function which returnoriginal charge distribution
173 //this function is just normalised for fKnorm
175 return fkNorm*fGRF->Eval(xin)/fInteg;
181 void AliTPCRF1D::SetParam( TF1 * GRF,Float_t padwidth,
182 Float_t kNorm, Float_t sigma)
184 //adjust parameters of the original charge distribution
185 //and pad size parameters
186 fpadWidth = padwidth;
189 if (sigma==0) sigma= fpadWidth/TMath::Sqrt(12.);
191 fDSTEPM1 = 10/TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
192 sprintf(fType,"User");
197 void AliTPCRF1D::SetGauss(Float_t sigma, Float_t padWidth,
201 // set parameters for Gauss generic charge distribution
203 fpadWidth = padWidth;
205 if (fGRF !=0 ) fGRF->Delete();
206 fGRF = new TF1("funGauss",funGauss,-5,5,1);
209 fGRF->SetParameters(funParam);
210 fDSTEPM1 = 10./TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
211 //by default I set the step as one tenth of sigma
212 sprintf(fType,"Gauss");
215 void AliTPCRF1D::SetCosh(Float_t sigma, Float_t padWidth,
219 // set parameters for Cosh generic charge distribution
221 fpadWidth = padWidth;
223 if (fGRF !=0 ) fGRF->Delete();
224 fGRF = new TF1("funCosh", funCosh, -5.,5.,2);
226 fGRF->SetParameters(funParam);
228 fDSTEPM1 = 10./TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
229 //by default I set the step as one tenth of sigma
230 sprintf(fType,"Cosh");
233 void AliTPCRF1D::SetGati(Float_t K3, Float_t padDistance, Float_t padWidth,
237 // set parameters for Gati generic charge distribution
239 fpadWidth = padWidth;
241 if (fGRF !=0 ) fGRF->Delete();
242 fGRF = new TF1("funGati", funGati, -5.,5.,2);
243 funParam[0]=padDistance;
245 fGRF->SetParameters(funParam);
246 forigsigma=padDistance;
247 fDSTEPM1 = 10./TMath::Sqrt(padDistance*padDistance+fpadWidth*fpadWidth/12);
248 //by default I set the step as one tenth of sigma
249 sprintf(fType,"Gati");
254 void AliTPCRF1D::DrawRF(Float_t x1,Float_t x2,Int_t N)
257 //Draw prf in selected region <x1,x2> with nuber of diviision = n
260 TCanvas * c1 = new TCanvas("canRF","Pad response function",700,900);
262 TPad * pad1 = new TPad("pad1RF","",0.05,0.55,0.95,0.95,21);
264 TPad * pad2 = new TPad("pad2RF","",0.05,0.05,0.95,0.45,21);
267 sprintf(s,"RF response function for %1.2f cm pad width",
270 TH1F * hRFo = new TH1F("hRFo","Original charge distribution",N+1,x1,x2);
272 gStyle->SetOptFit(1);
273 gStyle->SetOptStat(0);
274 TH1F * hRFc = new TH1F("hRFc",s,N+1,x1,x2);
279 for (Float_t i = 0;i<N+1;i++)
281 x+=(x2-x1)/Float_t(N);
293 void AliTPCRF1D::Update()
296 //update fields with interpolated values for
299 //at the begining initialize to 0
300 for (Int_t i =0; i<fNRF;i++) fcharge[i] = 0;
301 if ( fGRF == 0 ) return;
302 fInteg = fGRF->Integral(-5*forigsigma,5*forigsigma,funParam,0.00001);
303 if ( fInteg == 0 ) fInteg = 1;
304 if (fDirect==kFALSE){
305 //integrate charge over pad for different distance of pad
306 for (Int_t i =0; i<fNRF;i++)
307 { //x in cm fpadWidth in cm
308 Float_t x = (Float_t)(i-fNRF/2)/fDSTEPM1;
309 Float_t x1=TMath::Max(x-fpadWidth/2,-5*forigsigma);
310 Float_t x2=TMath::Min(x+fpadWidth/2,5*forigsigma);
312 fkNorm*fGRF->Integral(x1,x2,funParam,0.0001)/fInteg;
316 for (Int_t i =0; i<fNRF;i++)
317 { //x in cm fpadWidth in cm
318 Float_t x = (Float_t)(i-fNRF/2)/fDSTEPM1;
319 fcharge[i] = fkNorm*fGRF->Eval(x);
325 for (Float_t x =-fNRF/fDSTEPM1; x<fNRF/fDSTEPM1;x+=1/fDSTEPM1)
326 { //x in cm fpadWidth in cm
327 Float_t weight = GetRF(x+fOffset);
334 fSigma = TMath::Sqrt(fSigma/sum-mean*mean);
339 void AliTPCRF1D::Streamer(TBuffer &R__b)
341 // Stream an object of class AliTPCRF1D.
342 if (R__b.IsReading()) {
343 AliTPCRF1D::Class()->ReadBuffer(R__b, this);
346 if (strncmp(fType,"Gauss",3)==0) {delete fGRF; fGRF = new TF1("funGauss",funGauss,-5.,5.,4);}
347 if (strncmp(fType,"Cosh",3)==0) {delete fGRF; fGRF = new TF1("funCosh",funCosh,-5.,5.,4);}
348 if (strncmp(fType,"Gati",3)==0) {delete fGRF; fGRF = new TF1("funGati",funGati,-5.,5.,4);}
349 if (fGRF) fGRF->SetParameters(funParam);
352 AliTPCRF1D::Class()->WriteBuffer(R__b, this);
357 Double_t AliTPCRF1D::Gamma4(Double_t x, Double_t p0, Double_t p1){
359 // Gamma 4 Time response function of ALTRO
362 Double_t g1 = TMath::Exp(-4.*x/p1);
363 Double_t g2 = TMath::Power(x/p1,4);