]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/TPCbase/AliTPCRF1D.cxx
CMake: Retrieve Git information
[u/mrichter/AliRoot.git] / TPC / TPCbase / AliTPCRF1D.cxx
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
34 #include "AliTPCRF1D.h"
35
36 extern TStyle * gStyle;
37
38 Int_t   AliTPCRF1D::fgNRF=100;  ///< default  number of interpolation points
39 Float_t AliTPCRF1D::fgRFDSTEP=0.01; ///< default step in cm
40
41 static Double_t funGauss(Double_t *x, Double_t * par)
42 {
43   /// Gauss function  -needde by the generic function object
44
45   return TMath::Exp(-(x[0]*x[0])/(2*par[0]*par[0]));
46 }
47
48 static Double_t funCosh(Double_t *x, Double_t * par)
49 {
50   /// Cosh function  -needde by the generic function object
51
52   return 1/TMath::CosH(3.14159*x[0]/(2*par[0]));
53 }
54
55 static Double_t funGati(Double_t *x, Double_t * par)
56 {
57   /// Gati function  -needde by the generic function object
58
59   Float_t k3=par[1];
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);
65   tan2*=tan2;
66   Float_t res = k1*(1-tan2)/(1+k3*tan2);
67   return res;
68 }
69
70 ///////////////////////////////////////////////////////////////////////////
71 ///////////////////////////////////////////////////////////////////////////
72
73 /// \cond CLASSIMP
74 ClassImp(AliTPCRF1D)
75 /// \endcond
76
77
78 AliTPCRF1D::AliTPCRF1D(Bool_t direct,Int_t np,Float_t step)
79            :TObject(),
80             fNRF(0),
81             fDSTEPM1(0.),
82             fcharge(0),
83             forigsigma(0.),
84             fpadWidth(3.5),
85             fkNorm(0.5),
86             fInteg(0.),
87             fGRF(0),
88             fSigma(0.),
89             fOffset(0.),
90             fDirect(kFALSE),
91             fPadDistance(0.)
92 {
93   //default constructor for response function object
94   fDirect=direct;
95   if (np!=0)fNRF = np;
96   else (fNRF=fgNRF);
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++) {
101     funParam[i]=0.;
102     fType[i]=0;
103   }
104
105 }
106
107 AliTPCRF1D::AliTPCRF1D(const AliTPCRF1D &prf)
108            :TObject(prf),
109             fNRF(prf.fNRF),
110             fDSTEPM1(prf.fDSTEPM1),
111             fcharge(0),
112             forigsigma(prf.forigsigma),
113             fpadWidth(prf.fpadWidth),
114             fkNorm(prf.fkNorm),
115             fInteg(prf.fInteg),
116             fGRF(new TF1(*(prf.fGRF))),
117             fSigma(prf.fSigma),
118             fOffset(prf.fOffset),
119             fDirect(prf.fDirect),
120             fPadDistance(prf.fPadDistance)
121 {
122   //
123   //
124   for(Int_t i=0;i<5;i++) {
125     funParam[i]=0.;
126     fType[i]=0;
127   }
128   fcharge = new Float_t[fNRF];
129   memcpy(fcharge,prf.fcharge, fNRF*sizeof(Float_t));
130
131   //PH Change the name (add 0 to the end)
132   TString s(fGRF->GetName());
133   s+="0";
134   fGRF->SetName(s.Data());
135 }
136
137 AliTPCRF1D & AliTPCRF1D::operator = (const AliTPCRF1D &prf)
138 {
139   if(this!=&prf) {
140     TObject::operator=(prf);
141     fNRF=prf.fNRF;
142     fDSTEPM1=prf.fDSTEPM1;
143     delete [] fcharge;
144     fcharge = new Float_t[fNRF];
145     memcpy(fcharge,prf.fcharge, fNRF*sizeof(Float_t));
146     forigsigma=prf.forigsigma;
147     fpadWidth=prf.fpadWidth;
148     fkNorm=prf.fkNorm;
149     fInteg=prf.fInteg;
150     delete fGRF;
151     fGRF=new TF1(*(prf.fGRF));
152    //PH Change the name (add 0 to the end)
153     TString s(fGRF->GetName());
154     s+="0";
155     fGRF->SetName(s.Data());
156     fSigma=prf.fSigma;
157     fOffset=prf.fOffset;
158     fDirect=prf.fDirect;
159     fPadDistance=prf.fPadDistance;
160   }
161   return *this;
162 }
163
164
165
166 AliTPCRF1D::~AliTPCRF1D()
167 {
168   //
169   delete [] fcharge;
170   delete fGRF;
171 }
172
173 Float_t AliTPCRF1D::GetRF(Float_t xin)
174 {
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;
179   Int_t i1=Int_t(x);
180   if (x<0) i1-=1;
181   Float_t res=0;
182   if (i1+1<fNRF &&i1>0)
183     res = fcharge[i1]*(Float_t(i1+1)-x)+fcharge[i1+1]*(x-Float_t(i1));
184   return res;
185 }
186
187 Float_t  AliTPCRF1D::GetGRF(Float_t xin)
188 {
189   //function which returnoriginal charge distribution
190   //this function is just normalised for fKnorm
191   if (fGRF != 0 )
192     return fkNorm*fGRF->Eval(xin)/fInteg;
193       else
194     return 0.;
195 }
196
197
198 void AliTPCRF1D::SetParam( TF1 * GRF,Float_t padwidth,
199                        Float_t kNorm, Float_t sigma)
200 {
201   //adjust parameters of the original charge distribution
202   //and pad size parameters
203    fpadWidth = padwidth;
204    fGRF = GRF;
205    fkNorm = kNorm;
206    if (sigma==0) sigma= fpadWidth/TMath::Sqrt(12.);
207    forigsigma=sigma;
208    fDSTEPM1 = 10/TMath::Sqrt(sigma*sigma+fpadWidth*fpadWidth/12);
209    //sprintf(fType,"User");
210    snprintf(fType,5,"User");
211    //   Update();
212 }
213
214
215 void AliTPCRF1D::SetGauss(Float_t sigma, Float_t padWidth,
216                       Float_t kNorm)
217 {
218   //
219   // set parameters for Gauss generic charge distribution
220   //
221   fpadWidth = padWidth;
222   fkNorm = kNorm;
223   if (fGRF !=0 ) fGRF->Delete();
224   fGRF = new TF1("funGauss",funGauss,-5,5,1);
225   funParam[0]=sigma;
226   forigsigma=sigma;
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");
232 }
233
234 void AliTPCRF1D::SetCosh(Float_t sigma, Float_t padWidth,
235                      Float_t kNorm)
236 {
237   //
238   // set parameters for Cosh generic charge distribution
239   //
240   fpadWidth = padWidth;
241   fkNorm = kNorm;
242   if (fGRF !=0 ) fGRF->Delete();
243   fGRF = new TF1("funCosh",     funCosh, -5.,5.,2);
244   funParam[0]=sigma;
245   fGRF->SetParameters(funParam);
246   forigsigma=sigma;
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");
251 }
252
253 void AliTPCRF1D::SetGati(Float_t K3, Float_t padDistance, Float_t padWidth,
254                      Float_t kNorm)
255 {
256   //
257   // set parameters for Gati generic charge distribution
258   //
259   fpadWidth = padWidth;
260   fkNorm = kNorm;
261   if (fGRF !=0 ) fGRF->Delete();
262   fGRF = new TF1("funGati",     funGati, -5.,5.,2);
263   funParam[0]=padDistance;
264   funParam[1]=K3;
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");
271 }
272
273
274
275 void AliTPCRF1D::DrawRF(Float_t x1,Float_t x2,Int_t N)
276 {
277   //
278   //Draw prf in selected region <x1,x2> with nuber of diviision = n
279   //
280   char s[100];
281   TCanvas  * c1 = new TCanvas("canRF","Pad response function",700,900);
282   c1->cd();
283   TPad * pad1 = new TPad("pad1RF","",0.05,0.55,0.95,0.95,21);
284   pad1->Draw();
285   TPad * pad2 = new TPad("pad2RF","",0.05,0.05,0.95,0.45,21);
286   pad2->Draw();
287
288   //sprintf(s,"RF response function for %1.2f cm pad width",
289   //      fpadWidth);
290   snprintf(s,60,"RF response function for %1.2f cm pad width",fpadWidth);
291   pad1->cd();
292   TH1F * hRFo = new TH1F("hRFo","Original charge distribution",N+1,x1,x2);
293   pad2->cd();
294    gStyle->SetOptFit(1);
295    gStyle->SetOptStat(0);
296   TH1F * hRFc = new TH1F("hRFc",s,N+1,x1,x2);
297   Float_t x=x1;
298   Float_t y1;
299   Float_t y2;
300
301   for (Float_t i = 0;i<N+1;i++)
302     {
303       x+=(x2-x1)/Float_t(N);
304       y1 = GetRF(x);
305       hRFc->Fill(x,y1);
306       y2 = GetGRF(x);
307       hRFo->Fill(x,y2);
308     };
309   pad1->cd();
310   hRFo->Fit("gaus");
311   pad2->cd();
312   hRFc->Fit("gaus");
313 }
314
315 void AliTPCRF1D::Update()
316 {
317   //
318   //update fields  with interpolated values for
319   //PRF calculation
320
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);
327 #else
328   TArrayD savParam(fGRF->GetNpar(), fGRF->GetParameters());
329   fGRF->SetParameters(funParam);
330   fInteg  = fGRF->Integral(-5*forigsigma,5*forigsigma,0.00001);
331 #endif
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;
342 #else
343       fcharge[i] = fkNorm*fGRF->Integral(x1,x2,0.0001)/fInteg;
344 #endif
345     };
346   }
347   else{
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);
352       };
353   }
354   fSigma = 0;
355   Float_t sum =0;
356   Float_t mean=0;
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);
360       fSigma+=x*x*weight;
361       mean+=x*weight;
362       sum+=weight;
363     };
364   if (sum>0){
365     mean/=sum;
366     fSigma = TMath::Sqrt(fSigma/sum-mean*mean);
367   }
368   else fSigma=0;
369 #if ROOT_VERSION_CODE >= ROOT_VERSION(5,99,0)
370   fGRF->SetParameters(savParam.GetArray());
371 #endif
372 }
373
374 void AliTPCRF1D::Streamer(TBuffer &R__b)
375 {
376    // Stream an object of class AliTPCRF1D.
377    if (R__b.IsReading()) {
378       AliTPCRF1D::Class()->ReadBuffer(R__b, this);
379       //read functions
380
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);
385
386    } else {
387       AliTPCRF1D::Class()->WriteBuffer(R__b, this);
388    }
389 }
390
391
392 Double_t  AliTPCRF1D::Gamma4(Double_t x, Double_t p0, Double_t p1){
393   //
394   // Gamma 4 Time response function of ALTRO
395   //
396   if (x<0) return 0;
397   Double_t g1 = TMath::Exp(-4.*x/p1);
398   Double_t g2 = TMath::Power(x/p1,4);
399   return p0*g1*g2;
400 }
401