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 /// \brief Declaration of class AliTPCRF1D
21 /// \author Marian Ivanov, Uni. of Bratislava, ivanov@fmph.uniba.sk
24 #include <Riostream.h>
34 #include "AliTPCRF1D.h"
36 extern TStyle * gStyle;
38 Int_t AliTPCRF1D::fgNRF=100; ///< default number of interpolation points
39 Float_t AliTPCRF1D::fgRFDSTEP=0.01; ///< default step in cm
41 static Double_t funGauss(Double_t *x, Double_t * par)
43 /// Gauss function -needde by the generic function object
45 return TMath::Exp(-(x[0]*x[0])/(2*par[0]*par[0]));
48 static Double_t funCosh(Double_t *x, Double_t * par)
50 /// Cosh function -needde by the generic function object
52 return 1/TMath::CosH(3.14159*x[0]/(2*par[0]));
55 static Double_t funGati(Double_t *x, Double_t * par)
57 /// Gati function -needde by the generic function object
60 Float_t k3R=TMath::Sqrt(k3);
61 Float_t k2=(TMath::Pi()/2)*(1-k3R/2.);
62 Float_t k1=k2*k3R/(4*TMath::ATan(k3R));
63 Float_t l=x[0]/par[0];
64 Float_t tan2=TMath::TanH(k2*l);
66 Float_t res = k1*(1-tan2)/(1+k3*tan2);
70 ///////////////////////////////////////////////////////////////////////////
71 ///////////////////////////////////////////////////////////////////////////
78 AliTPCRF1D::AliTPCRF1D(Bool_t direct,Int_t np,Float_t step)
93 //default constructor for response function object
97 fcharge = new Float_t[fNRF];
98 if (step>0) fDSTEPM1=1./step;
99 else fDSTEPM1 = 1./fgRFDSTEP;
100 for(Int_t i=0;i<5;i++) {
107 AliTPCRF1D::AliTPCRF1D(const AliTPCRF1D &prf)
110 fDSTEPM1(prf.fDSTEPM1),
112 forigsigma(prf.forigsigma),
113 fpadWidth(prf.fpadWidth),
116 fGRF(new TF1(*(prf.fGRF))),
118 fOffset(prf.fOffset),
119 fDirect(prf.fDirect),
120 fPadDistance(prf.fPadDistance)
124 for(Int_t i=0;i<5;i++) {
128 fcharge = new Float_t[fNRF];
129 memcpy(fcharge,prf.fcharge, fNRF*sizeof(Float_t));
131 //PH Change the name (add 0 to the end)
132 TString s(fGRF->GetName());
134 fGRF->SetName(s.Data());
137 AliTPCRF1D & AliTPCRF1D::operator = (const AliTPCRF1D &prf)
140 TObject::operator=(prf);
142 fDSTEPM1=prf.fDSTEPM1;
144 fcharge = new Float_t[fNRF];
145 memcpy(fcharge,prf.fcharge, fNRF*sizeof(Float_t));
146 forigsigma=prf.forigsigma;
147 fpadWidth=prf.fpadWidth;
151 fGRF=new TF1(*(prf.fGRF));
152 //PH Change the name (add 0 to the end)
153 TString s(fGRF->GetName());
155 fGRF->SetName(s.Data());
159 fPadDistance=prf.fPadDistance;
166 AliTPCRF1D::~AliTPCRF1D()
173 Float_t AliTPCRF1D::GetRF(Float_t xin)
175 //function which return response
176 //for the charge in distance xin
177 //return linear aproximation of RF
178 Float_t x = (xin-fOffset)*fDSTEPM1+fNRF/2;
182 if (i1+1<fNRF &&i1>0)
183 res = fcharge[i1]*(Float_t(i1+1)-x)+fcharge[i1+1]*(x-Float_t(i1));
187 Float_t AliTPCRF1D::GetGRF(Float_t xin)
189 //function which returnoriginal charge distribution
190 //this function is just normalised for fKnorm
192 return fkNorm*fGRF->Eval(xin)/fInteg;
198 void AliTPCRF1D::SetParam( TF1 * GRF,Float_t padwidth,
199 Float_t kNorm, Float_t sigma)
201 //adjust parameters of the original charge distribution
202 //and pad size parameters
203 fpadWidth = padwidth;
206 if (sigma==0) sigma= fpadWidth/TMath::Sqrt(12.);
208 fDSTEPM1 = 10/TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
209 //sprintf(fType,"User");
210 snprintf(fType,5,"User");
215 void AliTPCRF1D::SetGauss(Float_t sigma, Float_t padWidth,
219 // set parameters for Gauss generic charge distribution
221 fpadWidth = padWidth;
223 if (fGRF !=0 ) fGRF->Delete();
224 fGRF = new TF1("funGauss",funGauss,-5,5,1);
227 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,"Gauss");
231 snprintf(fType,5,"Gauss");
234 void AliTPCRF1D::SetCosh(Float_t sigma, Float_t padWidth,
238 // set parameters for Cosh generic charge distribution
240 fpadWidth = padWidth;
242 if (fGRF !=0 ) fGRF->Delete();
243 fGRF = new TF1("funCosh", funCosh, -5.,5.,2);
245 fGRF->SetParameters(funParam);
247 fDSTEPM1 = 10./TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
248 //by default I set the step as one tenth of sigma
249 //sprintf(fType,"Cosh");
250 snprintf(fType,5,"Cosh");
253 void AliTPCRF1D::SetGati(Float_t K3, Float_t padDistance, Float_t padWidth,
257 // set parameters for Gati generic charge distribution
259 fpadWidth = padWidth;
261 if (fGRF !=0 ) fGRF->Delete();
262 fGRF = new TF1("funGati", funGati, -5.,5.,2);
263 funParam[0]=padDistance;
265 fGRF->SetParameters(funParam);
266 forigsigma=padDistance;
267 fDSTEPM1 = 10./TMath::Sqrt(padDistance*padDistance+fpadWidth*fpadWidth/12);
268 //by default I set the step as one tenth of sigma
269 //sprintf(fType,"Gati");
270 snprintf(fType,5,"Gati");
275 void AliTPCRF1D::DrawRF(Float_t x1,Float_t x2,Int_t N)
278 //Draw prf in selected region <x1,x2> with nuber of diviision = n
281 TCanvas * c1 = new TCanvas("canRF","Pad response function",700,900);
283 TPad * pad1 = new TPad("pad1RF","",0.05,0.55,0.95,0.95,21);
285 TPad * pad2 = new TPad("pad2RF","",0.05,0.05,0.95,0.45,21);
288 //sprintf(s,"RF response function for %1.2f cm pad width",
290 snprintf(s,60,"RF response function for %1.2f cm pad width",fpadWidth);
292 TH1F * hRFo = new TH1F("hRFo","Original charge distribution",N+1,x1,x2);
294 gStyle->SetOptFit(1);
295 gStyle->SetOptStat(0);
296 TH1F * hRFc = new TH1F("hRFc",s,N+1,x1,x2);
301 for (Float_t i = 0;i<N+1;i++)
303 x+=(x2-x1)/Float_t(N);
315 void AliTPCRF1D::Update()
318 //update fields with interpolated values for
321 //at the begining initialize to 0
322 for (Int_t i =0; i<fNRF;i++) fcharge[i] = 0;
323 if ( fGRF == 0 ) return;
324 // This form is no longer available
325 #if ROOT_VERSION_CODE < ROOT_VERSION(5,99,0)
326 fInteg = fGRF->Integral(-5*forigsigma,5*forigsigma,funParam,0.00001);
328 TArrayD savParam(fGRF->GetNpar(), fGRF->GetParameters());
329 fGRF->SetParameters(funParam);
330 fInteg = fGRF->Integral(-5*forigsigma,5*forigsigma,0.00001);
332 if ( fInteg == 0 ) fInteg = 1;
333 if (fDirect==kFALSE){
334 //integrate charge over pad for different distance of pad
335 for (Int_t i =0; i<fNRF;i++)
336 { //x in cm fpadWidth in cm
337 Float_t x = (Float_t)(i-fNRF/2)/fDSTEPM1;
338 Float_t x1=TMath::Max(x-fpadWidth/2,-5*forigsigma);
339 Float_t x2=TMath::Min(x+fpadWidth/2,5*forigsigma);
340 #if ROOT_VERSION_CODE < ROOT_VERSION(5,99,0)
341 fcharge[i] = fkNorm*fGRF->Integral(x1,x2,funParam,0.0001)/fInteg;
343 fcharge[i] = fkNorm*fGRF->Integral(x1,x2,0.0001)/fInteg;
348 for (Int_t i =0; i<fNRF;i++)
349 { //x in cm fpadWidth in cm
350 Float_t x = (Float_t)(i-fNRF/2)/fDSTEPM1;
351 fcharge[i] = fkNorm*fGRF->Eval(x);
357 for (Float_t x =-fNRF/fDSTEPM1; x<fNRF/fDSTEPM1;x+=1/fDSTEPM1)
358 { //x in cm fpadWidth in cm
359 Float_t weight = GetRF(x+fOffset);
366 fSigma = TMath::Sqrt(fSigma/sum-mean*mean);
369 #if ROOT_VERSION_CODE >= ROOT_VERSION(5,99,0)
370 fGRF->SetParameters(savParam.GetArray());
374 void AliTPCRF1D::Streamer(TBuffer &R__b)
376 // Stream an object of class AliTPCRF1D.
377 if (R__b.IsReading()) {
378 AliTPCRF1D::Class()->ReadBuffer(R__b, this);
381 if (strncmp(fType,"Gauss",3)==0) {delete fGRF; fGRF = new TF1("funGauss",funGauss,-5.,5.,4);}
382 if (strncmp(fType,"Cosh",3)==0) {delete fGRF; fGRF = new TF1("funCosh",funCosh,-5.,5.,4);}
383 if (strncmp(fType,"Gati",3)==0) {delete fGRF; fGRF = new TF1("funGati",funGati,-5.,5.,4);}
384 if (fGRF) fGRF->SetParameters(funParam);
387 AliTPCRF1D::Class()->WriteBuffer(R__b, this);
392 Double_t AliTPCRF1D::Gamma4(Double_t x, Double_t p0, Double_t p1){
394 // Gamma 4 Time response function of ALTRO
397 Double_t g1 = TMath::Exp(-4.*x/p1);
398 Double_t g2 = TMath::Power(x/p1,4);