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 //-----------------------------------------------------------------------------
32 #include <Riostream.h>
42 #include "AliTPCRF1D.h"
44 extern TStyle * gStyle;
46 Int_t AliTPCRF1D::fgNRF=100; //default number of interpolation points
47 Float_t AliTPCRF1D::fgRFDSTEP=0.01; //default step in cm
49 static Double_t funGauss(Double_t *x, Double_t * par)
51 //Gauss function -needde by the generic function object
52 return TMath::Exp(-(x[0]*x[0])/(2*par[0]*par[0]));
55 static Double_t funCosh(Double_t *x, Double_t * par)
57 //Cosh function -needde by the generic function object
58 return 1/TMath::CosH(3.14159*x[0]/(2*par[0]));
61 static Double_t funGati(Double_t *x, Double_t * par)
63 //Gati function -needde by the generic function object
65 Float_t k3R=TMath::Sqrt(k3);
66 Float_t k2=(TMath::Pi()/2)*(1-k3R/2.);
67 Float_t k1=k2*k3R/(4*TMath::ATan(k3R));
68 Float_t l=x[0]/par[0];
69 Float_t tan2=TMath::TanH(k2*l);
71 Float_t res = k1*(1-tan2)/(1+k3*tan2);
75 ///////////////////////////////////////////////////////////////////////////
76 ///////////////////////////////////////////////////////////////////////////
81 AliTPCRF1D::AliTPCRF1D(Bool_t direct,Int_t np,Float_t step)
96 //default constructor for response function object
100 fcharge = new Float_t[fNRF];
101 if (step>0) fDSTEPM1=1./step;
102 else fDSTEPM1 = 1./fgRFDSTEP;
103 for(Int_t i=0;i<5;i++) {
110 AliTPCRF1D::AliTPCRF1D(const AliTPCRF1D &prf)
113 fDSTEPM1(prf.fDSTEPM1),
115 forigsigma(prf.forigsigma),
116 fpadWidth(prf.fpadWidth),
119 fGRF(new TF1(*(prf.fGRF))),
121 fOffset(prf.fOffset),
122 fDirect(prf.fDirect),
123 fPadDistance(prf.fPadDistance)
127 for(Int_t i=0;i<5;i++) {
131 fcharge = new Float_t[fNRF];
132 memcpy(fcharge,prf.fcharge, fNRF*sizeof(Float_t));
134 //PH Change the name (add 0 to the end)
135 TString s(fGRF->GetName());
137 fGRF->SetName(s.Data());
140 AliTPCRF1D & AliTPCRF1D::operator = (const AliTPCRF1D &prf)
143 TObject::operator=(prf);
145 fDSTEPM1=prf.fDSTEPM1;
147 fcharge = new Float_t[fNRF];
148 memcpy(fcharge,prf.fcharge, fNRF*sizeof(Float_t));
149 forigsigma=prf.forigsigma;
150 fpadWidth=prf.fpadWidth;
154 fGRF=new TF1(*(prf.fGRF));
155 //PH Change the name (add 0 to the end)
156 TString s(fGRF->GetName());
158 fGRF->SetName(s.Data());
162 fPadDistance=prf.fPadDistance;
169 AliTPCRF1D::~AliTPCRF1D()
176 Float_t AliTPCRF1D::GetRF(Float_t xin)
178 //function which return response
179 //for the charge in distance xin
180 //return linear aproximation of RF
181 Float_t x = (xin-fOffset)*fDSTEPM1+fNRF/2;
185 if (i1+1<fNRF &&i1>0)
186 res = fcharge[i1]*(Float_t(i1+1)-x)+fcharge[i1+1]*(x-Float_t(i1));
190 Float_t AliTPCRF1D::GetGRF(Float_t xin)
192 //function which returnoriginal charge distribution
193 //this function is just normalised for fKnorm
195 return fkNorm*fGRF->Eval(xin)/fInteg;
201 void AliTPCRF1D::SetParam( TF1 * GRF,Float_t padwidth,
202 Float_t kNorm, Float_t sigma)
204 //adjust parameters of the original charge distribution
205 //and pad size parameters
206 fpadWidth = padwidth;
209 if (sigma==0) sigma= fpadWidth/TMath::Sqrt(12.);
211 fDSTEPM1 = 10/TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
212 //sprintf(fType,"User");
213 snprintf(fType,5,"User");
218 void AliTPCRF1D::SetGauss(Float_t sigma, Float_t padWidth,
222 // set parameters for Gauss generic charge distribution
224 fpadWidth = padWidth;
226 if (fGRF !=0 ) fGRF->Delete();
227 fGRF = new TF1("funGauss",funGauss,-5,5,1);
230 fGRF->SetParameters(funParam);
231 fDSTEPM1 = 10./TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
232 //by default I set the step as one tenth of sigma
233 //sprintf(fType,"Gauss");
234 snprintf(fType,5,"Gauss");
237 void AliTPCRF1D::SetCosh(Float_t sigma, Float_t padWidth,
241 // set parameters for Cosh generic charge distribution
243 fpadWidth = padWidth;
245 if (fGRF !=0 ) fGRF->Delete();
246 fGRF = new TF1("funCosh", funCosh, -5.,5.,2);
248 fGRF->SetParameters(funParam);
250 fDSTEPM1 = 10./TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
251 //by default I set the step as one tenth of sigma
252 //sprintf(fType,"Cosh");
253 snprintf(fType,5,"Cosh");
256 void AliTPCRF1D::SetGati(Float_t K3, Float_t padDistance, Float_t padWidth,
260 // set parameters for Gati generic charge distribution
262 fpadWidth = padWidth;
264 if (fGRF !=0 ) fGRF->Delete();
265 fGRF = new TF1("funGati", funGati, -5.,5.,2);
266 funParam[0]=padDistance;
268 fGRF->SetParameters(funParam);
269 forigsigma=padDistance;
270 fDSTEPM1 = 10./TMath::Sqrt(padDistance*padDistance+fpadWidth*fpadWidth/12);
271 //by default I set the step as one tenth of sigma
272 //sprintf(fType,"Gati");
273 snprintf(fType,5,"Gati");
278 void AliTPCRF1D::DrawRF(Float_t x1,Float_t x2,Int_t N)
281 //Draw prf in selected region <x1,x2> with nuber of diviision = n
284 TCanvas * c1 = new TCanvas("canRF","Pad response function",700,900);
286 TPad * pad1 = new TPad("pad1RF","",0.05,0.55,0.95,0.95,21);
288 TPad * pad2 = new TPad("pad2RF","",0.05,0.05,0.95,0.45,21);
291 //sprintf(s,"RF response function for %1.2f cm pad width",
293 snprintf(s,60,"RF response function for %1.2f cm pad width",fpadWidth);
295 TH1F * hRFo = new TH1F("hRFo","Original charge distribution",N+1,x1,x2);
297 gStyle->SetOptFit(1);
298 gStyle->SetOptStat(0);
299 TH1F * hRFc = new TH1F("hRFc",s,N+1,x1,x2);
304 for (Float_t i = 0;i<N+1;i++)
306 x+=(x2-x1)/Float_t(N);
318 void AliTPCRF1D::Update()
321 //update fields with interpolated values for
324 //at the begining initialize to 0
325 for (Int_t i =0; i<fNRF;i++) fcharge[i] = 0;
326 if ( fGRF == 0 ) return;
327 // This form is no longer available
328 #if ROOT_VERSION_CODE < ROOT_VERSION(5,99,0)
329 fInteg = fGRF->Integral(-5*forigsigma,5*forigsigma,funParam,0.00001);
331 TArrayD savParam(fGRF->GetNpar(), fGRF->GetParameters());
332 fGRF->SetParameters(funParam);
333 fInteg = fGRF->Integral(-5*forigsigma,5*forigsigma,0.00001);
335 if ( fInteg == 0 ) fInteg = 1;
336 if (fDirect==kFALSE){
337 //integrate charge over pad for different distance of pad
338 for (Int_t i =0; i<fNRF;i++)
339 { //x in cm fpadWidth in cm
340 Float_t x = (Float_t)(i-fNRF/2)/fDSTEPM1;
341 Float_t x1=TMath::Max(x-fpadWidth/2,-5*forigsigma);
342 Float_t x2=TMath::Min(x+fpadWidth/2,5*forigsigma);
343 #if ROOT_VERSION_CODE < ROOT_VERSION(5,99,0)
344 fcharge[i] = fkNorm*fGRF->Integral(x1,x2,funParam,0.0001)/fInteg;
346 fcharge[i] = fkNorm*fGRF->Integral(x1,x2,0.0001)/fInteg;
351 for (Int_t i =0; i<fNRF;i++)
352 { //x in cm fpadWidth in cm
353 Float_t x = (Float_t)(i-fNRF/2)/fDSTEPM1;
354 fcharge[i] = fkNorm*fGRF->Eval(x);
360 for (Float_t x =-fNRF/fDSTEPM1; x<fNRF/fDSTEPM1;x+=1/fDSTEPM1)
361 { //x in cm fpadWidth in cm
362 Float_t weight = GetRF(x+fOffset);
369 fSigma = TMath::Sqrt(fSigma/sum-mean*mean);
372 #if ROOT_VERSION_CODE >= ROOT_VERSION(5,99,0)
373 fGRF->SetParameters(savParam.GetArray());
377 void AliTPCRF1D::Streamer(TBuffer &R__b)
379 // Stream an object of class AliTPCRF1D.
380 if (R__b.IsReading()) {
381 AliTPCRF1D::Class()->ReadBuffer(R__b, this);
384 if (strncmp(fType,"Gauss",3)==0) {delete fGRF; fGRF = new TF1("funGauss",funGauss,-5.,5.,4);}
385 if (strncmp(fType,"Cosh",3)==0) {delete fGRF; fGRF = new TF1("funCosh",funCosh,-5.,5.,4);}
386 if (strncmp(fType,"Gati",3)==0) {delete fGRF; fGRF = new TF1("funGati",funGati,-5.,5.,4);}
387 if (fGRF) fGRF->SetParameters(funParam);
390 AliTPCRF1D::Class()->WriteBuffer(R__b, this);
395 Double_t AliTPCRF1D::Gamma4(Double_t x, Double_t p0, Double_t p1){
397 // Gamma 4 Time response function of ALTRO
400 Double_t g1 = TMath::Exp(-4.*x/p1);
401 Double_t g2 = TMath::Power(x/p1,4);