]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/TRD/AliTRDcheckTRK.cxx
full resolution monitoring for clusters, tracklets and trackIN
[u/mrichter/AliRoot.git] / PWG1 / TRD / AliTRDcheckTRK.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-commercialf 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 //                                                                        //
19 //  TRD tracker systematic                                                //
20 //
21 //
22 //  Authors:                                                              //
23 //    Alexandru Bercuci <A.Bercuci@gsi.de>                                //
24 //                                                                        //
25 ////////////////////////////////////////////////////////////////////////////
26
27 #include "TROOT.h"
28 #include "TAxis.h"
29 #include "TH1.h"
30 #include "TH2.h"
31 #include "TH3.h"
32 #include "TObjArray.h"
33 #include "THnSparse.h"
34 #include <TVectorT.h>
35
36 #include "AliLog.h"
37
38 #include "AliTRDcluster.h"
39 #include "AliTRDseedV1.h"
40 #include "AliTRDtrackV1.h"
41 #include "AliTRDtrackerV1.h"
42 #include "AliTRDtransform.h"
43
44 #include "AliTRDcheckTRK.h"
45
46 ClassImp(AliTRDcheckTRK)
47
48 Bool_t  AliTRDcheckTRK::fgKalmanUpdate = kTRUE;
49 Bool_t  AliTRDcheckTRK::fgClRecalibrate = kFALSE;
50 Float_t AliTRDcheckTRK::fgKalmanStep = 2.;
51 Float_t AliTRDcheckTRK::fgCalib[540][2];
52 //__________________________________________________________________________
53 AliTRDcheckTRK::AliTRDcheckTRK()
54   : AliTRDresolution()
55 {
56 // Default constructor
57   SetNameTitle("TRDtrackerSys", "TRD Tracker Systematic");
58   memset(AliTRDcheckTRK::fgCalib, 0, 540*2*sizeof(Float_t));
59 }
60
61 //__________________________________________________________________________
62 AliTRDcheckTRK::AliTRDcheckTRK(char* name)
63   : AliTRDresolution(name, kFALSE)
64 {
65 // User constructor
66   SetTitle("TRD Tracker Systematic");
67   memset(AliTRDcheckTRK::fgCalib, 0, 540*2*sizeof(Float_t));
68   InitFunctorList();
69 }
70
71 //__________________________________________________________________________
72 AliTRDcheckTRK::~AliTRDcheckTRK()
73 {
74 // Destructor
75 }
76
77 //__________________________________________________________________________
78 Int_t AliTRDcheckTRK::GetSpeciesByMass(Float_t m)
79 {
80 // Find particle index by mass
81 // 0 electron
82 // 1 muon
83 // 2 pion
84 // 3 kaon
85 // 4 proton
86
87   for(Int_t is(0); is<AliPID::kSPECIES; is++) if(TMath::Abs(m-AliPID::ParticleMass(is))<1.e-4) return is;
88   return -1;
89 }
90
91
92
93 //__________________________________________________________________________
94 TH1* AliTRDcheckTRK::PlotTrack(const AliTRDtrackV1 *track)
95 {
96 // comment needed
97
98   if(track) fkTrack = track;
99   if(!fkTrack){
100     AliDebug(4, "No Track defined.");
101     return NULL;
102   }
103   // make a local copy of current track
104   AliTRDtrackV1 lt(*fkTrack);
105   if(!PropagateKalman(lt)) return NULL;
106   PlotCluster(&lt);
107   PlotTracklet(&lt);
108   PlotTrackIn(&lt);
109   return NULL;
110 }
111
112
113 //___________________________________________________
114 Bool_t AliTRDcheckTRK::PropagateKalman(AliTRDtrackV1 &t)
115 {
116 // Propagate Back Kalman from the TPC input parameter to the last tracklet attached to track.
117 // On the propagation recalibration of clusters, tracklet refit and material budget are recalculated (on demand)
118 // On output the track is updated with the new info
119 //
120 // A.Bercuci@gsi.de
121
122 //  printf("PropagateKalman()\n");
123   Int_t ntracklets(t.GetNumberOfTracklets());
124   if(!ntracklets){
125     printf("E - AliTRDcheckTRK::PropagateKalman :: No tracklets attached to track.\n");
126     return kFALSE;
127   }
128
129   AliTRDseedV1 *tr(NULL);
130   AliExternalTrackParam *ref(NULL);
131   if(!(ref = t.GetTrackIn())){
132     printf("E - AliTRDcheckTRK::PropagateKalman :: Track did not entered TRD fiducial volume.\n");
133     return NULL;
134   }
135   if(ref->Pt()<1.e-3){printf("small refpt\n"); return kFALSE;}
136
137
138   // Initialize TRD track to the reference
139   AliTRDtrackV1 tt;
140   tt.Set(ref->GetX(), ref->GetAlpha(), ref->GetParameter(), ref->GetCovariance());
141   tt.SetMass(t.GetMass());
142   tt.SetTrackIn();tt.SetTrackOut(t.GetTrackOut());
143
144   for(Int_t ily(0); ily<AliTRDgeometry::kNlayer; ily++){
145     if(!(tr = t.GetTracklet(ily))) continue;
146     Int_t det(tr->GetDetector());
147     Float_t *calib = GetCalib(det);
148     if(fgClRecalibrate && calib[0]>0.){
149       AliTRDtransform trans(det);
150       AliTRDcluster *c(NULL);
151       Float_t exb, vd, t0, s2, dl, dt; tr->GetCalibParam(exb, vd, t0, s2, dl, dt);
152       tr->ResetClusterIter(kFALSE);
153       while((c = tr->PrevCluster())){
154         if(!trans.Transform(c/*, GetCalib(det)*/)){
155           printf("W - AliTRDcheckTRK::PropagateKalman :: Transform() failed for Det[%03d]\n", det);
156           break;
157         }
158       }
159       if(!tr->FitRobust()) printf("W - AliTRDcheckTRK::PropagateKalman :: FitRobust() failed for Det[%03d]\n", det);
160     }
161     if(!AliTRDtrackerV1::PropagateToX(tt, tr->GetX0(), fgKalmanStep)) continue;
162     tr->Update(&tt);
163     if(fgKalmanUpdate){
164       Double_t x(tr->GetX0()),
165                p[2] = { tr->GetYfit(0), tr->GetZfit(0)},
166                covTrklt[3];
167       tr->GetCovAt(x, covTrklt);
168       if(!((AliExternalTrackParam&)tt).Update(p, covTrklt)) continue;
169       //tr->Update(&tt);
170       tt.SetTracklet(tr, 0);
171       tt.SetNumberOfClusters();
172       tt.UpdateChi2(((AliExternalTrackParam)tt).GetPredictedChi2(p, covTrklt));
173     }
174   }
175   //tt.Print("a");
176   t.~AliTRDtrackV1();
177   new(&t) AliTRDtrackV1(tt);
178   return kTRUE;
179 }
180
181