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