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