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