]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/Cal/AliTRDCalTrkAttach.cxx
88dd3610abd79ee97c3c582d8d417ab22ae6f954
[u/mrichter/AliRoot.git] / TRD / Cal / AliTRDCalTrkAttach.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 //                                                                      //
18 // Container for calibration parameters for AliTRDseedV1::AttachClusters()
19 // .... Longer description of content ...
20 //
21 // For calibration procedure check AliTRDtrackleOflHelper::CalibrateAttach()
22 // .... some reference to the calibration procedure .........
23 //                                                                      //
24 // Authors:                                                             //
25 //   Alex Bercuci <a.bercuci@gsi.de>                                    //
26 //                                                                      //
27 //////////////////////////////////////////////////////////////////////////
28
29 #include <TFile.h>
30 #include <TROOT.h>
31 #include <TMath.h>
32 #include "TGraphErrors.h"
33 #include "TObjArray.h"
34
35 #include "AliLog.h"
36
37 #include "AliTRDCalTrkAttach.h"
38 #include "AliTRDcalibDB.h"
39
40 ClassImp(AliTRDCalTrkAttach)
41
42 //______________________________________________________________________________
43 AliTRDCalTrkAttach::AliTRDCalTrkAttach()
44   :TNamed("AliTRDCalTrkAttach", "Calibration of AliTRDseedV1::AttachClusters")
45   ,fRClikeLimit(0.65)
46   ,fScaleCov(2.)
47   ,fLike(NULL)
48 {
49 // Default constructor
50
51   fNsgmDy[0] = 5; fNsgmDy[1] = 7;
52   fLikeMinRelDecrease[0] = .2; fLikeMinRelDecrease[1] = .3;
53 }
54
55 //______________________________________________________________________________
56 AliTRDCalTrkAttach::~AliTRDCalTrkAttach()
57 {
58 // Destructor
59   if(fLike) delete fLike;
60 }
61
62
63 //______________________________________________________________________________
64 Double_t AliTRDCalTrkAttach::CookLikelihood(Bool_t chg, Int_t ly, Float_t pt, Float_t phiTrk, Int_t n, Double_t dyr, Double_t dphi, Double_t sr) const
65 {
66 // Calculate likelihood for a segment to belong to a tracklet
67 // Based on calibrated values
68
69   if(n<4){
70     AliDebug(4, Form("Failed basic cut[s] : n[%d] ...", n));
71     return 0.;
72   }
73   //check likelihood array  
74   if (!fLike || fLike->GetEntries() != 6*2*kNcalib) {
75     AliError("No usable AttachClusters calib object.");
76     return 0.;
77   }
78   
79   TGraphErrors *g(NULL);
80   if(!(g = (TGraphErrors*)fLike->At(ly*8+2*Int_t(kResPos)+Int_t(chg)))){
81     AliError("Failed retrieving AttachClusters graph.");
82     return 0.;
83   }
84   // Interpolate p_t
85   Int_t npts(g->GetN()), ip(0), jp(-1);
86   Double_t x0, y0, x1, y1, dd(0.), invdx(0.), f[4]={0., 0., 0., 0.};
87   for(Int_t kp(0); kp<npts; kp++){
88     g->GetPoint(kp, x1, y1);
89     if(x1>=pt){jp=kp; break;}
90   }
91   Bool_t boundary(kFALSE);
92   if(jp<0){
93     jp = npts-1; 
94     g->GetPoint(jp, x1, y1);
95     ip = npts-1;
96     boundary = kTRUE;
97   }else if(jp==0){ 
98     ip = 0;
99     boundary = kTRUE;
100   }else{ 
101     ip = jp-1;
102     g->GetPoint(ip, x0, y0);
103     invdx = 1./(x0-x1);
104   }
105   // process pt dependences
106   // +++ process dy
107   Double_t dym = boundary?y1:((pt*(y0-y1) + (x0*y1-x1*y0))*invdx),
108            sym = 0.5*(g->GetErrorY(ip)+g->GetErrorY(jp));
109   dd      = (dyr - dym)/sym; dd*=dd;
110   f[0] = TMath::Exp(-0.5*dd);
111   // +++ process dphi
112   if(!(g = (TGraphErrors*)fLike->At(ly*8+2*Int_t(kResAng)+Int_t(chg)))){
113     AliError("Failed retrieving AttachClusters graph.");
114     return 0.;
115   }
116   g->GetPoint(ip, x0, y0);g->GetPoint(jp, x1, y1);
117   Double_t dpm = boundary?y1:((pt*(y0-y1) + (x0*y1-x1*y0))*invdx),
118            spm = 0.5*(g->GetErrorY(ip)+g->GetErrorY(jp));
119   dd      = (dphi - dpm)/spm; dd*=dd;
120   f[1] = TMath::Exp(-0.5*dd);
121   // +++ process no of clusters
122   if(!(g = (TGraphErrors*)fLike->At(ly*8+2*Int_t(kNclMean)+Int_t(chg)))){
123     AliError("Failed retrieving AttachClusters graph.");
124     return 0.;
125   }
126   g->GetPoint(ip, x0, y0);g->GetPoint(jp, x1, y1);
127   Double_t nm = boundary?y1:((pt*(y0-y1) + (x0*y1-x1*y0))*invdx);
128   f[2] = (nm-TMath::Abs(n-nm))/nm;
129  
130   // process phi dependences
131   // +++ process <s>/s
132   if(!(g = (TGraphErrors*)fLike->At(ly*8+2*Int_t(kSigma)+Int_t(chg)))){
133     AliError("Failed retrieving AttachClusters graph.");
134     return 0.;
135   }
136   // Interpolate phi [deg]
137   npts=g->GetN(); jp=-1;
138   for(Int_t kp(0); kp<npts; kp++){
139     g->GetPoint(kp, x1, y1);
140     if(x1>=phiTrk){jp=kp; break;}
141   }
142   if(jp<0){
143     jp = npts-1; 
144     g->GetPoint(jp, x1, y1);
145     ip = jp;
146     boundary = kTRUE;
147   }else if(jp==0){ 
148     ip = jp;
149     boundary = kTRUE;
150   }else{ 
151     ip = jp-1;
152     g->GetPoint(ip, x0, y0);
153     invdx = 1./(x0-x1);
154     boundary = kFALSE;
155   }
156   Double_t sm = boundary?y1:((phiTrk*(y0-y1) + (x0*y1-x1*y0))*invdx),
157            ssm = 0.5*(g->GetErrorY(ip)+g->GetErrorY(jp));
158   dd      = (sr - sm)/ssm; dd*=dd;
159   f[3] = TMath::Exp(-0.5*dd);
160
161   // Calculate likelihood
162   Double_t length = f[0]*f[0]+f[1]*f[1]+f[2]*f[2]+f[3]*f[3];
163   length = TMath::Sqrt(length);
164   Double_t cosTht = f[0]+f[1]+f[2]+f[3];
165   cosTht /= (4.*length);
166   AliDebug(2, Form("Like[%5.3f] ThtLike[%6.2f](deg)\n"
167     "    F_dy (%+6.2f)=%4.2f\n" 
168     "    F_phi(%+6.2f)=%4.2f\n"
169     "    F_ncl(%+6d)=%4.2f\n"
170     "    F_<s>(%+6.2f)=%4.2f", 
171     length, TMath::ACos(cosTht)*TMath::RadToDeg(),
172     dyr, f[0], dphi, f[1], n, f[2], sr, f[3]));
173
174   return length;
175 }
176
177 //______________________________________________________________________________
178 Bool_t AliTRDCalTrkAttach::LoadReferences(const Char_t *file)
179 {
180 // Load calibration data from file
181
182   if(!file || !TFile::Open(file)){
183     AliError("Parametrization file missing or unreadable.");
184     return kFALSE;
185   }
186   TGraphErrors *g(NULL);
187   Char_t co[kNcalib] = {'y', 'p', 's', 'n'},
188          cs[2] = {'n', 'p'};
189
190   if(fLike) fLike->Clear();
191   else fLike = new TObjArray(6*2*kNcalib);
192   for(Int_t ily(0); ily<6; ily++){
193     for(Int_t icalib(0); icalib<kNcalib; icalib++){
194       for(Int_t isgn(0); isgn<2; isgn++){
195         if(!(g = (TGraphErrors*)gFile->Get(Form("%c%c%d", co[icalib], cs[isgn], ily)))) return kFALSE;
196         fLike->AddAt(g, ily*8+2*icalib+isgn);
197       }
198     }
199   }
200   return kTRUE;
201 }