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