Belle II Software development
V0Fitter Class Reference

V0Fitter class to create V0 mdst's from reconstructed tracks. More...

#include <V0Fitter.h>

Public Member Functions

 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.
 
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.
 
void setFitterMode (int fitterMode)
 set V0 fitter mode.
 
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.
 
std::pair< Const::ParticleType, Const::ParticleTypegetTrackHypotheses (const Const::ParticleType &v0Hypothesis) const
 Get track hypotheses for a given v0 hypothesis.
 

Private Member Functions

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.
 
RecoTrackcopyRecoTrack (RecoTrack *origRecoTrack)
 Create a copy of RecoTrack.
 
RecoTrackcopyRecoTrackAndFit (RecoTrack *origRecoTrack, const int trackPDG)
 Create a copy of RecoTrack and fit the Track.
 
bool removeInnerHits (RecoTrack *prevRecoTrack, RecoTrack *recoTrack, const int trackPDG, const ROOT::Math::XYZVector &vertexPosition)
 Remove inner hits from RecoTrack at once.
 
int checkSharedInnermostCluster (const RecoTrack *recoTrackPlus, const RecoTrack *recoTrackMinus)
 Compare innermost hits of daughter pairs to check if they are the same (shared) or not.
 
bool fitGFRaveVertex (genfit::Track &trackPlus, genfit::Track &trackMinus, genfit::GFRaveVertex &vertex)
 Fit the V0 vertex.
 
bool extrapolateToVertex (genfit::MeasuredStateOnPlane &stPlus, genfit::MeasuredStateOnPlane &stMinus, const ROOT::Math::XYZVector &vertexPosition)
 Extrapolate the fit results to the perigee to the vertex.
 
bool extrapolateToVertex (genfit::MeasuredStateOnPlane &stPlus, genfit::MeasuredStateOnPlane &stMinus, const ROOT::Math::XYZVector &vertexPosition, unsigned int &hasInnerHitStatus)
 Extrapolate the fit results to the perigee to the vertex.
 
TrackFitResultbuildTrackFitResult (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.
 

Private Attributes

bool m_validation
 Validation flag.
 
std::string m_recoTracksName
 RecoTrackColName (input).
 
StoreArray< RecoTrackm_recoTracks
 RecoTrack (input)
 
StoreArray< TrackFitResultm_trackFitResults
 TrackFitResult (output).
 
StoreArray< V0m_v0s
 V0 (output).
 
StoreArray< V0ValidationVertexm_validationV0s
 V0ValidationVertex (output, optional).
 
StoreArray< RecoTrackm_copiedRecoTracks
 RecoTrack used to refit tracks (output)
 
double m_beamPipeRadius
 Radius where inside/outside beampipe is defined.
 
double m_vertexChi2CutOutside
 Chi2 cut outside beampipe.
 
std::tuple< double, double > m_invMassRangeKshort
 invariant mass cut for Kshort.
 
std::tuple< double, double > m_invMassRangeLambda
 invariant mass cut for Lambda.
 
std::tuple< double, double > m_invMassRangePhoton
 invariant mass cut for Photon.
 
int m_v0FitterMode
 0: store V0 at the first vertex fit, regardless of inner hits, 1: remove hits inside the V0 vertex position, 2: mode 1 + don't use SVD hits if there is only one available SVD hit-pair (default)
 
bool m_forcestore
 true only if the V0Fitter mode is 1
 
bool m_useOnlyOneSVDHitPair
 false only if the V0Fitter mode is 3
 

Friends

class V0FitterTest_GetTrackHypotheses_Test
 
class V0FitterTest_EnableValidation_Test
 
class V0FitterTest_InitializeCuts_Test
 

Detailed Description

V0Fitter class to create V0 mdst's from reconstructed tracks.

To use this class, give the V0Fitter a positive and a negative charged track and call the fitAndStore function.

v0Fitter.fitAndStore(B2TrackPositive, B2TrackNegative, V0HypothesisParticleType);

Definition at line 39 of file V0Fitter.h.

Constructor & Destructor Documentation

◆ V0Fitter()

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.

m_recoTracks

register m_copiedRecoTracks

relation : m_recoTracks <--> m_copiedRecoTracks

Definition at line 39 of file V0Fitter.cc.

42 : m_validation(enableValidation), m_recoTracksName(recoTracksName), m_v0FitterMode(1), m_forcestore(false),
44{
45 m_trackFitResults.isRequired(trackFitResultsName);
47 //Relation to RecoTracks from Tracks is already tested at the module level.
48
49 if (m_validation) {
50 B2DEBUG(24, "Register DataStore for validation.");
51 m_validationV0s.registerInDataStore(v0ValidationVerticesName, DataStore::c_ErrorIfAlreadyRegistered);
52 m_v0s.registerRelationTo(m_validationV0s);
53 }
54
57
59 m_copiedRecoTracks.registerInDataStore(copiedRecoTracksName,
62
64 m_copiedRecoTracks.registerRelationTo(m_recoTracks);
65
66
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());
71}
@ 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
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
bool m_useOnlyOneSVDHitPair
false only if the V0Fitter mode is 3
Definition V0Fitter.h:170
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
bool m_validation
Validation flag.
Definition V0Fitter.h:155
StoreArray< V0 > m_v0s
V0 (output).
Definition V0Fitter.h:159
StoreArray< TrackFitResult > m_trackFitResults
TrackFitResult (output).
Definition V0Fitter.h:158
std::string m_recoTracksName
RecoTrackColName (input).
Definition V0Fitter.h:156
StoreArray< RecoTrack > m_recoTracks
RecoTrack (input)
Definition V0Fitter.h:157
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

Member Function Documentation

◆ buildTrackFitResult()

TrackFitResult * buildTrackFitResult ( const genfit::Track & track,
const RecoTrack * recoTrack,
const genfit::MeasuredStateOnPlane & msop,
const double Bz,
const Const::ParticleType & trackHypothesis,
const int sharedInnermostCluster )
private

Build TrackFitResult of V0 Track and set relation to genfit Track.

Definition at line 161 of file V0Fitter.cc.

165{
166 const uint64_t hitPatternCDCInitializer = TrackBuilder::getHitPatternCDCInitializer(*recoTrack);
167 uint32_t hitPatternVXDInitializer = TrackBuilder::getHitPatternVXDInitializer(*recoTrack);
168
169 // If the innermost hit is shared among V0 daughters, assign flag in the infoLayer.
170 if (sharedInnermostCluster > 0 && sharedInnermostCluster < 4) {
171 HitPatternVXD hitPatternVXD_forflag = HitPatternVXD(hitPatternVXDInitializer);
172 hitPatternVXD_forflag.setInnermostHitShareStatus(sharedInnermostCluster);
173 hitPatternVXDInitializer = hitPatternVXD_forflag.getInteger();
174 }
175
176 TrackFitResult* v0TrackFitResult
177 = m_trackFitResults.appendNew(ROOT::Math::XYZVector(msop.getPos()), ROOT::Math::XYZVector(msop.getMom()),
178 msop.get6DCov(), msop.getCharge(),
179 trackHypothesis,
180 track.getFitStatus()->getPVal(),
181 Bz, hitPatternCDCInitializer, hitPatternVXDInitializer, track.getFitStatus()->getNdf());
182 return v0TrackFitResult;
183}
unsigned int getInteger() const
Getter for the underlying integer.
void setInnermostHitShareStatus(const unsigned short innermostHitShareStatus)
Set the innermost hit share flags for V0 daughters.
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.

◆ checkSharedInnermostCluster()

int checkSharedInnermostCluster ( const RecoTrack * recoTrackPlus,
const RecoTrack * recoTrackMinus )
private

Compare innermost hits of daughter pairs to check if they are the same (shared) or not.

For SVD hits, compare U- and V- hit pair.

Parameters
recoTrackPlus,recoTrackMinusinput RecoTrack pair
Returns
If 1D- or 2D-hits are shared as the innermost hits among V0 daughters. 0x1(0x2) bit represents V/z(U/r-phi)-hit share. -1 for exception.

If 1D- or 2D-hits are shared as the innermost hits among V0 daughters. 0x1(0x2) bit in flag represents V/z(U/r-phi)-hit share.

if the innermost hit is a SVD U-hit, check the next hit in addition note: for the SVD pair hits, U-hit should be first and the V-hit the next in the sorted RecoHitInformation array.

Definition at line 621 of file V0Fitter.cc.

622{
625 int flag = 0; // -1 for exception
626
627 // get the innermost hit for plus/minus-daughter
628 const std::vector<RecoHitInformation*>& recoHitInformationsPlus = recoTrackPlus->getRecoHitInformations(
629 true); // true for sorted info.
630 const std::vector<RecoHitInformation*>& recoHitInformationsMinus = recoTrackMinus->getRecoHitInformations(
631 true);// true for sorted info.
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!");
639 return -1;
640 }
641 const auto& recoHitInfoPlus = recoHitInformationsPlus[iInnermostHitPlus];
642 const auto& recoHitInfoMinus = recoHitInformationsMinus[iInnermostHitMinus];
643
644 if (recoHitInfoPlus->getTrackingDetector() == recoHitInfoMinus->getTrackingDetector()) {
645 if (recoHitInfoPlus->getTrackingDetector() == RecoHitInformation::c_PXD) {
646 const PXDCluster* clusterPlus = recoHitInfoPlus->getRelatedTo<PXDCluster>();
647 const PXDCluster* clusterMinus = recoHitInfoMinus->getRelatedTo<PXDCluster>();
648 if (clusterPlus == clusterMinus) { // if they share a same PXDCluster, set the flag
649 flag = 3; // PXD cluster is a 2D-hit
650 }
651 } else if (recoHitInfoPlus->getTrackingDetector() == RecoHitInformation::c_SVD) {
654 const SVDCluster* clusterPlus = recoHitInfoPlus->getRelatedTo<SVDCluster>();
655 const SVDCluster* clusterMinus = recoHitInfoMinus->getRelatedTo<SVDCluster>();
656 if (clusterPlus->isUCluster() && clusterMinus->isUCluster()) {
657 if (recoHitInformationsPlus.size() > iInnermostHitPlus + 1
658 && recoHitInformationsMinus.size() > iInnermostHitMinus + 1) { // not to access an array out of boundary
659 const auto& recoHitInfoNextPlus = recoHitInformationsPlus[iInnermostHitPlus + 1];
660 const auto& recoHitInfoNextMinus = recoHitInformationsMinus[iInnermostHitMinus + 1];
661 // sanity check to access next hits
662 if (recoHitInfoNextPlus->useInFit() && recoHitInfoNextMinus->useInFit() // this should be satisfied by default
663 && recoHitInfoNextPlus->getTrackingDetector() == RecoHitInformation::c_SVD
664 && recoHitInfoNextMinus->getTrackingDetector() == RecoHitInformation::c_SVD) {
665 const SVDCluster* clusterNextPlus = recoHitInfoNextPlus->getRelatedTo<SVDCluster>();
666 const SVDCluster* clusterNextMinus = recoHitInfoNextMinus->getRelatedTo<SVDCluster>();
667 if (!(clusterNextPlus->isUCluster()) && !(clusterNextMinus->isUCluster())
668 && clusterPlus->getSensorID() == clusterNextPlus->getSensorID()
669 && clusterMinus->getSensorID() == clusterNextMinus->getSensorID()) {
670 if (clusterPlus == clusterMinus)
671 flag += 2; // SVD U-cluster is shared
672 if (clusterNextPlus == clusterNextMinus)
673 flag += 1; // SVD V-cluster is shared
674 } else {
675 B2WARNING("SVD cluster to be paired is not on V-side, or not on the same sensor.");
676 return -1;
677 }
678 } else {
679 B2WARNING("No SVD cluster to be paired.");
680 return -1;
681 }
682 } else {
683 B2WARNING("Innermost SVD U-cluster is the only hit in a daughter track. This should not happen.");
684 return -1;
685 }
686 } else {
687 B2WARNING("No SVD U-cluster in the innermost cluster.");
688 return -1;
689 }
690 }
691 }
692 return flag;
693
694}
std::vector< RecoHitInformation * > getRecoHitInformations(bool getSorted=false) const
Return a list of all RecoHitInformations associated with the RecoTrack.
Definition RecoTrack.cc:557
TO * getRelatedTo(const std::string &name="", const std::string &namedRelation="") const
Get the object to which this object has a relation.
VxdID getSensorID() const
Get the sensor ID.
Definition SVDCluster.h:102
bool isUCluster() const
Get the direction of strips.
Definition SVDCluster.h:110

◆ copyRecoTrack()

RecoTrack * copyRecoTrack ( RecoTrack * origRecoTrack)
private

Create a copy of RecoTrack.

Track fit should be executed in removeInnerHits function.

Parameters
origRecoTrackoriginal RecoTrack
Returns
copied RecoTrack stored in the m_copiedRecoTracks, nullptr if track fit fails (this should not happen)

Definition at line 483 of file V0Fitter.cc.

484{
485 RecoTrack* newRecoTrack = origRecoTrack->copyToStoreArray(m_copiedRecoTracks);
486 newRecoTrack->addHitsFromRecoTrack(origRecoTrack);
487 newRecoTrack->addRelationTo(origRecoTrack);
488 return newRecoTrack;
489}
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
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
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).

◆ copyRecoTrackAndFit()

RecoTrack * copyRecoTrackAndFit ( RecoTrack * origRecoTrack,
const int trackPDG )
private

Create a copy of RecoTrack and fit the Track.

Parameters
origRecoTrackoriginal RecoTrack
trackPDGsigned PDG used for the track fit hypothesis
Returns
copied RecoTrack stored in the m_copiedRecoTracks, nullptr if track fit fails (this should not happen)

original track information

only a positive PDG number is allowed for the input

fit newRecoTrack

check fit status of original track

Definition at line 491 of file V0Fitter.cc.

492{
494 Const::ChargedStable particleUsedForFitting(std::abs(trackPDG));
495 const genfit::AbsTrackRep* origTrackRep = origRecoTrack->getTrackRepresentationForPDG(std::abs(
496 trackPDG));
497
498 RecoTrack* newRecoTrack = copyRecoTrack(origRecoTrack);
499
501 TrackFitter fitter;
502 if (not fitter.fit(*newRecoTrack, particleUsedForFitting)) {
503 // This is not expected, but happens sometimes.
504 B2DEBUG(20, "track fit failed for copied RecoTrack.");
506 if (not origRecoTrack->wasFitSuccessful(origTrackRep))
507 B2DEBUG(20, "\t original track fit was also failed.");
508 return nullptr;
509 }
510
511 return newRecoTrack;
512}
bool wasFitSuccessful(const genfit::AbsTrackRep *representation=nullptr) const
Returns true if the last fit with the given representation was successful.
Definition RecoTrack.cc:336
genfit::AbsTrackRep * getTrackRepresentationForPDG(int pdgCode) const
Return an already created track representation of the given reco track for the PDG.
Definition RecoTrack.cc:475
RecoTrack * copyRecoTrack(RecoTrack *origRecoTrack)
Create a copy of RecoTrack.
Definition V0Fitter.cc:483

◆ extrapolateToVertex()

bool extrapolateToVertex ( genfit::MeasuredStateOnPlane & stPlus,
genfit::MeasuredStateOnPlane & stMinus,
const ROOT::Math::XYZVector & vertexPosition,
unsigned int & hasInnerHitStatus )
private

Extrapolate the fit results to the perigee to the vertex.

used in the fitAndStore function.

If the daughter tracks have hits inside the V0 vertex, bits in the hasInnerHiStatus variable are set.

if tracks have hits inside the V0 vertex position, bits in hasInnerHitStatus are set.

initialize

extrapolate the first (innermost) hit to the V0 vertex position the value will be positive (negative) if the direction of the extrapolation is (counter)momentum-wise

plus track has hits inside the V0 vertex.

minus track has hits inside the V0 vertex.

This shouldn't ever happen, but I can see the extrapolation code trying several windings before giving up, so this happens occasionally. Something more stable would perhaps be desirable.

Definition at line 137 of file V0Fitter.cc.

139{
141 hasInnerHitStatus = 0;
142
143 try {
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);
151 } catch (...) {
155 B2DEBUG(22, "Could not extrapolate track to vertex.");
156 return false;
157 }
158 return true;
159}
static constexpr auto XYZToTVector
Helper function to convert XYZVector to TVector3.
Definition VectorUtil.h:26

◆ fitAndStore()

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.

remove hits inside the V0 vertex position

Initialize status flag for counting

Existence of corresponding RecoTrack already checked at the module level;

hasInnerHitStatus: 0x1: plus track has hits inside the V0 vertex. 0x2: minus track has hits inside the V0 vertex.

fitted vertex position

Try V0 vertex fit. If the fit fails, return false. If hasInnerHitStatus evaluated in the vertexFitWithRecoTracks function results to be 0 (no hits inside the V0 vertex position in both tracks), TrackFitResult and V0 objects are build and stored, and then terminate this function with returning true. If hasInnerHitStatus!=0, objects are not stored and the tracks are refitted with removing inner hit later.

If one or two tracks have hits inside the V0 vertex. (0x1(0x2) bit in hasInnerHitStatus represents plus(minus)-track), redo vertex fit after refitting the track(s) with removing innermost hit in the track(s). Repeat until all the inner hits are removed or the vertex fit fails.

PDG code (positive number) used as the track hypothesis in the track fitting.

positive number

positive number

If the track has a hit inside the V0 vertex position, use refitted RecoTrack with removing inner hits in the next V0 vertex fit. Else, use the original RecoTrack.

for plus-charged track

if the track refit fails, break out of this loop and revert back to the original vertex fit with the original tracks.

for minus-charged track

if the track refit fails, break out of this loop and revert back to the original vertex fit with the original tracks.

V0 vertex fit

if vertex fit fails, break out of this loop and revert back to the original vertex fit with the original tracks.

end of the while loop

if failflag==true, revert back to the original vertex fit with the original tracks.

Definition at line 201 of file V0Fitter.cc.

204{
206 isForceStored = false;
207 isHitRemoved = false;
208
210 RecoTrack* recoTrackPlus = trackPlus->getRelated<RecoTrack>(m_recoTracksName);
211 RecoTrack* recoTrackMinus = trackMinus->getRelated<RecoTrack>(m_recoTracksName);
212
216 unsigned int hasInnerHitStatus = 0;
217
219 ROOT::Math::XYZVector vertexPos(0, 0, 0);
220
226 if (not vertexFitWithRecoTracks(trackPlus, trackMinus, recoTrackPlus, recoTrackMinus, v0Hypothesis,
227 hasInnerHitStatus, vertexPos, m_forcestore))
228 return false;
229
230 if ((hasInnerHitStatus == 0) or m_forcestore)
231 return true;
232
237 bool failflag = false;
239 const auto trackHypotheses = getTrackHypotheses(v0Hypothesis);
240 const int pdg_trackPlus = trackPlus->getTrackFitResultWithClosestMass(
241 trackHypotheses.first)->getParticleType().getPDGCode();
242 const int pdg_trackMinus = trackMinus->getTrackFitResultWithClosestMass(
243 trackHypotheses.second)->getParticleType().getPDGCode();
244
245 RecoTrack* recoTrackPlus_forRefit = nullptr;
246 RecoTrack* recoTrackMinus_forRefit = nullptr;
247 RecoTrack* cache_recoTrackPlus = nullptr;
248 RecoTrack* cache_recoTrackMinus = nullptr;
249
250 unsigned int count_removeInnerHits = 0;
251 while (hasInnerHitStatus != 0) {
252 ++count_removeInnerHits;
255
257 if (hasInnerHitStatus & 0x1) {
258 // create a copy of the original RecoTrack w/o track fit
259 recoTrackPlus_forRefit = copyRecoTrack(recoTrackPlus);
260 if (recoTrackPlus_forRefit == nullptr)
261 return false;
264 if (!cache_recoTrackPlus) {
265 if (not removeInnerHits(recoTrackPlus, recoTrackPlus_forRefit, pdg_trackPlus, vertexPos)) {
266 failflag = true;
267 break;
268 }
269 } else {
270 if (not removeInnerHits(cache_recoTrackPlus, recoTrackPlus_forRefit, pdg_trackPlus, vertexPos)) {
271 failflag = true;
272 break;
273 }
274 }
275 cache_recoTrackPlus = recoTrackPlus_forRefit;
276 } else if (recoTrackPlus_forRefit == nullptr) {
277 // create a copy of the original RecoTrack w/ track fit (only once)
278 recoTrackPlus_forRefit = copyRecoTrackAndFit(recoTrackPlus, pdg_trackPlus);
279 if (recoTrackPlus_forRefit == nullptr)
280 return false;
281 }
282
284 if (hasInnerHitStatus & 0x2) {
285 // create a copy of the original RecoTrack w/o track fit
286 recoTrackMinus_forRefit = copyRecoTrack(recoTrackMinus);
287 if (recoTrackMinus_forRefit == nullptr)
288 return false;
291 if (!cache_recoTrackMinus) {
292 if (not removeInnerHits(recoTrackMinus, recoTrackMinus_forRefit, pdg_trackMinus, vertexPos)) {
293 failflag = true;
294 break;
295 }
296 } else {
297 if (not removeInnerHits(cache_recoTrackMinus, recoTrackMinus_forRefit, pdg_trackMinus, vertexPos)) {
298 failflag = true;
299 break;
300 }
301 }
302 cache_recoTrackMinus = recoTrackMinus_forRefit;
303 } else if (recoTrackMinus_forRefit == nullptr) {
304 // create a copy of the original RecoTrack w/ track fit (only once)
305 recoTrackMinus_forRefit = copyRecoTrackAndFit(recoTrackMinus, pdg_trackMinus);
306 if (recoTrackMinus_forRefit == nullptr)
307 return false;
308 }
309
311 hasInnerHitStatus = 0;
314 if (not vertexFitWithRecoTracks(trackPlus, trackMinus, recoTrackPlus_forRefit, recoTrackMinus_forRefit,
315 v0Hypothesis, hasInnerHitStatus, vertexPos, false)) {
316 B2DEBUG(22, "Vertex refit failed, or rejected by invariant mass cut.");
317 failflag = true;
318 break;
319 } else if (hasInnerHitStatus == 0)
320 isHitRemoved = true;
321 if (count_removeInnerHits >= 5) {
322 B2WARNING("Inner hits remained after " << count_removeInnerHits << " times of removing inner hits!");
323 failflag = true;
324 break;
325 }
326 }
327
329 if (failflag) {
330 bool forcestore = true;
331 if (not vertexFitWithRecoTracks(trackPlus, trackMinus, recoTrackPlus, recoTrackMinus, v0Hypothesis,
332 hasInnerHitStatus, vertexPos, forcestore)) {
333 B2DEBUG(22, "Original vertex fit fails. Possibly rejected by invariant mass cut.");
334 return false;
335 } else
336 isForceStored = true;
337 }
338
339 return true;
340}
int getPDGCode() const
PDG code.
Definition Const.h:473
T * getRelated(const std::string &name="", const std::string &namedRelation="") const
Get the object to or from which this object has a relation.
Const::ParticleType getParticleType() const
Getter for ParticleType of the mass hypothesis of the track fit.
const TrackFitResult * getTrackFitResultWithClosestMass(const Const::ChargedStable &requestedType) const
Return the track fit for the fit hypothesis with the closest mass.
Definition Track.cc:104
RecoTrack * copyRecoTrackAndFit(RecoTrack *origRecoTrack, const int trackPDG)
Create a copy of RecoTrack and fit the Track.
Definition V0Fitter.cc:491
std::pair< Const::ParticleType, Const::ParticleType > getTrackHypotheses(const Const::ParticleType &v0Hypothesis) const
Get track hypotheses for a given v0 hypothesis.
Definition V0Fitter.cc:185
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:344
bool removeInnerHits(RecoTrack *prevRecoTrack, RecoTrack *recoTrack, const int trackPDG, const ROOT::Math::XYZVector &vertexPosition)
Remove inner hits from RecoTrack at once.
Definition V0Fitter.cc:514

◆ fitGFRaveVertex()

bool fitGFRaveVertex ( genfit::Track & trackPlus,
genfit::Track & trackMinus,
genfit::GFRaveVertex & vertex )
private

Fit the V0 vertex.

Parameters
trackPlus
trackMinus
vertexResult of the fit is returned via reference.
Returns

Definition at line 104 of file V0Fitter.cc.

105{
106 VertexVector vertexVector;
107 std::vector<genfit::Track*> trackPair {&trackPlus, &trackMinus};
108
109 try {
110 IOIntercept::OutputToLogMessages
111 logCapture("V0Fitter GFRaveVertexFactory", LogConfig::c_Debug, LogConfig::c_Debug);
112 logCapture.start();
113
114 genfit::GFRaveVertexFactory vertexFactory;
115 vertexFactory.findVertices(&vertexVector.v, trackPair);
116 } catch (...) {
117 B2ERROR("Exception during vertex fit.");
118 return false;
119 }
120
121 if (vertexVector.size() != 1) {
122 B2DEBUG(21, "Vertex fit failed. Size of vertexVector not 1, but: " << vertexVector.size());
123 return false;
124 }
125
126 if ((*vertexVector[0]).getNTracks() != 2) {
127 B2DEBUG(20, "Wrong number of tracks in vertex.");
128 return false;
129 }
130
131 vertex = *vertexVector[0];
132 return true;
133}
@ c_Debug
Debug: for code development.
Definition LogConfig.h:26
std::vector< genfit::GFRaveVertex * > v
Fitted vertices.
size_t size() const noexcept
Return size of vertex vector.

◆ getTrackHypotheses()

std::pair< Const::ParticleType, Const::ParticleType > getTrackHypotheses ( const Const::ParticleType & v0Hypothesis) const

Get track hypotheses for a given v0 hypothesis.

Definition at line 185 of file V0Fitter.cc.

186{
187 if (v0Hypothesis == Const::Kshort) {
188 return std::make_pair(Const::pion, Const::pion);
189 } else if (v0Hypothesis == Const::photon) {
190 return std::make_pair(Const::electron, Const::electron);
191 } else if (v0Hypothesis == Const::Lambda) {
192 return std::make_pair(Const::proton, Const::pion);
193 } else if (v0Hypothesis == Const::antiLambda) {
194 return std::make_pair(Const::pion, Const::proton);
195 }
196 B2FATAL("Given V0Hypothesis not available.");
197 return std::make_pair(Const::invalidParticle, Const::invalidParticle); // return something to avoid triggering cppcheck
198}
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

◆ initializeCuts()

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 at line 90 of file V0Fitter.cc.

95{
96 m_beamPipeRadius = beamPipeRadius;
97 m_vertexChi2CutOutside = vertexChi2CutOutside;
98 m_invMassRangeKshort = invMassRangeKshort;
99 m_invMassRangeLambda = invMassRangeLambda;
100 m_invMassRangePhoton = invMassRangePhoton;
101}
double m_beamPipeRadius
Radius where inside/outside beampipe is defined.
Definition V0Fitter.h:163
std::tuple< double, double > m_invMassRangePhoton
invariant mass cut for Photon.
Definition V0Fitter.h:167
double m_vertexChi2CutOutside
Chi2 cut outside beampipe.
Definition V0Fitter.h:164
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

◆ removeInnerHits()

bool removeInnerHits ( RecoTrack * prevRecoTrack,
RecoTrack * recoTrack,
const int trackPDG,
const ROOT::Math::XYZVector & vertexPosition )
private

Remove inner hits from RecoTrack at once.

Hits are removed from the minus-end of the momentum direction. For SVD hits, remove U- and V- hit pair at once. Input RecoTrack is fitted in the function. If track fit fails, return false.

Parameters
prevRecoTrackoriginal RecoTrack
recoTrackinput RecoTrack, updated in this function
trackPDGsigned PDG used for the track fit hypothesis
vertexPositionV0 vertex position
Returns

original track information

only a positive PDG number is allowed for the input

disable inner hits

true for sorted info.

true for sorted info.

check recohits one by one whether the vertex is outside/inside them

make a clone, not use the reference so that the genfit::MeasuredStateOnPlane will not be altered.

extrapolate the hit to the V0 vertex position the value will be positive (negative) if the direction of the extrapolation is (counter)momentum-wise

This shouldn't ever happen, but I can see the extrapolation code trying several windings before giving up, so this happens occasionally. Something more stable would perhaps be desirable.

N removed hits should not reach N hits in the track

if the last removed hit is a SVD U-hit, remove the next hit (SVD V-hit as the pair) in addition note: for the SVD pair hits, U-hit should be first and the V-hit the next in the sorted RecoHitInformation array.

if N of remaining SVD hit-pair is only one, don't use the SVD hits

SVD U-hit

SVD V-hit

not SVD hit (CDC)

fit recoTrack

check fit status of original track

Definition at line 514 of file V0Fitter.cc.

516{
517 if (!prevRecoTrack || !recoTrack) {
518 B2ERROR("Input recotrack is nullptr!");
519 return false;
520 }
522 Const::ChargedStable particleUsedForFitting(std::abs(trackPDG));
523 const genfit::AbsTrackRep* prevTrackRep = prevRecoTrack->getTrackRepresentationForPDG(std::abs(
524 trackPDG));
525
527 const std::vector<RecoHitInformation*>& recoHitInformations = recoTrack->getRecoHitInformations(true);
528 const std::vector<RecoHitInformation*>& prevRecoHitInformations = prevRecoTrack->getRecoHitInformations(
529 true);
530 unsigned int nRemoveHits = 0;
531 if (recoHitInformations.size() != prevRecoHitInformations.size()) {
532 B2WARNING("Copied RecoTrack has different number of hits from its original RecoTrack!");
533 return false;
534 }
536 for (nRemoveHits = 0; nRemoveHits < recoHitInformations.size(); ++nRemoveHits) {
537 if (!prevRecoHitInformations[nRemoveHits]->useInFit()) {
538 recoHitInformations[nRemoveHits]->setUseInFit(false);
539 continue;
540 }
541 try {
543 genfit::MeasuredStateOnPlane stPrevRecoHit = prevRecoTrack->getMeasuredStateOnPlaneFromRecoHit(
544 prevRecoHitInformations[nRemoveHits]);
547 double extralength = stPrevRecoHit.extrapolateToPoint(XYZToTVector(vertexPosition));
548 if (extralength > 0) {
549 recoHitInformations[nRemoveHits]->setUseInFit(false);
550 } else
551 break;
552 } catch (NoTrackFitResult()) {
553 B2WARNING("Exception: no FitterInfo assigned for TrackPoint created from this RecoHit.");
554 recoHitInformations[nRemoveHits]->setUseInFit(false);
555 continue;
556 } catch (...) {
560 B2DEBUG(22, "Could not extrapolate track to vertex when removing inner hits, aborting.");
561 return false;
562 }
563 }
564
565 if (nRemoveHits == 0) {
566 // This is not expected, but can happen if the track fit is different between copied and original RecoTrack.
567 B2DEBUG(20, "No hits removed in removeInnerHits, aborted. Switching to use the original RecoTrack.");
568 return false;
569 }
570
571 if (recoHitInformations.size() <= nRemoveHits) {
572 B2DEBUG(20, "Removed all the RecoHits in the RecoTrack, aborted. Switching to use the original RecoTrack.");
573 return false;
574 }
575
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>();
581 const SVDCluster* nextSVDHit = recoHitInformations[nRemoveHits] ->getRelatedTo<SVDCluster>();
582 if (!lastRemovedSVDHit || !nextSVDHit) B2ERROR("Last/Next SVD hit is null!");
583 else {
584 if (lastRemovedSVDHit->getSensorID() == nextSVDHit->getSensorID() &&
585 lastRemovedSVDHit->isUCluster() && !(nextSVDHit->isUCluster())) {
586 recoHitInformations[nRemoveHits]->setUseInFit(false);
587 ++nRemoveHits;
588 }
589 }
590 }
591 }
592
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);
602 nRemoveHits += 2;
603 }
604 }
605
606 B2DEBUG(22, nRemoveHits << " inner hits removed.");
607
609 TrackFitter fitter;
610 if (not fitter.fit(*recoTrack, particleUsedForFitting)) {
611 B2DEBUG(20, "track fit failed after removing inner hits.");
613 if (not prevRecoTrack->wasFitSuccessful(prevTrackRep))
614 B2DEBUG(20, "\t previous track fit was also failed.");
615 return false;
616 }
617
618 return true;
619}
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

◆ setFitterMode()

void setFitterMode ( int fitterMode)

set V0 fitter mode.

switch the mode of fitAndStore function. 0: store V0 at the first vertex fit, regardless of inner hits 1: remove hits inside the V0 vertex position 2: mode 1 + don't use SVD hits if there is only one available SVD hit-pair (default)

Definition at line 73 of file V0Fitter.cc.

74{
75 if (not(0 <= fitterMode && fitterMode <= 2)) {
76 B2FATAL("Invalid fitter mode!");
77 } else {
78 m_v0FitterMode = fitterMode;
79 if (fitterMode == 0)
80 m_forcestore = true;
81 else
82 m_forcestore = false;
83 if (fitterMode == 2)
85 else
87 }
88}

◆ vertexFitWithRecoTracks()

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 )
private

fit V0 vertex using RecoTrack's as inputs.

Fit V0 vertex using RecoTracks as inputs.

Return true (false) if the vertex fit has done well (failed). If RecoTracks have hits inside the fitted V0 vertex position, bits in hasInnerHitStatus are set. If there are no inside hits, store the V0 to the DataStore.

Parameters
trackPlusTrack of positively-charged daughter
trackMinusTrack of negatively-charged daughter
recoTrackPlusRecoTrack of positively-charged daughter
recoTrackMinusRecoTrack of negatively-charged daughter
v0HypothesisParticleType used in vertex fitting
hasInnerHitStatusstore a result of this function. if the plus(minus) track has hits inside the V0 vertex position, 0x1(0x2) bit is set.
vertexPosstore a result of this function. the fitted vertex position is stored.
forceStoreif true, store the fitted V0 to the DataStore even if there are some inside hits.
Returns

Store V0 and TrackFitResults if hasInnerHitStatus==0 or forceStore==true.

If existing, pass to the genfit::Track the correct cardinal representation

The two vectors should always have the same size, and this never happen

make a clone, not use the reference so that the genfit::MeasuredStateOnPlane and its TrackReps will not be altered.

Apply cuts. We have one set of cuts inside the beam pipe, the other outside.

store V0

Before storing, apply invariant mass cut.

Reconstruct invariant mass.

To build the trackFitResult, use the magnetic field at the origin = (0, 0, 0); the helix is extrapolated to the IP in a constant magnetic field and material effects are neglected so that the vertexing tool executed on the MDST object will find again this vertex position

Definition at line 344 of file V0Fitter.cc.

349{
350 const auto trackHypotheses = getTrackHypotheses(v0Hypothesis);
351
352 // make a clone, not use the reference so that the genfit::Track and its TrackReps will not be altered.
353 genfit::Track gfTrackPlus = RecoTrackGenfitAccess::getGenfitTrack(*recoTrackPlus);
354 const int pdgTrackPlus = trackPlus->getTrackFitResultWithClosestMass(trackHypotheses.first)->getParticleType().getPDGCode();
355 genfit::AbsTrackRep* plusRepresentation = recoTrackPlus->getTrackRepresentationForPDG(pdgTrackPlus);
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.");
358 return false;
359 }
360
361 // make a clone, not use the reference so that the genfit::Track and its TrackReps will not be altered.
362 genfit::Track gfTrackMinus = RecoTrackGenfitAccess::getGenfitTrack(*recoTrackMinus);
363 const int pdgTrackMinus = trackMinus->getTrackFitResultWithClosestMass(trackHypotheses.second)->getParticleType().getPDGCode();
364 genfit::AbsTrackRep* minusRepresentation = recoTrackMinus->getTrackRepresentationForPDG(pdgTrackMinus);
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.");
367 return false;
368 }
369
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);
379 }
380 }
382 else {
383 for (unsigned int id = 0; id < repsPlus.size(); id++) {
384 if (abs(repsPlus[id]->getPDG()) == pdgTrackPlus)
385 gfTrackPlus.setCardinalRep(id);
386 }
387 for (unsigned int id = 0; id < repsMinus.size(); id++) {
388 if (abs(repsMinus[id]->getPDG()) == pdgTrackMinus)
389 gfTrackMinus.setCardinalRep(id);
390 }
391 }
392
394 genfit::MeasuredStateOnPlane stPlus = recoTrackPlus->getMeasuredStateOnPlaneFromFirstHit(plusRepresentation);
395 genfit::MeasuredStateOnPlane stMinus = recoTrackMinus->getMeasuredStateOnPlaneFromFirstHit(minusRepresentation);
396
397 genfit::GFRaveVertex vert;
398 if (not fitGFRaveVertex(gfTrackPlus, gfTrackMinus, vert)) {
399 return false;
400 }
401
402 const ROOT::Math::XYZVector& posVert = ROOT::Math::XYZVector(vert.getPos());
403 vertexPos = posVert;
404
407 if (posVert.Rho() < m_beamPipeRadius) {
408 return false;
409 } else {
410 if (vert.getChi2() > m_vertexChi2CutOutside) {
411 B2DEBUG(22, "Vertex outside beam pipe, chi^2 too large.");
412 return false;
413 }
414 }
415
416 B2DEBUG(22, "Vertex accepted.");
417
418 if (not extrapolateToVertex(stPlus, stMinus, posVert, hasInnerHitStatus)) {
419 return false;
420 }
421
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();
431 if (v0Hypothesis == Const::Kshort) {
432 if (v0InvMass < std::get<0>(m_invMassRangeKshort) || v0InvMass > std::get<1>(m_invMassRangeKshort)) {
433 B2DEBUG(22, "Kshort vertex rejected, invariant mass out of range.");
434 return false;
435 }
436 } else if (v0Hypothesis == Const::Lambda || v0Hypothesis == Const::antiLambda) {
437 if (v0InvMass < std::get<0>(m_invMassRangeLambda) || v0InvMass > std::get<1>(m_invMassRangeLambda)) {
438 B2DEBUG(22, "Lambda vertex rejected, invariant mass out of range.");
439 return false;
440 }
441 } else if (v0Hypothesis == Const::photon) {
442 if (v0InvMass < std::get<0>(m_invMassRangePhoton) || v0InvMass > std::get<1>(m_invMassRangePhoton)) {
443 B2DEBUG(22, "Photon vertex rejected, invariant mass out of range.");
444 return false;
445 }
446 }
447
448
452 const double Bz = BFieldManager::getFieldInTesla({0, 0, 0}).Z();
453
454 int sharedInnermostCluster = checkSharedInnermostCluster(recoTrackPlus, recoTrackMinus);
455
456 TrackFitResult* tfrPlusVtx = buildTrackFitResult(gfTrackPlus, recoTrackPlus, stPlus, Bz, trackHypotheses.first,
457 sharedInnermostCluster);
458 TrackFitResult* tfrMinusVtx = buildTrackFitResult(gfTrackMinus, recoTrackMinus, stMinus, Bz, trackHypotheses.second,
459 sharedInnermostCluster);
460
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());
465
466 if (m_validation) {
467 B2DEBUG(24, "Create StoreArray and Output for validation.");
468 auto validationV0 = m_validationV0s.appendNew(
469 std::make_pair(trackPlus, tfrPlusVtx),
470 std::make_pair(trackMinus, tfrMinusVtx),
471 ROOT::Math::XYZVector(vert.getPos()),
472 vert.getCov(),
473 (lv0 + lv1).P(),
474 (lv0 + lv1).M(),
475 vert.getChi2()
476 );
477 v0->addRelationTo(validationV0);
478 }
479 }
480 return true;
481}
static ROOT::Math::XYZVector getFieldInTesla(const ROOT::Math::XYZVector &pos)
return the magnetic field at a given position in Tesla.
static genfit::Track & getGenfitTrack(RecoTrack &recoTrack)
Give access to the RecoTrack's genfit::Track.
Definition RecoTrack.cc:404
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
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:621
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:161
bool extrapolateToVertex(genfit::MeasuredStateOnPlane &stPlus, genfit::MeasuredStateOnPlane &stMinus, const ROOT::Math::XYZVector &vertexPosition)
Extrapolate the fit results to the perigee to the vertex.
bool fitGFRaveVertex(genfit::Track &trackPlus, genfit::Track &trackMinus, genfit::GFRaveVertex &vertex)
Fit the V0 vertex.
Definition V0Fitter.cc:104

Friends And Related Symbol Documentation

◆ V0FitterTest_EnableValidation_Test

friend class V0FitterTest_EnableValidation_Test
friend

Definition at line 41 of file V0Fitter.h.

◆ V0FitterTest_GetTrackHypotheses_Test

friend class V0FitterTest_GetTrackHypotheses_Test
friend

Definition at line 40 of file V0Fitter.h.

◆ V0FitterTest_InitializeCuts_Test

friend class V0FitterTest_InitializeCuts_Test
friend

Definition at line 42 of file V0Fitter.h.

Member Data Documentation

◆ m_beamPipeRadius

double m_beamPipeRadius
private

Radius where inside/outside beampipe is defined.

Definition at line 163 of file V0Fitter.h.

◆ m_copiedRecoTracks

StoreArray<RecoTrack> m_copiedRecoTracks
private

RecoTrack used to refit tracks (output)

Definition at line 161 of file V0Fitter.h.

◆ m_forcestore

bool m_forcestore
private

true only if the V0Fitter mode is 1

Definition at line 169 of file V0Fitter.h.

◆ m_invMassRangeKshort

std::tuple<double, double> m_invMassRangeKshort
private

invariant mass cut for Kshort.

Definition at line 165 of file V0Fitter.h.

◆ m_invMassRangeLambda

std::tuple<double, double> m_invMassRangeLambda
private

invariant mass cut for Lambda.

Definition at line 166 of file V0Fitter.h.

◆ m_invMassRangePhoton

std::tuple<double, double> m_invMassRangePhoton
private

invariant mass cut for Photon.

Definition at line 167 of file V0Fitter.h.

◆ m_recoTracks

StoreArray<RecoTrack> m_recoTracks
private

RecoTrack (input)

Definition at line 157 of file V0Fitter.h.

◆ m_recoTracksName

std::string m_recoTracksName
private

RecoTrackColName (input).

Definition at line 156 of file V0Fitter.h.

◆ m_trackFitResults

StoreArray<TrackFitResult> m_trackFitResults
private

TrackFitResult (output).

Definition at line 158 of file V0Fitter.h.

◆ m_useOnlyOneSVDHitPair

bool m_useOnlyOneSVDHitPair
private

false only if the V0Fitter mode is 3

Definition at line 170 of file V0Fitter.h.

◆ m_v0FitterMode

int m_v0FitterMode
private

0: store V0 at the first vertex fit, regardless of inner hits, 1: remove hits inside the V0 vertex position, 2: mode 1 + don't use SVD hits if there is only one available SVD hit-pair (default)

Definition at line 168 of file V0Fitter.h.

◆ m_v0s

StoreArray<V0> m_v0s
private

V0 (output).

Definition at line 159 of file V0Fitter.h.

◆ m_validation

bool m_validation
private

Validation flag.

Definition at line 155 of file V0Fitter.h.

◆ m_validationV0s

StoreArray<V0ValidationVertex> m_validationV0s
private

V0ValidationVertex (output, optional).

Definition at line 160 of file V0Fitter.h.

◆ m_vertexChi2CutOutside

double m_vertexChi2CutOutside
private

Chi2 cut outside beampipe.

Definition at line 164 of file V0Fitter.h.


The documentation for this class was generated from the following files: