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