]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGPP/TRD/AliTRDcheckTRK.cxx
add protection for backward data compatibility
[u/mrichter/AliRoot.git] / PWGPP / TRD / AliTRDcheckTRK.cxx
CommitLineData
3ceb45ae 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"
3ceb45ae 34
35#include "AliLog.h"
36
37#include "AliTRDcluster.h"
38#include "AliTRDseedV1.h"
566c3d46 39#include "AliTRDtrackletOflHelper.h"
3ceb45ae 40#include "AliTRDtrackV1.h"
41#include "AliTRDtrackerV1.h"
42#include "AliTRDtransform.h"
43
44#include "AliTRDcheckTRK.h"
45
46ClassImp(AliTRDcheckTRK)
47
eb05d549 48UChar_t AliTRDcheckTRK::fgSteer= 0;
3ceb45ae 49Float_t AliTRDcheckTRK::fgKalmanStep = 2.;
50//__________________________________________________________________________
51AliTRDcheckTRK::AliTRDcheckTRK()
3ed01fbe 52 : AliTRDresolution()
3ceb45ae 53{
54// Default constructor
55 SetNameTitle("TRDtrackerSys", "TRD Tracker Systematic");
3ceb45ae 56}
57
58//__________________________________________________________________________
59AliTRDcheckTRK::AliTRDcheckTRK(char* name)
3ed01fbe 60 : AliTRDresolution(name, kFALSE)
3ceb45ae 61{
62// User constructor
3ed01fbe 63 SetTitle("TRD Tracker Systematic");
566c3d46 64 MakePtCalib();
3ceb45ae 65 InitFunctorList();
3ceb45ae 66}
67
68//__________________________________________________________________________
69AliTRDcheckTRK::~AliTRDcheckTRK()
70{
71// Destructor
72}
73
566c3d46 74
3ceb45ae 75//__________________________________________________________________________
566c3d46 76TObjArray* AliTRDcheckTRK::Histos()
3ceb45ae 77{
566c3d46 78// Build extra calibration plots
79 if(!(fContainer = AliTRDresolution::Histos())) return NULL;
80 //fContainer->Expand(AliTRDresolution::kNclasses+1);
81
82 THnSparse *H(NULL);
b9058d72 83 if(!(H = (THnSparseI*)gROOT->FindObject("hRoads"))){
566c3d46 84 const Char_t *title[kNdim] = {"layer", "charge", fgkTitle[kPt], fgkTitle[kYrez], fgkTitle[kPrez], "#sigma^{*}/<#sigma_{y}> [a.u.]", "n_{cl}"};
85 const Int_t nbins[kNdim] = {AliTRDgeometry::kNlayer, 2, kNptBins, fgkNbins[kYrez], fgkNbins[kPrez], kNSigmaBins, kNclusters};
86 const Double_t min[kNdim] = {-0.5, -0.5, -0.5, -1., -5., 0., 8.5},
87 max[kNdim] = {AliTRDgeometry::kNlayer-0.5, 1.5, kNptBins-0.5, 1., 5., 5., min[6]+kNclusters};
88 TString st("Tracking Roads Calib;");
89 // define minimum info to be saved in non debug mode
90 for(Int_t idim(0); idim<kNdim; idim++){ st += title[idim]; st+=";";}
b9058d72 91 H = new THnSparseI("hRoads", st.Data(), kNdim, nbins, min, max);
566c3d46 92 } else H->Reset();
93 fContainer->AddAt(H, fContainer->GetEntries()/*AliTRDresolution::kNclasses*/);
94 return fContainer;
3ceb45ae 95}
96
3ceb45ae 97//__________________________________________________________________________
3ed01fbe 98TH1* AliTRDcheckTRK::PlotTrack(const AliTRDtrackV1 *track)
3ceb45ae 99{
100// comment needed
101
102 if(track) fkTrack = track;
103 if(!fkTrack){
104 AliDebug(4, "No Track defined.");
105 return NULL;
106 }
3ed01fbe 107 // make a local copy of current track
108 AliTRDtrackV1 lt(*fkTrack);
b9058d72 109 if(!PropagateKalman(lt, UseITS()?fkESD->GetITSoutParam():fkESD->GetTPCoutParam())) return NULL;
3ed01fbe 110 PlotCluster(&lt);
111 PlotTracklet(&lt);
112 PlotTrackIn(&lt);
566c3d46 113 DoRoads(&lt);
3ceb45ae 114 return NULL;
115}
116
566c3d46 117//________________________________________________________
118void AliTRDcheckTRK::MakePtCalib(Float_t pt0, Float_t dpt)
119{
120// Build pt segments
121 for(Int_t j(0); j<=kNptBins; j++){
122 pt0+=(TMath::Exp(j*j*dpt)-1.);
123 fPtBinCalib[j]=pt0;
124 }
125}
126
127//__________________________________________________________________________
128Int_t AliTRDcheckTRK::GetPtBinCalib(Float_t pt)
129{
130// Find pt bin according to local pt segmentation
131 Int_t ipt(-1);
132 while(ipt<kNptBins){
133 if(pt<fPtBinCalib[ipt+1]) break;
134 ipt++;
135 }
136 return ipt;
137}
138
139
140//__________________________________________________________________________
141TH1* AliTRDcheckTRK::DoRoads(const AliTRDtrackV1 *track)
142{
143// comment needed
144 if(track) fkTrack = track;
145 if(!fkTrack){
146 AliDebug(4, "No Track defined.");
147 return NULL;
148 }
149 if(TMath::Abs(fkESD->GetTOFbc())>1){
150 AliDebug(4, Form("Track with BC_index[%d] not used.", fkESD->GetTOFbc()));
151 return NULL;
152 }
153 THnSparse *H(NULL);
154 if(!fContainer || !(H = (THnSparse*)fContainer->At(3))){
155 AliWarning("No output container defined.");
156 return NULL;
157 }
158// return NULL;
159 Double_t val[kNdim];
160 AliTRDseedV1 *fTracklet(NULL);
161 AliTRDtrackletOflHelper helper; TObjArray cl(AliTRDseedV1::kNclusters); cl.SetOwner(kFALSE);
162 for(Int_t il(0); il<AliTRDgeometry::kNlayer; il++){
163 if(!(fTracklet = fkTrack->GetTracklet(il))) continue;
164 if(!fTracklet->IsOK() || !fTracklet->IsChmbGood()) continue;
165 Int_t ipt(GetPtBinCalib(fTracklet->GetPt()));
166 if(ipt<0) continue;
167 val[0] = il;
168 val[1] = fkTrack->Charge()<0?0:1;
169 val[2] = ipt;
170 Double_t dyt(fTracklet->GetYfit(0) - fTracklet->GetYref(0)),
171 dzt(fTracklet->GetZfit(0) - fTracklet->GetZref(0)),
172 dydx(fTracklet->GetYfit(1)),
173 tilt(fTracklet->GetTilt());
174 // correct for tilt rotation
175 val[3] = dyt - dzt*tilt;
176 dydx+= tilt*fTracklet->GetZref(1);
177 val[4] = TMath::ATan((dydx - fTracklet->GetYref(1))/(1.+ fTracklet->GetYref(1)*dydx)) * TMath::RadToDeg();
178 fTracklet->ResetClusterIter(kTRUE); AliTRDcluster *c(NULL);
179 while((c = fTracklet->NextCluster())) cl.AddLast(c);
180 helper.Init(AliTRDtransform::Geometry().GetPadPlane(fTracklet->GetDetector()), &cl);
181 Double_t r, y, s; helper.GetRMS(r, y, s, fTracklet->GetX0());
182 val[5] = s/helper.GetSyMean();
183 val[6] = fTracklet->GetN();
184 H->Fill(val);
185 }
186 return NULL;
187}
3ceb45ae 188
189//___________________________________________________
35983729 190Bool_t AliTRDcheckTRK::PropagateKalman(AliTRDtrackV1 &t, AliExternalTrackParam *ref)
3ceb45ae 191{
3ed01fbe 192// Propagate Back Kalman from the TPC input parameter to the last tracklet attached to track.
193// On the propagation recalibration of clusters, tracklet refit and material budget are recalculated (on demand)
194// On output the track is updated with the new info
3ceb45ae 195//
3ed01fbe 196// A.Bercuci@gsi.de
3ceb45ae 197
3ed01fbe 198// printf("PropagateKalman()\n");
199 Int_t ntracklets(t.GetNumberOfTracklets());
200 if(!ntracklets){
201 printf("E - AliTRDcheckTRK::PropagateKalman :: No tracklets attached to track.\n");
202 return kFALSE;
203 }
3ceb45ae 204
205 AliTRDseedV1 *tr(NULL);
01ccc21a 206 if(!t.GetTrackIn()){
3ed01fbe 207 printf("E - AliTRDcheckTRK::PropagateKalman :: Track did not entered TRD fiducial volume.\n");
943761d3 208 return kFALSE;
3ceb45ae 209 }
35983729 210 if(!ref){
b9058d72 211 //printf("W - AliTRDcheckTRK::PropagateKalman :: Missing starting param.\n");
35983729 212 return kFALSE;
213 }
214 if(ref->Pt()<1.e-3) return kFALSE;
3ed01fbe 215
3ceb45ae 216
217 // Initialize TRD track to the reference
218 AliTRDtrackV1 tt;
219 tt.Set(ref->GetX(), ref->GetAlpha(), ref->GetParameter(), ref->GetCovariance());
3ed01fbe 220 tt.SetMass(t.GetMass());
35983729 221 tt.SetTrackOut(t.GetTrackOut());
b9058d72 222 if(UseITS()){
223 if(!tt.Rotate((Int_t(ref->GetAlpha()/AliTRDgeometry::GetAlpha()) +(ref->GetAlpha()>0.?1:-1)* 0.5)*AliTRDgeometry::GetAlpha()-ref->GetAlpha())) return kFALSE;
224 }
3ceb45ae 225 for(Int_t ily(0); ily<AliTRDgeometry::kNlayer; ily++){
3ed01fbe 226 if(!(tr = t.GetTracklet(ily))) continue;
227 Int_t det(tr->GetDetector());
b9058d72 228 if(HasClRecalibrate()){
3ed01fbe 229 AliTRDtransform trans(det);
3ceb45ae 230 AliTRDcluster *c(NULL);
3ed01fbe 231 Float_t exb, vd, t0, s2, dl, dt; tr->GetCalibParam(exb, vd, t0, s2, dl, dt);
3ceb45ae 232 tr->ResetClusterIter(kFALSE);
3ed01fbe 233 while((c = tr->PrevCluster())){
234 if(!trans.Transform(c/*, GetCalib(det)*/)){
b9058d72 235 printf("W - AliTRDcheckTRK::PropagateKalman :: Transform() failed for Det[%03d %02d_%d_%d]\n", det,
236 AliTRDgeometry::GetSector(det), AliTRDgeometry::GetStack(det), AliTRDgeometry::GetLayer(det));
3ed01fbe 237 break;
238 }
239 }
01ccc21a 240 }
b9058d72 241 if(HasTrkltRefit()){
01ccc21a 242 if(!tr->FitRobust(tt.Charge()>0.)) printf("W - AliTRDcheckTRK::PropagateKalman :: FitRobust() failed for Det[%03d]\n", det);
3ceb45ae 243 }
244 if(!AliTRDtrackerV1::PropagateToX(tt, tr->GetX0(), fgKalmanStep)) continue;
35983729 245 if(!tt.GetTrackIn()) tt.SetTrackIn();
3ed01fbe 246 tr->Update(&tt);
b9058d72 247 if(HasKalmanUpdate()){
3ceb45ae 248 Double_t x(tr->GetX0()),
249 p[2] = { tr->GetYfit(0), tr->GetZfit(0)},
250 covTrklt[3];
251 tr->GetCovAt(x, covTrklt);
252 if(!((AliExternalTrackParam&)tt).Update(p, covTrklt)) continue;
3ed01fbe 253 tt.SetTracklet(tr, 0);
254 tt.SetNumberOfClusters();
255 tt.UpdateChi2(((AliExternalTrackParam)tt).GetPredictedChi2(p, covTrklt));
3ceb45ae 256 }
3ceb45ae 257 }
3ed01fbe 258 //tt.Print("a");
259 t.~AliTRDtrackV1();
260 new(&t) AliTRDtrackV1(tt);
3ceb45ae 261 return kTRUE;
262}
263
264