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