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