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