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