BmnRoot
Loading...
Searching...
No Matches
BmnStsKFTrackFitter.cxx
Go to the documentation of this file.
1// -------------------------------------------------------------------------
2// ----- BmnStsKFTrackFitter source file -----
3// ----- Created 27/03/21 by A.Zinchenko, D.Zinchenko -----
4// ----- (from CbmStsKFTrackFitter.cxx) -----
5// -------------------------------------------------------------------------
6
7// AZ #include "CbmStsKFTrackFitter.h"
9
10#include "CbmKF.h"
11#include "CbmKFMath.h"
12#include "CbmKFStsHit.h"
13#include "CbmKFTrack.h"
14#include "CbmKFVertex.h"
15#include "CbmMvdHit.h"
16#include "CbmStsTrack.h"
17#include "CbmVertex.h"
18#include "FairRootManager.h"
19#include "TClonesArray.h"
20#include "TMath.h"
21#include "math.h"
22
23#include <TFile.h>
24#include <TKey.h>
25#include <TProfile2D.h>
26#include <iostream>
27
28using std::cout;
29using std::endl;
30
31// AZ CbmStsKFTrackFitter::CbmStsKFTrackFitter():
33 : fHits()
34 , fMvdHitsArray(0)
35 , fStsHitsArray(0)
37 , fMatBudgetFile(0)
38{}
39
41{
42 /* AZ-060825 - this break macro processing at the end
43 for (auto mp : fMatHistos)
44 delete mp.second;
45 if (fMatBudgetFile)
46 delete fMatBudgetFile; */
47}
48
49// -------------------------------------------------------------------------
50
51// AZ void CbmStsKFTrackFitter::Init(){
53{
54 // Initialisation
55 FairRootManager* rootMgr = FairRootManager::Instance();
56 if (NULL == rootMgr) {
57 // AZ cout << "-E- CbmStsKFTrackFitter::Init(): "
58 cout << "-E- BmnStsKFTrackFitter::Init(): "
59 << "ROOT manager is not instantiated!" << endl;
60 return;
61 }
62 fStsHitsArray = reinterpret_cast<TClonesArray*>(rootMgr->GetObject("StsHit"));
63 if (!fStsHitsArray) {
64 // AZ cout << "-W- CbmStsKFTrackFitter::Init: "
65 cout << "-W- BmnStsKFTrackFitter::Init: "
66 << "no STS hits array" << endl;
67 // return;
68 }
69 // AZ-040623 fMvdHitsArray = reinterpret_cast<TClonesArray*>(rootMgr->GetObject("MvdHit"));
70 fMvdHitsArray = nullptr; // AZ-040623
71 if (!fMvdHitsArray) {
72 // AZ cout << "-W- CbmStsKFTrackFitter::Init: "
73 // cout << "-W- BmnStsKFTrackFitter::Init: "
74 // << "no MVD hits array" << endl;
75 // return;
76 }
77 fIsInitialised = 1;
78};
79
80// -------------------------------------------------------------------------
81
82// AZ void CbmStsKFTrackFitter::SetKFHits(CbmKFTrack &T, CbmStsTrack* track){
84{
85
86 T.fHits.clear();
87
88 if (!fIsInitialised)
89 Init();
90
91 Int_t NStsHits = (fStsHitsArray) ? track->GetNStsHits() : 0;
92 Int_t NMvdHits = (fMvdHitsArray) ? track->GetNMvdHits() : 0;
93
94 fHits.resize(NMvdHits + NStsHits);
95 if (NMvdHits > 0) {
96 for (Int_t i = 0; i < NMvdHits; i++) {
97 Int_t j = track->GetMvdHitIndex(i);
98 fHits[i].Create(reinterpret_cast<CbmMvdHit*>(fMvdHitsArray->At(j)));
99 T.fHits.push_back(&(fHits[i]));
100 }
101 }
102 if (NStsHits > 0 && fStsHitsArray) {
103 for (Int_t i = 0; i < NStsHits; i++) {
104 Int_t j = track->GetStsHitIndex(i);
105 fHits[NMvdHits + i].Create(reinterpret_cast<CbmStsHit*>(fStsHitsArray->At(j)));
106 T.fHits.push_back(&(fHits[NMvdHits + i]));
107 }
108 }
109}
110
111// -------------------------------------------------------------------------
112
113// AZ Int_t CbmStsKFTrackFitter::DoFit( CbmStsTrack* track, Int_t pidHypo )
114Int_t BmnStsKFTrackFitter::DoFit(CbmStsTrack* track, Int_t pidHypo)
115{
116 track->SetPidHypo(pidHypo);
117
118 CbmKFTrack T;
119 T.SetPID(pidHypo);
120 SetKFHits(T, track);
121 for (Int_t i = 0; i < 6; i++)
122 T.GetTrack()[i] = 0.; // no guess taken
123 T.Fit(1); // fit downstream
124 CheckTrack(T);
125 T.Fit(0); // fit upstream
126 CheckTrack(T);
127 Int_t err = T.Fit(1); // fit downstream
128 Bool_t ok = (!err) && CheckTrack(T);
129 if (ok) {
130 T.GetTrackParam(*track->GetParamLast()); // store fitted track & cov.matrix
131 err = T.Fit(0); // fit upstream
132 ok = ok && (!err) && CheckTrack(T);
133 if (ok)
134 T.GetStsTrack(*track, 1); // store fitted track & cov.matrix & chi2 & NDF
135 }
136 if (!ok) {
137 Double_t* t = T.GetTrack();
138 Double_t* c = T.GetCovMatrix();
139 for (int i = 0; i < 6; i++)
140 t[i] = 0;
141 for (int i = 0; i < 15; i++)
142 c[i] = 0;
143 c[0] = c[2] = c[5] = c[9] = c[14] = 100.;
144 T.GetRefChi2() = 100.;
145 T.GetRefNDF() = 0;
146 T.GetStsTrack(*track, 0);
147 T.GetStsTrack(*track, 1);
148 track->SetFlag(1);
149 } else {
150 track->SetFlag(0);
151 }
152 // cout << " zzz " << track->GetParamFirst()->GetZ() << " " << track->GetParamLast()->GetZ() << " " <<
153 // track->GetNStsHits() << endl; //AZ
154 return !ok;
155}
156
157// -------------------------------------------------------------------------
158
160{
161 // Get material budget histograms
162
163 cout << "-Info- BmnStsKFTrackFitter::ReadMatBudget: "
164 << "reading mat. budget histograms from file " << fileName << endl;
165
166 if (fMatBudgetFile)
167 delete fMatBudgetFile;
168 fMatBudgetFile = TFile::Open(fileName, "read");
169 const TList* keys = fMatBudgetFile->GetListOfKeys();
170 TIter next(keys);
171 TObject* obj;
172
173 while ((obj = next())) {
174 TKey* key = (TKey*)obj;
175 // key->Print();
176 if (!(TString(key->GetName()).Contains("Station")))
177 continue;
178 TProfile2D* prof = (TProfile2D*)key->ReadObj();
179 // cout << prof->GetName() << endl;
180 TString title = prof->GetTitle();
181 int pos = title.Last('_');
182 TString zstr = title(pos + 1, title.Length() - pos);
183 Double_t zpos = zstr.Atof() / 10;
184 // cout << zstr << " " << zpos << endl;
185 fMatHistos[zpos] = prof;
186 // cout << prof->GetNbinsX() << " " << prof->GetNbinsY() << endl;
187 // file.Get(key->GetName())->Write(key->GetName(),TObject::kSingleKey);
188 }
189 // fitter.SetMatHistos(&fMatHistos);
190}
191
192// -------------------------------------------------------------------------
193
194Int_t BmnStsKFTrackFitter::Fit(CbmStsTrack* track, Int_t pidHypo)
195{
196 // Track fit with material budget
197
198 track->SetPidHypo(pidHypo);
199
200 CbmKFTrack T;
201 T.SetPID(pidHypo);
202 SetKFHits(T, track);
203 for (Int_t i = 0; i < 6; i++)
204 T.GetTrack()[i] = 0.; // no guess taken
205 // AZ T.Fit( 1 ); // fit downstream
206 FitWithMat(T, 1);
207 CheckTrack(T);
208 // AZ T.Fit( 0 ); // fit upstream
209 FitWithMat(T, 0);
210 CheckTrack(T);
211 // AZ Int_t err = T.Fit( 1 ); // fit downstream
212 Int_t err = FitWithMat(T, 1); // fit downstream
213 Bool_t ok = (!err) && CheckTrack(T);
214 if (ok) {
215 T.GetTrackParam(*track->GetParamLast()); // store fitted track & cov.matrix
216 // AZ err = T.Fit( 0 ); // fit upstream
217 err = FitWithMat(T, 0); // fit upstream
218 ok = ok && (!err) && CheckTrack(T);
219 if (ok)
220 T.GetStsTrack(*track, 1); // store fitted track & cov.matrix & chi2 & NDF
221 // cout << " zzz " << 1/track->GetParamLast()->GetQp() << " " << 1/track->GetParamFirst()->GetQp() << " "
222 // << 1/track->GetParamFirst()->GetQp()-1/track->GetParamLast()->GetQp() << endl; //AZ
223 }
224 if (!ok) {
225 Double_t* t = T.GetTrack();
226 Double_t* c = T.GetCovMatrix();
227 for (int i = 0; i < 6; i++)
228 t[i] = 0;
229 for (int i = 0; i < 15; i++)
230 c[i] = 0;
231 c[0] = c[2] = c[5] = c[9] = c[14] = 100.;
232 T.GetRefChi2() = 100.;
233 T.GetRefNDF() = 0;
234 T.GetStsTrack(*track, 0);
235 T.GetStsTrack(*track, 1);
236 track->SetFlag(1);
237 } else {
238 track->SetFlag(0);
239 }
240 return !ok;
241}
242
243// -------------------------------------------------------------------------
244
245Int_t BmnStsKFTrackFitter::FitWithMat(CbmKFTrack& track, Int_t downstream)
246{
247 // Track fitting procedure with material budget
248
249 CbmKF* KF = CbmKF::Instance();
250 Double_t* T = track.GetTrack();
251 Double_t* C = track.GetCovMatrix();
252 Int_t NHits = track.GetNOfHits();
253
254 // if(line) KF->SetMethod(0);
255
256 Bool_t err = 0;
257 if (NHits == 0)
258 return 1;
259
260 // use initial momentum
261 // this fixed value will be used during fit instead of T[4]
262
263 Double_t qp0 = T[4];
264
265 const Double_t INF = 10000.;
266
267 track.GetRefChi2() = 0;
268 track.GetRefNDF() = 0;
269
270 // initialize covariance matrix
271
272 C[0] = INF;
273 C[1] = 0.;
274 C[2] = INF;
275 C[3] = 0.;
276 C[4] = 0.;
277 C[5] = INF;
278 C[6] = 0.;
279 C[7] = 0.;
280 C[8] = 0.;
281 C[9] = INF;
282 C[10] = 0.;
283 C[11] = 0.;
284 C[12] = 0.;
285 C[13] = 0.;
286 C[14] = INF;
287
288 CbmKFMaterial* mat = nullptr;
289 CbmKFMaterial matEval;
290
291 try {
292
293 if (downstream) {
294 CbmKFHit* h = track.GetHit(0);
295 // AZ err = h->Filter( *this, downstream, qp0 );
296 err = h->Filter(track, downstream, qp0);
297 Int_t istold = h->MaterialIndex;
298 for (Int_t i = 1; i < NHits; i++) {
299 h = track.GetHit(i);
300 Int_t ist = h->MaterialIndex;
301 // cout << " xxxxxxxxxxxxxxxxxxx " << istold << " " << ist << endl;
302 // AZ for (Int_t j = istold+1; j < ist; j++) {
303 for (Int_t j = istold + 0; j < ist; j++) { // AZ
304 // AZ err = err || KF->vMaterial[j]->Pass( *this, downstream, qp0 );
305 // err = err || KF->vMaterial[j]->Pass( track, downstream, qp0 ); //AZ
306 if (fMatHistos.size()) {
307 EvalMaterial(track, i, downstream, matEval); // AZ
308 mat = &matEval;
309 } else
310 mat = KF->vMaterial[j];
311 // cout << " Before: " << qp0 << endl;
312 err = err || mat->Pass(track, downstream, qp0); // AZ
313 // cout << " After: " << qp0 << endl;
314 }
315 // AZ err = err || h->Filter( *this, downstream, qp0 );
316 err = err || h->Filter(track, downstream, qp0);
317 istold = ist;
318 }
319 } else {
320 CbmKFHit* h = track.GetHit(NHits - 1);
321 err = h->Filter(track, downstream, qp0);
322 Int_t istold = h->MaterialIndex;
323 for (Int_t i = NHits - 2; i >= 0; i--) {
324 h = track.GetHit(i);
325 Int_t ist = h->MaterialIndex;
326 // cout << " yyyyyyyyyyyyyyyyyyy " << istold << " " << ist << endl;
327 // AZ for(Int_t j=istold-1; j>ist; j--) {
328 for (Int_t j = istold - 0; j > ist; j--) { // AZ
329 // AZ err = err || KF->vMaterial[j]->Pass( track, downstream, qp0 );
330 if (fMatHistos.size()) {
331 EvalMaterial(track, i, downstream, matEval); // AZ
332 mat = &matEval;
333 } else
334 mat = KF->vMaterial[j];
335 err = err || mat->Pass(track, downstream, qp0); // AZ
336 }
337 err = err || h->Filter(track, downstream, qp0);
338 istold = ist;
339 }
340 }
341
342 // correct NDF value to number of fitted track parameters (straight line(4) or helix(5) )
343
344 track.GetRefNDF() -= (KF->GetMethod() == 0) ? 4 : 5;
345 } catch (...) {
346 track.GetRefChi2() = 0;
347 track.GetRefNDF() = 0;
348 C[0] = INF;
349 C[1] = 0.;
350 C[2] = INF;
351 C[3] = 0.;
352 C[4] = 0.;
353 C[5] = INF;
354 C[6] = 0.;
355 C[7] = 0.;
356 C[8] = 0.;
357 C[9] = INF;
358 C[10] = 0.;
359 C[11] = 0.;
360 C[12] = 0.;
361 C[13] = 0.;
362 C[14] = INF;
363 T[0] = T[1] = T[2] = T[3] = T[4] = T[5] = 0;
364 return 1;
365 }
366 if (err) {
367 track.GetRefChi2() = 0;
368 track.GetRefNDF() = 0;
369 C[0] = INF;
370 C[1] = 0.;
371 C[2] = INF;
372 C[3] = 0.;
373 C[4] = 0.;
374 C[5] = INF;
375 C[6] = 0.;
376 C[7] = 0.;
377 C[8] = 0.;
378 C[9] = INF;
379 C[10] = 0.;
380 C[11] = 0.;
381 C[12] = 0.;
382 C[13] = 0.;
383 C[14] = INF;
384 T[0] = T[1] = T[2] = T[3] = T[4] = T[5] = 0;
385 }
386 return err;
387}
388
389// -------------------------------------------------------------------------
390
391void BmnStsKFTrackFitter::EvalMaterial(CbmKFTrack& track, int ihit, Int_t downstream, CbmKFMaterial& mat)
392{
393 // AZ Evaluate material between hits (stations)
394
395 if ((ihit == 0 && downstream) || (ihit == track.GetNOfHits() - 1 && !downstream)) {
396 cout << "!!! BmnStsKFTrackFitter::EvalMaterial - this should not happen " << ihit << " " << downstream << " "
397 << track.GetNOfHits() << endl;
398 exit(0);
399 }
400
401 // CbmKFHit *h = track.GetHit (ihit);
402 // CbmKFHit *hprev = track.GetHit (ihit-1);
403 BmnKFStsHit* h = (BmnKFStsHit*)track.GetHit(ihit);
404 int iprev = (downstream) ? ihit - 1 : ihit + 1;
405 BmnKFStsHit* hprev = (BmnKFStsHit*)track.GetHit(iprev);
406 if (!downstream) {
407 BmnKFStsHit* htmp = h;
408 h = hprev;
409 hprev = htmp;
410 }
411 Double_t zOffset =
412 (h->FitPoint[0].sigma2 < 1.e-4) ? 0.5 : 1.5; // !!! Si or GEM - Attention for sigma2 (based on strip pitch)
413 // TProfile2D *histo = nullptr;
414 TProfile2D* histo = fMatHistos.lower_bound(h->GetZ() - zOffset)->second;
415 TAxis* axis[2] = {histo->GetXaxis(), histo->GetYaxis()};
416 int ixy[2] = {axis[0]->FindBin(h->GetX()), axis[1]->FindBin(h->GetY())};
417 int ixyprev[2] = {axis[0]->FindBin(hprev->GetX()), axis[1]->FindBin(hprev->GetY())};
418 Double_t thick = histo->GetBinContent(ixyprev[0], ixyprev[1]);
419 int nbins = 1;
420 if ((ixy[0] != ixyprev[0]) || (ixy[1] != ixyprev[1])) {
421 int idxy[2], ilong = 0, ishor = 1;
422 Double_t dxy[2];
423 for (int j = 0; j < 2; ++j)
424 idxy[j] = ixy[j] - ixyprev[j];
425 for (int j = 0; j < 2; ++j)
426 dxy[j] = idxy[j] * axis[j]->GetBinWidth(1);
427 if (TMath::Abs(idxy[0]) < TMath::Abs(idxy[1])) {
428 ilong = 1;
429 ishor = 0;
430 }
431 Double_t slope = dxy[ishor] / dxy[ilong], xyshort = axis[ishor]->GetBinCenter(ixyprev[ishor]);
432 int istep = TMath::Sign(1, idxy[ilong]), i = ixyprev[ilong];
433 while (1) {
434 i += istep;
435 ++nbins;
436 xyshort += slope;
437 int j = axis[ishor]->FindBin(xyshort);
438 if (ilong == 0)
439 thick += histo->GetBinContent(i, j);
440 else
441 thick += histo->GetBinContent(j, i);
442 // if (ilong == 0) thick = TMath::Max (histo->GetBinContent(i,j), thick); // max. value
443 // else thick = TMath::Max (histo->GetBinContent(j,i), thick); // max. value
444 if (i == ixy[ilong])
445 break;
446 }
447 }
448 thick /= (nbins * 100); // conversion from %
449 // thick *= 10; //!!! - test
450 mat.ID = 77777;
451 mat.ZReference = (h->GetZ() + hprev->GetZ()) / 2;
452 // mat.ZReference = (downstream) ? hprev->GetZ() + 0.1 : h->GetZ() - 0.1;
453 // mat.ZThickness = TMath::Abs (h->GetZ() - hprev->GetZ());
454 mat.ZThickness = TMath::Abs(h->GetZ() - hprev->GetZ()) - h->tube->dz;
455 mat.RadLength = mat.ZThickness / thick;
456 // mat.Fe = 0.10; // energy loss correction factor - what does it mean?
457 // cout << nbins << " " << thick << " " << mat.RadLength << " " << mat.ZThickness << " " << ixy[0] << " " <<
458 // ixyprev[0]
459 // << " " << ixy[1] << " " << ixyprev[1] << " " << h->GetZ() << " " << hprev->GetZ() << " " << histo->GetTitle()
460 // << endl;
461}
462
463// -------------------------------------------------------------------------
464
465// AZ void CbmStsKFTrackFitter::Extrapolate( FairTrackParam* track, Double_t z, FairTrackParam* e_track )
466void BmnStsKFTrackFitter::Extrapolate(FairTrackParam* track, Double_t z, FairTrackParam* e_track)
467{
468 if (!track)
469 return;
470 CbmKFTrack T;
471 T.SetTrackParam(*track);
472 T.Extrapolate(z);
473 if (e_track)
474 T.GetTrackParam(*e_track);
475}
476
477// -------------------------------------------------------------------------
478
479// AZ void CbmStsKFTrackFitter::Extrapolate( CbmStsTrack* track, Double_t z, FairTrackParam* e_track )
480void BmnStsKFTrackFitter::Extrapolate(CbmStsTrack* track, Double_t z, FairTrackParam* e_track)
481{
482 if (!track)
483 return;
484 CbmKFTrack T;
485 T.SetPID(track->GetPidHypo());
486 FairTrackParam *fpar = track->GetParamFirst(), *lpar = track->GetParamLast();
487
488 if (z <= fpar->GetZ()) { // extrapolate first parameters
489 T.SetTrackParam(*fpar);
490 T.Extrapolate(z);
491 } else if (z < fpar->GetZ() + 0.1) { // extrapolate first parameters
492 T.SetTrackParam(*fpar);
493 T.Propagate(z);
494 } else if (lpar->GetZ() <= z) { // extrapolate last parameters
495 T.SetTrackParam(*lpar);
496 T.Extrapolate(z);
497 } else if (lpar->GetZ() - 0.1 < z) { // extrapolate last parameters
498 T.SetTrackParam(*lpar);
499 T.Propagate(z);
500 } else { // refit with smoother
501 SetKFHits(T, track);
502 T.SetTrackParam(*fpar);
503 T.Smooth(z);
504 }
505 if (e_track)
506 T.GetTrackParam(*e_track);
507}
508
509// -------------------------------------------------------------------------
510
511// AZ Double_t CbmStsKFTrackFitter::GetChiToVertex( CbmStsTrack* track, CbmVertex *vtx )
513{
514 if (!vtx) {
515 FairRootManager* fManger = FairRootManager::Instance();
516 vtx = reinterpret_cast<CbmVertex*>(fManger->GetObject("PrimaryVertex"));
517 if (!vtx) {
518 cout << "-W- CbmStsKFTrackFitter::GetChiToVertex: No Primary Vertex found!" << endl;
519 return 100.;
520 }
521 }
522 CbmKFTrack T;
523 T.SetStsTrack(*track, 1);
524 T.Extrapolate(vtx->GetZ());
525
526 TMatrixFSym tmp(3);
527 vtx->CovMatrix(tmp);
528 Double_t Cv[3] = {tmp(0, 0), tmp(0, 1), tmp(1, 1)};
529
530 return CbmKFMath::getDeviation(T.GetTrack()[0], T.GetTrack()[1], T.GetCovMatrix(), vtx->GetX(), vtx->GetY(), Cv);
531}
532
533// -------------------------------------------------------------------------
534
535// AZ Double_t CbmStsKFTrackFitter::FitToVertex( CbmStsTrack* track, CbmVertex *vtx, FairTrackParam *v_track )
536Double_t BmnStsKFTrackFitter::FitToVertex(CbmStsTrack* track, CbmVertex* vtx, FairTrackParam* v_track)
537{
538 Double_t ret = 100.;
539 if (!track || !vtx || !v_track)
540 return ret;
541 CbmKFTrack T(*track);
542 CbmKFVertex V(*vtx);
543 T.Fit2Vertex(V);
544 T.GetTrackParam(*v_track);
545 if (T.GetRefNDF() > 0 && T.GetRefChi2() >= 0) {
546 ret = T.GetRefChi2() / T.GetRefNDF();
547 if (finite(ret))
548 ret = sqrt(ret);
549 }
550 return ret;
551}
552
553// -------------------------------------------------------------------------
554
555// AZ Bool_t CbmStsKFTrackFitter::CheckTrack( CbmKFTrack &T )
556Bool_t BmnStsKFTrackFitter::CheckTrack(CbmKFTrack& T)
557{
558 Bool_t ok = 1;
559 Double_t* t = T.GetTrack();
560 Double_t* c = T.GetCovMatrix();
561 for (int i = 0; i < 6; i++)
562 ok = ok && finite(t[i]) && TMath::Abs(t[i]) < 1.e5;
563 for (int i = 0; i < 15; i++)
564 ok = ok && finite(c[i]);
565 ok = ok && finite(T.GetMass()) && finite(T.GetRefChi2());
566 if (ok) {
567 ok = ok && (c[0] > 0);
568 ok = ok && (c[2] > 0);
569 ok = ok && (c[5] > 0);
570 ok = ok && (c[9] > 0);
571 ok = ok && (c[14] > 0);
572 }
573 if (ok) { // correct the cov matrix
574 Double_t c00 = TMath::Sqrt(c[0]);
575 Double_t c11 = TMath::Sqrt(c[2]);
576 Double_t c22 = TMath::Sqrt(c[5]);
577 Double_t c33 = TMath::Sqrt(c[9]);
578 Double_t c44 = TMath::Sqrt(c[14]);
579 Double_t a = c11 * c00;
580 if (c[1] > a)
581 c[1] = a;
582 if (c[1] < -a)
583 c[1] = -a;
584 a = c22 * c00;
585 if (c[3] > a)
586 c[3] = a;
587 if (c[3] < -a)
588 c[3] = -a;
589 a = c22 * c11;
590 if (c[4] > a)
591 c[4] = a;
592 if (c[4] < -a)
593 c[4] = -a;
594 a = c33 * c00;
595 if (c[6] > a)
596 c[6] = a;
597 if (c[6] < -a)
598 c[6] = -a;
599 a = c33 * c11;
600 if (c[7] > a)
601 c[7] = a;
602 if (c[7] < -a)
603 c[7] = -a;
604 a = c33 * c22;
605 if (c[8] > a)
606 c[8] = a;
607 if (c[8] < -a)
608 c[8] = -a;
609 a = c44 * c00;
610 if (c[10] > a)
611 c[10] = a;
612 if (c[10] < -a)
613 c[10] = -a;
614 a = c44 * c11;
615 if (c[11] > a)
616 c[11] = a;
617 if (c[11] < -a)
618 c[11] = -a;
619 a = c44 * c22;
620 if (c[12] > a)
621 c[12] = a;
622 if (c[12] < -a)
623 c[12] = -a;
624 a = c44 * c33;
625 if (c[13] > a)
626 c[13] = a;
627 if (c[13] < -a)
628 c[13] = -a;
629 }
630 if (!ok) {
631 for (int i = 0; i < 6; i++)
632 t[i] = 0;
633 for (int i = 0; i < 15; i++)
634 c[i] = 0;
635 c[0] = c[2] = c[5] = c[9] = c[14] = 100.;
636 T.GetRefChi2() = 100.;
637 T.GetRefNDF() = 0;
638 }
639 return ok;
640}
Bool_t fIsInitialised
kTRUE if Init() was called
friend F32vec4 sqrt(const F32vec4 &a)
Definition P4_F32vec4.h:34
int i
Definition P4_F32vec4.h:22
Double_t GetY() const
Definition BmnKFStsHit.h:32
CbmKFTube * tube
Definition BmnKFStsHit.h:27
Double_t GetX() const
Definition BmnKFStsHit.h:31
CbmKFUMeasurement FitPoint[2]
Definition BmnKFStsHit.h:26
Double_t GetZ() const
Definition BmnKFStsHit.h:33
void Extrapolate(CbmStsTrack *track, Double_t z, FairTrackParam *e_track)
Int_t FitWithMat(CbmKFTrack &track, Int_t downstream)
void EvalMaterial(CbmKFTrack &track, int ihit, Int_t downstream, CbmKFMaterial &mat)
void ReadMatBudget(TString &matBudgetFileName)
Double_t GetChiToVertex(CbmStsTrack *track, CbmVertex *vtx=0)
void SetKFHits(CbmKFTrack &T, CbmStsTrack *track)
Int_t DoFit(CbmStsTrack *track, Int_t pidHypo=211)
Int_t Fit(CbmStsTrack *track, Int_t pidHypo=211)
Double_t FitToVertex(CbmStsTrack *track, CbmVertex *vtx, FairTrackParam *v_track)
Int_t MaterialIndex
Definition CbmKFHit.h:23
virtual Int_t Filter(CbmKFTrackInterface &track, Bool_t downstream, Double_t &QP0)=0
Double_t ZReference
virtual Int_t Pass(Double_t ZCross, Double_t ZThick, CbmKFTrackInterface &track, Bool_t downstream, Double_t &QP0)
Double_t RadLength
Double_t ZThickness
static Double_t getDeviation(Double_t x, Double_t y, Double_t C[], Double_t vx, Double_t vy, Double_t Cv[]=0)
Int_t Propagate(Double_t z_out, Double_t QP0, Bool_t line=false)
Int_t Extrapolate(Double_t z, Double_t *QP0=0, Bool_t line=false)
Access to i-th hit.
Int_t Fit(Bool_t downstream=1, Bool_t line=false)
void Fit2Vertex(CbmKFVertexInterface &vtx)
void GetStsTrack(CbmStsTrack &track, bool first=1)
std::vector< CbmKFHit * > fHits
Definition CbmKFTrack.h:31
Double_t GetMass()
Definition CbmKFTrack.h:55
void SetStsTrack(CbmStsTrack &track, bool first=1)
CbmKFHit * GetHit(Int_t i)
Number of hits.
Definition CbmKFTrack.h:58
void SetTrackParam(FairTrackParam &track)
void GetTrackParam(FairTrackParam &track)
Double_t * GetTrack()
Is it electron.
Definition CbmKFTrack.h:51
Double_t & GetRefChi2()
array[15] of covariance matrix
Definition CbmKFTrack.h:53
Int_t GetNOfHits()
Number of Degrees of Freedom after fit.
Definition CbmKFTrack.h:57
Int_t & GetRefNDF()
Chi^2 after fit.
Definition CbmKFTrack.h:54
Double_t * GetCovMatrix()
array[6] of track parameters(x,y,tx,ty,qp,z)
Definition CbmKFTrack.h:52
void SetPID(Int_t pidHypo)
Double_t dz
Definition CbmKF.h:28
std::vector< CbmKFMaterial * > vMaterial
Definition CbmKF.h:58
static CbmKF * Instance()
Definition CbmKF.h:35
Int_t GetMethod()
Definition CbmKF.h:79
FairTrackParam * GetParamLast()
Definition CbmStsTrack.h:70
Int_t GetNMvdHits() const
Definition CbmStsTrack.h:61
Int_t GetPidHypo() const
Definition CbmStsTrack.h:64
Int_t GetMvdHitIndex(Int_t iHit) const
Definition CbmStsTrack.h:63
void SetFlag(Int_t flag)
Definition CbmStsTrack.h:85
Int_t GetNStsHits() const
Definition CbmStsTrack.h:60
Int_t GetStsHitIndex(Int_t iHit) const
Definition CbmStsTrack.h:62
FairTrackParam * GetParamFirst()
Definition CbmStsTrack.h:69
void SetPidHypo(Int_t pid)
Definition CbmStsTrack.h:82
Double_t GetZ() const
Definition CbmVertex.h:60
void CovMatrix(TMatrixFSym &covMat) const
Double_t GetX() const
Definition CbmVertex.h:58
Double_t GetY() const
Definition CbmVertex.h:59