Belle II Software development
V0Fitter.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8#include <tracking/v0Finding/fitter/V0Fitter.h>
9
10#include <framework/logging/Logger.h>
11#include <framework/datastore/StoreArray.h>
12#include <framework/geometry/VectorUtil.h>
13#include <framework/geometry/BFieldManager.h>
14#include <tracking/v0Finding/dataobjects/VertexVector.h>
15#include <tracking/dataobjects/RecoTrack.h>
16#include <tracking/trackFitting/fitter/base/TrackFitter.h>
17#include <tracking/trackFitting/trackBuilder/factories/TrackBuilder.h>
18
19#include <mdst/dataobjects/HitPatternVXD.h>
20#include <mdst/dataobjects/HitPatternCDC.h>
21#include <mdst/dataobjects/Track.h>
22#include <mdst/dataobjects/TrackFitResult.h>
23#include <cdc/dataobjects/CDCRecoHit.h>
24#include <pxd/reconstruction/PXDRecoHit.h>
25#include <svd/reconstruction/SVDRecoHit.h>
26#include <svd/reconstruction/SVDRecoHit2D.h>
27
28#include <genfit/Track.h>
29#include <genfit/TrackPoint.h>
30#include <genfit/MeasuredStateOnPlane.h>
31#include <genfit/GFRaveVertexFactory.h>
32#include <genfit/GFRaveVertex.h>
33#include <genfit/FieldManager.h>
34#include <genfit/MaterialEffects.h>
35#include <genfit/FitStatus.h>
36#include <genfit/KalmanFitStatus.h>
37#include <genfit/KalmanFitterInfo.h>
38
39#include <framework/utilities/IOIntercept.h>
40
41using namespace Belle2;
42
43V0Fitter::V0Fitter(const std::string& trackFitResultsName, const std::string& v0sName,
44 const std::string& v0ValidationVerticesName, const std::string& recoTracksName,
45 const std::string& copiedRecoTracksName, bool enableValidation)
46 : m_validation(enableValidation), m_recoTracksName(recoTracksName), m_v0FitterMode(1), m_forcestore(false),
47 m_useOnlyOneSVDHitPair(true)
48{
49 m_trackFitResults.isRequired(trackFitResultsName);
51 //Relation to RecoTracks from Tracks is already tested at the module level.
52
53 if (m_validation) {
54 B2DEBUG(24, "Register DataStore for validation.");
55 m_validationV0s.registerInDataStore(v0ValidationVerticesName, DataStore::c_ErrorIfAlreadyRegistered);
56 m_v0s.registerRelationTo(m_validationV0s);
57 }
58
61
63 m_copiedRecoTracks.registerInDataStore(copiedRecoTracksName,
66
69
70
71 B2ASSERT("Material effects not set up. Please use SetupGenfitExtrapolationModule.",
72 genfit::MaterialEffects::getInstance()->isInitialized());
73 B2ASSERT("Magnetic field not set up. Please use SetupGenfitExtrapolationModule.",
74 genfit::FieldManager::getInstance()->isInitialized());
75}
76
77void V0Fitter::setFitterMode(int fitterMode)
78{
79 if (not(0 <= fitterMode && fitterMode <= 2)) {
80 B2FATAL("Invalid fitter mode!");
81 } else {
82 m_v0FitterMode = fitterMode;
83 if (fitterMode == 0)
84 m_forcestore = true;
85 else
86 m_forcestore = false;
87 if (fitterMode == 2)
89 else
91 }
92}
93
94void V0Fitter::initializeCuts(double beamPipeRadius,
95 double vertexChi2CutOutside,
96 std::tuple<double, double> invMassRangeKshort,
97 std::tuple<double, double> invMassRangeLambda,
98 std::tuple<double, double> invMassRangePhoton)
99{
100 m_beamPipeRadius = beamPipeRadius;
101 m_vertexChi2CutOutside = vertexChi2CutOutside;
102 m_invMassRangeKshort = invMassRangeKshort;
103 m_invMassRangeLambda = invMassRangeLambda;
104 m_invMassRangePhoton = invMassRangePhoton;
105}
106
107
108bool V0Fitter::fitGFRaveVertex(genfit::Track& trackPlus, genfit::Track& trackMinus, genfit::GFRaveVertex& vertex)
109{
110 VertexVector vertexVector;
111 std::vector<genfit::Track*> trackPair {&trackPlus, &trackMinus};
112
113 try {
115 logCapture("V0Fitter GFRaveVertexFactory", LogConfig::c_Debug, LogConfig::c_Debug);
116 logCapture.start();
117
118 genfit::GFRaveVertexFactory vertexFactory;
119 vertexFactory.findVertices(&vertexVector.v, trackPair);
120 } catch (...) {
121 B2ERROR("Exception during vertex fit.");
122 return false;
123 }
124
125 if (vertexVector.size() != 1) {
126 B2DEBUG(21, "Vertex fit failed. Size of vertexVector not 1, but: " << vertexVector.size());
127 return false;
128 }
129
130 if ((*vertexVector[0]).getNTracks() != 2) {
131 B2DEBUG(20, "Wrong number of tracks in vertex.");
132 return false;
133 }
134
135 vertex = *vertexVector[0];
136 return true;
137}
138
141bool V0Fitter::extrapolateToVertex(genfit::MeasuredStateOnPlane& stPlus, genfit::MeasuredStateOnPlane& stMinus,
142 const ROOT::Math::XYZVector& vertexPosition, unsigned int& hasInnerHitStatus)
143{
145 hasInnerHitStatus = 0;
146
147 try {
150 double extralengthPlus = stPlus.extrapolateToPoint(XYZToTVector(vertexPosition));
151 double extralengthMinus = stMinus.extrapolateToPoint(XYZToTVector(vertexPosition));
152 if (extralengthPlus > 0) hasInnerHitStatus |= 0x1;
153 if (extralengthMinus > 0) hasInnerHitStatus |= 0x2;
154 B2DEBUG(22, "extralengthPlus=" << extralengthPlus << ", extralengthMinus=" << extralengthMinus);
155 } catch (...) {
159 B2DEBUG(22, "Could not extrapolate track to vertex.");
160 return false;
161 }
162 return true;
163}
164
165TrackFitResult* V0Fitter::buildTrackFitResult(const genfit::Track& track, const RecoTrack* recoTrack,
166 const genfit::MeasuredStateOnPlane& msop, const double Bz,
167 const Const::ParticleType& trackHypothesis,
168 const int sharedInnermostCluster)
169{
170 const uint64_t hitPatternCDCInitializer = TrackBuilder::getHitPatternCDCInitializer(*recoTrack);
171 uint32_t hitPatternVXDInitializer = TrackBuilder::getHitPatternVXDInitializer(*recoTrack);
172
173 // If the innermost hit is shared among V0 daughters, assign flag in the infoLayer.
174 if (sharedInnermostCluster > 0 && sharedInnermostCluster < 4) {
175 HitPatternVXD hitPatternVXD_forflag = HitPatternVXD(hitPatternVXDInitializer);
176 hitPatternVXD_forflag.setInnermostHitShareStatus(sharedInnermostCluster);
177 hitPatternVXDInitializer = hitPatternVXD_forflag.getInteger();
178 }
179
180 TrackFitResult* v0TrackFitResult
181 = m_trackFitResults.appendNew(ROOT::Math::XYZVector(msop.getPos()), ROOT::Math::XYZVector(msop.getMom()),
182 msop.get6DCov(), msop.getCharge(),
183 trackHypothesis,
184 track.getFitStatus()->getPVal(),
185 Bz, hitPatternCDCInitializer, hitPatternVXDInitializer, track.getFitStatus()->getNdf());
186 return v0TrackFitResult;
187}
188
189std::pair<Const::ParticleType, Const::ParticleType> V0Fitter::getTrackHypotheses(const Const::ParticleType& v0Hypothesis) const
190{
191 if (v0Hypothesis == Const::Kshort) {
192 return std::make_pair(Const::pion, Const::pion);
193 } else if (v0Hypothesis == Const::photon) {
194 return std::make_pair(Const::electron, Const::electron);
195 } else if (v0Hypothesis == Const::Lambda) {
196 return std::make_pair(Const::proton, Const::pion);
197 } else if (v0Hypothesis == Const::antiLambda) {
198 return std::make_pair(Const::pion, Const::proton);
199 }
200 B2FATAL("Given V0Hypothesis not available.");
201 return std::make_pair(Const::invalidParticle, Const::invalidParticle); // return something to avoid triggering cppcheck
202}
203
205bool V0Fitter::fitAndStore(const Track* trackPlus, const Track* trackMinus,
206 const Const::ParticleType& v0Hypothesis,
207 bool& isForceStored, bool& isHitRemoved)
208{
210 isForceStored = false;
211 isHitRemoved = false;
212
214 RecoTrack* recoTrackPlus = trackPlus->getRelated<RecoTrack>(m_recoTracksName);
215 RecoTrack* recoTrackMinus = trackMinus->getRelated<RecoTrack>(m_recoTracksName);
216
220 unsigned int hasInnerHitStatus = 0;
221
223 ROOT::Math::XYZVector vertexPos(0, 0, 0);
224
230 if (not vertexFitWithRecoTracks(trackPlus, trackMinus, recoTrackPlus, recoTrackMinus, v0Hypothesis,
231 hasInnerHitStatus, vertexPos, m_forcestore))
232 return false;
233
234 if ((hasInnerHitStatus == 0) or m_forcestore)
235 return true;
236
241 bool failflag = false;
243 const auto trackHypotheses = getTrackHypotheses(v0Hypothesis);
244 const int pdg_trackPlus = trackPlus->getTrackFitResultWithClosestMass(
245 trackHypotheses.first)->getParticleType().getPDGCode();
246 const int pdg_trackMinus = trackMinus->getTrackFitResultWithClosestMass(
247 trackHypotheses.second)->getParticleType().getPDGCode();
248
249 RecoTrack* recoTrackPlus_forRefit = nullptr;
250 RecoTrack* recoTrackMinus_forRefit = nullptr;
251 RecoTrack* cache_recoTrackPlus = nullptr;
252 RecoTrack* cache_recoTrackMinus = nullptr;
253
254 unsigned int count_removeInnerHits = 0;
255 while (hasInnerHitStatus != 0) {
256 ++count_removeInnerHits;
259
261 if (hasInnerHitStatus & 0x1) {
262 // create a copy of the original RecoTrack w/o track fit
263 recoTrackPlus_forRefit = copyRecoTrack(recoTrackPlus);
264 if (recoTrackPlus_forRefit == nullptr)
265 return false;
268 if (!cache_recoTrackPlus) {
269 if (not removeInnerHits(recoTrackPlus, recoTrackPlus_forRefit, pdg_trackPlus, vertexPos)) {
270 failflag = true;
271 break;
272 }
273 } else {
274 if (not removeInnerHits(cache_recoTrackPlus, recoTrackPlus_forRefit, pdg_trackPlus, vertexPos)) {
275 failflag = true;
276 break;
277 }
278 }
279 cache_recoTrackPlus = recoTrackPlus_forRefit;
280 } else if (recoTrackPlus_forRefit == nullptr) {
281 // create a copy of the original RecoTrack w/ track fit (only once)
282 recoTrackPlus_forRefit = copyRecoTrackAndFit(recoTrackPlus, pdg_trackPlus);
283 if (recoTrackPlus_forRefit == nullptr)
284 return false;
285 }
286
288 if (hasInnerHitStatus & 0x2) {
289 // create a copy of the original RecoTrack w/o track fit
290 recoTrackMinus_forRefit = copyRecoTrack(recoTrackMinus);
291 if (recoTrackMinus_forRefit == nullptr)
292 return false;
295 if (!cache_recoTrackMinus) {
296 if (not removeInnerHits(recoTrackMinus, recoTrackMinus_forRefit, pdg_trackMinus, vertexPos)) {
297 failflag = true;
298 break;
299 }
300 } else {
301 if (not removeInnerHits(cache_recoTrackMinus, recoTrackMinus_forRefit, pdg_trackMinus, vertexPos)) {
302 failflag = true;
303 break;
304 }
305 }
306 cache_recoTrackMinus = recoTrackMinus_forRefit;
307 } else if (recoTrackMinus_forRefit == nullptr) {
308 // create a copy of the original RecoTrack w/ track fit (only once)
309 recoTrackMinus_forRefit = copyRecoTrackAndFit(recoTrackMinus, pdg_trackMinus);
310 if (recoTrackMinus_forRefit == nullptr)
311 return false;
312 }
313
315 hasInnerHitStatus = 0;
318 if (not vertexFitWithRecoTracks(trackPlus, trackMinus, recoTrackPlus_forRefit, recoTrackMinus_forRefit,
319 v0Hypothesis, hasInnerHitStatus, vertexPos, false)) {
320 B2DEBUG(22, "Vertex refit failed, or rejected by invariant mass cut.");
321 failflag = true;
322 break;
323 } else if (hasInnerHitStatus == 0)
324 isHitRemoved = true;
325 if (count_removeInnerHits >= 5) {
326 B2WARNING("Inner hits remained after " << count_removeInnerHits << " times of removing inner hits!");
327 failflag = true;
328 break;
329 }
330 }
331
333 if (failflag) {
334 bool forcestore = true;
335 if (not vertexFitWithRecoTracks(trackPlus, trackMinus, recoTrackPlus, recoTrackMinus, v0Hypothesis,
336 hasInnerHitStatus, vertexPos, forcestore)) {
337 B2DEBUG(22, "Original vertex fit fails. Possibly rejected by invariant mass cut.");
338 return false;
339 } else
340 isForceStored = true;
341 }
342
343 return true;
344}
345
348bool V0Fitter::vertexFitWithRecoTracks(const Track* trackPlus, const Track* trackMinus,
349 RecoTrack* recoTrackPlus, RecoTrack* recoTrackMinus,
350 const Const::ParticleType& v0Hypothesis,
351 unsigned int& hasInnerHitStatus, ROOT::Math::XYZVector& vertexPos,
352 const bool forceStore)
353{
354 const auto trackHypotheses = getTrackHypotheses(v0Hypothesis);
355
356 // make a clone, not use the reference so that the genfit::Track and its TrackReps will not be altered.
357 genfit::Track gfTrackPlus = RecoTrackGenfitAccess::getGenfitTrack(*recoTrackPlus);
358 const int pdgTrackPlus = trackPlus->getTrackFitResultWithClosestMass(trackHypotheses.first)->getParticleType().getPDGCode();
359 genfit::AbsTrackRep* plusRepresentation = recoTrackPlus->getTrackRepresentationForPDG(pdgTrackPlus);
360 if ((plusRepresentation == nullptr) or (not recoTrackPlus->wasFitSuccessful(plusRepresentation))) {
361 B2ERROR("Track hypothesis with closest mass not available. Should never happen, but I can continue savely anyway.");
362 return false;
363 }
364
365 // make a clone, not use the reference so that the genfit::Track and its TrackReps will not be altered.
366 genfit::Track gfTrackMinus = RecoTrackGenfitAccess::getGenfitTrack(*recoTrackMinus);
367 const int pdgTrackMinus = trackMinus->getTrackFitResultWithClosestMass(trackHypotheses.second)->getParticleType().getPDGCode();
368 genfit::AbsTrackRep* minusRepresentation = recoTrackMinus->getTrackRepresentationForPDG(pdgTrackMinus);
369 if ((minusRepresentation == nullptr) or (not recoTrackMinus->wasFitSuccessful(minusRepresentation))) {
370 B2ERROR("Track hypothesis with closest mass not available. Should never happen, but I can continue savely anyway.");
371 return false;
372 }
373
375 std::vector<genfit::AbsTrackRep*> repsPlus = gfTrackPlus.getTrackReps();
376 std::vector<genfit::AbsTrackRep*> repsMinus = gfTrackMinus.getTrackReps();
377 if (repsPlus.size() == repsMinus.size()) {
378 for (unsigned int id = 0; id < repsPlus.size(); id++) {
379 if (abs(repsPlus[id]->getPDG()) == pdgTrackPlus)
380 gfTrackPlus.setCardinalRep(id);
381 if (abs(repsMinus[id]->getPDG()) == pdgTrackMinus)
382 gfTrackMinus.setCardinalRep(id);
383 }
384 }
386 else {
387 for (unsigned int id = 0; id < repsPlus.size(); id++) {
388 if (abs(repsPlus[id]->getPDG()) == pdgTrackPlus)
389 gfTrackPlus.setCardinalRep(id);
390 }
391 for (unsigned int id = 0; id < repsMinus.size(); id++) {
392 if (abs(repsMinus[id]->getPDG()) == pdgTrackMinus)
393 gfTrackMinus.setCardinalRep(id);
394 }
395 }
396
398 genfit::MeasuredStateOnPlane stPlus = recoTrackPlus->getMeasuredStateOnPlaneFromFirstHit(plusRepresentation);
399 genfit::MeasuredStateOnPlane stMinus = recoTrackMinus->getMeasuredStateOnPlaneFromFirstHit(minusRepresentation);
400
401 genfit::GFRaveVertex vert;
402 if (not fitGFRaveVertex(gfTrackPlus, gfTrackMinus, vert)) {
403 return false;
404 }
405
406 const ROOT::Math::XYZVector& posVert = ROOT::Math::XYZVector(vert.getPos());
407 vertexPos = posVert;
408
411 if (posVert.Rho() < m_beamPipeRadius) {
412 return false;
413 } else {
414 if (vert.getChi2() > m_vertexChi2CutOutside) {
415 B2DEBUG(22, "Vertex outside beam pipe, chi^2 too large.");
416 return false;
417 }
418 }
419
420 B2DEBUG(22, "Vertex accepted.");
421
422 if (not extrapolateToVertex(stPlus, stMinus, posVert, hasInnerHitStatus)) {
423 return false;
424 }
425
427 if (forceStore || hasInnerHitStatus == 0) {
429 const genfit::GFRaveTrackParameters* tr0 = vert.getParameters(0);
430 const genfit::GFRaveTrackParameters* tr1 = vert.getParameters(1);
431 ROOT::Math::PxPyPzMVector lv0(tr0->getMom().Px(), tr0->getMom().Py(), tr0->getMom().Pz(), trackHypotheses.first.getMass());
432 ROOT::Math::PxPyPzMVector lv1(tr1->getMom().Px(), tr1->getMom().Py(), tr1->getMom().Pz(), trackHypotheses.second.getMass());
434 double v0InvMass = (lv0 + lv1).M();
435 if (v0Hypothesis == Const::Kshort) {
436 if (v0InvMass < std::get<0>(m_invMassRangeKshort) || v0InvMass > std::get<1>(m_invMassRangeKshort)) {
437 B2DEBUG(22, "Kshort vertex rejected, invariant mass out of range.");
438 return false;
439 }
440 } else if (v0Hypothesis == Const::Lambda || v0Hypothesis == Const::antiLambda) {
441 if (v0InvMass < std::get<0>(m_invMassRangeLambda) || v0InvMass > std::get<1>(m_invMassRangeLambda)) {
442 B2DEBUG(22, "Lambda vertex rejected, invariant mass out of range.");
443 return false;
444 }
445 } else if (v0Hypothesis == Const::photon) {
446 if (v0InvMass < std::get<0>(m_invMassRangePhoton) || v0InvMass > std::get<1>(m_invMassRangePhoton)) {
447 B2DEBUG(22, "Photon vertex rejected, invariant mass out of range.");
448 return false;
449 }
450 }
451
452
456 const double Bz = BFieldManager::getFieldInTesla({0, 0, 0}).Z();
457
458 int sharedInnermostCluster = checkSharedInnermostCluster(recoTrackPlus, recoTrackMinus);
459
460 TrackFitResult* tfrPlusVtx = buildTrackFitResult(gfTrackPlus, recoTrackPlus, stPlus, Bz, trackHypotheses.first,
461 sharedInnermostCluster);
462 TrackFitResult* tfrMinusVtx = buildTrackFitResult(gfTrackMinus, recoTrackMinus, stMinus, Bz, trackHypotheses.second,
463 sharedInnermostCluster);
464
465 B2DEBUG(20, "Creating new V0.");
466 auto v0 = m_v0s.appendNew(std::make_pair(trackPlus, tfrPlusVtx),
467 std::make_pair(trackMinus, tfrMinusVtx),
468 posVert.X(), posVert.Y(), posVert.Z());
469
470 if (m_validation) {
471 B2DEBUG(24, "Create StoreArray and Output for validation.");
472 auto validationV0 = m_validationV0s.appendNew(
473 std::make_pair(trackPlus, tfrPlusVtx),
474 std::make_pair(trackMinus, tfrMinusVtx),
475 ROOT::Math::XYZVector(vert.getPos()),
476 vert.getCov(),
477 (lv0 + lv1).P(),
478 (lv0 + lv1).M(),
479 vert.getChi2()
480 );
481 v0->addRelationTo(validationV0);
482 }
483 }
484 return true;
485}
486
488{
489 RecoTrack* newRecoTrack = origRecoTrack->copyToStoreArray(m_copiedRecoTracks);
490 newRecoTrack->addHitsFromRecoTrack(origRecoTrack);
491 newRecoTrack->addRelationTo(origRecoTrack);
492 return newRecoTrack;
493}
494
495RecoTrack* V0Fitter::copyRecoTrackAndFit(RecoTrack* origRecoTrack, const int trackPDG)
496{
498 Const::ChargedStable particleUsedForFitting(std::abs(trackPDG));
499 const genfit::AbsTrackRep* origTrackRep = origRecoTrack->getTrackRepresentationForPDG(std::abs(
500 trackPDG));
501
502 RecoTrack* newRecoTrack = copyRecoTrack(origRecoTrack);
503
505 TrackFitter fitter;
506 if (not fitter.fit(*newRecoTrack, particleUsedForFitting)) {
507 // This is not expected, but happens sometimes.
508 B2DEBUG(20, "track fit failed for copied RecoTrack.");
510 if (not origRecoTrack->wasFitSuccessful(origTrackRep))
511 B2DEBUG(20, "\t original track fit was also failed.");
512 return nullptr;
513 }
514
515 return newRecoTrack;
516}
517
518bool V0Fitter::removeInnerHits(RecoTrack* prevRecoTrack, RecoTrack* recoTrack,
519 const int trackPDG, const ROOT::Math::XYZVector& vertexPosition)
520{
521 if (!prevRecoTrack || !recoTrack) {
522 B2ERROR("Input recotrack is nullptr!");
523 return false;
524 }
526 Const::ChargedStable particleUsedForFitting(std::abs(trackPDG));
527 const genfit::AbsTrackRep* prevTrackRep = prevRecoTrack->getTrackRepresentationForPDG(std::abs(
528 trackPDG));
529
531 const std::vector<RecoHitInformation*>& recoHitInformations = recoTrack->getRecoHitInformations(true);
532 const std::vector<RecoHitInformation*>& prevRecoHitInformations = prevRecoTrack->getRecoHitInformations(
533 true);
534 unsigned int nRemoveHits = 0;
535 if (recoHitInformations.size() != prevRecoHitInformations.size()) {
536 B2WARNING("Copied RecoTrack has different number of hits from its original RecoTrack!");
537 return false;
538 }
540 for (nRemoveHits = 0; nRemoveHits < recoHitInformations.size(); ++nRemoveHits) {
541 if (!prevRecoHitInformations[nRemoveHits]->useInFit()) {
542 recoHitInformations[nRemoveHits]->setUseInFit(false);
543 continue;
544 }
545 try {
547 genfit::MeasuredStateOnPlane stPrevRecoHit = prevRecoTrack->getMeasuredStateOnPlaneFromRecoHit(
548 prevRecoHitInformations[nRemoveHits]);
551 double extralength = stPrevRecoHit.extrapolateToPoint(XYZToTVector(vertexPosition));
552 if (extralength > 0) {
553 recoHitInformations[nRemoveHits]->setUseInFit(false);
554 } else
555 break;
556 } catch (NoTrackFitResult()) {
557 B2WARNING("Exception: no FitterInfo assigned for TrackPoint created from this RecoHit.");
558 recoHitInformations[nRemoveHits]->setUseInFit(false);
559 continue;
560 } catch (...) {
564 B2DEBUG(22, "Could not extrapolate track to vertex when removing inner hits, aborting.");
565 return false;
566 }
567 }
568
569 if (nRemoveHits == 0) {
570 // This is not expected, but can happen if the track fit is different between copied and original RecoTrack.
571 B2DEBUG(20, "No hits removed in removeInnerHits, aborted. Switching to use the original RecoTrack.");
572 return false;
573 }
574
575 if (recoHitInformations.size() <= nRemoveHits) {
576 B2DEBUG(20, "Removed all the RecoHits in the RecoTrack, aborted. Switching to use the original RecoTrack.");
577 return false;
578 }
579
582 if (recoHitInformations[nRemoveHits - 1]->getTrackingDetector() == RecoHitInformation::RecoHitDetector::c_SVD) {
583 if (recoHitInformations[nRemoveHits]->getTrackingDetector() == RecoHitInformation::RecoHitDetector::c_SVD) {
584 const SVDCluster* lastRemovedSVDHit = recoHitInformations[nRemoveHits - 1]->getRelatedTo<SVDCluster>();
585 const SVDCluster* nextSVDHit = recoHitInformations[nRemoveHits] ->getRelatedTo<SVDCluster>();
586 if (!lastRemovedSVDHit || !nextSVDHit) B2ERROR("Last/Next SVD hit is null!");
587 else {
588 if (lastRemovedSVDHit->getSensorID() == nextSVDHit->getSensorID() &&
589 lastRemovedSVDHit->isUCluster() && !(nextSVDHit->isUCluster())) {
590 recoHitInformations[nRemoveHits]->setUseInFit(false);
591 ++nRemoveHits;
592 }
593 }
594 }
595 }
596
599 recoHitInformations[nRemoveHits - 1]->getTrackingDetector() == RecoHitInformation::RecoHitDetector::c_SVD &&
600 recoHitInformations.size() > nRemoveHits + 2) {
601 if (recoHitInformations[nRemoveHits ]->getTrackingDetector() == RecoHitInformation::RecoHitDetector::c_SVD &&
602 recoHitInformations[nRemoveHits + 1]->getTrackingDetector() == RecoHitInformation::RecoHitDetector::c_SVD &&
603 recoHitInformations[nRemoveHits + 2]->getTrackingDetector() != RecoHitInformation::RecoHitDetector::c_SVD) {
604 recoHitInformations[nRemoveHits ]->setUseInFit(false);
605 recoHitInformations[nRemoveHits + 1]->setUseInFit(false);
606 nRemoveHits += 2;
607 }
608 }
609
610 B2DEBUG(22, nRemoveHits << " inner hits removed.");
611
613 TrackFitter fitter;
614 if (not fitter.fit(*recoTrack, particleUsedForFitting)) {
615 B2DEBUG(20, "track fit failed after removing inner hits.");
617 if (not prevRecoTrack->wasFitSuccessful(prevTrackRep))
618 B2DEBUG(20, "\t previous track fit was also failed.");
619 return false;
620 }
621
622 return true;
623}
624
625int V0Fitter::checkSharedInnermostCluster(const RecoTrack* recoTrackPlus, const RecoTrack* recoTrackMinus)
626{
629 int flag = 0; // -1 for exception
630
631 // get the innermost hit for plus/minus-daughter
632 const std::vector<RecoHitInformation*>& recoHitInformationsPlus = recoTrackPlus->getRecoHitInformations(
633 true); // true for sorted info.
634 const std::vector<RecoHitInformation*>& recoHitInformationsMinus = recoTrackMinus->getRecoHitInformations(
635 true);// true for sorted info.
636 unsigned int iInnermostHitPlus, iInnermostHitMinus;
637 for (iInnermostHitPlus = 0 ; iInnermostHitPlus < recoHitInformationsPlus.size() ; ++iInnermostHitPlus)
638 if (recoHitInformationsPlus[iInnermostHitPlus]->useInFit()) break;
639 for (iInnermostHitMinus = 0 ; iInnermostHitMinus < recoHitInformationsMinus.size() ; ++iInnermostHitMinus)
640 if (recoHitInformationsMinus[iInnermostHitMinus]->useInFit()) break;
641 if (iInnermostHitPlus == recoHitInformationsPlus.size() || iInnermostHitMinus == recoHitInformationsMinus.size()) {
642 B2WARNING("checkSharedInnermostCluster function called for recoTrack including no hit used for fit! This should not happen!");
643 return -1;
644 }
645 const auto& recoHitInfoPlus = recoHitInformationsPlus[iInnermostHitPlus];
646 const auto& recoHitInfoMinus = recoHitInformationsMinus[iInnermostHitMinus];
647
648 if (recoHitInfoPlus->getTrackingDetector() == recoHitInfoMinus->getTrackingDetector()) {
649 if (recoHitInfoPlus->getTrackingDetector() == RecoHitInformation::c_PXD) {
650 const PXDCluster* clusterPlus = recoHitInfoPlus->getRelatedTo<PXDCluster>();
651 const PXDCluster* clusterMinus = recoHitInfoMinus->getRelatedTo<PXDCluster>();
652 if (clusterPlus == clusterMinus) { // if they share a same PXDCluster, set the flag
653 flag = 3; // PXD cluster is a 2D-hit
654 }
655 } else if (recoHitInfoPlus->getTrackingDetector() == RecoHitInformation::c_SVD) {
658 const SVDCluster* clusterPlus = recoHitInfoPlus->getRelatedTo<SVDCluster>();
659 const SVDCluster* clusterMinus = recoHitInfoMinus->getRelatedTo<SVDCluster>();
660 if (clusterPlus->isUCluster() && clusterMinus->isUCluster()) {
661 if (recoHitInformationsPlus.size() > iInnermostHitPlus + 1
662 && recoHitInformationsMinus.size() > iInnermostHitMinus + 1) { // not to access an array out of boundary
663 const auto& recoHitInfoNextPlus = recoHitInformationsPlus[iInnermostHitPlus + 1];
664 const auto& recoHitInfoNextMinus = recoHitInformationsMinus[iInnermostHitMinus + 1];
665 // sanity check to access next hits
666 if (recoHitInfoNextPlus->useInFit() && recoHitInfoNextMinus->useInFit() // this should be satisfied by default
667 && recoHitInfoNextPlus->getTrackingDetector() == RecoHitInformation::c_SVD
668 && recoHitInfoNextMinus->getTrackingDetector() == RecoHitInformation::c_SVD) {
669 const SVDCluster* clusterNextPlus = recoHitInfoNextPlus->getRelatedTo<SVDCluster>();
670 const SVDCluster* clusterNextMinus = recoHitInfoNextMinus->getRelatedTo<SVDCluster>();
671 if (!(clusterNextPlus->isUCluster()) && !(clusterNextMinus->isUCluster())
672 && clusterPlus->getSensorID() == clusterNextPlus->getSensorID()
673 && clusterMinus->getSensorID() == clusterNextMinus->getSensorID()) {
674 if (clusterPlus == clusterMinus)
675 flag += 2; // SVD U-cluster is shared
676 if (clusterNextPlus == clusterNextMinus)
677 flag += 1; // SVD V-cluster is shared
678 } else {
679 B2WARNING("SVD cluster to be paired is not on V-side, or not on the same sensor.");
680 return -1;
681 }
682 } else {
683 B2WARNING("No SVD cluster to be paired.");
684 return -1;
685 }
686 } else {
687 B2WARNING("Innermost SVD U-cluster is the only hit in a daughter track. This should not happen.");
688 return -1;
689 }
690 } else {
691 B2WARNING("No SVD U-cluster in the innermost cluster.");
692 return -1;
693 }
694 }
695 }
696 return flag;
697
698}
static ROOT::Math::XYZVector getFieldInTesla(const ROOT::Math::XYZVector &pos)
return the magnetic field at a given position in Tesla.
Definition: BFieldManager.h:61
Provides a type-safe way to pass members of the chargedStableSet set.
Definition: Const.h:589
The ParticleType class for identifying different particle types.
Definition: Const.h:408
int getPDGCode() const
PDG code.
Definition: Const.h:473
static const ParticleType Lambda
Lambda particle.
Definition: Const.h:679
static const ChargedStable pion
charged pion particle
Definition: Const.h:661
static const ParticleType antiLambda
Anti-Lambda particle.
Definition: Const.h:680
static const ChargedStable proton
proton particle
Definition: Const.h:663
static const ParticleType invalidParticle
Invalid particle, used internally.
Definition: Const.h:681
static const ParticleType Kshort
K^0_S particle.
Definition: Const.h:677
static const ParticleType photon
photon particle
Definition: Const.h:673
static const ChargedStable electron
electron particle
Definition: Const.h:659
@ c_WriteOut
Object/array should be saved by output modules.
Definition: DataStore.h:70
@ c_ErrorIfAlreadyRegistered
If the object/array was already registered, produce an error (aborting initialisation).
Definition: DataStore.h:72
Hit pattern of the VXD within a track.
Definition: HitPatternVXD.h:37
unsigned int getInteger() const
Getter for the underlying integer.
Definition: HitPatternVXD.h:58
void setInnermostHitShareStatus(const unsigned short innermostHitShareStatus)
Set the innermost hit share flags for V0 daughters.
bool start()
Start intercepting the output.
Definition: IOIntercept.h:155
Capture stdout and stderr and convert into log messages.
Definition: IOIntercept.h:226
@ c_Debug
Debug: for code development.
Definition: LogConfig.h:26
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
Definition: PXDCluster.h:30
static genfit::Track & getGenfitTrack(RecoTrack &recoTrack)
Give access to the RecoTrack's genfit::Track.
Definition: RecoTrack.cc:404
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
size_t addHitsFromRecoTrack(const RecoTrack *recoTrack, unsigned int sortingParameterOffset=0, bool reversed=false, std::optional< double > optionalMinimalWeight=std::nullopt)
Add all hits from another RecoTrack to this RecoTrack.
Definition: RecoTrack.cc:240
bool wasFitSuccessful(const genfit::AbsTrackRep *representation=nullptr) const
Returns true if the last fit with the given representation was successful.
Definition: RecoTrack.cc:336
RecoTrack * copyToStoreArray(StoreArray< RecoTrack > &storeArray) const
Append a new RecoTrack to the given store array and copy its general properties, but not the hits the...
Definition: RecoTrack.cc:529
std::vector< RecoHitInformation * > getRecoHitInformations(bool getSorted=false) const
Return a list of all RecoHitInformations associated with the RecoTrack.
Definition: RecoTrack.cc:557
static void registerRequiredRelations(StoreArray< RecoTrack > &recoTracks, std::string const &pxdHitsStoreArrayName="", std::string const &svdHitsStoreArrayName="", std::string const &cdcHitsStoreArrayName="", std::string const &bklmHitsStoreArrayName="", std::string const &eklmHitsStoreArrayName="", std::string const &recoHitInformationStoreArrayName="")
Convenience method which registers all relations required to fully use a RecoTrack.
Definition: RecoTrack.cc:53
const genfit::MeasuredStateOnPlane & getMeasuredStateOnPlaneFromFirstHit(const genfit::AbsTrackRep *representation=nullptr) const
Return genfit's MeasuredStateOnPlane for the first hit in a fit useful for extrapolation of measureme...
Definition: RecoTrack.cc:605
genfit::AbsTrackRep * getTrackRepresentationForPDG(int pdgCode) const
Return an already created track representation of the given reco track for the PDG.
Definition: RecoTrack.cc:475
const genfit::MeasuredStateOnPlane & getMeasuredStateOnPlaneFromRecoHit(const RecoHitInformation *recoHitInfo, const genfit::AbsTrackRep *representation=nullptr) const
Return genfit's MeasuredStateOnPlane on plane for associated with one RecoHitInformation.
Definition: RecoTrack.cc:579
void addRelationTo(const RelationsInterface< BASE > *object, float weight=1.0, const std::string &namedRelation="") const
Add a relation from this object to another object (with caching).
TO * getRelatedTo(const std::string &name="", const std::string &namedRelation="") const
Get the object to which this object has a relation.
T * getRelated(const std::string &name="", const std::string &namedRelation="") const
Get the object to or from which this object has a relation.
The SVD Cluster class This class stores all information about reconstructed SVD clusters.
Definition: SVDCluster.h:29
VxdID getSensorID() const
Get the sensor ID.
Definition: SVDCluster.h:102
bool isUCluster() const
Get the direction of strips.
Definition: SVDCluster.h:110
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.
bool registerRelationTo(const StoreArray< TO > &toArray, DataStore::EDurability durability=DataStore::c_Event, DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut, const std::string &namedRelation="") const
Register a relation to the given StoreArray.
Definition: StoreArray.h:140
static uint32_t getHitPatternVXDInitializer(const RecoTrack &recoTrack, const genfit::AbsTrackRep *representation=nullptr)
Get the HitPattern in the VXD.
static uint64_t getHitPatternCDCInitializer(const RecoTrack &recoTrack, const genfit::AbsTrackRep *representation=nullptr)
Get the HitPattern in the CDC.
Values of the result of a track fit with a given particle hypothesis.
Const::ParticleType getParticleType() const
Getter for ParticleType of the mass hypothesis of the track fit.
Algorithm class to handle the fitting of RecoTrack objects.
Definition: TrackFitter.h:121
Class that bundles various TrackFitResults.
Definition: Track.h:25
const TrackFitResult * getTrackFitResultWithClosestMass(const Const::ChargedStable &requestedType) const
Return the track fit for a fit hypothesis with the closest mass.
Definition: Track.cc:104
bool m_useOnlyOneSVDHitPair
false only if the V0Fitter mode is 3
Definition: V0Fitter.h:170
void setFitterMode(int fitterMode)
set V0 fitter mode.
Definition: V0Fitter.cc:77
RecoTrack * copyRecoTrackAndFit(RecoTrack *origRecoTrack, const int trackPDG)
Create a copy of RecoTrack and fit the Track.
Definition: V0Fitter.cc:495
StoreArray< V0ValidationVertex > m_validationV0s
V0ValidationVertex (output, optional).
Definition: V0Fitter.h:160
StoreArray< RecoTrack > m_copiedRecoTracks
RecoTrack used to refit tracks (output)
Definition: V0Fitter.h:161
bool m_forcestore
true only if the V0Fitter mode is 1
Definition: V0Fitter.h:169
void initializeCuts(double beamPipeRadius, double vertexChi2CutOutside, std::tuple< double, double > invMassRangeKshort, std::tuple< double, double > invMassRangeLambda, std::tuple< double, double > invMassRangePhoton)
Initialize the cuts which will be applied during the fit and store process.
Definition: V0Fitter.cc:94
RecoTrack * copyRecoTrack(RecoTrack *origRecoTrack)
Create a copy of RecoTrack.
Definition: V0Fitter.cc:487
bool m_validation
Validation flag.
Definition: V0Fitter.h:155
std::pair< Const::ParticleType, Const::ParticleType > getTrackHypotheses(const Const::ParticleType &v0Hypothesis) const
Get track hypotheses for a given v0 hypothesis.
Definition: V0Fitter.cc:189
StoreArray< V0 > m_v0s
V0 (output).
Definition: V0Fitter.h:159
int checkSharedInnermostCluster(const RecoTrack *recoTrackPlus, const RecoTrack *recoTrackMinus)
Compare innermost hits of daughter pairs to check if they are the same (shared) or not.
Definition: V0Fitter.cc:625
TrackFitResult * buildTrackFitResult(const genfit::Track &track, const RecoTrack *recoTrack, const genfit::MeasuredStateOnPlane &msop, const double Bz, const Const::ParticleType &trackHypothesis, const int sharedInnermostCluster)
Build TrackFitResult of V0 Track and set relation to genfit Track.
Definition: V0Fitter.cc:165
V0Fitter(const std::string &trackFitResultsName="", const std::string &v0sName="", const std::string &v0ValidationVerticesName="", const std::string &recoTracksName="", const std::string &copiedRecoTracksName="CopiedRecoTracks", bool enableValidation=false)
Constructor for the V0Fitter.
Definition: V0Fitter.cc:43
StoreArray< TrackFitResult > m_trackFitResults
TrackFitResult (output).
Definition: V0Fitter.h:158
bool fitAndStore(const Track *trackPlus, const Track *trackMinus, const Const::ParticleType &v0Hypothesis, bool &isForceStored, bool &isHitRemoved)
Fit V0 with given hypothesis and store if fit was successful.
Definition: V0Fitter.cc:205
double m_beamPipeRadius
Radius where inside/outside beampipe is defined.
Definition: V0Fitter.h:163
bool vertexFitWithRecoTracks(const Track *trackPlus, const Track *trackMinus, RecoTrack *recoTrackPlus, RecoTrack *recoTrackMinus, const Const::ParticleType &v0Hypothesis, unsigned int &hasInnerHitStatus, ROOT::Math::XYZVector &vertexPos, const bool forceStore)
fit V0 vertex using RecoTrack's as inputs.
Definition: V0Fitter.cc:348
std::tuple< double, double > m_invMassRangePhoton
invariant mass cut for Photon.
Definition: V0Fitter.h:167
bool removeInnerHits(RecoTrack *prevRecoTrack, RecoTrack *recoTrack, const int trackPDG, const ROOT::Math::XYZVector &vertexPosition)
Remove inner hits from RecoTrack at once.
Definition: V0Fitter.cc:518
std::string m_recoTracksName
RecoTrackColName (input).
Definition: V0Fitter.h:156
bool extrapolateToVertex(genfit::MeasuredStateOnPlane &stPlus, genfit::MeasuredStateOnPlane &stMinus, const ROOT::Math::XYZVector &vertexPosition)
Extrapolate the fit results to the perigee to the vertex.
double m_vertexChi2CutOutside
Chi2 cut outside beampipe.
Definition: V0Fitter.h:164
StoreArray< RecoTrack > m_recoTracks
RecoTrack (input)
Definition: V0Fitter.h:157
bool fitGFRaveVertex(genfit::Track &trackPlus, genfit::Track &trackMinus, genfit::GFRaveVertex &vertex)
Fit the V0 vertex.
Definition: V0Fitter.cc:108
int m_v0FitterMode
0: store V0 at the first vertex fit, regardless of inner hits, 1: remove hits inside the V0 vertex po...
Definition: V0Fitter.h:168
std::tuple< double, double > m_invMassRangeKshort
invariant mass cut for Kshort.
Definition: V0Fitter.h:165
std::tuple< double, double > m_invMassRangeLambda
invariant mass cut for Lambda.
Definition: V0Fitter.h:166
Need this container for exception-safe cleanup, GFRave's interface isn't exception-safe as is.
Definition: VertexVector.h:24
std::vector< genfit::GFRaveVertex * > v
Fitted vertices.
Definition: VertexVector.h:42
size_t size() const noexcept
Return size of vertex vector.
Definition: VertexVector.h:36
static constexpr auto XYZToTVector
Helper function to convert XYZVector to TVector3.
Definition: VectorUtil.h:24
Abstract base class for different kinds of events.