Line data Source code
1 : /**************************************************************************
2 : * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 : * *
4 : * Author: The ALICE Off-line Project. *
5 : * Contributors are mentioned in the code where appropriate. *
6 : * *
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 : **************************************************************************/
15 :
16 :
17 :
18 : /// \class AliTPCRF1D
19 : /// \brief Declaration of class AliTPCRF1D
20 : ///
21 : /// \author Marian Ivanov, Uni. of Bratislava, ivanov@fmph.uniba.sk
22 :
23 : #include <RVersion.h>
24 : #include <Riostream.h>
25 : #include <TCanvas.h>
26 : #include <TClass.h>
27 : #include <TF2.h>
28 : #include <TH1.h>
29 : #include <TMath.h>
30 : #include <TPad.h>
31 : #include <TString.h>
32 : #include <TStyle.h>
33 : #include "TBuffer.h"
34 :
35 : #include "AliTPCRF1D.h"
36 :
37 : extern TStyle * gStyle;
38 :
39 : Int_t AliTPCRF1D::fgNRF=100; ///< default number of interpolation points
40 : Float_t AliTPCRF1D::fgRFDSTEP=0.01; ///< default step in cm
41 :
42 : static Double_t funGauss(Double_t *x, Double_t * par)
43 : {
44 : /// Gauss function -needde by the generic function object
45 :
46 0 : return TMath::Exp(-(x[0]*x[0])/(2*par[0]*par[0]));
47 : }
48 :
49 : static Double_t funCosh(Double_t *x, Double_t * par)
50 : {
51 : /// Cosh function -needde by the generic function object
52 :
53 0 : return 1/TMath::CosH(3.14159*x[0]/(2*par[0]));
54 : }
55 :
56 : static Double_t funGati(Double_t *x, Double_t * par)
57 : {
58 : /// Gati function -needde by the generic function object
59 :
60 0 : Float_t k3=par[1];
61 0 : Float_t k3R=TMath::Sqrt(k3);
62 0 : Float_t k2=(TMath::Pi()/2)*(1-k3R/2.);
63 0 : Float_t k1=k2*k3R/(4*TMath::ATan(k3R));
64 0 : Float_t l=x[0]/par[0];
65 0 : Float_t tan2=TMath::TanH(k2*l);
66 0 : tan2*=tan2;
67 0 : Float_t res = k1*(1-tan2)/(1+k3*tan2);
68 0 : return res;
69 : }
70 :
71 : ///////////////////////////////////////////////////////////////////////////
72 : ///////////////////////////////////////////////////////////////////////////
73 :
74 : /// \cond CLASSIMP
75 24 : ClassImp(AliTPCRF1D)
76 : /// \endcond
77 :
78 :
79 : AliTPCRF1D::AliTPCRF1D(Bool_t direct,Int_t np,Float_t step)
80 1 : :TObject(),
81 1 : fNRF(0),
82 1 : fDSTEPM1(0.),
83 1 : fcharge(0),
84 1 : forigsigma(0.),
85 1 : fpadWidth(3.5),
86 1 : fkNorm(0.5),
87 1 : fInteg(0.),
88 1 : fGRF(0),
89 1 : fSigma(0.),
90 1 : fOffset(0.),
91 1 : fDirect(kFALSE),
92 1 : fPadDistance(0.)
93 3 : {
94 : //default constructor for response function object
95 1 : fDirect=direct;
96 2 : if (np!=0)fNRF = np;
97 0 : else (fNRF=fgNRF);
98 2 : fcharge = new Float_t[fNRF];
99 1 : if (step>0) fDSTEPM1=1./step;
100 1 : else fDSTEPM1 = 1./fgRFDSTEP;
101 12 : for(Int_t i=0;i<5;i++) {
102 5 : funParam[i]=0.;
103 5 : fType[i]=0;
104 : }
105 :
106 2 : }
107 :
108 : AliTPCRF1D::AliTPCRF1D(const AliTPCRF1D &prf)
109 0 : :TObject(prf),
110 0 : fNRF(prf.fNRF),
111 0 : fDSTEPM1(prf.fDSTEPM1),
112 0 : fcharge(0),
113 0 : forigsigma(prf.forigsigma),
114 0 : fpadWidth(prf.fpadWidth),
115 0 : fkNorm(prf.fkNorm),
116 0 : fInteg(prf.fInteg),
117 0 : fGRF(new TF1(*(prf.fGRF))),
118 0 : fSigma(prf.fSigma),
119 0 : fOffset(prf.fOffset),
120 0 : fDirect(prf.fDirect),
121 0 : fPadDistance(prf.fPadDistance)
122 0 : {
123 : //
124 : //
125 0 : for(Int_t i=0;i<5;i++) {
126 0 : funParam[i]=0.;
127 0 : fType[i]=0;
128 : }
129 0 : fcharge = new Float_t[fNRF];
130 0 : memcpy(fcharge,prf.fcharge, fNRF*sizeof(Float_t));
131 :
132 : //PH Change the name (add 0 to the end)
133 0 : TString s(fGRF->GetName());
134 0 : s+="0";
135 0 : fGRF->SetName(s.Data());
136 0 : }
137 :
138 : AliTPCRF1D & AliTPCRF1D::operator = (const AliTPCRF1D &prf)
139 : {
140 0 : if(this!=&prf) {
141 0 : TObject::operator=(prf);
142 0 : fNRF=prf.fNRF;
143 0 : fDSTEPM1=prf.fDSTEPM1;
144 0 : delete [] fcharge;
145 0 : fcharge = new Float_t[fNRF];
146 0 : memcpy(fcharge,prf.fcharge, fNRF*sizeof(Float_t));
147 0 : forigsigma=prf.forigsigma;
148 0 : fpadWidth=prf.fpadWidth;
149 0 : fkNorm=prf.fkNorm;
150 0 : fInteg=prf.fInteg;
151 0 : delete fGRF;
152 0 : fGRF=new TF1(*(prf.fGRF));
153 : //PH Change the name (add 0 to the end)
154 0 : TString s(fGRF->GetName());
155 0 : s+="0";
156 0 : fGRF->SetName(s.Data());
157 0 : fSigma=prf.fSigma;
158 0 : fOffset=prf.fOffset;
159 0 : fDirect=prf.fDirect;
160 0 : fPadDistance=prf.fPadDistance;
161 0 : }
162 0 : return *this;
163 0 : }
164 :
165 :
166 :
167 : AliTPCRF1D::~AliTPCRF1D()
168 6 : {
169 : //
170 2 : delete [] fcharge;
171 2 : delete fGRF;
172 3 : }
173 :
174 : Float_t AliTPCRF1D::GetRF(Float_t xin)
175 : {
176 : //function which return response
177 : //for the charge in distance xin
178 : //return linear aproximation of RF
179 9002 : Float_t x = (xin-fOffset)*fDSTEPM1+fNRF/2;
180 4501 : Int_t i1=Int_t(x);
181 5002 : if (x<0) i1-=1;
182 : Float_t res=0;
183 4501 : if (i1+1<fNRF &&i1>0)
184 3498 : res = fcharge[i1]*(Float_t(i1+1)-x)+fcharge[i1+1]*(x-Float_t(i1));
185 4501 : return res;
186 : }
187 :
188 : Float_t AliTPCRF1D::GetGRF(Float_t xin)
189 : {
190 : //function which returnoriginal charge distribution
191 : //this function is just normalised for fKnorm
192 0 : if (fGRF != 0 )
193 0 : return fkNorm*fGRF->Eval(xin)/fInteg;
194 : else
195 0 : return 0.;
196 0 : }
197 :
198 :
199 : void AliTPCRF1D::SetParam( TF1 * GRF,Float_t padwidth,
200 : Float_t kNorm, Float_t sigma)
201 : {
202 : //adjust parameters of the original charge distribution
203 : //and pad size parameters
204 2 : fpadWidth = padwidth;
205 1 : fGRF = GRF;
206 1 : fkNorm = kNorm;
207 1 : if (sigma==0) sigma= fpadWidth/TMath::Sqrt(12.);
208 1 : forigsigma=sigma;
209 1 : fDSTEPM1 = 10/TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
210 : //sprintf(fType,"User");
211 1 : snprintf(fType,5,"User");
212 : // Update();
213 1 : }
214 :
215 :
216 : void AliTPCRF1D::SetGauss(Float_t sigma, Float_t padWidth,
217 : Float_t kNorm)
218 : {
219 : //
220 : // set parameters for Gauss generic charge distribution
221 : //
222 0 : fpadWidth = padWidth;
223 0 : fkNorm = kNorm;
224 0 : if (fGRF !=0 ) fGRF->Delete();
225 0 : fGRF = new TF1("funGauss",funGauss,-5,5,1);
226 0 : funParam[0]=sigma;
227 0 : forigsigma=sigma;
228 0 : fGRF->SetParameters(funParam);
229 0 : fDSTEPM1 = 10./TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
230 : //by default I set the step as one tenth of sigma
231 : //sprintf(fType,"Gauss");
232 0 : snprintf(fType,5,"Gauss");
233 0 : }
234 :
235 : void AliTPCRF1D::SetCosh(Float_t sigma, Float_t padWidth,
236 : Float_t kNorm)
237 : {
238 : //
239 : // set parameters for Cosh generic charge distribution
240 : //
241 0 : fpadWidth = padWidth;
242 0 : fkNorm = kNorm;
243 0 : if (fGRF !=0 ) fGRF->Delete();
244 0 : fGRF = new TF1("funCosh", funCosh, -5.,5.,2);
245 0 : funParam[0]=sigma;
246 0 : fGRF->SetParameters(funParam);
247 0 : forigsigma=sigma;
248 0 : fDSTEPM1 = 10./TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
249 : //by default I set the step as one tenth of sigma
250 : //sprintf(fType,"Cosh");
251 0 : snprintf(fType,5,"Cosh");
252 0 : }
253 :
254 : void AliTPCRF1D::SetGati(Float_t K3, Float_t padDistance, Float_t padWidth,
255 : Float_t kNorm)
256 : {
257 : //
258 : // set parameters for Gati generic charge distribution
259 : //
260 0 : fpadWidth = padWidth;
261 0 : fkNorm = kNorm;
262 0 : if (fGRF !=0 ) fGRF->Delete();
263 0 : fGRF = new TF1("funGati", funGati, -5.,5.,2);
264 0 : funParam[0]=padDistance;
265 0 : funParam[1]=K3;
266 0 : fGRF->SetParameters(funParam);
267 0 : forigsigma=padDistance;
268 0 : fDSTEPM1 = 10./TMath::Sqrt(padDistance*padDistance+fpadWidth*fpadWidth/12);
269 : //by default I set the step as one tenth of sigma
270 : //sprintf(fType,"Gati");
271 0 : snprintf(fType,5,"Gati");
272 0 : }
273 :
274 :
275 :
276 : void AliTPCRF1D::DrawRF(Float_t x1,Float_t x2,Int_t N)
277 : {
278 : //
279 : //Draw prf in selected region <x1,x2> with nuber of diviision = n
280 : //
281 0 : char s[100];
282 0 : TCanvas * c1 = new TCanvas("canRF","Pad response function",700,900);
283 0 : c1->cd();
284 0 : TPad * pad1 = new TPad("pad1RF","",0.05,0.55,0.95,0.95,21);
285 0 : pad1->Draw();
286 0 : TPad * pad2 = new TPad("pad2RF","",0.05,0.05,0.95,0.45,21);
287 0 : pad2->Draw();
288 :
289 : //sprintf(s,"RF response function for %1.2f cm pad width",
290 : // fpadWidth);
291 0 : snprintf(s,60,"RF response function for %1.2f cm pad width",fpadWidth);
292 0 : pad1->cd();
293 0 : TH1F * hRFo = new TH1F("hRFo","Original charge distribution",N+1,x1,x2);
294 0 : pad2->cd();
295 0 : gStyle->SetOptFit(1);
296 0 : gStyle->SetOptStat(0);
297 0 : TH1F * hRFc = new TH1F("hRFc",s,N+1,x1,x2);
298 : Float_t x=x1;
299 : Float_t y1;
300 : Float_t y2;
301 :
302 0 : for (Float_t i = 0;i<N+1;i++)
303 : {
304 0 : x+=(x2-x1)/Float_t(N);
305 0 : y1 = GetRF(x);
306 0 : hRFc->Fill(x,y1);
307 0 : y2 = GetGRF(x);
308 0 : hRFo->Fill(x,y2);
309 : };
310 0 : pad1->cd();
311 0 : hRFo->Fit("gaus");
312 0 : pad2->cd();
313 0 : hRFc->Fit("gaus");
314 0 : }
315 :
316 : void AliTPCRF1D::Update()
317 : {
318 : //
319 : //update fields with interpolated values for
320 : //PRF calculation
321 :
322 : //at the begining initialize to 0
323 2003 : for (Int_t i =0; i<fNRF;i++) fcharge[i] = 0;
324 1 : if ( fGRF == 0 ) return;
325 : // This form is no longer available
326 : #if ROOT_VERSION_CODE < ROOT_VERSION(5,99,0)
327 1 : fInteg = fGRF->Integral(-5*forigsigma,5*forigsigma,funParam,0.00001);
328 : #else
329 : TArrayD savParam(fGRF->GetNpar(), fGRF->GetParameters());
330 : fGRF->SetParameters(funParam);
331 : fInteg = fGRF->Integral(-5*forigsigma,5*forigsigma,0.00001);
332 : #endif
333 1 : if ( fInteg == 0 ) fInteg = 1;
334 1 : if (fDirect==kFALSE){
335 : //integrate charge over pad for different distance of pad
336 0 : for (Int_t i =0; i<fNRF;i++)
337 : { //x in cm fpadWidth in cm
338 0 : Float_t x = (Float_t)(i-fNRF/2)/fDSTEPM1;
339 0 : Float_t x1=TMath::Max(x-fpadWidth/2,-5*forigsigma);
340 0 : Float_t x2=TMath::Min(x+fpadWidth/2,5*forigsigma);
341 : #if ROOT_VERSION_CODE < ROOT_VERSION(5,99,0)
342 0 : fcharge[i] = fkNorm*fGRF->Integral(x1,x2,funParam,0.0001)/fInteg;
343 : #else
344 : fcharge[i] = fkNorm*fGRF->Integral(x1,x2,0.0001)/fInteg;
345 : #endif
346 : };
347 0 : }
348 : else{
349 2002 : for (Int_t i =0; i<fNRF;i++)
350 : { //x in cm fpadWidth in cm
351 1000 : Float_t x = (Float_t)(i-fNRF/2)/fDSTEPM1;
352 1000 : fcharge[i] = fkNorm*fGRF->Eval(x);
353 : };
354 : }
355 1 : fSigma = 0;
356 : Float_t sum =0;
357 : Float_t mean=0;
358 4004 : for (Float_t x =-fNRF/fDSTEPM1; x<fNRF/fDSTEPM1;x+=1/fDSTEPM1)
359 : { //x in cm fpadWidth in cm
360 2001 : Float_t weight = GetRF(x+fOffset);
361 2001 : fSigma+=x*x*weight;
362 2001 : mean+=x*weight;
363 2001 : sum+=weight;
364 : };
365 1 : if (sum>0){
366 1 : mean/=sum;
367 1 : fSigma = TMath::Sqrt(fSigma/sum-mean*mean);
368 1 : }
369 0 : else fSigma=0;
370 : #if ROOT_VERSION_CODE >= ROOT_VERSION(5,99,0)
371 : fGRF->SetParameters(savParam.GetArray());
372 : #endif
373 2 : }
374 :
375 : void AliTPCRF1D::Streamer(TBuffer &R__b)
376 : {
377 : // Stream an object of class AliTPCRF1D.
378 0 : if (R__b.IsReading()) {
379 0 : AliTPCRF1D::Class()->ReadBuffer(R__b, this);
380 : //read functions
381 :
382 0 : if (strncmp(fType,"Gauss",3)==0) {delete fGRF; fGRF = new TF1("funGauss",funGauss,-5.,5.,4);}
383 0 : if (strncmp(fType,"Cosh",3)==0) {delete fGRF; fGRF = new TF1("funCosh",funCosh,-5.,5.,4);}
384 0 : if (strncmp(fType,"Gati",3)==0) {delete fGRF; fGRF = new TF1("funGati",funGati,-5.,5.,4);}
385 0 : if (fGRF) fGRF->SetParameters(funParam);
386 :
387 : } else {
388 0 : AliTPCRF1D::Class()->WriteBuffer(R__b, this);
389 : }
390 0 : }
391 :
392 :
393 : Double_t AliTPCRF1D::Gamma4(Double_t x, Double_t p0, Double_t p1){
394 : //
395 : // Gamma 4 Time response function of ALTRO
396 : //
397 2808 : if (x<0) return 0;
398 624 : Double_t g1 = TMath::Exp(-4.*x/p1);
399 624 : Double_t g2 = TMath::Power(x/p1,4);
400 624 : return p0*g1*g2;
401 1144 : }
402 :
|