]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCcalibCosmic.cxx
AliTPCcalibAlign.cxx - Protection again 0 pointer
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibCosmic.cxx
CommitLineData
f7f33dec 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
54b76c13 16/*
108953e9 17 Comments to be written here:
54b76c13 18 1. What do we calibrate.
19 2. How to interpret results
20 3. Simple example
21 4. Analysis using debug streamers.
22
23
24
25 3.Simple example
26 // To make cosmic scan the user interaction neccessary
27 //
28 .x ~/UliStyle.C
29 gSystem->Load("libANALYSIS");
30 gSystem->Load("libTPCcalib");
31 TFile fcalib("CalibObjects.root");
32 TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
33 AliTPCcalibCosmic * cosmic = ( AliTPCcalibCosmic *)array->FindObject("cosmicTPC");
34
35
3e55050f 36 //
37 //
38 gSystem->AddIncludePath("-I$ALICE_ROOT/TPC/macros");
39 gROOT->LoadMacro("$ALICE_ROOT/TPC/macros/AliXRDPROOFtoolkit.cxx+")
40 AliXRDPROOFtoolkit tool;
b9908d0b 41 TChain * chainCosmic = tool.MakeChainRandom("cosmicF.txt","Track0",0,10000);
a390f11f 42 //
b9908d0b 43 TCut cutITSN="min(Orig0.fITSncls,Orig1.fITSncls)>2";
44 TCut cutTPCN="min(Orig0.fTPCncls,Orig1.fTPCncls)>120";
45
46 chainCosmic->Draw(">>listITS",cutITSN+cutTPCN,"entryList");
47 TEntryList *elistITS = (TEntryList*)gDirectory->Get("listITS");
48 chainCosmic->SetEntryList(elistITS);
a390f11f 49
b9908d0b 50 */
54b76c13 51
52
53
f7f33dec 54#include "Riostream.h"
55#include "TChain.h"
56#include "TTree.h"
57#include "TH1F.h"
58#include "TH2F.h"
59#include "TList.h"
60#include "TMath.h"
61#include "TCanvas.h"
62#include "TFile.h"
54b76c13 63#include "TF1.h"
91fd44c9 64#include "THnSparse.h"
9963b5e2 65#include "TDatabasePDG.h"
f7f33dec 66
54b76c13 67#include "AliTPCclusterMI.h"
f7f33dec 68#include "AliTPCseed.h"
69#include "AliESDVertex.h"
70#include "AliESDEvent.h"
71#include "AliESDfriend.h"
72#include "AliESDInputHandler.h"
54b76c13 73#include "AliAnalysisManager.h"
f7f33dec 74
75#include "AliTracker.h"
f7a1cc68 76#include "AliMagF.h"
54b76c13 77#include "AliTPCCalROC.h"
f7f33dec 78
79#include "AliLog.h"
80
81#include "AliTPCcalibCosmic.h"
f7f33dec 82#include "TTreeStream.h"
83#include "AliTPCTracklet.h"
3326b323 84//#include "AliESDcosmic.h"
15e48021 85
f7f33dec 86
87ClassImp(AliTPCcalibCosmic)
88
89
90AliTPCcalibCosmic::AliTPCcalibCosmic()
91 :AliTPCcalibBase(),
9b27d39b 92 fHistNTracks(0),
93 fClusters(0),
94 fModules(0),
95 fHistPt(0),
9b27d39b 96 fDeDx(0),
54b76c13 97 fDeDxMIP(0),
98 fMIPvalue(1),
9b27d39b 99 fCutMaxD(5), // maximal distance in rfi ditection
a6dc0cf6 100 fCutMaxDz(40), // maximal distance in z ditection
9b27d39b 101 fCutTheta(0.03), // maximal distan theta
102 fCutMinDir(-0.99) // direction vector products
f7f33dec 103{
54b76c13 104 AliInfo("Default Constructor");
91fd44c9 105 for (Int_t ihis=0; ihis<6;ihis++){
106 fHistoDelta[ihis]=0;
107 fHistoPull[ihis]=0;
8a92e133 108 }
109 for (Int_t ihis=0; ihis<4;ihis++){
110 fHistodEdxMax[ihis] =0;
111 fHistodEdxTot[ihis] =0;
91fd44c9 112 }
f7f33dec 113}
114
115
116AliTPCcalibCosmic::AliTPCcalibCosmic(const Text_t *name, const Text_t *title)
117 :AliTPCcalibBase(),
9b27d39b 118 fHistNTracks(0),
119 fClusters(0),
120 fModules(0),
121 fHistPt(0),
9b27d39b 122 fDeDx(0),
54b76c13 123 fDeDxMIP(0),
124 fMIPvalue(1),
a6dc0cf6 125 fCutMaxD(5), // maximal distance in rfi ditection
126 fCutMaxDz(40), // maximal distance in z ditection
9b27d39b 127 fCutTheta(0.03), // maximal distan theta
128 fCutMinDir(-0.99) // direction vector products
f7f33dec 129{
130 SetName(name);
131 SetTitle(title);
54b76c13 132
133 fHistNTracks = new TH1F("ntracks","Number of Tracks per Event; number of tracks per event; number of tracks",501,-0.5,500.5);
134 fClusters = new TH1F("signal","Number of Clusters per track; number of clusters per track n_{cl}; counts",160,0,160);
135 fModules = new TH2F("sector","Acorde hits; z (cm); x(cm)",1200,-650,650,600,-700,700);
136 fHistPt = new TH1F("Pt","Pt distribution; p_{T} (GeV); counts",2000,0,50);
137 fDeDx = new TH2F("DeDx","dEdx; momentum p (GeV); TPC signal (a.u.)",500,0.01,100.,500,2.,1000);
f7f33dec 138 BinLogX(fDeDx);
54b76c13 139 fDeDxMIP = new TH1F("DeDxMIP","MIP region; TPC signal (a.u.);counts ",500,2.,1000);
91fd44c9 140 Init();
9b27d39b 141 AliInfo("Non Default Constructor");
bca20570 142 //
f7f33dec 143}
144
145AliTPCcalibCosmic::~AliTPCcalibCosmic(){
146 //
147 //
148 //
91fd44c9 149 for (Int_t ihis=0; ihis<6;ihis++){
150 delete fHistoDelta[ihis];
151 delete fHistoPull[ihis];
91fd44c9 152 }
8a92e133 153 for (Int_t ihis=0; ihis<4;ihis++){
154 delete fHistodEdxTot[ihis];
155 delete fHistodEdxMax[ihis];
156 }
157
158 delete fHistNTracks; // histogram showing number of ESD tracks per event
159 delete fClusters; // histogram showing the number of clusters per track
160 delete fModules; // 2d histogram of tracks which are propagated to the ACORDE scintillator array
161 delete fHistPt; // Pt histogram of reconstructed tracks
162 delete fDeDx; // dEdx spectrum showing the different particle types
163 delete fDeDxMIP; // TPC signal close to the MIP region of muons 0.4 < p < 0.45 GeV
f7f33dec 164}
165
166
91fd44c9 167void AliTPCcalibCosmic::Init(){
168 //
169 // init component
170 // Make performance histograms
171 //
172
173 // tracking performance bins
174 // 0 - delta of interest
175 // 1 - min (track0, track1) number of clusters
176 // 2 - R - vertex radius
177 // 3 - P1 - mean z
178 // 4 - P2 - snp(phi) at inner wall of TPC
179 // 5 - P3 - tan(theta) at inner wall of TPC
180 // 6 - P4 - 1/pt mean
181 // 7 - pt - pt mean
182 // 8 - alpha
183
184 Double_t xminTrack[9], xmaxTrack[9];
185 Int_t binsTrack[9];
186 TString axisName[9];
187 //
188 binsTrack[0] =100;
189 axisName[0] ="#Delta";
190 //
191 binsTrack[1] =8;
192 xminTrack[1] =80; xmaxTrack[1]=160;
193 axisName[1] ="N_{cl}";
194 //
195 binsTrack[2] =10;
196 xminTrack[2] =0; xmaxTrack[2]=90; //
197 axisName[2] ="dca_{r} (cm)";
198 //
199 binsTrack[3] =25;
200 xminTrack[3] =-250; xmaxTrack[3]=250; //
201 axisName[3] ="z (cm)";
202 //
203 binsTrack[4] =10;
204 xminTrack[4] =-0.8; xmaxTrack[4]=0.8; //
205 axisName[4] ="sin(#phi)";
206 //
207 binsTrack[5] =10;
208 xminTrack[5] =-1; xmaxTrack[5]=1; //
8a92e133 209 axisName[5] ="tan(#theta)";
91fd44c9 210 //
a390f11f 211 binsTrack[6] =40;
212 xminTrack[6] =-2; xmaxTrack[6]=2; //
91fd44c9 213 axisName[6] ="1/pt (1/GeV)";
214 //
8a92e133 215 binsTrack[7] =40;
91fd44c9 216 xminTrack[7] =0.2; xmaxTrack[7]=50; //
8a92e133 217 axisName[7] ="pt (GeV)";
91fd44c9 218 //
a390f11f 219 binsTrack[8] =18;
91fd44c9 220 xminTrack[8] =0; xmaxTrack[8]=TMath::Pi(); //
221 axisName[8] ="alpha";
222 //
223 // delta y
224 xminTrack[0] =-1; xmaxTrack[0]=1; //
225 fHistoDelta[0] = new THnSparseS("#Delta_{Y} (cm)","#Delta_{Y} (cm)", 9, binsTrack,xminTrack, xmaxTrack);
226 xminTrack[0] =-5; xmaxTrack[0]=5; //
227 fHistoPull[0] = new THnSparseS("#Delta_{Y} (unit)","#Delta_{Y} (unit)", 9, binsTrack,xminTrack, xmaxTrack);
228 //
229 // delta z
230 xminTrack[0] =-1; xmaxTrack[0]=1; //
231 fHistoDelta[1] = new THnSparseS("#Delta_{Z} (cm)","#Delta_{Z} (cm)", 9, binsTrack,xminTrack, xmaxTrack);
232 xminTrack[0] =-5; xmaxTrack[0]=5; //
233 fHistoPull[1] = new THnSparseS("#Delta_{Z} (unit)","#Delta_{Z} (unit)", 9, binsTrack,xminTrack, xmaxTrack);
234 //
235 // delta P2
236 xminTrack[0] =-10; xmaxTrack[0]=10; //
237 fHistoDelta[2] = new THnSparseS("#Delta_{#phi} (mrad)","#Delta_{#phi} (mrad)", 9, binsTrack,xminTrack, xmaxTrack);
238 xminTrack[0] =-5; xmaxTrack[0]=5; //
239 fHistoPull[2] = new THnSparseS("#Delta_{#phi} (unit)","#Delta_{#phi} (unit)", 9, binsTrack,xminTrack, xmaxTrack);
240 //
241 // delta P3
242 xminTrack[0] =-10; xmaxTrack[0]=10; //
243 fHistoDelta[3] = new THnSparseS("#Delta_{#theta} (mrad)","#Delta_{#theta} (mrad)", 9, binsTrack,xminTrack, xmaxTrack);
244 xminTrack[0] =-5; xmaxTrack[0]=5; //
245 fHistoPull[3] = new THnSparseS("#Delta_{#theta} (unit)","#Delta_{#theta} (unit)", 9, binsTrack,xminTrack, xmaxTrack);
246 //
247 // delta P4
248 xminTrack[0] =-0.2; xmaxTrack[0]=0.2; //
249 fHistoDelta[4] = new THnSparseS("#Delta_{1/pt} (1/GeV)","#Delta_{1/pt} (1/GeV)", 9, binsTrack,xminTrack, xmaxTrack);
250 xminTrack[0] =-5; xmaxTrack[0]=5; //
251 fHistoPull[4] = new THnSparseS("#Delta_{1/pt} (unit)","#Delta_{1/pt} (unit)", 9, binsTrack,xminTrack, xmaxTrack);
252
253 //
254 // delta Pt
255 xminTrack[0] =-0.5; xmaxTrack[0]=0.5; //
256 fHistoDelta[5] = new THnSparseS("#Delta_{pt}/p_{t}","#Delta_{pt}/p_{t}", 9, binsTrack,xminTrack, xmaxTrack);
257 xminTrack[0] =-5; xmaxTrack[0]=5; //
258 fHistoPull[5] = new THnSparseS("#Delta_{pt}/p_{t} (unit)","#Delta_{pt}/p_{t} (unit)", 9, binsTrack,xminTrack, xmaxTrack);
259 //
8a92e133 260
261 for (Int_t idedx=0;idedx<4;idedx++){
262 xminTrack[0] =0.5; xmaxTrack[0]=1.5; //
263 binsTrack[1] =40;
264 xminTrack[1] =10; xmaxTrack[1]=160;
265
266 fHistodEdxMax[idedx] = new THnSparseS(Form("dEdx_{MaxUp}/dEdx_{MaxDown}_Pad%d",idedx),
267 Form("dEdx_{MaxUp}/dEdx_{MaxDown}_Pad%d",idedx),
268 9, binsTrack,xminTrack, xmaxTrack);
269 fHistodEdxTot[idedx] = new THnSparseS(Form("dEdx_{TotUp}/dEdx_{TotDown}_Pad%d",idedx),
270 Form("dEdx_{TotUp}/dEdx_{TotDown}_Pad%d",idedx),
271 9, binsTrack,xminTrack, xmaxTrack);
272 }
273
274
275
91fd44c9 276 for (Int_t ivar=0;ivar<6;ivar++){
8a92e133 277 for (Int_t ivar2=0;ivar2<9;ivar2++){
91fd44c9 278 fHistoDelta[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
279 fHistoDelta[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
280 fHistoPull[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
281 fHistoPull[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
282 BinLogX(fHistoDelta[ivar],7);
283 BinLogX(fHistoPull[ivar],7);
8a92e133 284 if (ivar<4){
285 fHistodEdxMax[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
286 fHistodEdxMax[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
287 fHistodEdxTot[ivar]->GetAxis(ivar2)->SetName(axisName[ivar2].Data());
288 fHistodEdxTot[ivar]->GetAxis(ivar2)->SetTitle(axisName[ivar2].Data());
289 BinLogX(fHistodEdxMax[ivar],7);
290 BinLogX(fHistodEdxTot[ivar],7);
291 }
91fd44c9 292 }
293 }
294}
295
296
297void AliTPCcalibCosmic::Add(const AliTPCcalibCosmic* cosmic){
298 //
299 //
300 //
301 for (Int_t ivar=0; ivar<6;ivar++){
302 if (fHistoDelta[ivar] && cosmic->fHistoDelta[ivar]){
303 fHistoDelta[ivar]->Add(cosmic->fHistoDelta[ivar]);
304 }
305 if (fHistoPull[ivar] && cosmic->fHistoPull[ivar]){
306 fHistoPull[ivar]->Add(cosmic->fHistoPull[ivar]);
307 }
308 }
8a92e133 309 for (Int_t ivar=0; ivar<4;ivar++){
310 if (fHistodEdxMax[ivar] && cosmic->fHistodEdxMax[ivar]){
311 fHistodEdxMax[ivar]->Add(cosmic->fHistodEdxMax[ivar]);
312 }
313 if (fHistodEdxTot[ivar] && cosmic->fHistodEdxTot[ivar]){
314 fHistodEdxTot[ivar]->Add(cosmic->fHistodEdxTot[ivar]);
315 }
316 }
91fd44c9 317}
318
9b27d39b 319
320
321
f7f33dec 322void AliTPCcalibCosmic::Process(AliESDEvent *event) {
9b27d39b 323 //
324 //
325 //
f7f33dec 326 if (!event) {
327 Printf("ERROR: ESD not available");
328 return;
9b27d39b 329 }
76273318 330 AliESDfriend *esdFriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
331 if (!esdFriend) {
332 Printf("ERROR: esdFriend not available");
f7f33dec 333 return;
334 }
2acad464 335
108953e9 336
54b76c13 337 FindPairs(event); // nearly everything takes place in find pairs...
338
2acad464 339 if (GetDebugLevel()>20) printf("Hallo world: Im here and processing an event\n");
9b27d39b 340 Int_t ntracks=event->GetNumberOfTracks();
341 fHistNTracks->Fill(ntracks);
9b27d39b 342 if (ntracks==0) return;
cc65e4f5 343 // AliESDcosmic cosmicESD;
344// TTreeSRedirector * cstream = GetDebugStreamer();
345// cosmicESD.SetDebugStreamer(cstream);
346// cosmicESD.ProcessEvent(event);
347// if (cstream) cosmicESD.DumpToTree();
15e48021 348
349
54b76c13 350}
9b27d39b 351
f7f33dec 352
3326b323 353void AliTPCcalibCosmic::FillHistoPerformance(const AliExternalTrackParam *par0, const AliExternalTrackParam *par1, const AliExternalTrackParam *inner0, const AliExternalTrackParam */*inner1*/, AliTPCseed *seed0, AliTPCseed *seed1, const AliExternalTrackParam *param0Combined ){
91fd44c9 354 //
a390f11f 355 // par0,par1 - parameter of tracks at DCA 0
356 // inner0,inner1 - parameter of tracks at the TPC entrance
357 // seed0, seed1 - detailed track information
358 // param0Combined - Use combined track parameters for binning
91fd44c9 359 //
8a92e133 360 Int_t kMinCldEdx =20;
361 Int_t ncl0 = seed0->GetNumberOfClusters();
362 Int_t ncl1 = seed1->GetNumberOfClusters();
363
91fd44c9 364 const Double_t kpullCut = 10;
365 Double_t x[9];
366 Double_t xyz0[3];
367 Double_t xyz1[3];
368 par0->GetXYZ(xyz0);
369 par1->GetXYZ(xyz1);
370 Double_t radius0 = TMath::Sqrt(xyz0[0]*xyz0[0]+xyz0[1]*xyz0[1]);
371 Double_t radius1 = TMath::Sqrt(xyz1[0]*xyz1[0]+xyz1[1]*xyz1[1]);
372 inner0->GetXYZ(xyz0);
373 Double_t alpha = TMath::ATan2(xyz0[1],xyz0[0]);
374 // bin parameters
375 x[1] = TMath::Min(ncl0,ncl1);
376 x[2] = (radius0+radius1)*0.5;
a390f11f 377 x[3] = param0Combined->GetZ();
378 x[4] = inner0->GetSnp();
379 x[5] = param0Combined->GetTgl();
380 x[6] = param0Combined->GetSigned1Pt();
381 x[7] = param0Combined->Pt();
91fd44c9 382 x[8] = alpha;
383 // deltas
384 Double_t delta[6];
385 Double_t sigma[6];
386 delta[0] = (par0->GetY()+par1->GetY());
387 delta[1] = (par0->GetZ()-par1->GetZ());
388 delta[2] = (par0->GetAlpha()-par1->GetAlpha()-TMath::Pi());
389 delta[3] = (par0->GetTgl()+par1->GetTgl());
390 delta[4] = (par0->GetParameter()[4]+par1->GetParameter()[4]);
391 delta[5] = (par0->Pt()-par1->Pt())/((par0->Pt()+par1->Pt())*0.5);
392 //
393 sigma[0] = TMath::Sqrt(par0->GetSigmaY2()+par1->GetSigmaY2());
394 sigma[1] = TMath::Sqrt(par0->GetSigmaZ2()+par1->GetSigmaZ2());
395 sigma[2] = TMath::Sqrt(par0->GetSigmaSnp2()+par1->GetSigmaSnp2());
396 sigma[3] = TMath::Sqrt(par0->GetSigmaTgl2()+par1->GetSigmaTgl2());
397 sigma[4] = TMath::Sqrt(par0->GetSigma1Pt2()+par1->GetSigma1Pt2());
398 sigma[5] = sigma[4]*((par0->Pt()+par1->Pt())*0.5);
399 //
400 Bool_t isOK = kTRUE;
401 for (Int_t ivar=0;ivar<6;ivar++){
402 if (sigma[ivar]==0) isOK=kFALSE;
403 x[0]= delta[ivar]/sigma[ivar];
404 if (TMath::Abs(x[0])>kpullCut) isOK = kFALSE;
405 }
406 //
407
408 if (isOK) for (Int_t ivar=0;ivar<6;ivar++){
409 x[0]= delta[ivar]/TMath::Sqrt(2);
410 if (ivar==2 || ivar ==3) x[0]*=1000;
411 fHistoDelta[ivar]->Fill(x);
412 if (sigma[ivar]>0){
413 x[0]= delta[ivar]/sigma[ivar];
414 fHistoPull[ivar]->Fill(x);
415 }
416 }
8a92e133 417
418 //
419 // Fill dedx performance
420 //
421 for (Int_t ipad=0; ipad<4;ipad++){
422 //
423 //
424 //
425 Int_t row0=0;
426 Int_t row1=160;
427 if (ipad==0) row1=63;
428 if (ipad==1) {row0=63; row1=63+64;}
429 if (ipad==2) {row0=128;}
430 Int_t nclUp = TMath::Nint(seed0->CookdEdxAnalytical(0.01,0.7,0,row0,row1,2));
431 Int_t nclDown = TMath::Nint(seed1->CookdEdxAnalytical(0.01,0.7,0,row0,row1,2));
432 Int_t minCl = TMath::Min(nclUp,nclDown);
433 if (minCl<kMinCldEdx) continue;
434 x[1] = minCl;
435 //
436 Float_t dEdxTotUp = seed0->CookdEdxAnalytical(0.01,0.7,0,row0,row1);
437 Float_t dEdxTotDown = seed1->CookdEdxAnalytical(0.01,0.7,0,row0,row1);
438 Float_t dEdxMaxUp = seed0->CookdEdxAnalytical(0.01,0.7,1,row0,row1);
439 Float_t dEdxMaxDown = seed1->CookdEdxAnalytical(0.01,0.7,1,row0,row1);
440 //
441 if (dEdxTotDown<=0) continue;
442 if (dEdxMaxDown<=0) continue;
443 x[0]=dEdxTotUp/dEdxTotDown;
444 fHistodEdxTot[ipad]->Fill(x);
445 x[0]=dEdxMaxUp/dEdxMaxDown;
446 fHistodEdxMax[ipad]->Fill(x);
447 }
448
449
91fd44c9 450
451}
452
3326b323 453void AliTPCcalibCosmic::MaterialBudgetDump(AliExternalTrackParam *const par0, AliExternalTrackParam *const par1, const AliExternalTrackParam *inner0, const AliExternalTrackParam *inner1, AliTPCseed *const seed0, AliTPCseed *const seed1, AliExternalTrackParam *const param0Combined, AliExternalTrackParam *const param1Combined){
a390f11f 454 //
455 // matrial budget AOD dump
456 //
457 // par0,par1 - parameter of tracks at DCA 0
458 // inner0,inner1 - parameter of tracks at the TPC entrance
459 // seed0, seed1 - detailed track information
460 // param0Combined - Use combined track parameters for binning
461 // param1Combined -
462 Double_t p0In = inner0->GetP();
463 Double_t p1In = inner1->GetP();
464 Double_t p0V = par0->GetP();
465 Double_t p1V = par1->GetP();
466 //
467 Double_t pt0In = inner0->Pt();
468 Double_t pt1In = inner1->Pt();
469 Double_t pt0V = par0->Pt();
470 Double_t pt1V = par1->Pt();
471 Int_t ncl0 = seed0->GetNumberOfClusters();
472 Int_t ncl1 = seed1->GetNumberOfClusters();
473 Int_t nclmin=TMath::Min(ncl0,ncl1);
474 Double_t sign = (param0Combined->GetSigned1Pt()>0) ? 1:-1.;
475 //
476 TTreeSRedirector * pcstream = GetDebugStreamer();
477 if (pcstream){
478 (*pcstream)<<"material"<<
479 "run="<<fRun<< // run number
480 "event="<<fEvent<< // event number
481 "time="<<fTime<< // time stamp of event
482 "trigger="<<fTrigger<< // trigger
483 "triggerClass="<<&fTriggerClass<< // trigger
484 "mag="<<fMagF<< // magnetic field
485 "sign="<<sign<< // sign of the track
486 //
487 "ncl0="<<ncl0<<
488 "ncl1="<<ncl1<<
489 "nclmin="<<nclmin<<
490 //
491 "p0In="<<p0In<<
492 "p1In="<<p1In<<
493 "p0V="<<p0V<<
494 "p1V="<<p1V<<
495 "pt0In="<<pt0In<<
496 "pt1In="<<pt1In<<
497 "pt0V="<<pt0V<<
498 "pt1V="<<pt1V<<
499 "p0.="<<par0<<
500 "p1.="<<par1<<
501 "up0.="<<param0Combined<<
502 "up1.="<<param1Combined<<
503 "\n";
504 }
505
506}
f7f33dec 507
54b76c13 508void AliTPCcalibCosmic::Analyze() {
509
510 fMIPvalue = CalculateMIPvalue(fDeDxMIP);
511
512 return;
513
514}
515
f7f33dec 516
f7f33dec 517
9b27d39b 518void AliTPCcalibCosmic::FindPairs(AliESDEvent *event) {
519 //
520 // Find cosmic pairs
04087794 521 //
522 // Track0 is choosen in upper TPC part
523 // Track1 is choosen in lower TPC part
9b27d39b 524 //
2acad464 525 if (GetDebugLevel()>20) printf("Hallo world: Im here\n");
76273318 526 AliESDfriend *esdFriend=static_cast<AliESDfriend*>(event->FindListObject("AliESDfriend"));
9b27d39b 527 Int_t ntracks=event->GetNumberOfTracks();
9b27d39b 528 TObjArray tpcSeeds(ntracks);
529 if (ntracks==0) return;
530 Double_t vtxx[3]={0,0,0};
531 Double_t svtxx[3]={0.000001,0.000001,100.};
532 AliESDVertex vtx(vtxx,svtxx);
533 //
534 //track loop
535 //
54b76c13 536 for (Int_t i=0;i<ntracks;++i) {
537 AliESDtrack *track = event->GetTrack(i);
538 fClusters->Fill(track->GetTPCNcls());
539
540 const AliExternalTrackParam * trackIn = track->GetInnerParam();
541 const AliExternalTrackParam * trackOut = track->GetOuterParam();
542 if (!trackIn) continue;
543 if (!trackOut) continue;
860b3d93 544 if (ntracks>4 && TMath::Abs(trackIn->GetTgl())<0.0015) continue; // filter laser
545
546
76273318 547 AliESDfriendTrack *friendTrack = esdFriend->GetTrack(i);
b9908d0b 548 if (!friendTrack) continue;
9b27d39b 549 TObject *calibObject;
b9908d0b 550 AliTPCseed *seed = 0;
9b27d39b 551 for (Int_t l=0;(calibObject=friendTrack->GetCalibObject(l));++l) {
552 if ((seed=dynamic_cast<AliTPCseed*>(calibObject))) break;
553 }
554 if (seed) tpcSeeds.AddAt(seed,i);
54b76c13 555
556 Double_t meanP = 0.5*(trackIn->GetP() + trackOut->GetP());
557 if (seed && track->GetTPCNcls() > 80 + 60/(1+TMath::Exp(-meanP+5))) {
15e48021 558 fDeDx->Fill(meanP, seed->CookdEdxNorm(0.0,0.45,0,0,159));
54b76c13 559 //
15e48021 560 if (meanP > 0.4 && meanP < 0.45) fDeDxMIP->Fill(seed->CookdEdxNorm(0.0,0.45,0,0,159));
54b76c13 561 //
82628455 562 // if (GetDebugLevel()>0&&meanP>0.2&&seed->CookdEdxNorm(0.0,0.45,0,0,159)>300) {
563// //TFile *curfile = AliAnalysisManager::GetAnalysisManager()->GetTree()->GetCurrentFile();
564// //if (curfile) printf(">>> p+ in file: %s \t event: %i \t Number of ESD tracks: %i \n", curfile->GetName(), (int)event->GetEventNumberInFile(), (int)ntracks);
565// // if (track->GetOuterParam()->GetAlpha()<0) cout << " Polartiy: " << track->GetSign() << endl;
566// }
54b76c13 567
568 }
569
9b27d39b 570 }
54b76c13 571
9b27d39b 572 if (ntracks<2) return;
573 //
574 // Find pairs
575 //
576 for (Int_t i=0;i<ntracks;++i) {
577 AliESDtrack *track0 = event->GetTrack(i);
04087794 578 // track0 - choosen upper part
579 if (!track0) continue;
580 if (!track0->GetOuterParam()) continue;
581 if (track0->GetOuterParam()->GetAlpha()<0) continue;
e9362f9d 582 Double_t dir0[3];
583 track0->GetDirection(dir0);
04087794 584 for (Int_t j=0;j<ntracks;++j) {
585 if (i==j) continue;
586 AliESDtrack *track1 = event->GetTrack(j);
587 //track 1 lower part
588 if (!track1) continue;
589 if (!track1->GetOuterParam()) continue;
590 if (track1->GetOuterParam()->GetAlpha()>0) continue;
591 //
e9362f9d 592 Double_t dir1[3];
593 track1->GetDirection(dir1);
54b76c13 594
9b27d39b 595 AliTPCseed * seed0 = (AliTPCseed*) tpcSeeds.At(i);
596 AliTPCseed * seed1 = (AliTPCseed*) tpcSeeds.At(j);
597 if (! seed0) continue;
598 if (! seed1) continue;
15e48021 599 Float_t dedx0 = seed0->CookdEdxNorm(0.05,0.55,0,0,159);
600 Float_t dedx1 = seed1->CookdEdxNorm(0.05,0.55,0,0,159);
bca20570 601 //
15e48021 602 Float_t dedx0I = seed0->CookdEdxNorm(0.05,0.55,0,0,63);
603 Float_t dedx1I = seed1->CookdEdxNorm(0.05,0.55,0,0,63);
bca20570 604 //
15e48021 605 Float_t dedx0O = seed0->CookdEdxNorm(0.05,0.55,0,64,159);
606 Float_t dedx1O = seed1->CookdEdxNorm(0.05,0.55,0,64,159);
bca20570 607 //
e9362f9d 608 Float_t dir = (dir0[0]*dir1[0] + dir0[1]*dir1[1] + dir0[2]*dir1[2]);
9b27d39b 609 Float_t d0 = track0->GetLinearD(0,0);
610 Float_t d1 = track1->GetLinearD(0,0);
611 //
612 // conservative cuts - convergence to be guarantied
613 // applying before track propagation
614 if (TMath::Abs(d0+d1)>fCutMaxD) continue; // distance to the 0,0
615 if (dir>fCutMinDir) continue; // direction vector product
616 Float_t bz = AliTracker::GetBz();
617 Float_t dvertex0[2]; //distance to 0,0
618 Float_t dvertex1[2]; //distance to 0,0
619 track0->GetDZ(0,0,0,bz,dvertex0);
620 track1->GetDZ(0,0,0,bz,dvertex1);
621 if (TMath::Abs(dvertex0[1])>250) continue;
622 if (TMath::Abs(dvertex1[1])>250) continue;
623 //
624 //
625 //
626 Float_t dmax = TMath::Max(TMath::Abs(d0),TMath::Abs(d1));
627 AliExternalTrackParam param0(*track0);
628 AliExternalTrackParam param1(*track1);
629 //
630 // Propagate using Magnetic field and correct fo material budget
631 //
a390f11f 632 Double_t sign0=-1;
633 Double_t sign1=1;
634 Double_t maxsnp=0.90;
9963b5e2 635 AliTracker::PropagateTrackToBxByBz(&param0,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE,maxsnp,sign0);
636 AliTracker::PropagateTrackToBxByBz(&param1,dmax+1,TDatabasePDG::Instance()->GetParticle("e-")->Mass(),3,kTRUE,maxsnp,sign1);
9b27d39b 637 //
638 // Propagate rest to the 0,0 DCA - z should be ignored
639 //
640 Bool_t b0 = param0.PropagateToDCA(&vtx,bz,1000);
641 Bool_t b1 = param1.PropagateToDCA(&vtx,bz,1000);
642 //
643 param0.GetDZ(0,0,0,bz,dvertex0);
644 param1.GetDZ(0,0,0,bz,dvertex1);
a6dc0cf6 645 if (TMath::Abs(param0.GetZ()-param1.GetZ())>fCutMaxDz) continue;
9b27d39b 646 //
647 Double_t xyz0[3];//,pxyz0[3];
648 Double_t xyz1[3];//,pxyz1[3];
649 param0.GetXYZ(xyz0);
650 param1.GetXYZ(xyz1);
651 Bool_t isPair = IsPair(&param0,&param1);
652 //
54b76c13 653 if (isPair) FillAcordeHist(track0);
15e48021 654 //
655 // combined track params
656 //
657 AliExternalTrackParam *par0U=MakeCombinedTrack(&param0,&param1);
658 AliExternalTrackParam *par1U=MakeCombinedTrack(&param1,&param0);
659
660
54b76c13 661 //
9b27d39b 662 if (fStreamLevel>0){
663 TTreeSRedirector * cstream = GetDebugStreamer();
108953e9 664 //printf("My stream=%p\n",(void*)cstream);
e9362f9d 665 AliExternalTrackParam *ip0 = (AliExternalTrackParam *)track0->GetInnerParam();
666 AliExternalTrackParam *ip1 = (AliExternalTrackParam *)track1->GetInnerParam();
667 AliExternalTrackParam *op0 = (AliExternalTrackParam *)track0->GetOuterParam();
668 AliExternalTrackParam *op1 = (AliExternalTrackParam *)track1->GetOuterParam();
669 Bool_t isCrossI = ip0->GetZ()*ip1->GetZ()<0;
670 Bool_t isCrossO = op0->GetZ()*op1->GetZ()<0;
671 Double_t alpha0 = TMath::ATan2(dir0[1],dir0[0]);
672 Double_t alpha1 = TMath::ATan2(dir1[1],dir1[0]);
91fd44c9 673 //
674 //
675 //
a390f11f 676 FillHistoPerformance(&param0, &param1, ip0, ip1, seed0, seed1,par0U);
677 MaterialBudgetDump(&param0, &param1, ip0, ip1, seed0, seed1,par0U,par1U);
9b27d39b 678 if (cstream) {
679 (*cstream) << "Track0" <<
108953e9 680 "run="<<fRun<< // run number
681 "event="<<fEvent<< // event number
682 "time="<<fTime<< // time stamp of event
683 "trigger="<<fTrigger<< // trigger
15e48021 684 "triggerClass="<<&fTriggerClass<< // trigger
108953e9 685 "mag="<<fMagF<< // magnetic field
9b27d39b 686 "dir="<<dir<< // direction
54b76c13 687 "OK="<<isPair<< // will be accepted
9b27d39b 688 "b0="<<b0<< // propagate status
689 "b1="<<b1<< // propagate status
e9362f9d 690 "crossI="<<isCrossI<< // cross inner
691 "crossO="<<isCrossO<< // cross outer
692 //
9b27d39b 693 "Orig0.=" << track0 << // original track 0
694 "Orig1.=" << track1 << // original track 1
695 "Tr0.="<<&param0<< // track propagated to the DCA 0,0
696 "Tr1.="<<&param1<< // track propagated to the DCA 0,0
e9362f9d 697 "Ip0.="<<ip0<< // inner param - upper
698 "Ip1.="<<ip1<< // inner param - lower
699 "Op0.="<<op0<< // outer param - upper
700 "Op1.="<<op1<< // outer param - lower
15e48021 701 "Up0.="<<par0U<< // combined track 0
702 "Up1.="<<par1U<< // combined track 1
e9362f9d 703 //
9b27d39b 704 "v00="<<dvertex0[0]<< // distance using kalman
705 "v01="<<dvertex0[1]<< //
706 "v10="<<dvertex1[0]<< //
707 "v11="<<dvertex1[1]<< //
708 "d0="<<d0<< // linear distance to 0,0
709 "d1="<<d1<< // linear distance to 0,0
710 //
e9362f9d 711 //
712 //
713 "x00="<<xyz0[0]<< // global position close to vertex
9b27d39b 714 "x01="<<xyz0[1]<<
715 "x02="<<xyz0[2]<<
716 //
e9362f9d 717 "x10="<<xyz1[0]<< // global position close to vertex
9b27d39b 718 "x11="<<xyz1[1]<<
719 "x12="<<xyz1[2]<<
720 //
e9362f9d 721 "alpha0="<<alpha0<<
722 "alpha1="<<alpha1<<
723 "dir00="<<dir0[0]<< // direction upper
724 "dir01="<<dir0[1]<<
725 "dir02="<<dir0[2]<<
726 //
727 "dir10="<<dir1[0]<< // direction lower
728 "dir11="<<dir1[1]<<
729 "dir12="<<dir1[2]<<
730 //
731 //
54b76c13 732 "Seed0.=" << seed0 << // original seed 0
733 "Seed1.=" << seed1 << // original seed 1
bca20570 734 //
735 "dedx0="<<dedx0<< // dedx0 - all
736 "dedx1="<<dedx1<< // dedx1 - all
737 //
54b76c13 738 "dedx0I="<<dedx0I<< // dedx0 - inner ROC
739 "dedx1I="<<dedx1I<< // dedx1 - inner ROC
bca20570 740 //
54b76c13 741 "dedx0O="<<dedx0O<< // dedx0 - outer ROC
742 "dedx1O="<<dedx1O<< // dedx1 - outer ROC
9b27d39b 743 "\n";
744 }
15e48021 745 }
746 delete par0U;
747 delete par1U;
9b27d39b 748 }
749 }
750}
751
9b27d39b 752
bca20570 753
754
54b76c13 755void AliTPCcalibCosmic::FillAcordeHist(AliESDtrack *upperTrack) {
756
757 // Pt cut to select straight tracks which can be easily propagated to ACORDE which is outside the magnetic field
758 if (upperTrack->Pt() < 10 || upperTrack->GetTPCNcls() < 80) return;
759
76273318 760 const Double_t acordePlane = 850.; // distance of the central Acorde detectors to the beam line at y =0
54b76c13 761 const Double_t roof = 210.5; // distance from x =0 to end of magnet roof
762
763 Double_t r[3];
764 upperTrack->GetXYZ(r);
765 Double_t d[3];
766 upperTrack->GetDirection(d);
767 Double_t x,z;
76273318 768 z = r[2] + (d[2]/d[1])*(acordePlane - r[1]);
769 x = r[0] + (d[0]/d[1])*(acordePlane - r[1]);
54b76c13 770
771 if (x > roof) {
76273318 772 x = r[0] + (d[0]/(d[0]+d[1]))*(acordePlane+roof-r[0]-r[1]);
773 z = r[2] + (d[2]/(d[0]+d[1]))*(acordePlane+roof-r[0]-r[1]);
54b76c13 774 }
775 if (x < -roof) {
76273318 776 x = r[0] + (d[0]/(d[1]-d[0]))*(acordePlane+roof+r[0]-r[1]);
777 z = r[2] + (d[2]/(d[1]-d[0]))*(acordePlane+roof+r[0]-r[1]);
54b76c13 778 }
779
780 fModules->Fill(z, x);
781
782}
783
784
785
3326b323 786Long64_t AliTPCcalibCosmic::Merge(TCollection *const li) {
bca20570 787
788 TIterator* iter = li->MakeIterator();
789 AliTPCcalibCosmic* cal = 0;
790
791 while ((cal = (AliTPCcalibCosmic*)iter->Next())) {
792 if (!cal->InheritsFrom(AliTPCcalibCosmic::Class())) {
860b3d93 793 //Error("Merge","Attempt to add object of class %s to a %s", cal->ClassName(), this->ClassName());
bca20570 794 return -1;
795 }
796
54b76c13 797 fHistNTracks->Add(cal->GetHistNTracks());
798 fClusters->Add(cal-> GetHistClusters());
799 fModules->Add(cal->GetHistAcorde());
800 fHistPt->Add(cal->GetHistPt());
801 fDeDx->Add(cal->GetHistDeDx());
802 fDeDxMIP->Add(cal->GetHistMIP());
91fd44c9 803 Add(cal);
bca20570 804 }
bca20570 805 return 0;
9b27d39b 806
f7f33dec 807}
808
54b76c13 809
9b27d39b 810Bool_t AliTPCcalibCosmic::IsPair(AliExternalTrackParam *tr0, AliExternalTrackParam *tr1){
811 //
812 //
813 /*
814 // 0. Same direction - OPOSITE - cutDir +cutT
815 TCut cutDir("cutDir","dir<-0.99")
816 // 1.
817 TCut cutT("cutT","abs(Tr1.fP[3]+Tr0.fP[3])<0.03")
818 //
819 // 2. The same rphi
820 TCut cutD("cutD","abs(Tr0.fP[0]+Tr1.fP[0])<5")
821 //
822 //
823 //
824 TCut cutPt("cutPt","abs(Tr1.fP[4]+Tr0.fP[4])<1&&abs(Tr0.fP[4])+abs(Tr1.fP[4])<10");
825 // 1/Pt diff cut
826 */
827 const Double_t *p0 = tr0->GetParameter();
828 const Double_t *p1 = tr1->GetParameter();
829 if (TMath::Abs(p0[3]+p1[3])>fCutTheta) return kFALSE;
a6dc0cf6 830 if (TMath::Abs(p0[1]-p1[1])>fCutMaxDz) return kFALSE;
9b27d39b 831 if (TMath::Abs(p0[0]+p1[0])>fCutMaxD) return kFALSE;
a6dc0cf6 832
9b27d39b 833 Double_t d0[3], d1[3];
834 tr0->GetDirection(d0);
835 tr1->GetDirection(d1);
836 if (d0[0]*d1[0] + d0[1]*d1[1] + d0[2]*d1[2] >fCutMinDir) return kFALSE;
837 //
838 return kTRUE;
839}
54b76c13 840
841
842
843Double_t AliTPCcalibCosmic::CalculateMIPvalue(TH1F * hist) {
844
845 TF1 * funcDoubleGaus = new TF1("funcDoubleGaus", "gaus(0)+gaus(3)",0,1000);
846 funcDoubleGaus->SetParameters(hist->GetEntries()*0.75,hist->GetMean()/1.3,hist->GetMean()*0.10,
847 hist->GetEntries()*0.25,hist->GetMean()*1.3,hist->GetMean()*0.10);
848 hist->Fit(funcDoubleGaus);
3326b323 849 Double_t aMIPvalue = TMath::Min(funcDoubleGaus->GetParameter(1),funcDoubleGaus->GetParameter(4));
54b76c13 850
851 delete funcDoubleGaus;
852
3326b323 853 return aMIPvalue;
54b76c13 854
855}
9b27d39b 856
857
f7f33dec 858
54b76c13 859
57dc06f2 860void AliTPCcalibCosmic::CalculateBetheParams(TH2F */*hist*/, Double_t * /*initialParam*/) {
861 //
862 // Not implemented yet
863 //
54b76c13 864 return;
865
866}
867
868
3326b323 869void AliTPCcalibCosmic::BinLogX(THnSparse *const h, Int_t axisDim) {
91fd44c9 870
871 // Method for the correct logarithmic binning of histograms
872
873 TAxis *axis = h->GetAxis(axisDim);
874 int bins = axis->GetNbins();
875
876 Double_t from = axis->GetXmin();
877 Double_t to = axis->GetXmax();
3326b323 878 Double_t *newBins = new Double_t[bins + 1];
91fd44c9 879
3326b323 880 newBins[0] = from;
91fd44c9 881 Double_t factor = pow(to/from, 1./bins);
882
883 for (int i = 1; i <= bins; i++) {
3326b323 884 newBins[i] = factor * newBins[i-1];
91fd44c9 885 }
3326b323 886 axis->Set(bins, newBins);
887 delete newBins;
91fd44c9 888
889}
890
891
3326b323 892void AliTPCcalibCosmic::BinLogX(TH1 *const h) {
f7f33dec 893
894 // Method for the correct logarithmic binning of histograms
895
896 TAxis *axis = h->GetXaxis();
897 int bins = axis->GetNbins();
898
899 Double_t from = axis->GetXmin();
900 Double_t to = axis->GetXmax();
3326b323 901 Double_t *newBins = new Double_t[bins + 1];
f7f33dec 902
3326b323 903 newBins[0] = from;
f7f33dec 904 Double_t factor = pow(to/from, 1./bins);
905
906 for (int i = 1; i <= bins; i++) {
3326b323 907 newBins[i] = factor * newBins[i-1];
f7f33dec 908 }
3326b323 909 axis->Set(bins, newBins);
910 delete newBins;
f7f33dec 911
912}
913
7b18d067 914
860b3d93 915AliExternalTrackParam *AliTPCcalibCosmic::MakeTrack(const AliExternalTrackParam *track0, const AliExternalTrackParam *track1){
916 //
917 //
918 //
919 AliExternalTrackParam *par1R= new AliExternalTrackParam(*track1);
920 par1R->Rotate(track0->GetAlpha());
15e48021 921 par1R->PropagateTo(track0->GetX(),AliTracker::GetBz());
860b3d93 922 //
923 //
924 Double_t * param = (Double_t*)par1R->GetParameter();
925 Double_t * covar = (Double_t*)par1R->GetCovariance();
15e48021 926
860b3d93 927 param[0]*=1; //OK
928 param[1]*=1; //OK
929 param[2]*=1; //?
930 param[3]*=-1; //OK
931 param[4]*=-1; //OK
932 //
933 covar[6] *=-1.; covar[7] *=-1.; covar[8] *=-1.;
934 //covar[10]*=-1.; covar[11]*=-1.; covar[12]*=-1.;
935 covar[13]*=-1.;
860b3d93 936 return par1R;
937}
938
15e48021 939AliExternalTrackParam *AliTPCcalibCosmic::MakeCombinedTrack(const AliExternalTrackParam *track0, const AliExternalTrackParam *track1){
940 //
941 // Make combined track
942 //
943 //
944 AliExternalTrackParam * par1T = MakeTrack(track0,track1);
945 AliExternalTrackParam * par0U = new AliExternalTrackParam(*track0);
946 //
947 UpdateTrack(*par0U,*par1T);
948 delete par1T;
949 return par0U;
950}
951
952
860b3d93 953void AliTPCcalibCosmic::UpdateTrack(AliExternalTrackParam &track1, const AliExternalTrackParam &track2){
954 //
955 // Update track 1 with track 2
956 //
957 //
958 //
959 TMatrixD vecXk(5,1); // X vector
960 TMatrixD covXk(5,5); // X covariance
961 TMatrixD matHk(5,5); // vector to mesurement
962 TMatrixD measR(5,5); // measurement error
963 TMatrixD vecZk(5,1); // measurement
964 //
965 TMatrixD vecYk(5,1); // Innovation or measurement residual
966 TMatrixD matHkT(5,5);
967 TMatrixD matSk(5,5); // Innovation (or residual) covariance
968 TMatrixD matKk(5,5); // Optimal Kalman gain
969 TMatrixD mat1(5,5); // update covariance matrix
970 TMatrixD covXk2(5,5); //
971 TMatrixD covOut(5,5);
972 //
973 Double_t *param1=(Double_t*) track1.GetParameter();
974 Double_t *covar1=(Double_t*) track1.GetCovariance();
975 Double_t *param2=(Double_t*) track2.GetParameter();
976 Double_t *covar2=(Double_t*) track2.GetCovariance();
977 //
978 // copy data to the matrix
979 for (Int_t ipar=0; ipar<5; ipar++){
860b3d93 980 for (Int_t jpar=0; jpar<5; jpar++){
981 covXk(ipar,jpar) = covar1[track1.GetIndex(ipar, jpar)];
982 measR(ipar,jpar) = covar2[track2.GetIndex(ipar, jpar)];
15e48021 983 matHk(ipar,jpar)=0;
984 mat1(ipar,jpar)=0;
860b3d93 985 }
15e48021 986 vecXk(ipar,0) = param1[ipar];
987 vecZk(ipar,0) = param2[ipar];
988 matHk(ipar,ipar)=1;
989 mat1(ipar,ipar)=0;
860b3d93 990 }
991 //
992 //
993 //
994 //
860b3d93 995 //
996 vecYk = vecZk-matHk*vecXk; // Innovation or measurement residual
997 matHkT=matHk.T(); matHk.T();
998 matSk = (matHk*(covXk*matHkT))+measR; // Innovation (or residual) covariance
999 matSk.Invert();
1000 matKk = (covXk*matHkT)*matSk; // Optimal Kalman gain
1001 vecXk += matKk*vecYk; // updated vector
860b3d93 1002 covXk2 = (mat1-(matKk*matHk));
1003 covOut = covXk2*covXk;
1004 //
1005 //
1006 //
1007 // copy from matrix to parameters
1008 if (0) {
1009 vecXk.Print();
1010 vecZk.Print();
1011 //
1012 measR.Print();
1013 covXk.Print();
1014 covOut.Print();
1015 //
1016 track1.Print();
1017 track2.Print();
1018 }
1019
1020 for (Int_t ipar=0; ipar<5; ipar++){
1021 param1[ipar]= vecXk(ipar,0) ;
1022 for (Int_t jpar=0; jpar<5; jpar++){
1023 covar1[track1.GetIndex(ipar, jpar)]=covOut(ipar,jpar);
1024 }
1025 }
1026}
1027
860b3d93 1028
1029
1030