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