Belle II Software development
NewV0Fitter.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/NewV0Fitter.h>
9
10#include <framework/logging/Logger.h>
11#include <framework/geometry/BFieldManager.h>
12#include <tracking/v0Finding/dataobjects/VertexVector.h>
13#include <tracking/trackFitting/fitter/base/TrackFitter.h>
14#include <tracking/trackFitting/trackBuilder/factories/TrackBuilder.h>
15
16#include <mdst/dataobjects/HitPatternVXD.h>
17#include <mdst/dataobjects/HitPatternCDC.h>
18
19#include <genfit/GFRaveVertexFactory.h>
20#include <genfit/FieldManager.h>
21#include <genfit/MaterialEffects.h>
22
23#include <framework/utilities/IOIntercept.h>
24
25using namespace Belle2;
26
27NewV0Fitter::NewV0Fitter(const std::string& trackFitResultsName, const std::string& v0sName,
28 const std::string& v0ValidationVerticesName, const std::string& recoTracksName,
29 const std::string& copiedRecoTracksName, bool enableValidation)
30 : m_recoTracksName(recoTracksName), m_validation(enableValidation)
31{
32 B2ASSERT("V0Fitter: material effects not set up. Please use SetupGenfitExtrapolationModule.",
33 genfit::MaterialEffects::getInstance()->isInitialized());
34 B2ASSERT("V0Fitter: magnetic field not set up. Please use SetupGenfitExtrapolationModule.",
35 genfit::FieldManager::getInstance()->isInitialized());
36
38 m_trackFitResults.isRequired(trackFitResultsName);
39
44
45 if (m_validation) {
46 m_validationV0s.registerInDataStore(v0ValidationVerticesName, DataStore::c_ErrorIfAlreadyRegistered);
47 m_v0s.registerRelationTo(m_validationV0s);
48 }
49
50}
51
52
53void NewV0Fitter::initializeCuts(double vertexDistanceCut,
54 double vertexChi2Cut,
55 const std::tuple<double, double>& invMassRangeKshort,
56 const std::tuple<double, double>& invMassRangeLambda,
57 const std::tuple<double, double>& invMassRangePhoton)
58{
59 m_vertexDistanceCut = vertexDistanceCut;
60 m_vertexChi2Cut = vertexChi2Cut;
61 m_invMassCuts[22] = std::make_pair(std::get<0>(invMassRangePhoton), std::get<1>(invMassRangePhoton));
62 m_invMassCuts[310] = std::make_pair(std::get<0>(invMassRangeKshort), std::get<1>(invMassRangeKshort));
63 m_invMassCuts[3122] = std::make_pair(std::get<0>(invMassRangeLambda), std::get<1>(invMassRangeLambda));
64}
65
66
67std::pair<Const::ParticleType, Const::ParticleType> NewV0Fitter::getTrackHypotheses(const Const::ParticleType& v0Hypothesis)
68{
69 if (v0Hypothesis == Const::Kshort) {
70 return std::make_pair(Const::pion, Const::pion);
71 } else if (v0Hypothesis == Const::photon) {
72 return std::make_pair(Const::electron, Const::electron);
73 } else if (v0Hypothesis == Const::Lambda) {
74 return std::make_pair(Const::proton, Const::pion);
75 } else if (v0Hypothesis == Const::antiLambda) {
76 return std::make_pair(Const::pion, Const::proton);
77 }
78 B2FATAL("V0Fitter: given V0 hypothesis not available.");
79 return std::make_pair(Const::invalidParticle, Const::invalidParticle); // return something to avoid triggering cppcheck
80}
81
82
83bool NewV0Fitter::fitAndStore(const Track* trackPlus, const Track* trackMinus, const Const::ParticleType& v0Hypothesis,
84 bool& isForceStored, bool& isHitRemoved)
85{
87 isForceStored = false;
88 isHitRemoved = false;
89 m_trkPlus.recoTrack = nullptr;
90 m_trkMinus.recoTrack = nullptr;
91
93 const RecoTrack* recoTrackPlus = trackPlus->getRelated<RecoTrack>(m_recoTracksName);
94 if (not recoTrackPlus) return false;
95 const RecoTrack* recoTrackMinus = trackMinus->getRelated<RecoTrack>(m_recoTracksName);
96 if (not recoTrackMinus) return false;
97
99 auto trackHypotheses = getTrackHypotheses(v0Hypothesis);
100 auto ptypeTrackPlus = trackPlus->getTrackFitResultWithClosestMass(trackHypotheses.first)->getParticleType();
101 auto ptypeTrackMinus = trackMinus->getTrackFitResultWithClosestMass(trackHypotheses.second)->getParticleType();
102 int pdgTrackPlus = std::abs(ptypeTrackPlus.getPDGCode());
103 int pdgTrackMinus = std::abs(ptypeTrackMinus.getPDGCode());
104
106 int status = vertexFit(recoTrackPlus, recoTrackMinus, pdgTrackPlus, pdgTrackMinus, v0Hypothesis);
107 if (status < 0) return false;
108
109 if (m_fitterMode > 0) {
111 const RecoTrack* origRecoTrackPlus = recoTrackPlus;
112 const RecoTrack* origRecoTrackMinus = recoTrackMinus;
113 int counter = 0;
114 while (status != 0 and counter < 5) {
115 counter++;
118 const RecoTrack* trkPlus = recoTrackPlus;
119 if (status & c_BitTrackPlus) {
120 trkPlus = removeHitsAndRefit(origRecoTrackPlus, recoTrackPlus, ptypeTrackPlus);
121 if (not trkPlus) return false;
122 isHitRemoved = true;
123 }
124 const RecoTrack* trkMinus = recoTrackMinus;
125 if (status & c_BitTrackMinus) {
126 trkMinus = removeHitsAndRefit(origRecoTrackMinus, recoTrackMinus, ptypeTrackMinus);
127 if (not trkMinus) return false;
128 isHitRemoved = true;
129 }
130 if (trkPlus == recoTrackPlus and trkMinus == recoTrackMinus) break; // vertex fit already done, so exit the loop
131 status = vertexFit(trkPlus, trkMinus, pdgTrackPlus, pdgTrackMinus, v0Hypothesis);
132 if (status < 0) break; // save the results of the last successful iteration
133 recoTrackPlus = trkPlus;
134 recoTrackMinus = trkMinus;
135 }
136
137 isForceStored = (status != 0);
138 }
139
141 int sharedCluster = isInnermostClusterShared(recoTrackPlus, recoTrackMinus);
142 const auto* fitPlus = saveTrackFitResult(m_trkPlus, sharedCluster);
143 if (not fitPlus) return false;
144 const auto* fitMinus = saveTrackFitResult(m_trkMinus, sharedCluster);
145 if (not fitMinus) return false;
146
147 auto* v0 = m_v0s.appendNew(std::make_pair(trackPlus, fitPlus),
148 std::make_pair(trackMinus, fitMinus),
149 m_fittedVertex.getPos().X(), m_fittedVertex.getPos().Y(), m_fittedVertex.getPos().Z());
150
151 if (m_validation) {
152 auto* validationV0 = m_validationV0s.appendNew(std::make_pair(trackPlus, fitPlus),
153 std::make_pair(trackMinus, fitMinus),
154 ROOT::Math::XYZVector(m_fittedVertex.getPos()),
155 m_fittedVertex.getCov(),
157 m_fittedVertex.getChi2());
158 v0->addRelationTo(validationV0);
159 }
160
161 return true;
162}
163
164
165int NewV0Fitter::vertexFit(const RecoTrack* recoTrackPlus, const RecoTrack* recoTrackMinus,
166 int pdgTrackPlus, int pdgTrackMinus, const Const::ParticleType& v0Hypothesis)
167{
168
169 // get track representations for given PDG codes and check their existence
170
171 const auto* plusRepresentation = getTrackRepresentation(recoTrackPlus, pdgTrackPlus);
172 if (not plusRepresentation) return c_NoTrackRepresentation;
173
174 const auto* minusRepresentation = getTrackRepresentation(recoTrackMinus, pdgTrackMinus);
175 if (not minusRepresentation) return c_NoTrackRepresentation;
176
177 // make copies of genfit tracks which will be passed to vertex fit and set the cardinal representations
178
179 auto gfTrackPlus = recoTrackPlus->getGenfitTrack(); // a copy of
180 if (not setCardinalRep(gfTrackPlus, pdgTrackPlus)) return c_NoTrackRepresentation;
181
182 auto gfTrackMinus = recoTrackMinus->getGenfitTrack(); // a copy of
183 if (not setCardinalRep(gfTrackMinus, pdgTrackMinus)) return c_NoTrackRepresentation;
184
185 // fit vertex
186
187 genfit::GFRaveVertex vert;
188 if (not fitGFRaveVertex(gfTrackPlus, gfTrackMinus, vert)) return c_VertexFitFailed;
189 auto vertexPos = ROOT::Math::XYZVector(vert.getPos());
190
191 // apply cuts on the vertex
192
193 if (vertexPos.Rho() < m_vertexDistanceCut) return c_NotSelected;
194 if (vert.getChi2() > m_vertexChi2Cut) return c_NotSelected;
195
196 // apply cut on the invariant mass
197
198 const auto& p1 = vert.getParameters(0)->getMom();
199 const auto& p2 = vert.getParameters(1)->getMom();
200 auto trackHypotheses = getTrackHypotheses(v0Hypothesis);
201 double mass1 = trackHypotheses.first.getMass();
202 double mass2 = trackHypotheses.second.getMass();
203 ROOT::Math::PxPyPzMVector fourVec1(p1.X(), p1.Y(), p1.Z(), mass1);
204 ROOT::Math::PxPyPzMVector fourVec2(p2.X(), p2.Y(), p2.Z(), mass2);
205 double invMass = (fourVec1 + fourVec2).M();
206 int pdgCode = abs(v0Hypothesis.getPDGCode());
207 const auto& cuts = m_invMassCuts[pdgCode];
208 if (invMass < cuts.first or invMass > cuts.second) return c_NotSelected;
209
210 // extrapolate tracks to fitted vertex; the return status, if positive, indicates whether there are inner hits
211
212 auto statePlus = recoTrackPlus->getMeasuredStateOnPlaneFromFirstHit(plusRepresentation); // a copy of
213 auto stateMinus = recoTrackMinus->getMeasuredStateOnPlaneFromFirstHit(minusRepresentation); // a copy of
214 int status = extrapolateToVertex(statePlus, stateMinus, vert);
215 if (status < 0) return c_ExtrapolationFailed;
216
217 // save fitted vertex and tracks
218
219 m_fittedVertex = vert;
220 m_momentum = (fourVec1 + fourVec2).P();
221 m_invMass = invMass;
222 m_trkPlus.set(recoTrackPlus, trackHypotheses.first, statePlus, plusRepresentation);
223 m_trkMinus.set(recoTrackMinus, trackHypotheses.second, stateMinus, minusRepresentation);
224
225 return status;
226}
227
228
229const genfit::AbsTrackRep* NewV0Fitter::getTrackRepresentation(const RecoTrack* recoTrack, int pdgCode)
230{
231 const auto* rep = recoTrack->getTrackRepresentationForPDG(pdgCode);
232 if (rep and recoTrack->wasFitSuccessful(rep)) return rep;
233
234 B2ERROR("V0Fitter: track hypothesis with closest mass not available. Should never happen!");
235 return nullptr;
236}
237
238
239bool NewV0Fitter::setCardinalRep(genfit::Track& gfTrack, int pdgCode)
240{
241 const auto& reps = gfTrack.getTrackReps();
242 for (unsigned id = 0; id < reps.size(); id++) {
243 if (abs(reps[id]->getPDG()) == pdgCode) {
244 gfTrack.setCardinalRep(id);
245 return true;
246 }
247 }
248
249 B2ERROR("V0Fitter: cannot set cardinal representation for PDG = " << pdgCode);
250 return false;
251}
252
253
254bool NewV0Fitter::fitGFRaveVertex(genfit::Track& trackPlus, genfit::Track& trackMinus, genfit::GFRaveVertex& vertex)
255{
256 VertexVector vertexVector;
257 std::vector<genfit::Track*> trackPair {&trackPlus, &trackMinus};
258
259 try {
261 logCapture("V0Fitter GFRaveVertexFactory", LogConfig::c_Debug, LogConfig::c_Debug);
262 logCapture.start();
263
264 genfit::GFRaveVertexFactory vertexFactory;
265 vertexFactory.findVertices(&vertexVector.v, trackPair);
266 } catch (...) {
267 B2ERROR("V0Fitter: exception during vertex fit.");
268 return false;
269 }
270
271 if (vertexVector.size() != 1) {
272 B2DEBUG(21, "Vertex fit failed. Size of vertexVector not 1, but: " << vertexVector.size());
273 return false;
274 }
275
276 if ((*vertexVector[0]).getNTracks() != 2) {
277 B2DEBUG(20, "Wrong number of tracks in vertex.");
278 return false;
279 }
280
281 vertex = *vertexVector[0];
282 return true;
283}
284
285
286int NewV0Fitter::extrapolateToVertex(genfit::MeasuredStateOnPlane& statePlus, genfit::MeasuredStateOnPlane& stateMinus,
287 const genfit::GFRaveVertex& vertex)
288{
289 try {
292 double extralengthPlus = statePlus.extrapolateToPoint(vertex.getPos());
293 double extralengthMinus = stateMinus.extrapolateToPoint(vertex.getPos());
294 unsigned status = 0;
295 if (extralengthPlus > 0) status |= c_BitTrackPlus;
296 if (extralengthMinus > 0) status |= c_BitTrackMinus;
297 return status;
298 } catch (...) {
302 B2DEBUG(22, "Could not extrapolate track to vertex.");
304 }
305}
306
307
308const RecoTrack* NewV0Fitter::removeHitsAndRefit(const RecoTrack* origRecoTrack, const RecoTrack* lastRecoTrack,
309 const Const::ParticleType& ptype)
310{
311 // make a copy of useInFit flags
312 std::vector<bool> useInFit;
313 const auto& recoHitInformations = origRecoTrack->getRecoHitInformations(true); // true to get sorted hits info
314 useInFit.reserve(recoHitInformations.size());
315 for (const auto& hitInfo : recoHitInformations) useInFit.push_back(hitInfo->useInFit());
316
317 // get track representation for a given particle
318 const auto* rep = getTrackRepresentation(origRecoTrack, abs(ptype.getPDGCode()));
319 if (not rep) return lastRecoTrack;
320
321 // remove inner hits by setting useInFit flags to false
322 int removedHits = 0;
323 unsigned firstHit = 0;
324 for (unsigned i = 0; i < recoHitInformations.size(); i++) {
325 const auto& hitInfo = recoHitInformations[i];
326 if (not hitInfo->useInFit()) continue;
327 try {
328 auto state = origRecoTrack->getMeasuredStateOnPlaneFromRecoHit(hitInfo, rep); // a copy of
329 double extraLength = state.extrapolateToPoint(m_fittedVertex.getPos());
330 if (extraLength > 0) {
331 useInFit[i] = false;
332 removedHits++;
333 } else {
334 firstHit = i;
335 break;
336 }
337 } catch (NoTrackFitResult()) {
338 B2WARNING("V0Fitter exception: no FitterInfo assigned for TrackPoint created from this RecoHit.");
339 useInFit[i] = false;
340 removedHits++;
341 continue;
342 } catch (...) {
346 B2DEBUG(22, "Could not extrapolate track to vertex when removing inner hits.");
347 return lastRecoTrack;
348 }
349 }
350
351 if (m_fitterMode == 2) {
352 // remove SVD hits if there is only one or two left just after the vertex
353 std::vector<unsigned> svdIndex;
354 for (unsigned i = 0; i < useInFit.size(); i++) {
355 if (not useInFit[i]) continue;
356 const auto& hitInfo = recoHitInformations[i];
357 if (hitInfo->getTrackingDetector() == RecoHitInformation::c_SVD) svdIndex.push_back(i);
358 else break;
359 }
360 if (not svdIndex.empty() and svdIndex.size() < 3) {
361 for (unsigned i : svdIndex) {
362 useInFit[i] = false;
363 removedHits++;
364 }
365 }
366 }
367
368 if (removedHits == 0) return origRecoTrack; // in this case track doesn't need to be refitted
369
370 // count remaining hits and return if there is no hope for the fit to succeed
371 int nHits = 0;
372 for (auto x : useInFit) if (x) nHits++;
373 if (nHits < 5) return lastRecoTrack;
374
375 // make a copy of recoTrack
376 const auto& state = origRecoTrack->getMeasuredStateOnPlaneFromRecoHit(recoHitInformations[firstHit], rep);
377 auto* recoTrack_copy = copyRecoTrack(origRecoTrack, state);
378
379 // set useInFit flags in a copy of recoTrack
380 const auto& recoHitInfos = recoTrack_copy->getRecoHitInformations(true); // true to get sorted hits info
381 if (recoHitInfos.size() != useInFit.size()) {
382 B2ERROR("V0Fitter: copied recoTrack has different number of hits than the original one");
383 return lastRecoTrack;
384 }
385 for (unsigned i = 0; i < recoHitInfos.size(); i++) recoHitInfos[i]->setUseInFit(useInFit[i]);
386
387 // fit a copy of recoTrack
388 TrackFitter fitter;
389 bool ok = fitter.fit(*recoTrack_copy, ptype);
390 if (not ok) return lastRecoTrack;
391
392 return recoTrack_copy;
393}
394
395
396RecoTrack* NewV0Fitter::copyRecoTrack(const RecoTrack* origRecoTrack, const genfit::MeasuredStateOnPlane& state)
397{
398 RecoTrack* newRecoTrack = origRecoTrack->copyToStoreArrayUsing(m_copiedRecoTracks,
399 ROOT::Math::XYZVector(state.getPos()),
400 ROOT::Math::XYZVector(state.getMom()),
401 state.getCharge(),
402 state.get6DCov(), state.getTime());
403 newRecoTrack->addHitsFromRecoTrack(origRecoTrack);
404 newRecoTrack->addRelationTo(origRecoTrack);
405 return newRecoTrack;
406}
407
408
409int NewV0Fitter::isInnermostClusterShared(const RecoTrack* recoTrackPlus, const RecoTrack* recoTrackMinus)
410{
411 const auto& recoHitInformationsPlus = recoTrackPlus->getRecoHitInformations(true); // true to get sorted hits info
412 std::vector<RecoHitInformation*> innerHitsPlus;
413 for (const auto& hitInfo : recoHitInformationsPlus) {
414 if (hitInfo->useInFit()) innerHitsPlus.push_back(hitInfo);
415 if (innerHitsPlus.size() == 2) break;
416 }
417 if (innerHitsPlus.empty()) return 0;
418
419 const auto& recoHitInformationsMinus = recoTrackMinus->getRecoHitInformations(true); // true to get sorted hits info
420 std::vector<RecoHitInformation*> innerHitsMinus;
421 for (const auto& hitInfo : recoHitInformationsMinus) {
422 if (hitInfo->useInFit()) innerHitsMinus.push_back(hitInfo);
423 if (innerHitsMinus.size() == 2) break;
424 }
425 if (innerHitsMinus.empty()) return 0;
426
427 if (innerHitsPlus.front()->getTrackingDetector() != innerHitsMinus.front()->getTrackingDetector()) return 0;
428
429 if (innerHitsPlus.front()->getTrackingDetector() == RecoHitInformation::c_PXD) {
430 const auto* clusterPlus = innerHitsPlus.front()->getRelatedTo<PXDCluster>();
431 const auto* clusterMinus = innerHitsMinus.front()->getRelatedTo<PXDCluster>();
432 if (clusterPlus and clusterPlus == clusterMinus) return 0x03; // PXD cluster the same: set both bits
433 return 0;
434 }
435
436 int flag = 0;
437 VxdID sensorID = 0;
438 if (innerHitsPlus.front()->getTrackingDetector() == RecoHitInformation::c_SVD) {
439 const auto* clusterPlus = innerHitsPlus.front()->getRelatedTo<SVDCluster>();
440 const auto* clusterMinus = innerHitsMinus.front()->getRelatedTo<SVDCluster>();
441 if (clusterPlus and clusterPlus == clusterMinus) {
442 sensorID = clusterPlus->getSensorID();
443 if (clusterPlus->isUCluster()) flag = 0x01; // SVD U-cluster the same: set first bit
444 else flag = 0x02; // SVD V-cluster the same: set second bit
445 }
446 }
447
448 if (innerHitsPlus.size() != 2 or innerHitsMinus.size() != 2) return flag;
449
450 if (innerHitsPlus.back()->getTrackingDetector() == RecoHitInformation::c_SVD) {
451 const auto* clusterPlus = innerHitsPlus.back()->getRelatedTo<SVDCluster>();
452 const auto* clusterMinus = innerHitsMinus.back()->getRelatedTo<SVDCluster>();
453 if (clusterPlus and clusterPlus == clusterMinus and clusterPlus->getSensorID() == sensorID) {
454 if (clusterPlus->isUCluster()) flag |= 0x01; // SVD U-cluster the same: set first bit
455 else flag |= 0x02; // SVD V-cluster the same: set second bit
456 }
457 }
458
459 return flag;
460}
461
462
463const TrackFitResult* NewV0Fitter::saveTrackFitResult(const FittedTrack& trk, int sharedInnermostCluster)
464{
465
466 const auto* recoTrack = trk.recoTrack;
467 if (not recoTrack) {
468 B2ERROR("V0Fitter: bug in saving track fit result, recoTrack is nullptr");
469 return nullptr;
470 }
471
472 auto hitPatternCDC = TrackBuilder::getHitPatternCDCInitializer(*recoTrack);
473 auto hitPatternVXD = TrackBuilder::getHitPatternVXDInitializer(*recoTrack);
474 if (sharedInnermostCluster > 0) {
475 auto pattern = HitPatternVXD(hitPatternVXD);
476 pattern.setInnermostHitShareStatus(sharedInnermostCluster);
477 hitPatternVXD = pattern.getInteger();
478 }
479 double Bz = BFieldManager::getFieldInTesla({0, 0, 0}).Z();
480 const auto& state = trk.state;
481
482 auto* fit = m_trackFitResults.appendNew(ROOT::Math::XYZVector(state.getPos()), ROOT::Math::XYZVector(state.getMom()),
483 state.get6DCov(), state.getCharge(), trk.ptype, trk.pValue,
484 Bz, hitPatternCDC, hitPatternVXD, trk.Ndf);
485 return fit;
486}
static ROOT::Math::XYZVector getFieldInTesla(const ROOT::Math::XYZVector &pos)
return the magnetic field at a given position in Tesla.
Definition: BFieldManager.h:61
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
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
StoreArray< V0ValidationVertex > m_validationV0s
V0ValidationVertex collection (optional)
Definition: NewV0Fitter.h:235
StoreArray< RecoTrack > m_copiedRecoTracks
copied RecoTracks collection
Definition: NewV0Fitter.h:236
FittedTrack m_trkPlus
positively charged track data of last successfully fitted vertex
Definition: NewV0Fitter.h:250
bool m_validation
validation flag
Definition: NewV0Fitter.h:243
StoreArray< V0 > m_v0s
V0s collection.
Definition: NewV0Fitter.h:234
double m_vertexDistanceCut
cut on the transverse radius
Definition: NewV0Fitter.h:238
double m_vertexChi2Cut
Chi2 cut.
Definition: NewV0Fitter.h:239
StoreArray< TrackFitResult > m_trackFitResults
TrackFitResults collection.
Definition: NewV0Fitter.h:233
bool fitAndStore(const Track *trackPlus, const Track *trackMinus, const Const::ParticleType &v0Hypothesis, bool &isForceStored, bool &isHitRemoved)
Fit V0 with given hypothesis and store results if fit is successful.
Definition: NewV0Fitter.cc:83
int extrapolateToVertex(genfit::MeasuredStateOnPlane &statePlus, genfit::MeasuredStateOnPlane &stateMinus, const genfit::GFRaveVertex &vertex)
Extrapolation of both tracks to the vertex.
Definition: NewV0Fitter.cc:286
double m_momentum
momentum of last successfully fitted vertex
Definition: NewV0Fitter.h:248
@ c_NoTrackRepresentation
no track representation for given PDG code
Definition: NewV0Fitter.h:129
@ c_NotSelected
fitted vertex not passing the cuts
Definition: NewV0Fitter.h:132
@ c_VertexFitFailed
vertex fit failed
Definition: NewV0Fitter.h:130
@ c_ExtrapolationFailed
track extrapolation failed
Definition: NewV0Fitter.h:131
bool setCardinalRep(genfit::Track &gfTrack, int pdgCode)
Sets cardinal representation of a given genfit track and PDG code.
Definition: NewV0Fitter.cc:239
std::string m_recoTracksName
name of the RecoTracks collection
Definition: NewV0Fitter.h:231
const TrackFitResult * saveTrackFitResult(const FittedTrack &trk, int sharedInnermostCluster)
Append track fit result to the collection.
Definition: NewV0Fitter.cc:463
int isInnermostClusterShared(const RecoTrack *recoTrackPlus, const RecoTrack *recoTrackMinus)
Returns bit flags indicating that the innermost cluster is shared between both tracks.
Definition: NewV0Fitter.cc:409
const genfit::AbsTrackRep * getTrackRepresentation(const RecoTrack *recoTrack, int pdgCode)
Returns track representation for a given PDG code.
Definition: NewV0Fitter.cc:229
std::map< int, std::pair< double, double > > m_invMassCuts
invariant mass cuts, key = abs(PDG)
Definition: NewV0Fitter.h:240
NewV0Fitter(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.
Definition: NewV0Fitter.cc:27
StoreArray< RecoTrack > m_recoTracks
RecoTracks collection.
Definition: NewV0Fitter.h:232
genfit::GFRaveVertex m_fittedVertex
last successfully fitted vertex
Definition: NewV0Fitter.h:247
static std::pair< Const::ParticleType, Const::ParticleType > getTrackHypotheses(const Const::ParticleType &v0Hypothesis)
Returns daughter particle types for a given V0 hypothesis.
Definition: NewV0Fitter.cc:67
bool fitGFRaveVertex(genfit::Track &trackPlus, genfit::Track &trackMinus, genfit::GFRaveVertex &vertex)
Genfit Rave vertex fit called by vertexFit method.
Definition: NewV0Fitter.cc:254
void initializeCuts(double vertexDistanceCut, double vertexChi2Cut, const std::tuple< double, double > &invMassRangeKshort, const std::tuple< double, double > &invMassRangeLambda, const std::tuple< double, double > &invMassRangePhoton)
Initialization of cuts applied during the fit and store process.
Definition: NewV0Fitter.cc:53
double m_invMass
invariant mass of last successfully fitted vertex
Definition: NewV0Fitter.h:249
int m_fitterMode
fitter mode
Definition: NewV0Fitter.h:242
FittedTrack m_trkMinus
negatively charged track data of last successfully fitted vertex
Definition: NewV0Fitter.h:251
@ c_BitTrackMinus
negative track has inner hits
Definition: NewV0Fitter.h:141
@ c_BitTrackPlus
positive track has inner hits
Definition: NewV0Fitter.h:140
int vertexFit(const RecoTrack *recoTrackPlus, const RecoTrack *recoTrackMinus, int pdgTrackPlus, int pdgTrackMinus, const Const::ParticleType &v0Hypothesis)
Performs a vertex fit.
Definition: NewV0Fitter.cc:165
RecoTrack * copyRecoTrack(const RecoTrack *origRecoTrack, const genfit::MeasuredStateOnPlane &state)
Make a copy of reco track.
Definition: NewV0Fitter.cc:396
const RecoTrack * removeHitsAndRefit(const RecoTrack *origRecoTrack, const RecoTrack *lastRecoTrack, const Const::ParticleType &ptype)
Remove track inner hits and refit the track.
Definition: NewV0Fitter.cc:308
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
Definition: PXDCluster.h:30
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 * copyToStoreArrayUsing(StoreArray< RecoTrack > &storeArray, const ROOT::Math::XYZVector &position, const ROOT::Math::XYZVector &momentum, short charge, const TMatrixDSym &covariance, double timeSeed) const
Append a new RecoTrack to the given store array and copy its general properties, but not the hits the...
Definition: RecoTrack.cc:510
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
const genfit::Track & getGenfitTrack() const
Returns genfit track.
Definition: RecoTrack.h:505
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 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
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
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
Abstract base class for different kinds of events.
Structure to save track data of the last successful iteration.
Definition: NewV0Fitter.h:37
void set(const RecoTrack *recoTrk, const Const::ParticleType &hypo, const genfit::MeasuredStateOnPlane &mSoP, const genfit::AbsTrackRep *rep)
Sets the data members.
Definition: NewV0Fitter.h:51
double pValue
p-value of track fit
Definition: NewV0Fitter.h:41
genfit::MeasuredStateOnPlane state
measured state at first hit, extrapolated to fitted vertex
Definition: NewV0Fitter.h:40
Const::ParticleType ptype
particle type of the V0 track
Definition: NewV0Fitter.h:39
const RecoTrack * recoTrack
reco track
Definition: NewV0Fitter.h:38
int Ndf
degrees-of-freedom of track fit
Definition: NewV0Fitter.h:42