CbmRoot
Loading...
Searching...
No Matches
CbmKFParticleFinder.cxx
Go to the documentation of this file.
1/* Copyright (C) 2014-2020 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Maksym Zyzak, Volker Friese [committer] */
4
5//-----------------------------------------------------------
6//-----------------------------------------------------------
7
8// Cbm Headers ----------------------
10
11#include "CbmEvent.h"
12#include "CbmKF.h"
14#include "CbmKFVertex.h"
15#include "CbmL1PFFitter.h"
16#include "CbmMCDataArray.h"
17#include "CbmMCDataManager.h"
18#include "CbmMCEventList.h"
19#include "CbmMCTrack.h"
20#include "CbmTrackMatchNew.h"
21#include "FairRunAna.h"
22
23//KF Particle headers
24#include "KFPTrackVector.h"
25#include "KFParticleTopoReconstructor.h"
26
27#include <Logger.h>
28
29//ROOT headers
30#include "TClonesArray.h" //to get arrays from the FairRootManager
31#include "TMath.h" //to calculate Prob function
32#include "TStopwatch.h" //to measure the time
33
34//c++ and std headers
35#include <cmath>
36#include <iostream>
37#include <vector>
38using std::vector;
39
40CbmKFParticleFinder::CbmKFParticleFinder(const char* name, Int_t iVerbose)
41 : FairTask(name, iVerbose)
42 , fStsTrackBranchName("StsTrack")
43 , fTrackArray(nullptr)
44 , fEvents(nullptr)
45 , fTopoReconstructor(nullptr)
46 , fPVFindMode(2)
47 , fPID(nullptr)
48 , fSuperEventAnalysis(0)
49 , fSETracks(0)
50 , fSEField(0)
51 , fSEpdg(0)
52 , fSETrackId(0)
53 , fSEChiPrim(0)
54{
55 fTopoReconstructor = new KFParticleTopoReconstructor;
56
57 // set default cuts
58 SetPrimaryProbCut(0.0001); // 0.01% to consider primary track as a secondary;
59}
60
67
69{
70 std::string prefix = std::string(GetName()) + "::Init: ";
71 //Get ROOT Manager
72 FairRootManager* ioman = FairRootManager::Instance();
73
74 if (ioman == nullptr) {
75 LOG(error) << prefix << "RootManager not instantiated!";
76 return kERROR;
77 }
78
79 // Get the event branch
80
81 if (!ioman->CheckBranch("CbmEvent")) {
82 LOG(error) << prefix << "No event branch found. Run the event builder first.";
83 return kERROR;
84 }
85
86 fEvents = (TClonesArray*) ioman->GetObject("CbmEvent");
87 if (fEvents == nullptr) {
88 LOG(fatal) << prefix << "No events available in event-by-event mode.";
89 return kERROR;
90 }
91
92 // Get input collection
93 fTrackArray = (TClonesArray*) ioman->GetObject(fStsTrackBranchName);
94 if (fTrackArray == nullptr) {
95 LOG(error) << "track-array not found!";
96 return kERROR;
97 }
98
99 //In case of reconstruction of pure signal no PV is defined. The MC PV is used.
100 if (fPVFindMode == 0) {
101
102 CbmMCDataManager* mcManager = (CbmMCDataManager*) ioman->GetObject("MCDataManager");
103 if (mcManager == nullptr) {
104 LOG(error) << prefix << "MC Data Manager not found!";
105 return kERROR;
106 }
107
108 fMcTrackArray = mcManager->InitBranch("MCTrack");
109
110 if (fMcTrackArray == nullptr) {
111 LOG(error) << prefix << "MC track array not found!";
112 return kERROR;
113 }
114
115 fMcEventList = (CbmMCEventList*) ioman->GetObject("MCEventList.");
116 if (fMcEventList == nullptr) {
117 LOG(error) << prefix << "MC Event List not found!";
118 return kERROR;
119 }
120
121 fTrackMatchArray = (TClonesArray*) ioman->GetObject("StsTrackMatch");
122 if (fTrackMatchArray == nullptr) {
123 LOG(error) << prefix << " Sts Track Match array not found!";
124 return kERROR;
125 }
126 } // if (fPVFindMode == 0)
127
128 auto& target = CbmKF::Instance()->vTargets[0];
129 const std::array<float, 3> targetXYZ{(float) target.x, (float) target.y, (float) target.z};
130 fTopoReconstructor->SetTarget(targetXYZ);
131
132 return kSUCCESS;
133}
134
135void CbmKFParticleFinder::Exec(Option_t* /*opt*/)
136{
137 fTopoReconstructor->Clear();
138
139 int nEvents = fEvents->GetEntriesFast();
140
141 vector<KFParticleTopoReconstructor> eventTopoReconstructor(nEvents);
142
143 for (int iEvent = 0, firstEventTrack = 0, nTracksEvent = 0; iEvent < nEvents;
144 iEvent++, firstEventTrack += nTracksEvent) {
145
146 CbmEvent* event = dynamic_cast<CbmEvent*>(fEvents->At(iEvent));
147
148 nTracksEvent = event->GetNofStsTracks();
149
150 eventTopoReconstructor[iEvent].SetTarget(fTopoReconstructor->GetTargetPosition());
151 eventTopoReconstructor[iEvent].SetChi2PrimaryCut(InversedChi2Prob(0.0001, 2));
152 eventTopoReconstructor[iEvent].CopyCuts(fTopoReconstructor);
153 eventTopoReconstructor[iEvent].GetKFParticleFinder()->SetReconstructionList(
154 fTopoReconstructor->GetKFParticleFinder()->GetReconstructionList());
155
156 Int_t ntracks = 0; //fTrackArray->GetEntriesFast();
157
158 //calculate number of d-He4 candidates
159 int nCandidatesDHe4 = 0;
160 if (fPID) {
161 if ((int) fPID->GetPID().size() >= firstEventTrack + nTracksEvent) {
162 for (int iTr = 0; iTr < nTracksEvent; iTr++) {
163 if (TMath::Abs(fPID->GetPID()[firstEventTrack + iTr]) == 1000010029) {
164 nCandidatesDHe4++;
165 }
166 }
167 }
168 else {
169 LOG(error) << "CbmKFParticleFinder::Event: PID task has a wrong number of tracks: " << fPID->GetPID().size()
170 << " of " << firstEventTrack + nTracksEvent;
171 }
172 }
173
174 vector<CbmStsTrack> vRTracks(nTracksEvent + nCandidatesDHe4);
175 vector<int> pdg(nTracksEvent + nCandidatesDHe4, -1);
176 vector<int> trackId(nTracksEvent + nCandidatesDHe4, -1);
177
178 for (int iTr = 0; iTr < nTracksEvent; iTr++) {
179 int stsTrackIndex = event->GetStsTrackIndex(iTr);
180 CbmStsTrack* stsTrack = dynamic_cast<CbmStsTrack*>(fTrackArray->At(stsTrackIndex));
181
182 const FairTrackParam* parameters = stsTrack->GetParamFirst();
183
184 Double_t V[15] = {0.f};
185 for (Int_t i = 0, iCov = 0; i < 5; i++) {
186 for (Int_t j = 0; j <= i; j++, iCov++) {
187 V[iCov] = parameters->GetCovariance(i, j);
188 }
189 }
190
191 if (stsTrack->GetTotalNofHits() < 3) {
192 continue;
193 }
194
195 bool ok = 1;
196 ok = ok && std::isfinite(parameters->GetX());
197 ok = ok && std::isfinite(parameters->GetY());
198 ok = ok && std::isfinite(parameters->GetZ());
199 ok = ok && std::isfinite(parameters->GetTx());
200 ok = ok && std::isfinite(parameters->GetTy());
201 ok = ok && std::isfinite(parameters->GetQp());
202
203 for (unsigned short iC = 0; iC < 15; iC++) {
204 ok = ok && std::isfinite(V[iC]);
205 }
206
207 ok = ok && (V[0] < 1. && V[0] > 0.) && (V[2] < 1. && V[2] > 0.) && (V[5] < 1. && V[5] > 0.)
208 && (V[9] < 1. && V[9] > 0.) && (V[14] < 1. && V[14] > 0.);
209 ok = ok && stsTrack->GetChiSq() < 10 * stsTrack->GetNDF();
210 if (!ok) {
211 continue;
212 }
213
214 if (fPID) {
215 if (fPID->GetPID()[stsTrackIndex] == -2) {
216 continue;
217 }
218
219 //not clear separation between d and He4
220 if (TMath::Abs(fPID->GetPID()[stsTrackIndex]) == 1000010029) {
221 int sgn = fPID->GetPID()[stsTrackIndex] / TMath::Abs(fPID->GetPID()[stsTrackIndex]);
222 pdg[ntracks] = sgn * 1000010020;
223 vRTracks[ntracks] = *stsTrack;
224 trackId[ntracks] = stsTrackIndex;
225 ntracks++;
226
227 pdg[ntracks] = sgn * 1000020040;
228 }
229 else {
230 pdg[ntracks] = fPID->GetPID()[stsTrackIndex];
231 }
232 }
233 else {
234 pdg[ntracks] = -1;
235 }
236
237 vRTracks[ntracks] = *stsTrack;
238 trackId[ntracks] = stsTrackIndex;
239
240 ntracks++;
241 } // track loop
242
243 vRTracks.resize(ntracks);
244 pdg.resize(ntracks);
245 trackId.resize(ntracks);
246
247 CbmL1PFFitter fitter;
248 vector<float> vChiToPrimVtx;
249
250 CbmKFVertex kfVertex;
251
252 if (fPVFindMode == 0) {
253
254 // find an MC event with that matches the reco event
255 int iMcEvent = -1;
256 {
257 int nMCEvents = fMcEventList->GetNofEvents();
258 std::vector<double> mcWeight(nMCEvents, 0.);
259 for (int trId : trackId) {
260 CbmTrackMatchNew* stsTrackMatch = dynamic_cast<CbmTrackMatchNew*>(fTrackMatchArray->At(trId));
261 for (int iLink = 0; iLink < stsTrackMatch->GetNofLinks(); iLink++) {
262 CbmLink link = stsTrackMatch->GetLink(iLink);
263 Int_t iTrackMcEvent = fMcEventList->GetEventIndex(link);
264 if (iTrackMcEvent < 0 || iTrackMcEvent >= nMCEvents) continue;
265 mcWeight[iTrackMcEvent] += link.GetWeight();
266 }
267 }
268
269 iMcEvent = std::distance(mcWeight.begin(), std::max_element(mcWeight.begin(), mcWeight.end()));
270
271 if (iMcEvent < 0 || iMcEvent >= nMCEvents) {
272 LOG(error) << "CbmKFParticleFinder::Event: No MC event found for event " << iEvent;
273 break;
274 }
275 }
276
277 bool isMCPVFound = false;
278
279 CbmLink mcEventLink = fMcEventList->GetEventLinkByIndex(iMcEvent);
280 int nMCTracks = fMcTrackArray->Size(mcEventLink);
281
282 for (Int_t iMC = 0; (iMC < nMCTracks) && (!isMCPVFound); iMC++) {
283 CbmLink mcTrackLink = mcEventLink;
284 mcTrackLink.SetIndex(iMC);
285 auto* cbmMCTrack = dynamic_cast<CbmMCTrack*>(fMcTrackArray->Get(mcTrackLink));
286 assert(cbmMCTrack);
287 if (cbmMCTrack->GetMotherId() < 0) {
288 kfVertex.GetRefX() = cbmMCTrack->GetStartX();
289 kfVertex.GetRefY() = cbmMCTrack->GetStartY();
290 kfVertex.GetRefZ() = cbmMCTrack->GetStartZ();
291 isMCPVFound = true;
292 }
293 }
294 if (!isMCPVFound) {
295 break;
296 }
297 } // if (fPVFindMode == 0)
298
299 if (fPVFindMode == 3) {
300 const CbmVertex* primVertex = event->GetVertex();
301 assert(primVertex);
302 kfVertex.GetRefX() = primVertex->GetX();
303 kfVertex.GetRefY() = primVertex->GetY();
304 kfVertex.GetRefZ() = primVertex->GetZ();
305 }
306
307 vector<CbmL1PFFitter::PFFieldRegion> vField, vFieldAtLastPoint;
308 fitter.Fit(vRTracks, pdg);
309 fitter.GetChiToVertex(vRTracks, vField, vChiToPrimVtx, kfVertex, 3);
310 fitter.CalculateFieldRegionAtLastPoint(vRTracks, vFieldAtLastPoint);
311 vector<KFFieldVector> vFieldVector(ntracks), vFieldVectorAtLastPoint(ntracks);
312 for (Int_t iTr = 0; iTr < ntracks; iTr++) {
313 for (int i = 0; i < 10; i++) {
314 vFieldVector[iTr].fField[i] = vField[iTr].fP[i];
315 }
316 }
317 for (Int_t iTr = 0; iTr < ntracks; iTr++) {
318 for (int i = 0; i < 10; i++) {
319 vFieldVectorAtLastPoint[iTr].fField[i] = vFieldAtLastPoint[iTr].fP[i];
320 }
321 }
322
323 if (!fSuperEventAnalysis) {
324 KFPTrackVector tracks;
325 FillKFPTrackVector(&tracks, vRTracks, vFieldVector, pdg, trackId, vChiToPrimVtx);
326 KFPTrackVector tracksAtLastPoint;
327 FillKFPTrackVector(&tracksAtLastPoint, vRTracks, vFieldVectorAtLastPoint, pdg, trackId, vChiToPrimVtx, 0);
328
329 TStopwatch timer;
330 timer.Start();
331
332 eventTopoReconstructor[iEvent].Init(tracks, tracksAtLastPoint);
333
334 if (fPVFindMode == 0 || fPVFindMode == 3) {
335 KFPVertex primVtx_tmp;
336 primVtx_tmp.SetXYZ(kfVertex.GetRefX(), kfVertex.GetRefY(), kfVertex.GetRefZ());
337 primVtx_tmp.SetCovarianceMatrix(0, 0, 0, 0, 0, 0);
338 primVtx_tmp.SetNContributors(0);
339 primVtx_tmp.SetChi2(-100);
340
341 vector<int> pvTrackIds;
342 KFVertex pv(primVtx_tmp);
343 eventTopoReconstructor[iEvent].AddPV(pv, pvTrackIds);
344 }
345 else if (fPVFindMode == 1) {
346 eventTopoReconstructor[iEvent].ReconstructPrimVertex();
347 }
348 else if (fPVFindMode == 2) {
349 eventTopoReconstructor[iEvent].ReconstructPrimVertex(0);
350 }
351
352 eventTopoReconstructor[iEvent].SortTracks();
353 eventTopoReconstructor[iEvent].ReconstructParticles();
354
355 timer.Stop();
356 eventTopoReconstructor[iEvent].SetTime(timer.RealTime());
357 }
358 else {
359 for (int iTr = 0; iTr < ntracks; iTr++) {
360 const FairTrackParam* parameters = vRTracks[iTr].GetParamFirst();
361 float a = parameters->GetTx(), b = parameters->GetTy(), qp = parameters->GetQp();
362 Int_t q = 0;
363 if (qp > 0.f) {
364 q = 1;
365 }
366 if (qp < 0.f) {
367 q = -1;
368 }
369 float c2 = 1.f / (1.f + a * a + b * b);
370 float pq = 1.f / qp * TMath::Abs(q);
371 float p2 = pq * pq;
372 float pz = sqrt(p2 * c2);
373 float px = a * pz;
374 float py = b * pz;
375 float pt = sqrt(px * px + py * py);
376
377 bool save = 0;
378
379 if (vChiToPrimVtx[iTr] < 3) {
380 if ((fabs(pdg[iTr]) == 11 && pt > 0.2f) || (fabs(pdg[iTr]) == 13) || (fabs(pdg[iTr]) == 19)) {
381 save = 1;
382 }
383 }
384
385 if (vChiToPrimVtx[iTr] > 3) {
386 if ((fabs(pdg[iTr]) == 211 || fabs(pdg[iTr]) == 321 || fabs(pdg[iTr]) == 2212 || pdg[iTr] == -1)
387 && pt > 0.2f) {
388 save = 1;
389 }
390 }
391
392 if (save) {
393 fSETracks.push_back(vRTracks[iTr]);
394 fSEField.push_back(vFieldVector[iTr]);
395 fSEpdg.push_back(pdg[iTr]);
396 fSETrackId.push_back(fSETrackId.size());
397 fSEChiPrim.push_back(vChiToPrimVtx[iTr]);
398 }
399 }
400 }
401 } // event loop
402
403
404 vector<int> PVTrackIndexArray;
405 int indexAdd = 0;
406 for (int iEvent = 0; iEvent < nEvents; iEvent++) {
407 const KFParticleTopoReconstructor* eventTR = &eventTopoReconstructor[iEvent];
408
409 for (int iPV = 0; iPV < eventTR->NPrimaryVertices(); iPV++) {
410 PVTrackIndexArray = eventTR->GetPVTrackIndexArray(iPV);
411 for (unsigned int iTr = 0; iTr < PVTrackIndexArray.size(); iTr++) {
412 PVTrackIndexArray[iTr] = PVTrackIndexArray[iTr] + indexAdd;
413 }
414 fTopoReconstructor->AddPV(eventTR->GetPrimKFVertex(iPV), PVTrackIndexArray);
415 PVTrackIndexArray.clear();
416 }
417
418 for (unsigned int iP = 0; iP < eventTR->GetParticles().size(); iP++) {
419 const KFParticle& particleEvent = eventTR->GetParticles()[iP];
420 KFParticle particle = particleEvent;
421 particle.CleanDaughtersId();
422 int idP = particleEvent.Id() + indexAdd;
423 int idDaughter = 0;
424 for (int nD = 0; nD < particleEvent.NDaughters(); nD++) {
425 if (particleEvent.NDaughters() == 1) {
426 idDaughter = particleEvent.DaughterIds()[nD];
427 }
428 if (particleEvent.NDaughters() > 1) {
429 idDaughter = particleEvent.DaughterIds()[nD] + indexAdd;
430 }
431 particle.AddDaughterId(idDaughter);
432 }
433 particle.SetId(idP);
434 fTopoReconstructor->AddParticle(particle);
435 }
436
437 indexAdd += eventTR->GetParticles().size();
438 }
439}
440
442{
444 KFPTrackVector tracks;
446 KFPTrackVector tracksAtLastPoint;
447
448 LOG(info) << "CbmKFParticleFinder: Start SE analysis";
449 TStopwatch timer;
450 timer.Start();
451
452 fTopoReconstructor->Init(tracks, tracksAtLastPoint);
453
454 KFPVertex primVtx_tmp;
455 primVtx_tmp.SetXYZ(0, 0, 0);
456 primVtx_tmp.SetCovarianceMatrix(0, 0, 0, 0, 0, 0);
457 primVtx_tmp.SetNContributors(0);
458 primVtx_tmp.SetChi2(-100);
459 vector<int> pvTrackIds;
460 KFVertex pv(primVtx_tmp);
461 fTopoReconstructor->AddPV(pv, pvTrackIds);
462
463 fTopoReconstructor->SortTracks();
464 fTopoReconstructor->ReconstructParticles();
465
466 timer.Stop();
467 fTopoReconstructor->SetTime(timer.RealTime());
468 LOG(info) << "CbmKFParticleFinder: Finish SE analysis" << timer.RealTime();
469 }
470}
471
472void CbmKFParticleFinder::FillKFPTrackVector(KFPTrackVector* tracks, const vector<CbmStsTrack>& vRTracks,
473 const vector<KFFieldVector>& vField, const vector<int>& pdg,
474 const vector<int>& trackId, const vector<float>& vChiToPrimVtx,
475 bool atFirstPoint) const
476{
477 int ntracks = vRTracks.size();
478 tracks->Resize(ntracks);
479 //fill vector with tracks
480 for (Int_t iTr = 0; iTr < ntracks; iTr++) {
481 const FairTrackParam* parameters;
482 if (atFirstPoint) {
483 parameters = vRTracks[iTr].GetParamFirst();
484 }
485 else {
486 parameters = vRTracks[iTr].GetParamLast();
487 }
488
489 double par[6] = {0.f};
490
491 double tx = parameters->GetTx(), ty = parameters->GetTy(), qp = parameters->GetQp();
492
493 Int_t q = 0;
494 if (qp > 0.f) {
495 q = 1;
496 }
497 if (qp < 0.f) {
498 q = -1;
499 }
500 if (TMath::Abs(pdg[iTr]) == 1000020030 || TMath::Abs(pdg[iTr]) == 1000020040) {
501 q *= 2;
502 }
503
504
505 double c2 = 1.f / (1.f + tx * tx + ty * ty);
506 double pq = 1.f / qp * TMath::Abs(q);
507 double p2 = pq * pq;
508 double pz = sqrt(p2 * c2);
509 double px = tx * pz;
510 double py = ty * pz;
511
512 par[0] = parameters->GetX();
513 par[1] = parameters->GetY();
514 par[2] = parameters->GetZ();
515 par[3] = px;
516 par[4] = py;
517 par[5] = pz;
518
519 //calculate covariance matrix
520 double t = sqrt(1.f + tx * tx + ty * ty);
521 double t3 = t * t * t;
522 double dpxdtx = q / qp * (1.f + ty * ty) / t3;
523 double dpxdty = -q / qp * tx * ty / t3;
524 double dpxdqp = -q / (qp * qp) * tx / t;
525 double dpydtx = -q / qp * tx * ty / t3;
526 double dpydty = q / qp * (1.f + tx * tx) / t3;
527 double dpydqp = -q / (qp * qp) * ty / t;
528 double dpzdtx = -q / qp * tx / t3;
529 double dpzdty = -q / qp * ty / t3;
530 double dpzdqp = -q / (qp * qp * t3);
531
532 double F[6][5] = {{1.f, 0.f, 0.f, 0.f, 0.f}, {0.f, 1.f, 0.f, 0.f, 0.f},
533 {0.f, 0.f, 0.f, 0.f, 0.f}, {0.f, 0.f, dpxdtx, dpxdty, dpxdqp},
534 {0.f, 0.f, dpydtx, dpydty, dpydqp}, {0.f, 0.f, dpzdtx, dpzdty, dpzdqp}};
535
536 double VFT[5][6];
537 for (int i = 0; i < 5; i++) {
538 for (int j = 0; j < 6; j++) {
539 VFT[i][j] = 0;
540 for (int k = 0; k < 5; k++) {
541 VFT[i][j] += parameters->GetCovariance(i, k) * F[j][k];
542 }
543 }
544 }
545
546 double cov[21];
547 for (int i = 0, l = 0; i < 6; i++) {
548 for (int j = 0; j <= i; j++, l++) {
549 cov[l] = 0;
550 for (int k = 0; k < 5; k++) {
551 cov[l] += F[i][k] * VFT[k][j];
552 }
553 }
554 }
555
556 for (Int_t iP = 0; iP < 6; iP++) {
557 tracks->SetParameter(par[iP], iP, iTr);
558 }
559 for (Int_t iC = 0; iC < 21; iC++) {
560 tracks->SetCovariance(cov[iC], iC, iTr);
561 }
562 for (Int_t iF = 0; iF < 10; iF++) {
563 tracks->SetFieldCoefficient(vField[iTr].fField[iF], iF, iTr);
564 }
565 tracks->SetId(trackId[iTr], iTr);
566 tracks->SetPDG(pdg[iTr], iTr);
567 tracks->SetQ(q, iTr);
568 tracks->SetNPixelHits(vRTracks[iTr].GetNofMvdHits(), iTr);
569
570 if (fPVFindMode == 0 || fPVFindMode == 3) {
571 if (vChiToPrimVtx[iTr] < 3) {
572 tracks->SetPVIndex(0, iTr);
573 }
574 else {
575 tracks->SetPVIndex(-1, iTr);
576 }
577 }
578 else {
579 tracks->SetPVIndex(-1, iTr);
580 }
581 }
582}
583
584double CbmKFParticleFinder::InversedChi2Prob(double p, int ndf) const
585{
586 double epsilon = 1.e-14;
587 double chi2Left = 0.f;
588 double chi2Right = 10000.f;
589
590 double probLeft = p - TMath::Prob(chi2Left, ndf);
591
592 double chi2Centr = (chi2Left + chi2Right) / 2.f;
593 double probCentr = p - TMath::Prob(chi2Centr, ndf);
594
595 while (TMath::Abs(chi2Right - chi2Centr) / chi2Centr > epsilon) {
596 if (probCentr * probLeft > 0.f) {
597 chi2Left = chi2Centr;
598 probLeft = probCentr;
599 }
600 else {
601 chi2Right = chi2Centr;
602 }
603
604 chi2Centr = (chi2Left + chi2Right) / 2.f;
605 probCentr = p - TMath::Prob(chi2Centr, ndf);
606 }
607
608 return chi2Centr;
609}
610
612{
613 fTopoReconstructor->SetChi2PrimaryCut(InversedChi2Prob(prob, 2));
614}
615
617{
619 fPVFindMode = 0;
620 fTopoReconstructor->SetMixedEventAnalysis();
621}
622
623
624void CbmKFParticleFinder::SetTarget(const std::array<float, 3>& target) { fTopoReconstructor->SetTarget(target); }
625KFParticleFinder* CbmKFParticleFinder::GetKFParticleFinder() { return fTopoReconstructor->GetKFParticleFinder(); }
627{
628 GetKFParticleFinder()->SetMaxDistanceBetweenParticlesCut(cut);
629}
630void CbmKFParticleFinder::SetLCut(float cut) { GetKFParticleFinder()->SetLCut(cut); }
631void CbmKFParticleFinder::SetChiPrimaryCut2D(float cut) { GetKFParticleFinder()->SetChiPrimaryCut2D(cut); }
632void CbmKFParticleFinder::SetChi2Cut2D(float cut) { GetKFParticleFinder()->SetChi2Cut2D(cut); }
633void CbmKFParticleFinder::SetLdLCut2D(float cut) { GetKFParticleFinder()->SetLdLCut2D(cut); }
634void CbmKFParticleFinder::SetLdLCutXiOmega(float cut) { GetKFParticleFinder()->SetLdLCutXiOmega(cut); }
635void CbmKFParticleFinder::SetChi2TopoCutXiOmega(float cut) { GetKFParticleFinder()->SetChi2TopoCutXiOmega(cut); }
636void CbmKFParticleFinder::SetChi2CutXiOmega(float cut) { GetKFParticleFinder()->SetChi2CutXiOmega(cut); }
637void CbmKFParticleFinder::SetChi2TopoCutResonances(float cut) { GetKFParticleFinder()->SetChi2TopoCutResonances(cut); }
638void CbmKFParticleFinder::SetChi2CutResonances(float cut) { GetKFParticleFinder()->SetChi2CutResonances(cut); }
639void CbmKFParticleFinder::SetPtCutLMVM(float cut) { GetKFParticleFinder()->SetPtCutLMVM(cut); }
640void CbmKFParticleFinder::SetPCutLMVM(float cut) { GetKFParticleFinder()->SetPCutLMVM(cut); }
641void CbmKFParticleFinder::SetPtCutJPsi(float cut) { GetKFParticleFinder()->SetPtCutJPsi(cut); }
642void CbmKFParticleFinder::SetPtCutCharm(float cut) { GetKFParticleFinder()->SetPtCutCharm(cut); }
643void CbmKFParticleFinder::SetChiPrimaryCutCharm(float cut) { GetKFParticleFinder()->SetChiPrimaryCutCharm(cut); }
645{
646 GetKFParticleFinder()->SetLdLCutCharmManybodyDecays(cut);
647}
649{
650 GetKFParticleFinder()->SetChi2TopoCutCharmManybodyDecays(cut);
651}
653{
654 GetKFParticleFinder()->SetChi2CutCharmManybodyDecays(cut);
655}
656void CbmKFParticleFinder::SetLdLCutCharm2D(float cut) { GetKFParticleFinder()->SetLdLCutCharm2D(cut); }
657void CbmKFParticleFinder::SetChi2TopoCutCharm2D(float cut) { GetKFParticleFinder()->SetChi2TopoCutCharm2D(cut); }
658void CbmKFParticleFinder::SetChi2CutCharm2D(float cut) { GetKFParticleFinder()->SetChi2CutCharm2D(cut); }
660{
661 GetKFParticleFinder()->AddDecayToReconstructionList(pdg);
662}
663
TClonesArray * tracks
Int_t nMCTracks
ClassImp(CbmKFParticleFinder)
friend fvec sqrt(const fvec &a)
friend fscal sgn(fscal x)
static constexpr size_t size()
Definition KfSimdPseudo.h:2
Class characterising one event by a collection of links (indices) to data objects,...
Definition CbmEvent.h:34
int32_t GetNofStsTracks() const
Definition CbmEvent.h:115
const std::vector< int > & GetPID() const
void SetChi2CutCharm2D(float cut)
void SetChi2TopoCutCharmManybodyDecays(float cut)
void SetLdLCutXiOmega(float cut)
TClonesArray * fTrackArray
Name of the input TCA with reco tracks.
CbmMCEventList * fMcEventList
void SetMaxDistanceBetweenParticlesCut(float cut)
CbmMCDataArray * fMcTrackArray
CbmKFParticleFinder(const char *name="CbmKFParticleFinder", Int_t iVerbose=0)
void AddDecayToReconstructionList(int pdg)
double InversedChi2Prob(double p, int ndf) const
virtual InitStatus Init()
KFParticleTopoReconstructor * fTopoReconstructor
void SetChi2CutCharmManybodyDecays(float cut)
void FillKFPTrackVector(KFPTrackVector *tracks, const std::vector< CbmStsTrack > &vRTracks, const std::vector< KFFieldVector > &vField, const std::vector< int > &pdg, const std::vector< int > &trackId, const std::vector< float > &vChiToPrimVtx, bool atFirstPoint=1) const
CbmKFParticleFinderPID * fPID
void SetChi2TopoCutXiOmega(float cut)
void SetLdLCutCharmManybodyDecays(float cut)
void SetChiPrimaryCut2D(float cut)
void SetChiPrimaryCutCharm(float cut)
KFParticleFinder * GetKFParticleFinder()
std::vector< KFFieldVector > fSEField
void SetChi2TopoCutResonances(float cut)
void SetChi2CutResonances(float cut)
virtual void Exec(Option_t *opt)
std::vector< int > fSEpdg
std::vector< CbmStsTrack > fSETracks
void SetChi2TopoCutCharm2D(float cut)
std::vector< int > fSETrackId
void SetLdLCutCharm2D(float cut)
void SetTarget(const std::array< float, 3 > &target)
TClonesArray * fTrackMatchArray
void SetPrimaryProbCut(float prob)
std::vector< float > fSEChiPrim
void SetChi2CutXiOmega(float cut)
Double_t & GetRefZ()
Definition CbmKFVertex.h:27
Double_t & GetRefX()
Definition CbmKFVertex.h:25
Double_t & GetRefY()
Definition CbmKFVertex.h:26
static CbmKF * Instance()
Definition CbmKF.h:40
std::vector< CbmKFTube > vTargets
Definition CbmKF.h:71
void GetChiToVertex(std::vector< CbmStsTrack > &Tracks, std::vector< PFFieldRegion > &field, std::vector< float > &chiToVtx, CbmKFVertex &primVtx, float chiPrim=-1)
void Fit(std::vector< CbmStsTrack > &Tracks, const std::vector< CbmMvdHit > &vMvdHits, const std::vector< CbmStsHit > &vStsHits, const std::vector< int > &pidHypo)
void CalculateFieldRegionAtLastPoint(std::vector< CbmStsTrack > &Tracks, std::vector< PFFieldRegion > &field)
TObject * Get(const CbmLink *lnk)
Int_t Size(Int_t fileNumber, Int_t eventNumber)
Task class creating and managing CbmMCDataArray objects.
CbmMCDataObject * GetObject(const char *name)
CbmMCDataArray * InitBranch(const char *name)
Container class for MC events with number, file and start time.
CbmLink GetEventLinkByIndex(uint32_t index)
Event file and event indices as CbmLink.
Int_t GetEventIndex(UInt_t event, UInt_t file)
Event index.
std::size_t GetNofEvents() const
Number of events in the list @value Number of events.
const CbmLink & GetLink(int32_t i) const
Definition CbmMatch.h:39
int32_t GetNofLinks() const
Definition CbmMatch.h:42
virtual int32_t GetTotalNofHits() const
Definition CbmStsTrack.h:81
int32_t GetNDF() const
Definition CbmTrack.h:64
const FairTrackParam * GetParamFirst() const
Definition CbmTrack.h:68
double GetChiSq() const
Definition CbmTrack.h:63
double GetZ() const
Definition CbmVertex.h:69
double GetY() const
Definition CbmVertex.h:68
double GetX() const
Definition CbmVertex.h:67