]> git.uio.no Git - u/mrichter/AliRoot.git/blame - FASTSIM/AliMUONFastTracking.cxx
Front absorber geometry as built.
[u/mrichter/AliRoot.git] / FASTSIM / AliMUONFastTracking.cxx
CommitLineData
6255180c 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
d53fc881 16/*
17$Log$
9a226b7e 18Revision 1.9 2004/02/03 16:45:20 morsch
19Unique name for functions (TF1).
20
c6ba3ff4 21Revision 1.8 2004/01/27 18:02:38 hristov
22Removing some warning (Sun)
23
15a060e4 24Revision 1.7 2003/11/13 14:21:57 morsch
25Coding Rule violation corrections.
26
041f7f97 27Revision 1.6 2003/08/12 15:16:25 morsch
28Saver initialisation of fFitp array. (Lenaic COUEDEL)
29
1e68a0c5 30Revision 1.5 2003/08/05 16:14:20 morsch
31Some problems with too big fluctuations corrected. (A. de Falco)
32
d53fc881 33Revision 1.2 2003/01/08 10:29:33 morsch
34Path to data file changed.
35
36Revision 1.1 2003/01/06 10:13:33 morsch
37First commit.
38
39*/
6255180c 40
041f7f97 41//-------------------------------------------------------------------------
42// Class AliMUONFastTracking
43//
44// Manager for the fast simulation of tracking in the muon spectrometer
45// This class reads the lookup tables containing the parameterization
46// of the deltap, deltatheta, deltaphi for different background levels
47// and provides the related smeared parameters.
48// Used by AliFastMuonTrackingEff, AliFastMuonTrackingAcc,
49// AliFastMuonTrackingRes.
50//-------------------------------------------------------------------------
51
6255180c 52#include "AliMUONFastTracking.h"
53#include "AliMUONFastTrackingEntry.h"
6255180c 54#include <TSpline.h>
55#include <TFile.h>
6255180c 56#include <TH3.h>
57#include <TF1.h>
58#include <TRandom.h>
59#include <stdlib.h>
60#include <stdio.h>
61#include <string.h>
62#include <Riostream.h>
63
64ClassImp(AliMUONFastTracking)
65
041f7f97 66
6255180c 67AliMUONFastTracking* AliMUONFastTracking::fgMUONFastTracking=NULL;
68
69static Double_t FitP(Double_t *x, Double_t *par){
70 Double_t dx = x[0] - par[0];
71 Double_t dx2 = x[0] - par[4];
72 Double_t sigma = par[1] * ( 1 + par[2] * dx);
73 if (sigma == 0) {
74
75 return 0.;
76 }
77 Double_t fasymm = TMath::Exp(-0.5 * dx * dx / (sigma * sigma));
78 Double_t sigma2 = par[1] * par[5];
79 Double_t fgauss = TMath::Exp(-0.5 * dx2 * dx2 / (sigma2 * sigma2));
d53fc881 80 Double_t value = fasymm + par[3] * fgauss;
9a226b7e 81 return TMath::Abs(value);
6255180c 82}
83
041f7f97 84AliMUONFastTracking::AliMUONFastTracking(const AliMUONFastTracking & ft):TObject()
85{
86// Copy constructor
87 ft.Copy(*this);
88}
89
90
6255180c 91AliMUONFastTracking* AliMUONFastTracking::Instance()
92{
93// Set random number generator
94 if (fgMUONFastTracking) {
95 return fgMUONFastTracking;
96 } else {
97 fgMUONFastTracking = new AliMUONFastTracking();
98 return fgMUONFastTracking;
99 }
100}
101
102AliMUONFastTracking::AliMUONFastTracking()
103{
041f7f97 104//
105// constructor
106//
1e68a0c5 107 for (Int_t i = 0; i<20;i++) {
108 for (Int_t j = 0; j<20; j++) {
109 for (Int_t k = 0; k<20; k++) {
110 fFitp[i][j][k] = 0x0;
111 }
112 }
113 }
114
d53fc881 115 fClusterFinder = kOld;
6255180c 116 fPrintLevel = 1;
117 // read binning; temporarily put by hand
118 Float_t pmin = 0, pmax = 200;
119 Int_t nbinp = 10;
120 Float_t thetamin = 2, thetamax = 9;
121 Int_t nbintheta=10;
122 Float_t phimin = -180, phimax =180;
123 Int_t nbinphi=10;
124 //--------------------------------------
125
126 fNbinp = nbinp;
127 fPmin = pmin;
128 fPmax = pmax;
129
130 fNbintheta = nbintheta;
131 fThetamin = thetamin;
132 fThetamax = thetamax;
133
134 fNbinphi = nbinphi;
135 fPhimin = phimin;
136 fPhimax = phimax;
137
138 fDeltaP = (fPmax-fPmin)/fNbinp;
139 fDeltaTheta = (fThetamax-fThetamin)/fNbintheta;
140 fDeltaPhi = (fPhimax-fPhimin)/fNbinphi;
141}
142
143void AliMUONFastTracking::Init(Float_t bkg)
144{
145 //
146 // Initialization
147 //
148 for (Int_t ip=0; ip< fNbinp; ip++){
149 for (Int_t itheta=0; itheta< fNbintheta; itheta++){
150 for (Int_t iphi=0; iphi< fNbinphi; iphi++){
151 fCurrentEntry[ip][itheta][iphi] = new AliMUONFastTrackingEntry;
152 for (Int_t ibkg=0; ibkg<4; ibkg++){
153 fEntry[ip][itheta][iphi][ibkg] = new AliMUONFastTrackingEntry;
154 }
155 }
156 }
157 }
158
159 char filename [100];
d53fc881 160 if (fClusterFinder==kOld) sprintf (filename,"$(ALICE_ROOT)/FASTSIM/data/MUONtrackLUT.root");
041f7f97 161 else sprintf (filename,"$(ALICE_ROOT)/FASTSIM/data/MUONtrackLUT-AZ.root");
d53fc881 162
6255180c 163 TFile *file = new TFile(filename);
164 ReadLUT(file);
165 SetBackground(bkg);
d53fc881 166 UseSpline(0);
6255180c 167}
168
169
170void AliMUONFastTracking::ReadLUT(TFile* file)
171{
041f7f97 172 //
173 // read the lookup tables from file
174 //
6255180c 175 TH3F *heff[5][3], *hacc[5][3], *hmeanp, *hsigmap, *hsigma1p, *hchi2p;
176 TH3F *hnormg2, *hmeang2, *hsigmag2, *hmeantheta, *hsigmatheta, *hchi2theta;
177 TH3F *hmeanphi, *hsigmaphi, *hchi2phi;
178 char tag[40], tag2[40];
179
d53fc881 180 printf ("Reading parameters from LUT file %s...\n",file->GetName());
181
041f7f97 182 const Float_t kBkg[4] = {0, 0.5, 1, 2};
6255180c 183 for (Int_t ibkg=0; ibkg<4; ibkg++) {
041f7f97 184 sprintf (tag,"BKG%g",kBkg[ibkg]);
6255180c 185 file->cd(tag);
186 for (Int_t isplp = 0; isplp<kSplitP; isplp++) {
187 for (Int_t ispltheta = 0; ispltheta<kSplitTheta; ispltheta++) {
188 sprintf (tag2,"heff[%d][%d]",isplp,ispltheta);
189 heff[isplp][ispltheta] = (TH3F*)gDirectory->Get(tag2);
190 sprintf (tag2,"hacc[%d][%d]",isplp,ispltheta);
191 hacc[isplp][ispltheta] = (TH3F*)gDirectory->Get(tag2);
192 }
193 }
194 hmeanp = (TH3F*)gDirectory->Get("hmeanp");
195 hsigmap = (TH3F*)gDirectory->Get("hsigmap");
196 hsigma1p = (TH3F*)gDirectory->Get("hsigma1p");
197 hchi2p = (TH3F*)gDirectory->Get("hchi2p");
198 hnormg2 = (TH3F*)gDirectory->Get("hnormg2");
199 hmeang2 = (TH3F*)gDirectory->Get("hmeang2");
200 hsigmag2 = (TH3F*)gDirectory->Get("hsigmag2");
201 hmeantheta = (TH3F*)gDirectory->Get("hmeantheta");
202 hsigmatheta = (TH3F*)gDirectory->Get("hsigmatheta");
203 hchi2theta = (TH3F*)gDirectory->Get("hchi2theta");
204 hmeanphi = (TH3F*)gDirectory->Get("hmeanphi");
205 hsigmaphi = (TH3F*)gDirectory->Get("hsigmaphi");
206 hchi2phi = (TH3F*)gDirectory->Get("hchi2phi");
207
6255180c 208 for (Int_t ip=0; ip<fNbinp ;ip++) {
209 for (Int_t itheta=0; itheta<fNbintheta ;itheta++) {
210 for (Int_t iphi=0; iphi<fNbinphi ;iphi++) {
211 Float_t p = fPmin + fDeltaP * (ip + 0.5);
212 Float_t theta = fThetamin + fDeltaTheta * (itheta + 0.5);
213 Float_t phi = fPhimin + fDeltaPhi * (iphi + 0.5);
214
215 fEntry[ip][itheta][iphi][ibkg]->fP = p;
216 fEntry[ip][itheta][iphi][ibkg]->fMeanp =
217 hmeanp->GetBinContent(ip+1,itheta+1,iphi+1);
218 fEntry[ip][itheta][iphi][ibkg]->fSigmap =
219 TMath::Abs(hsigmap->GetBinContent(ip+1,itheta+1,iphi+1));
220 fEntry[ip][itheta][iphi][ibkg]->fSigma1p =
221 hsigma1p->GetBinContent(ip+1,itheta+1,iphi+1);
222 fEntry[ip][itheta][iphi][ibkg]->fChi2p =
223 hchi2p->GetBinContent(ip+1,itheta+1,iphi+1);
224 fEntry[ip][itheta][iphi][ibkg]->fNormG2 =
225 hnormg2->GetBinContent(ip+1,itheta+1,iphi+1);
226 fEntry[ip][itheta][iphi][ibkg]->fMeanG2 =
227 hmeang2->GetBinContent(ip+1,itheta+1,iphi+1);
228 if (ibkg == 0) fEntry[ip][itheta][iphi][ibkg]->fSigmaG2 = 9999;
229 else fEntry[ip][itheta][iphi][ibkg]->fSigmaG2 =
230 hsigmag2->GetBinContent(ip+1,itheta+1,iphi+1);
231 fEntry[ip][itheta][iphi][ibkg]->fTheta = theta;
232 fEntry[ip][itheta][iphi][ibkg]->fMeantheta =
233 hmeantheta->GetBinContent(ip+1,itheta+1,iphi+1);
234 fEntry[ip][itheta][iphi][ibkg]->fSigmatheta =
235 TMath::Abs(hsigmatheta->GetBinContent(ip+1,itheta+1,iphi+1));
236 fEntry[ip][itheta][iphi][ibkg]->fChi2theta =
237 hchi2theta->GetBinContent(ip+1,itheta+1,iphi+1);
238 fEntry[ip][itheta][iphi][ibkg]->fPhi = phi;
239 fEntry[ip][itheta][iphi][ibkg]->fMeanphi =
240 hmeanphi->GetBinContent(ip+1,itheta+1,iphi+1);
241 fEntry[ip][itheta][iphi][ibkg]->fSigmaphi =
242 TMath::Abs(hsigmaphi->GetBinContent(ip+1,itheta+1,iphi+1));
243 fEntry[ip][itheta][iphi][ibkg]->fChi2phi =
244 hchi2phi->GetBinContent(ip+1,itheta+1,iphi+1);
245 for (Int_t i=0; i<kSplitP; i++) {
246 for (Int_t j=0; j<kSplitTheta; j++) {
247 fEntry[ip][itheta][iphi][ibkg]->fAcc[i][j] =
248 hacc[i][j]->GetBinContent(ip+1,itheta+1,iphi+1);
249 fEntry[ip][itheta][iphi][ibkg]->fEff[i][j] =
250 heff[i][j]->GetBinContent(ip+1,itheta+1,iphi+1);
251 }
252 }
253 } // iphi
254 } // itheta
255 } // ip
256 } // ibkg
257
258 TGraph *graph = new TGraph(3);
259 TF1 *f = new TF1("f","[0]+[1]*x");
260
261 for (Int_t ip=0; ip< fNbinp; ip++){
262 for (Int_t itheta=0; itheta< fNbintheta; itheta++){
263 for (Int_t iphi=0; iphi< fNbinphi; iphi++){
264 graph->SetPoint(0,0.5,fEntry[ip][itheta][iphi][1]->fSigmaG2);
265 graph->SetPoint(1,1,fEntry[ip][itheta][iphi][2]->fSigmaG2);
266 graph->SetPoint(2,2,fEntry[ip][itheta][iphi][3]->fSigmaG2);
267 graph->Fit("f","q");
268 fEntry[ip][itheta][iphi][0]->fSigmaG2 = f->Eval(0);
269 }
270 }
271 }
272 f->Delete();
273 graph->Delete();
274 printf ("parameters read. \n");
275}
276
277void AliMUONFastTracking::GetBinning(Int_t &nbinp, Float_t &pmin, Float_t &pmax,
278 Int_t &nbintheta, Float_t &thetamin,
279 Float_t &thetamax,
041f7f97 280 Int_t &nbinphi, Float_t &phimin, Float_t &phimax) const
6255180c 281{
041f7f97 282 //
283 // gets the binning for the discrete parametrizations in the lookup table
284 //
6255180c 285 nbinp = fNbinp;
286 pmin = fPmin;
287 pmax = fPmax;
288 nbintheta = fNbintheta;
289 thetamin = fThetamin;
290 thetamax = fThetamax;
291 nbinphi = fNbinphi;
292 phimin = fPhimin;
293 phimax = fPhimax;
294}
295
296
297void AliMUONFastTracking::GetIpIthetaIphi(Float_t p, Float_t theta, Float_t phi,
298 Int_t charge, Int_t &ip, Int_t &itheta,
041f7f97 299 Int_t &iphi) const
6255180c 300{
041f7f97 301 //
302 // gets the id of the cells in the LUT for a given (p,theta,phi, charge)
303 //
6255180c 304 if (charge < 0) phi = -phi;
305 ip = Int_t (( p - fPmin ) / fDeltaP);
306 itheta = Int_t (( theta - fThetamin ) / fDeltaTheta);
307 iphi = Int_t (( phi - fPhimin ) / fDeltaPhi);
308
1e68a0c5 309
d53fc881 310 if (ip< 0) ip = 0;
311 if (ip>= fNbinp) ip = fNbinp-1;
312 if (itheta< 0) itheta = 0;
313 if (itheta>= fNbintheta) itheta = fNbintheta-1;
6255180c 314
d53fc881 315 if (iphi< 0) iphi = 0;
316 if (iphi>= fNbinphi) iphi = fNbinphi-1;
6255180c 317}
318
319void AliMUONFastTracking::GetSplit(Int_t ip, Int_t itheta,
041f7f97 320 Int_t &nSplitP, Int_t &nSplitTheta) const
321{
322 //
323 // the first cell is splitted in more bins for theta and momentum
324 // parameterizations. Get the number of divisions for the splitted bins
325 //
6255180c 326 if (ip==0) nSplitP = 5;
327 else nSplitP = 2;
328 if (itheta==0) nSplitTheta = 3;
329 else nSplitTheta = 1;
330}
331
332Float_t AliMUONFastTracking::Efficiency(Float_t p, Float_t theta,
333 Float_t phi, Int_t charge){
041f7f97 334 //
335 // gets the tracking efficiency
336 //
6255180c 337 Int_t ip=0, itheta=0, iphi=0;
338 GetIpIthetaIphi(p,theta,phi,charge,ip,itheta,iphi);
339 Int_t nSplitP, nSplitTheta;
340 GetSplit(ip,itheta,nSplitP,nSplitTheta);
341
342 Float_t dp = p - fPmin;
343 Int_t ibinp = Int_t(nSplitP*(dp - fDeltaP * Int_t(dp / fDeltaP))/fDeltaP);
344 Float_t dtheta = theta - fThetamin;
345 Int_t ibintheta = Int_t(nSplitTheta*(dtheta - fDeltaTheta * Int_t(dtheta / fDeltaTheta))/fDeltaTheta);
346 Float_t eff = fCurrentEntry[ip][itheta][iphi]->fEff[ibinp][ibintheta];
347 return eff;
348}
349
350Float_t AliMUONFastTracking::Acceptance(Float_t p, Float_t theta,
351 Float_t phi, Int_t charge){
041f7f97 352 //
353 // gets the geometrical acceptance
354 //
6255180c 355 if (theta<fThetamin || theta>fThetamax) return 0;
356
357 Int_t ip=0, itheta=0, iphi=0;
358 GetIpIthetaIphi(p,theta,phi,charge,ip,itheta,iphi);
359 Int_t nSplitP, nSplitTheta;
360 GetSplit(ip,itheta,nSplitP,nSplitTheta);
361 // central value and corrections with spline
362
363 Float_t dp = p - fPmin;
364 Int_t ibinp = Int_t(nSplitP*(dp - fDeltaP * Int_t(dp / fDeltaP))/fDeltaP);
365 Float_t dtheta = theta - fThetamin;
366 Int_t ibintheta = Int_t(nSplitTheta*(dtheta - fDeltaTheta * Int_t(dtheta / fDeltaTheta))/fDeltaTheta);
367 Float_t acc = fCurrentEntry[ip][itheta][iphi]->fAcc[ibinp][ibintheta];
368 return acc;
369}
370
371Float_t AliMUONFastTracking::MeanP(Float_t p, Float_t theta,
041f7f97 372 Float_t phi, Int_t charge) const
6255180c 373{
041f7f97 374 //
375 // gets the mean value of the prec-pgen distribution
376 //
6255180c 377 Int_t ip=0, itheta=0, iphi=0;
378 GetIpIthetaIphi(p,theta,phi,charge,ip,itheta,iphi);
379 return fCurrentEntry[ip][itheta][iphi]->fMeanp;
380}
381
382Float_t AliMUONFastTracking::SigmaP(Float_t p, Float_t theta,
041f7f97 383 Float_t phi, Int_t charge) const
6255180c 384{
041f7f97 385 //
386 // gets the width of the prec-pgen distribution
387 //
6255180c 388 Int_t ip=0, itheta=0, iphi=0;
389 Int_t index;
390 GetIpIthetaIphi(p,theta,phi,charge,ip,itheta,iphi);
391 // central value and corrections with spline
392 Float_t sigmap = fCurrentEntry[ip][itheta][iphi]->fSigmap;
393 if (!fSpline) return sigmap;
394 // corrections vs p, theta, phi
395 index = iphi + fNbinphi * itheta;
396 Double_t xmin,ymin,xmax,ymax;
397 Float_t frac1 = fSplineSigmap[index][0]->Eval(p)/sigmap;
398
399 if (p>fPmax-fDeltaP/2.) {
400 Float_t s1 = fCurrentEntry[fNbinp-1][itheta][iphi]->fSigmap;
401 Float_t s2 = fCurrentEntry[fNbinp-2][itheta][iphi]->fSigmap;
402 Float_t s3 = fCurrentEntry[fNbinp-3][itheta][iphi]->fSigmap;
403 Float_t p1 = fDeltaP * (fNbinp - 1 + 0.5) + fPmin;
404 Float_t p2 = fDeltaP * (fNbinp - 2 + 0.5) + fPmin;
405 Float_t p3 = fDeltaP * (fNbinp - 3 + 0.5) + fPmin;
406 Float_t p12 = p1 * p1, p22 = p2 * p2, p32 = p3 * p3;
407 Float_t d = p12*p2 + p1*p32 + p22*p3 - p32*p2 - p3*p12 - p22*p1;
408 Float_t a = (s1*p2 + p1*s3 + s2*p3 - s3*p2 - p3*s1 - s2*p1) / d;
409 Float_t b = (p12*s2 + s1*p32 + p22*s3 - p32*s2 - s3*p12 - p22*s1)/d;
410 Float_t c = (p12*p2*s3 + p1*p32*s2 + p22*p3*s1
411 - p32*p2*s1 - p3*p12*s2 - p22*p1*s3) / d;
412 Float_t sigma = a * p * p + b * p + c;
413 frac1 = sigma/sigmap;
414 }
415 index = iphi + fNbinphi * ip;
416 fSplineEff[index][1]->GetKnot(0,xmin,ymin);
417 fSplineEff[index][1]->GetKnot(9,xmax,ymax);
418 if (theta>xmax) theta = xmax;
419 Float_t frac2 = fSplineSigmap[index][1]->Eval(theta)/sigmap;
420 index = itheta + fNbintheta * ip;
421 fSplineEff[index][2]->GetKnot(0,xmin,ymin);
422 fSplineEff[index][2]->GetKnot(9,xmax,ymax);
423 if (phi>xmax) phi = xmax;
424 Float_t frac3 = fSplineSigmap[index][2]->Eval(phi)/sigmap;
425 Float_t sigmatot = sigmap * frac1 * frac2 * frac3;
426 if (sigmatot<0) sigmatot = sigmap;
427 return sigmatot;
428}
429
430Float_t AliMUONFastTracking::Sigma1P(Float_t p, Float_t theta,
041f7f97 431 Float_t phi, Int_t charge) const
6255180c 432{
041f7f97 433 //
434 // gets the width correction of the prec-pgen distribution (see FitP)
435 //
6255180c 436 Int_t ip=0, itheta=0, iphi=0;
437 GetIpIthetaIphi(p,theta,phi,charge,ip,itheta,iphi);
438 if (p>fPmax) {
439 // linear extrapolation of sigmap for p out of range
440 Float_t s1 = fCurrentEntry[fNbinp-1][itheta][iphi]->fSigma1p;
441 Float_t s2 = fCurrentEntry[fNbinp-2][itheta][iphi]->fSigma1p;
442 Float_t p1 = fDeltaP * (fNbinp - 1 + 0.5) + fPmin;
443 Float_t p2 = fDeltaP * (fNbinp - 2 + 0.5) + fPmin;
444 Float_t sigma = 1./(p1-p2) * ( (s1-s2)*p + (s2-s1)*p1 + s1*(p1-p2) );
445 return sigma;
446 }
447 else return fCurrentEntry[ip][itheta][iphi]->fSigma1p;
448}
449
450Float_t AliMUONFastTracking::NormG2(Float_t p, Float_t theta,
041f7f97 451 Float_t phi, Int_t charge) const
6255180c 452{
041f7f97 453 //
454 // gets the relative normalization of the background
455 // (gaussian) component in the prec-pgen distribution
456 //
6255180c 457 Int_t ip=0, itheta=0, iphi=0;
458 GetIpIthetaIphi(p,theta,phi,charge,ip,itheta,iphi);
459 if (p>fPmax) {
460 // linear extrapolation of sigmap for p out of range
461 Float_t s1 = fCurrentEntry[fNbinp-1][itheta][iphi]->fNormG2;
462 Float_t s2 = fCurrentEntry[fNbinp-2][itheta][iphi]->fNormG2;
463 Float_t p1 = fDeltaP * (fNbinp - 1 + 0.5) + fPmin;
464 Float_t p2 = fDeltaP * (fNbinp - 2 + 0.5) + fPmin;
465 Float_t norm = 1./(p1-p2) * ( (s1-s2)*p + (s2-s1)*p1 + s1*(p1-p2) );
466 return norm;
467 }
468 else return fCurrentEntry[ip][itheta][iphi]->fNormG2;
469}
470
471Float_t AliMUONFastTracking::MeanG2(Float_t p, Float_t theta,
041f7f97 472 Float_t phi, Int_t charge) const
6255180c 473{
041f7f97 474 //
475 // gets the mean value of the background
476 // (gaussian) component in the prec-pgen distribution
477 //
6255180c 478 Int_t ip=0, itheta=0, iphi=0;
479 GetIpIthetaIphi(p,theta,phi,charge,ip,itheta,iphi);
480 if (p>fPmax) {
481 // linear extrapolation of sigmap for p out of range
482 Float_t s1 = fCurrentEntry[fNbinp-1][itheta][iphi]->fMeanG2;
483 Float_t s2 = fCurrentEntry[fNbinp-2][itheta][iphi]->fMeanG2;
484 Float_t p1 = fDeltaP * (fNbinp - 1 + 0.5) + fPmin;
485 Float_t p2 = fDeltaP * (fNbinp - 2 + 0.5) + fPmin;
486 Float_t norm = 1./(p1-p2) * ( (s1-s2)*p + (s2-s1)*p1 + s1*(p1-p2) );
487 return norm;
488 }
489 else return fCurrentEntry[ip][itheta][iphi]->fMeanG2;
490}
491
492Float_t AliMUONFastTracking::SigmaG2(Float_t p, Float_t theta,
041f7f97 493 Float_t phi, Int_t charge) const
6255180c 494{
041f7f97 495 //
496 // gets the width of the background
497 // (gaussian) component in the prec-pgen distribution
498 //
6255180c 499 Int_t ip=0, itheta=0, iphi=0;
500 GetIpIthetaIphi(p,theta,phi,charge,ip,itheta,iphi);
501 if (p>fPmax) {
502 // linear extrapolation of sigmap for p out of range
503 Float_t s1 = fCurrentEntry[fNbinp-1][itheta][iphi]->fSigmaG2;
504 Float_t s2 = fCurrentEntry[fNbinp-2][itheta][iphi]->fSigmaG2;
505 Float_t p1 = fDeltaP * (fNbinp - 1 + 0.5) + fPmin;
506 Float_t p2 = fDeltaP * (fNbinp - 2 + 0.5) + fPmin;
507 Float_t sigma = 1./(p1-p2) * ( (s1-s2)*p + (s2-s1)*p1 + s1*(p1-p2) );
508 return sigma;
509 }
510 else return fCurrentEntry[ip][itheta][iphi]->fSigmaG2;
511}
512
513
514Float_t AliMUONFastTracking::MeanTheta(Float_t p, Float_t theta,
041f7f97 515 Float_t phi, Int_t charge) const
6255180c 516{
041f7f97 517 //
518 // gets the mean value of the thetarec-thetagen distribution
519 //
6255180c 520 Int_t ip=0, itheta=0, iphi=0;
521 GetIpIthetaIphi(p,theta,phi,charge,ip,itheta,iphi);
522 return fCurrentEntry[ip][itheta][iphi]->fMeantheta;
523}
524
041f7f97 525Float_t AliMUONFastTracking::SigmaTheta(Float_t p, Float_t theta,
526 Float_t phi, Int_t charge) const
527{
528 //
529 // gets the width of the thetarec-thetagen distribution
530 //
6255180c 531 Int_t ip=0, itheta=0, iphi=0;
532 Int_t index;
533 GetIpIthetaIphi(p,theta,phi,charge,ip,itheta,iphi);
534 // central value and corrections with spline
535 Float_t sigmatheta = fCurrentEntry[ip][itheta][iphi]->fSigmatheta;
536 if (!fSpline) return sigmatheta;
537 // corrections vs p, theta, phi
538 index = iphi + fNbinphi * itheta;
539 Double_t xmin,ymin,xmax,ymax;
540 Float_t frac1 = fSplineSigmatheta[index][0]->Eval(p)/sigmatheta;
541 if (p>fPmax-fDeltaP/2.) {
542 // linear extrapolation of sigmap for p out of range
543 Float_t s1 = fCurrentEntry[fNbinp-1][itheta][iphi]->fSigmatheta;
544 Float_t s2 = fCurrentEntry[fNbinp-2][itheta][iphi]->fSigmatheta;
545 Float_t p1 = fDeltaP * (fNbinp - 1 + 0.5) + fPmin;
546 Float_t p2 = fDeltaP * (fNbinp - 2 + 0.5) + fPmin;
547 Float_t sigma = 1./(p1-p2) * ( (s1-s2)*p + (s2-s1)*p1 + s1*(p1-p2) );
548 frac1=sigma/sigmatheta;
549 }
550 index = iphi + fNbinphi * ip;
551 fSplineEff[index][1]->GetKnot(0,xmin,ymin);
552 fSplineEff[index][1]->GetKnot(9,xmax,ymax);
553 if (theta>xmax) theta = xmax;
554 Float_t frac2 = fSplineSigmatheta[index][1]->Eval(theta)/sigmatheta;
555 index = itheta + fNbintheta * ip;
556 fSplineEff[index][2]->GetKnot(0,xmin,ymin);
557 fSplineEff[index][2]->GetKnot(9,xmax,ymax);
558 if (phi>xmax) phi = xmax;
559 Float_t frac3 = fSplineSigmatheta[index][2]->Eval(phi)/sigmatheta;
560 return sigmatheta * frac1 * frac2 * frac3;
561}
562
563
564Float_t AliMUONFastTracking::MeanPhi(Float_t p, Float_t theta,
041f7f97 565 Float_t phi, Int_t charge) const
566{
567 //
568 // gets the mean value of the phirec-phigen distribution
569 //
6255180c 570 Int_t ip=0, itheta=0, iphi=0;
571 GetIpIthetaIphi(p,theta,phi,charge,ip,itheta,iphi);
572 return fCurrentEntry[ip][itheta][iphi]->fMeanphi;
573}
574
575Float_t AliMUONFastTracking::SigmaPhi(Float_t p, Float_t theta,
576 Float_t phi, Int_t charge){
041f7f97 577 //
578 // gets the width of the phirec-phigen distribution
579 //
6255180c 580 Int_t ip=0, itheta=0, iphi=0;
581 Int_t index;
582 GetIpIthetaIphi(p,theta,phi,charge,ip,itheta,iphi);
583 // central value and corrections with spline
584 Float_t sigmaphi = fCurrentEntry[ip][itheta][iphi]->fSigmaphi;
585 if (!fSpline) return sigmaphi;
586 // corrections vs p, theta, phi
587 index = iphi + fNbinphi * itheta;
588 Float_t frac1 = fSplineSigmaphi[index][0]->Eval(p)/sigmaphi;
589 Double_t xmin,ymin,xmax,ymax;
590 if (p>fPmax-fDeltaP/2.) {
591 Float_t s1 = fCurrentEntry[fNbinp-1][itheta][iphi]->fSigmaphi;
592 Float_t s2 = fCurrentEntry[fNbinp-2][itheta][iphi]->fSigmaphi;
593 Float_t p1 = fDeltaP * (fNbinp - 1 + 0.5) + fPmin;
594 Float_t p2 = fDeltaP * (fNbinp - 2 + 0.5) + fPmin;
595 Float_t sigma = 1./(p1-p2) * ( (s1-s2)*p + (s2-s1)*p1 + s1*(p1-p2) );
596 frac1 = sigma/sigmaphi;
597 }
598
599 index = iphi + fNbinphi * ip;
600 fSplineEff[index][1]->GetKnot(0,xmin,ymin);
601 fSplineEff[index][1]->GetKnot(9,xmax,ymax);
602 if (theta>xmax) theta = xmax;
603 Float_t frac2 = fSplineSigmaphi[index][1]->Eval(theta)/sigmaphi;
604 index = itheta + fNbintheta * ip;
605 fSplineEff[index][2]->GetKnot(0,xmin,ymin);
606 fSplineEff[index][2]->GetKnot(9,xmax,ymax);
607 if (phi>xmax) phi = xmax;
608 Float_t frac3 = fSplineSigmaphi[index][2]->Eval(phi)/sigmaphi;
609 return sigmaphi * frac1 * frac2 * frac3;
610}
611
612void AliMUONFastTracking::SetSpline(){
041f7f97 613 //
614 // sets the spline functions for a smooth behaviour of the parameters
615 // when going from one cell to another
616 //
6255180c 617 printf ("Setting spline functions...");
618 char splname[40];
619 Double_t x[20][3];
620 Double_t x2[50][3];
621 Int_t nbins[3] = {fNbinp, fNbintheta, fNbinphi};
622 Double_t xspl[20],yeff[50],ysigmap[20],ysigma1p[20];
623 Double_t yacc[50], ysigmatheta[20],ysigmaphi[20];
624 Double_t xsp2[50];
625 // let's calculate the x axis for p, theta, phi
626
627 Int_t i, ispline, ivar;
628 for (i=0; i< fNbinp; i++) x[i][0] = fPmin + fDeltaP * (i + 0.5);
629 for (i=0; i< fNbintheta; i++) x[i][1] = fThetamin + fDeltaTheta * (i + 0.5);
630 for (i=0; i< fNbinphi; i++) x[i][2] = fPhimin + fDeltaPhi * (i + 0.5);
631
632 for (i=0; i< 5 * fNbinp; i++) x2[i][0] = fPmin + fDeltaP * (i + 0.5)/5.;
633 for (i=0; i< 5 * fNbintheta; i++) x2[i][1] = fThetamin + fDeltaTheta * (i + 0.5)/5.;
634 for (i=0; i< 5 * fNbinphi; i++) x2[i][2] = fPhimin + fDeltaPhi * (i + 0.5)/5.;
635
636 // splines in p
637 ivar = 0;
638 for (i=0; i<nbins[ivar]; i++) xspl[i] = x[i][ivar];
639 for (i=0; i<5 * nbins[ivar]; i++) xsp2[i] = x2[i][ivar];
640 ispline=0;
641 for (Int_t itheta=0; itheta< fNbintheta; itheta++){
642 for (Int_t iphi=0; iphi< fNbinphi; iphi++){
643 for (Int_t ip=0; ip<fNbinp; ip++) {
644 // for (Int_t i=0; i<5; i++) {
645 // yeff[5 * ip + i] = fCurrentEntry[ip][itheta][iphi]->fEff[i];
646 // yacc[5 * ip + i] = fCurrentEntry[ip][itheta][iphi]->fAcc[i];
647 // }
648 ysigmap[ip] = fCurrentEntry[ip][itheta][iphi]->fSigmap;
649 ysigma1p[ip] = fCurrentEntry[ip][itheta][iphi]->fSigma1p;
650 ysigmatheta[ip] = fCurrentEntry[ip][itheta][iphi]->fSigmatheta;
651 ysigmaphi[ip] = fCurrentEntry[ip][itheta][iphi]->fSigmaphi;
652 }
653 if (fPrintLevel>3) cout << " creating new spline " << splname << endl;
654 sprintf (splname,"fSplineEff[%d][%d]",ispline,ivar);
655 fSplineEff[ispline][ivar] = new TSpline3(splname,xsp2,yeff,5 * nbins[ivar]);
656 sprintf (splname,"fSplineAcc[%d][%d]",ispline,ivar);
657 fSplineAcc[ispline][ivar] = new TSpline3(splname,xsp2,yacc,5 * nbins[ivar]);
658 sprintf (splname,"fSplineSigmap[%d][%d]",ispline,ivar);
659 fSplineSigmap[ispline][ivar] = new TSpline3(splname,xspl,ysigmap,nbins[ivar]);
660 sprintf (splname,"fSplineSigma1p[%d][%d]",ispline,ivar);
661 fSplineSigma1p[ispline][ivar] = new TSpline3(splname,xspl,ysigma1p,nbins[ivar]);
662 sprintf (splname,"fSplineSigmatheta[%d][%d]",ispline,ivar);
663 fSplineSigmatheta[ispline][ivar] = new TSpline3(splname,xspl,ysigmatheta,nbins[ivar]);
664 sprintf (splname,"fSplineSigmaphi[%d][%d]",ispline,ivar);
665 fSplineSigmaphi[ispline][ivar] = new TSpline3(splname,xspl,ysigmaphi,nbins[ivar]);
666 ispline++;
667 }
668 }
669
670 ivar = 1;
671 for (i=0; i<nbins[ivar]; i++) xspl[i] = x[i][ivar];
672 ispline=0;
673 for (Int_t ip=0; ip<fNbinp; ip++) {
674 for (Int_t iphi=0; iphi< fNbinphi; iphi++){
675 for (Int_t itheta=0; itheta< fNbintheta; itheta++){
676 // for efficiency and acceptance let's take the central value
677 // yeff[itheta] = fCurrentEntry[ip][itheta][iphi]->fEff[2];
678 // yacc[itheta] = fCurrentEntry[ip][itheta][iphi]->fAcc[2];
679 ysigmap[itheta] = fCurrentEntry[ip][itheta][iphi]->fSigmap;
680 ysigma1p[itheta] = fCurrentEntry[ip][itheta][iphi]->fSigma1p;
681 ysigmatheta[itheta] = fCurrentEntry[ip][itheta][iphi]->fSigmatheta;
682 ysigmaphi[itheta] = fCurrentEntry[ip][itheta][iphi]->fSigmaphi;
683 }
684 if (fPrintLevel>3) cout << " creating new spline " << splname << endl;
685 sprintf (splname,"fSplineEff[%d][%d]",ispline,ivar);
686 fSplineEff[ispline][ivar] = new TSpline3(splname,xspl,yeff, nbins[ivar]);
687 sprintf (splname,"fSplineAcc[%d][%d]",ispline,ivar);
688 fSplineAcc[ispline][ivar] = new TSpline3(splname,xspl,yacc, nbins[ivar]);
689 sprintf (splname,"fSplineSigmap[%d][%d]",ispline,ivar);
690 fSplineSigmap[ispline][ivar] = new TSpline3(splname,xspl,ysigmap,nbins[ivar]);
691 sprintf (splname,"fSplineSigma1p[%d][%d]",ispline,ivar);
692 fSplineSigma1p[ispline][ivar] = new TSpline3(splname,xspl,ysigma1p,nbins[ivar]);
693 sprintf (splname,"fSplineSigmatheta[%d][%d]",ispline,ivar);
694 fSplineSigmatheta[ispline][ivar] = new TSpline3(splname,xspl,ysigmatheta,nbins[ivar]);
695 sprintf (splname,"fSplineSigmaphi[%d][%d]",ispline,ivar);
696 fSplineSigmaphi[ispline][ivar] = new TSpline3(splname,xspl,ysigmaphi,nbins[ivar]);
697 ispline++;
698 }
699 }
700
701 ivar = 2;
702 for (i=0; i<nbins[ivar]; i++) xspl[i] = x[i][ivar];
703 ispline=0;
704 for (Int_t ip=0; ip<fNbinp; ip++) {
705 for (Int_t itheta=0; itheta< fNbintheta; itheta++){
706 for (Int_t iphi=0; iphi< fNbinphi; iphi++){
707 // for efficiency and acceptance let's take the central value
708 // yeff[iphi] = fCurrentEntry[ip][itheta][iphi]->fEff[2];
709 // yacc[iphi] = fCurrentEntry[ip][itheta][iphi]->fAcc[2];
710 ysigmap[iphi] = fCurrentEntry[ip][itheta][iphi]->fSigmap;
711 ysigma1p[iphi] = fCurrentEntry[ip][itheta][iphi]->fSigma1p;
712 ysigmatheta[iphi] = fCurrentEntry[ip][itheta][iphi]->fSigmatheta;
713 ysigmaphi[iphi] = fCurrentEntry[ip][itheta][iphi]->fSigmaphi;
714 }
715 if (fPrintLevel>3) cout << " creating new spline " << splname << endl;
716 sprintf (splname,"fSplineEff[%d][%d]",ispline,ivar);
717 fSplineEff[ispline][ivar] = new TSpline3(splname,xspl,yeff, nbins[ivar]);
718 sprintf (splname,"fSplineAcc[%d][%d]",ispline,ivar);
719 fSplineAcc[ispline][ivar] = new TSpline3(splname,xspl,yacc, nbins[ivar]);
720 sprintf (splname,"fSplineSigmap[%d][%d]",ispline,ivar);
721 fSplineSigmap[ispline][ivar] = new TSpline3(splname,xspl,ysigmap,nbins[ivar]);
722 sprintf (splname,"fSplineSigma1p[%d][%d]",ispline,ivar);
723 fSplineSigma1p[ispline][ivar] = new TSpline3(splname,xspl,ysigma1p,nbins[ivar]);
724 sprintf (splname,"fSplineSigmatheta[%d][%d]",ispline,ivar);
725 fSplineSigmatheta[ispline][ivar] = new TSpline3(splname,xspl,ysigmatheta,nbins[ivar]);
726 sprintf (splname,"fSplineSigmaphi[%d][%d]",ispline,ivar);
727 fSplineSigmaphi[ispline][ivar] = new TSpline3(splname,xspl,ysigmaphi,nbins[ivar]);
728 ispline++;
729 }
730 }
731 printf ("...done\n");
732}
733
6255180c 734void AliMUONFastTracking::SetBackground(Float_t bkg){
041f7f97 735 //
6255180c 736 // linear interpolation of the parameters in the LUT between 2 values where
737 // the background has been actually calculated
041f7f97 738 //
6255180c 739 if (bkg>2) printf ("WARNING: unsafe extrapolation!\n");
740 fBkg = bkg;
741
041f7f97 742 Float_t bkgLevel[4] = {0, 0.5, 1, 2}; // bkg values for which LUT is calculated
6255180c 743 Int_t ibkg;
041f7f97 744 for (ibkg=0; ibkg<4; ibkg++) if ( bkg < bkgLevel[ibkg]) break;
6255180c 745 if (ibkg == 4) ibkg--;
746 if (ibkg == 0) ibkg++;
747
041f7f97 748 Float_t x0 = bkgLevel[ibkg-1];
749 Float_t x1 = bkgLevel[ibkg];
6255180c 750 Float_t x = (bkg - x0) / (x1 - x0);
751
752 Float_t y0, y1;
753
754 for (Int_t ip=0; ip< fNbinp; ip++){
755 for (Int_t itheta=0; itheta< fNbintheta; itheta++){
756 for (Int_t iphi=0; iphi< fNbinphi; iphi++){
757 fCurrentEntry[ip][itheta][iphi]->fP = fEntry[ip][itheta][iphi][ibkg]->fP;
758 fCurrentEntry[ip][itheta][iphi]->fTheta = fEntry[ip][itheta][iphi][ibkg]->fTheta;
759 fCurrentEntry[ip][itheta][iphi]->fPhi = fEntry[ip][itheta][iphi][ibkg]->fPhi;
760 fCurrentEntry[ip][itheta][iphi]->fChi2p = -1;
761 fCurrentEntry[ip][itheta][iphi]->fChi2theta = -1;
762 fCurrentEntry[ip][itheta][iphi]->fChi2phi = -1;
763
764 y0 = fEntry[ip][itheta][iphi][ibkg-1]->fMeanp;
765 y1 = fEntry[ip][itheta][iphi][ibkg]->fMeanp;
766 fCurrentEntry[ip][itheta][iphi] ->fMeanp = (y1 - y0) * x + y0;
767 y0 = fEntry[ip][itheta][iphi][ibkg-1]->fMeantheta;
768 y1 = fEntry[ip][itheta][iphi][ibkg]->fMeantheta;
769 fCurrentEntry[ip][itheta][iphi] ->fMeantheta = (y1 - y0) * x + y0;
770 y0 = fEntry[ip][itheta][iphi][ibkg-1]->fMeanphi;
771 y1 = fEntry[ip][itheta][iphi][ibkg]->fMeanphi;
772 fCurrentEntry[ip][itheta][iphi] ->fMeanphi = (y1 - y0) * x + y0;
773 y0 = fEntry[ip][itheta][iphi][ibkg-1]->fSigmap;
774 y1 = fEntry[ip][itheta][iphi][ibkg]->fSigmap;
775 fCurrentEntry[ip][itheta][iphi] ->fSigmap = (y1 - y0) * x + y0;
776 y0 = fEntry[ip][itheta][iphi][ibkg-1]->fSigmatheta;
777 y1 = fEntry[ip][itheta][iphi][ibkg]->fSigmatheta;
778 fCurrentEntry[ip][itheta][iphi] ->fSigmatheta = (y1 - y0) * x + y0;
779 y0 = fEntry[ip][itheta][iphi][ibkg-1]->fSigmaphi;
780 y1 = fEntry[ip][itheta][iphi][ibkg]->fSigmaphi;
781 fCurrentEntry[ip][itheta][iphi] ->fSigmaphi = (y1 - y0) * x + y0;
782 y0 = fEntry[ip][itheta][iphi][ibkg-1]->fSigma1p;
783 y1 = fEntry[ip][itheta][iphi][ibkg]->fSigma1p;
784 fCurrentEntry[ip][itheta][iphi] ->fSigma1p = (y1 - y0) * x + y0;
785 y0 = fEntry[ip][itheta][iphi][ibkg-1]->fNormG2;
786 y1 = fEntry[ip][itheta][iphi][ibkg]->fNormG2;
787 fCurrentEntry[ip][itheta][iphi] ->fNormG2 = (y1 - y0) * x + y0;
788 y0 = fEntry[ip][itheta][iphi][ibkg-1]->fMeanG2;
789 y1 = fEntry[ip][itheta][iphi][ibkg]->fMeanG2;
790 fCurrentEntry[ip][itheta][iphi] ->fMeanG2 = (y1 - y0) * x + y0;
791
792 y0 = fEntry[ip][itheta][iphi][ibkg-1]->fSigmaG2;
793 y1 = fEntry[ip][itheta][iphi][ibkg]->fSigmaG2;
794 fCurrentEntry[ip][itheta][iphi] ->fSigmaG2 = (y1 - y0) * x + y0;
795 for (Int_t i=0; i<kSplitP; i++) {
796 for (Int_t j=0; j<kSplitTheta; j++) {
797 fCurrentEntry[ip][itheta][iphi]->fAcc[i][j] = fEntry[ip][itheta][iphi][ibkg]->fAcc[i][j];
798 y0 = fEntry[ip][itheta][iphi][ibkg-1]->fEff[i][j];
799 y1 = fEntry[ip][itheta][iphi][ibkg]->fEff[i][j];
800 fCurrentEntry[ip][itheta][iphi]->fEff[i][j] = (y1 - y0) * x + y0;
801 }
802 }
803 }
804 }
805 }
806 SetSpline();
807}
808
d53fc881 809TF1* AliMUONFastTracking::GetFitP(Int_t ip,Int_t itheta,Int_t iphi) {
041f7f97 810 // gets the correct prec-pgen distribution for a given LUT cell
c6ba3ff4 811 if (!fFitp[ip][itheta][iphi]) {
812 char name[256];
813 sprintf(name, "fit_%d_%d_%d", ip, itheta, iphi);
814 fFitp[ip][itheta][iphi] = new TF1(name ,FitP,-20.,20.,6);
d53fc881 815 fFitp[ip][itheta][iphi]->SetNpx(500);
1e68a0c5 816 fFitp[ip][itheta][iphi]->SetParameters(0.,0.,0.,0.,0.,0.);
d53fc881 817 }
818 return fFitp[ip][itheta][iphi];
819}
6255180c 820
041f7f97 821AliMUONFastTracking& AliMUONFastTracking::operator=(const AliMUONFastTracking& rhs)
822{
823// Assignment operator
824 rhs.Copy(*this);
825 return *this;
826}
827
15a060e4 828void AliMUONFastTracking::Copy(TObject&) const
041f7f97 829{
830 //
831 // Copy
832 //
833 Fatal("Copy","Not implemented!\n");
834}
6255180c 835
836
837
838
839
840
841
842