8#include <tracking/v0Finding/fitter/V0Fitter.h>
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>
19#include <mdst/dataobjects/HitPatternVXD.h>
20#include <mdst/dataobjects/HitPatternCDC.h>
21#include <mdst/dataobjects/Track.h>
22#include <mdst/dataobjects/TrackFitResult.h>
24#include <genfit/Track.h>
25#include <genfit/TrackPoint.h>
26#include <genfit/MeasuredStateOnPlane.h>
27#include <genfit/GFRaveVertexFactory.h>
28#include <genfit/GFRaveVertex.h>
29#include <genfit/FieldManager.h>
30#include <genfit/MaterialEffects.h>
31#include <genfit/FitStatus.h>
32#include <genfit/KalmanFitStatus.h>
33#include <genfit/KalmanFitterInfo.h>
35#include <framework/utilities/IOIntercept.h>
40 const std::string& v0ValidationVerticesName,
const std::string& recoTracksName,
41 const std::string& copiedRecoTracksName,
bool enableValidation)
50 B2DEBUG(24,
"Register DataStore for validation.");
67 B2ASSERT(
"Material effects not set up. Please use SetupGenfitExtrapolationModule.",
68 genfit::MaterialEffects::getInstance()->isInitialized());
69 B2ASSERT(
"Magnetic field not set up. Please use SetupGenfitExtrapolationModule.",
70 genfit::FieldManager::getInstance()->isInitialized());
75 if (not(0 <= fitterMode && fitterMode <= 2)) {
76 B2FATAL(
"Invalid fitter mode!");
91 double vertexChi2CutOutside,
92 std::tuple<double, double> invMassRangeKshort,
93 std::tuple<double, double> invMassRangeLambda,
94 std::tuple<double, double> invMassRangePhoton)
107 std::vector<genfit::Track*> trackPair {&trackPlus, &trackMinus};
114 genfit::GFRaveVertexFactory vertexFactory;
115 vertexFactory.findVertices(&vertexVector.
v, trackPair);
117 B2ERROR(
"Exception during vertex fit.");
121 if (vertexVector.
size() != 1) {
122 B2DEBUG(21,
"Vertex fit failed. Size of vertexVector not 1, but: " << vertexVector.
size());
126 if ((*vertexVector[0]).getNTracks() != 2) {
127 B2DEBUG(20,
"Wrong number of tracks in vertex.");
131 vertex = *vertexVector[0];
138 const ROOT::Math::XYZVector& vertexPosition,
unsigned int& hasInnerHitStatus)
141 hasInnerHitStatus = 0;
146 double extralengthPlus = stPlus.extrapolateToPoint(
XYZToTVector(vertexPosition));
147 double extralengthMinus = stMinus.extrapolateToPoint(
XYZToTVector(vertexPosition));
148 if (extralengthPlus > 0) hasInnerHitStatus |= 0x1;
149 if (extralengthMinus > 0) hasInnerHitStatus |= 0x2;
150 B2DEBUG(22,
"extralengthPlus=" << extralengthPlus <<
", extralengthMinus=" << extralengthMinus);
155 B2DEBUG(22,
"Could not extrapolate track to vertex.");
162 const genfit::MeasuredStateOnPlane& msop,
const double Bz,
164 const int sharedInnermostCluster)
170 if (sharedInnermostCluster > 0 && sharedInnermostCluster < 4) {
173 hitPatternVXDInitializer = hitPatternVXD_forflag.
getInteger();
177 =
m_trackFitResults.appendNew(ROOT::Math::XYZVector(msop.getPos()), ROOT::Math::XYZVector(msop.getMom()),
178 msop.get6DCov(), msop.getCharge(),
180 track.getFitStatus()->getPVal(),
181 Bz, hitPatternCDCInitializer, hitPatternVXDInitializer, track.getFitStatus()->getNdf());
182 return v0TrackFitResult;
196 B2FATAL(
"Given V0Hypothesis not available.");
203 bool& isForceStored,
bool& isHitRemoved)
206 isForceStored =
false;
207 isHitRemoved =
false;
216 unsigned int hasInnerHitStatus = 0;
219 ROOT::Math::XYZVector vertexPos(0, 0, 0);
237 bool failflag =
false;
245 RecoTrack* recoTrackPlus_forRefit =
nullptr;
246 RecoTrack* recoTrackMinus_forRefit =
nullptr;
247 RecoTrack* cache_recoTrackPlus =
nullptr;
248 RecoTrack* cache_recoTrackMinus =
nullptr;
250 unsigned int count_removeInnerHits = 0;
251 while (hasInnerHitStatus != 0) {
252 ++count_removeInnerHits;
257 if (hasInnerHitStatus & 0x1) {
260 if (recoTrackPlus_forRefit ==
nullptr)
264 if (!cache_recoTrackPlus) {
265 if (not
removeInnerHits(recoTrackPlus, recoTrackPlus_forRefit, pdg_trackPlus, vertexPos)) {
270 if (not
removeInnerHits(cache_recoTrackPlus, recoTrackPlus_forRefit, pdg_trackPlus, vertexPos)) {
275 cache_recoTrackPlus = recoTrackPlus_forRefit;
276 }
else if (recoTrackPlus_forRefit ==
nullptr) {
279 if (recoTrackPlus_forRefit ==
nullptr)
284 if (hasInnerHitStatus & 0x2) {
287 if (recoTrackMinus_forRefit ==
nullptr)
291 if (!cache_recoTrackMinus) {
292 if (not
removeInnerHits(recoTrackMinus, recoTrackMinus_forRefit, pdg_trackMinus, vertexPos)) {
297 if (not
removeInnerHits(cache_recoTrackMinus, recoTrackMinus_forRefit, pdg_trackMinus, vertexPos)) {
302 cache_recoTrackMinus = recoTrackMinus_forRefit;
303 }
else if (recoTrackMinus_forRefit ==
nullptr) {
306 if (recoTrackMinus_forRefit ==
nullptr)
311 hasInnerHitStatus = 0;
315 v0Hypothesis, hasInnerHitStatus, vertexPos,
false)) {
316 B2DEBUG(22,
"Vertex refit failed, or rejected by invariant mass cut.");
319 }
else if (hasInnerHitStatus == 0)
321 if (count_removeInnerHits >= 5) {
322 B2WARNING(
"Inner hits remained after " << count_removeInnerHits <<
" times of removing inner hits!");
330 bool forcestore =
true;
332 hasInnerHitStatus, vertexPos, forcestore)) {
333 B2DEBUG(22,
"Original vertex fit fails. Possibly rejected by invariant mass cut.");
336 isForceStored =
true;
347 unsigned int& hasInnerHitStatus, ROOT::Math::XYZVector& vertexPos,
348 const bool forceStore)
356 if ((plusRepresentation ==
nullptr) or (not recoTrackPlus->
wasFitSuccessful(plusRepresentation))) {
357 B2ERROR(
"Track hypothesis with closest mass not available. Should never happen, but I can continue safely anyway.");
365 if ((minusRepresentation ==
nullptr) or (not recoTrackMinus->
wasFitSuccessful(minusRepresentation))) {
366 B2ERROR(
"Track hypothesis with closest mass not available. Should never happen, but I can continue safely anyway.");
371 std::vector<genfit::AbsTrackRep*> repsPlus = gfTrackPlus.getTrackReps();
372 std::vector<genfit::AbsTrackRep*> repsMinus = gfTrackMinus.getTrackReps();
373 if (repsPlus.size() == repsMinus.size()) {
374 for (
unsigned int id = 0;
id < repsPlus.size();
id++) {
375 if (abs(repsPlus[
id]->getPDG()) == pdgTrackPlus)
376 gfTrackPlus.setCardinalRep(
id);
377 if (abs(repsMinus[
id]->getPDG()) == pdgTrackMinus)
378 gfTrackMinus.setCardinalRep(
id);
383 for (
unsigned int id = 0;
id < repsPlus.size();
id++) {
384 if (abs(repsPlus[
id]->getPDG()) == pdgTrackPlus)
385 gfTrackPlus.setCardinalRep(
id);
387 for (
unsigned int id = 0;
id < repsMinus.size();
id++) {
388 if (abs(repsMinus[
id]->getPDG()) == pdgTrackMinus)
389 gfTrackMinus.setCardinalRep(
id);
397 genfit::GFRaveVertex vert;
402 const ROOT::Math::XYZVector& posVert = ROOT::Math::XYZVector(vert.getPos());
411 B2DEBUG(22,
"Vertex outside beam pipe, chi^2 too large.");
416 B2DEBUG(22,
"Vertex accepted.");
423 if (forceStore || hasInnerHitStatus == 0) {
425 const genfit::GFRaveTrackParameters* tr0 = vert.getParameters(0);
426 const genfit::GFRaveTrackParameters* tr1 = vert.getParameters(1);
427 ROOT::Math::PxPyPzMVector lv0(tr0->getMom().Px(), tr0->getMom().Py(), tr0->getMom().Pz(), trackHypotheses.first.getMass());
428 ROOT::Math::PxPyPzMVector lv1(tr1->getMom().Px(), tr1->getMom().Py(), tr1->getMom().Pz(), trackHypotheses.second.getMass());
430 double v0InvMass = (lv0 + lv1).M();
433 B2DEBUG(22,
"Kshort vertex rejected, invariant mass out of range.");
438 B2DEBUG(22,
"Lambda vertex rejected, invariant mass out of range.");
443 B2DEBUG(22,
"Photon vertex rejected, invariant mass out of range.");
457 sharedInnermostCluster);
459 sharedInnermostCluster);
461 B2DEBUG(20,
"Creating new V0.");
462 auto v0 =
m_v0s.appendNew(std::make_pair(trackPlus, tfrPlusVtx),
463 std::make_pair(trackMinus, tfrMinusVtx),
464 posVert.X(), posVert.Y(), posVert.Z());
467 B2DEBUG(24,
"Create StoreArray and Output for validation.");
469 std::make_pair(trackPlus, tfrPlusVtx),
470 std::make_pair(trackMinus, tfrMinusVtx),
471 ROOT::Math::XYZVector(vert.getPos()),
477 v0->addRelationTo(validationV0);
502 if (not fitter.fit(*newRecoTrack, particleUsedForFitting)) {
504 B2DEBUG(20,
"track fit failed for copied RecoTrack.");
507 B2DEBUG(20,
"\t original track fit was also failed.");
515 const int trackPDG,
const ROOT::Math::XYZVector& vertexPosition)
517 if (!prevRecoTrack || !recoTrack) {
518 B2ERROR(
"Input recotrack is nullptr!");
530 unsigned int nRemoveHits = 0;
531 if (recoHitInformations.size() != prevRecoHitInformations.size()) {
532 B2WARNING(
"Copied RecoTrack has different number of hits from its original RecoTrack!");
536 for (nRemoveHits = 0; nRemoveHits < recoHitInformations.size(); ++nRemoveHits) {
537 if (!prevRecoHitInformations[nRemoveHits]->useInFit()) {
538 recoHitInformations[nRemoveHits]->setUseInFit(
false);
544 prevRecoHitInformations[nRemoveHits]);
547 double extralength = stPrevRecoHit.extrapolateToPoint(
XYZToTVector(vertexPosition));
548 if (extralength > 0) {
549 recoHitInformations[nRemoveHits]->setUseInFit(
false);
552 }
catch (NoTrackFitResult()) {
553 B2WARNING(
"Exception: no FitterInfo assigned for TrackPoint created from this RecoHit.");
554 recoHitInformations[nRemoveHits]->setUseInFit(
false);
560 B2DEBUG(22,
"Could not extrapolate track to vertex when removing inner hits, aborting.");
565 if (nRemoveHits == 0) {
567 B2DEBUG(20,
"No hits removed in removeInnerHits, aborted. Switching to use the original RecoTrack.");
571 if (recoHitInformations.size() <= nRemoveHits) {
572 B2DEBUG(20,
"Removed all the RecoHits in the RecoTrack, aborted. Switching to use the original RecoTrack.");
578 if (recoHitInformations[nRemoveHits - 1]->getTrackingDetector() == RecoHitInformation::RecoHitDetector::c_SVD) {
579 if (recoHitInformations[nRemoveHits]->getTrackingDetector() == RecoHitInformation::RecoHitDetector::c_SVD) {
580 const SVDCluster* lastRemovedSVDHit = recoHitInformations[nRemoveHits - 1]->getRelatedTo<
SVDCluster>();
582 if (!lastRemovedSVDHit || !nextSVDHit) B2ERROR(
"Last/Next SVD hit is null!");
584 if (lastRemovedSVDHit->
getSensorID() == nextSVDHit->getSensorID() &&
585 lastRemovedSVDHit->
isUCluster() && !(nextSVDHit->isUCluster())) {
586 recoHitInformations[nRemoveHits]->setUseInFit(
false);
595 recoHitInformations[nRemoveHits - 1]->getTrackingDetector() == RecoHitInformation::RecoHitDetector::c_SVD &&
596 recoHitInformations.size() > nRemoveHits + 2) {
597 if (recoHitInformations[nRemoveHits ]->getTrackingDetector() == RecoHitInformation::RecoHitDetector::c_SVD &&
598 recoHitInformations[nRemoveHits + 1]->getTrackingDetector() == RecoHitInformation::RecoHitDetector::c_SVD &&
599 recoHitInformations[nRemoveHits + 2]->getTrackingDetector() != RecoHitInformation::RecoHitDetector::c_SVD) {
600 recoHitInformations[nRemoveHits ]->setUseInFit(
false);
601 recoHitInformations[nRemoveHits + 1]->setUseInFit(
false);
606 B2DEBUG(22, nRemoveHits <<
" inner hits removed.");
610 if (not fitter.fit(*recoTrack, particleUsedForFitting)) {
611 B2DEBUG(20,
"track fit failed after removing inner hits.");
614 B2DEBUG(20,
"\t previous track fit was also failed.");
630 const std::vector<RecoHitInformation*>& recoHitInformationsMinus = recoTrackMinus->
getRecoHitInformations(
632 unsigned int iInnermostHitPlus, iInnermostHitMinus;
633 for (iInnermostHitPlus = 0 ; iInnermostHitPlus < recoHitInformationsPlus.size() ; ++iInnermostHitPlus)
634 if (recoHitInformationsPlus[iInnermostHitPlus]->useInFit())
break;
635 for (iInnermostHitMinus = 0 ; iInnermostHitMinus < recoHitInformationsMinus.size() ; ++iInnermostHitMinus)
636 if (recoHitInformationsMinus[iInnermostHitMinus]->useInFit())
break;
637 if (iInnermostHitPlus == recoHitInformationsPlus.size() || iInnermostHitMinus == recoHitInformationsMinus.size()) {
638 B2WARNING(
"checkSharedInnermostCluster function called for recoTrack including no hit used for fit! This should not happen!");
641 const auto& recoHitInfoPlus = recoHitInformationsPlus[iInnermostHitPlus];
642 const auto& recoHitInfoMinus = recoHitInformationsMinus[iInnermostHitMinus];
644 if (recoHitInfoPlus->getTrackingDetector() == recoHitInfoMinus->getTrackingDetector()) {
645 if (recoHitInfoPlus->getTrackingDetector() == RecoHitInformation::c_PXD) {
648 if (clusterPlus == clusterMinus) {
651 }
else if (recoHitInfoPlus->getTrackingDetector() == RecoHitInformation::c_SVD) {
656 if (clusterPlus->
isUCluster() && clusterMinus->isUCluster()) {
657 if (recoHitInformationsPlus.size() > iInnermostHitPlus + 1
658 && recoHitInformationsMinus.size() > iInnermostHitMinus + 1) {
659 const auto& recoHitInfoNextPlus = recoHitInformationsPlus[iInnermostHitPlus + 1];
660 const auto& recoHitInfoNextMinus = recoHitInformationsMinus[iInnermostHitMinus + 1];
662 if (recoHitInfoNextPlus->useInFit() && recoHitInfoNextMinus->useInFit()
663 && recoHitInfoNextPlus->getTrackingDetector() == RecoHitInformation::c_SVD
664 && recoHitInfoNextMinus->getTrackingDetector() == RecoHitInformation::c_SVD) {
667 if (!(clusterNextPlus->
isUCluster()) && !(clusterNextMinus->isUCluster())
669 && clusterMinus->getSensorID() == clusterNextMinus->getSensorID()) {
670 if (clusterPlus == clusterMinus)
672 if (clusterNextPlus == clusterNextMinus)
675 B2WARNING(
"SVD cluster to be paired is not on V-side, or not on the same sensor.");
679 B2WARNING(
"No SVD cluster to be paired.");
683 B2WARNING(
"Innermost SVD U-cluster is the only hit in a daughter track. This should not happen.");
687 B2WARNING(
"No SVD U-cluster in the innermost cluster.");
static ROOT::Math::XYZVector getFieldInTesla(const ROOT::Math::XYZVector &pos)
return the magnetic field at a given position in Tesla.
Provides a type-safe way to pass members of the chargedStableSet set.
The ParticleType class for identifying different particle types.
int getPDGCode() const
PDG code.
static const ParticleType Lambda
Lambda particle.
static const ChargedStable pion
charged pion particle
static const ParticleType antiLambda
Anti-Lambda particle.
static const ChargedStable proton
proton particle
static const ParticleType invalidParticle
Invalid particle, used internally.
static const ParticleType Kshort
K^0_S particle.
static const ParticleType photon
photon particle
static const ChargedStable electron
electron particle
@ c_WriteOut
Object/array should be saved by output modules.
@ c_ErrorIfAlreadyRegistered
If the object/array was already registered, produce an error (aborting initialisation).
Hit pattern of the VXD within a track.
unsigned int getInteger() const
Getter for the underlying integer.
void setInnermostHitShareStatus(const unsigned short innermostHitShareStatus)
Set the innermost hit share flags for V0 daughters.
bool start()
Start intercepting the output.
Capture stdout and stderr and convert into log messages.
@ c_Debug
Debug: for code development.
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
static genfit::Track & getGenfitTrack(RecoTrack &recoTrack)
Give access to the RecoTrack's genfit::Track.
This is the Reconstruction Event-Data Model Track.
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.
bool wasFitSuccessful(const genfit::AbsTrackRep *representation=nullptr) const
Returns true if the last fit with the given representation was successful.
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...
std::vector< RecoHitInformation * > getRecoHitInformations(bool getSorted=false) const
Return a list of all RecoHitInformations associated with the RecoTrack.
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.
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...
genfit::AbsTrackRep * getTrackRepresentationForPDG(int pdgCode) const
Return an already created track representation of the given reco track for the PDG.
const genfit::MeasuredStateOnPlane & getMeasuredStateOnPlaneFromRecoHit(const RecoHitInformation *recoHitInfo, const genfit::AbsTrackRep *representation=nullptr) const
Return genfit's MeasuredStateOnPlane on plane for associated with one RecoHitInformation.
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.
VxdID getSensorID() const
Get the sensor ID.
bool isUCluster() const
Get the direction of strips.
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.
Class that bundles various TrackFitResults.
const TrackFitResult * getTrackFitResultWithClosestMass(const Const::ChargedStable &requestedType) const
Return the track fit for the fit hypothesis with the closest mass.
bool m_useOnlyOneSVDHitPair
false only if the V0Fitter mode is 3
void setFitterMode(int fitterMode)
set V0 fitter mode.
RecoTrack * copyRecoTrackAndFit(RecoTrack *origRecoTrack, const int trackPDG)
Create a copy of RecoTrack and fit the Track.
StoreArray< V0ValidationVertex > m_validationV0s
V0ValidationVertex (output, optional).
StoreArray< RecoTrack > m_copiedRecoTracks
RecoTrack used to refit tracks (output)
bool m_forcestore
true only if the V0Fitter mode is 1
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.
RecoTrack * copyRecoTrack(RecoTrack *origRecoTrack)
Create a copy of RecoTrack.
bool m_validation
Validation flag.
std::pair< Const::ParticleType, Const::ParticleType > getTrackHypotheses(const Const::ParticleType &v0Hypothesis) const
Get track hypotheses for a given v0 hypothesis.
StoreArray< V0 > m_v0s
V0 (output).
int checkSharedInnermostCluster(const RecoTrack *recoTrackPlus, const RecoTrack *recoTrackMinus)
Compare innermost hits of daughter pairs to check if they are the same (shared) or not.
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.
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.
StoreArray< TrackFitResult > m_trackFitResults
TrackFitResult (output).
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.
double m_beamPipeRadius
Radius where inside/outside beampipe is defined.
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.
std::tuple< double, double > m_invMassRangePhoton
invariant mass cut for Photon.
bool removeInnerHits(RecoTrack *prevRecoTrack, RecoTrack *recoTrack, const int trackPDG, const ROOT::Math::XYZVector &vertexPosition)
Remove inner hits from RecoTrack at once.
std::string m_recoTracksName
RecoTrackColName (input).
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.
StoreArray< RecoTrack > m_recoTracks
RecoTrack (input)
bool fitGFRaveVertex(genfit::Track &trackPlus, genfit::Track &trackMinus, genfit::GFRaveVertex &vertex)
Fit the V0 vertex.
int m_v0FitterMode
0: store V0 at the first vertex fit, regardless of inner hits, 1: remove hits inside the V0 vertex po...
std::tuple< double, double > m_invMassRangeKshort
invariant mass cut for Kshort.
std::tuple< double, double > m_invMassRangeLambda
invariant mass cut for Lambda.
Need this container for exception-safe cleanup, GFRave's interface isn't exception-safe as is.
std::vector< genfit::GFRaveVertex * > v
Fitted vertices.
size_t size() const noexcept
Return size of vertex vector.
static constexpr auto XYZToTVector
Helper function to convert XYZVector to TVector3.
Abstract base class for different kinds of events.