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