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