Belle II Software  release-05-01-25
TrackSZFitter.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2016 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Oliver Frost *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 #include <tracking/trackFindingCDC/findlets/minimal/TrackSZFitter.h>
11 
12 #include <tracking/trackFindingCDC/eventdata/tracks/CDCTrack.h>
13 #include <tracking/trackFindingCDC/fitting/CDCSZFitter.h>
14 #include <tracking/trackFindingCDC/fitting/CDCAxialStereoFusion.h>
15 
16 #include <tracking/trackFindingCDC/eventdata/segments/CDCSegment2D.h>
17 
18 #include <tracking/trackFindingCDC/eventdata/hits/CDCRecoHit3D.h>
19 #include <tracking/trackFindingCDC/eventdata/hits/CDCRecoHit2D.h>
20 
21 #include <tracking/trackFindingCDC/eventdata/trajectories/CDCTrajectory3D.h>
22 #include <tracking/trackFindingCDC/eventdata/trajectories/CDCTrajectory2D.h>
23 #include <tracking/trackFindingCDC/eventdata/trajectories/CDCTrajectorySZ.h>
24 
25 #include <tracking/trackFindingCDC/geometry/UncertainHelix.h>
26 #include <tracking/trackFindingCDC/geometry/Vector2D.h>
27 
28 #include <vector>
29 
30 using namespace Belle2;
31 using namespace TrackFindingCDC;
32 
34 {
35  return "Use an SZFitter to create the 3D trajectory out of the 2D one.";
36 }
37 
38 void TrackSZFitter::apply(std::vector<CDCTrack>& tracks)
39 {
40  // Postprocess each track (=fit)
41  for (CDCTrack& track : tracks) {
42  const CDCSZFitter& szFitter = CDCSZFitter::getFitter();
43 
44  track.shiftToPositiveArcLengths2D();
45  track.sortByArcLength2D();
46 
47  CDCTrajectory2D originalTrajectory2D = track.getStartTrajectory3D().getTrajectory2D();
48  Vector2D localOrigin = originalTrajectory2D.getLocalOrigin();
49 
50  const CDCTrajectorySZ& szTrajectory = szFitter.fitWithStereoHits(track);
51  CDCTrajectory3D preliminaryTrajectory3D{originalTrajectory2D, szTrajectory};
52  track.setStartTrajectory3D(preliminaryTrajectory3D);
53 
54  // Estimate a better covariance matrix
55  CDCAxialStereoFusion axialStereoFusionFitter;
56  CDCSegment2D axialSegment2D;
57  CDCSegment2D stereoSegment2D;
58  for (const CDCRecoHit3D& recoHit3D : track) {
59  if (recoHit3D.isAxial()) {
60  axialSegment2D.push_back(recoHit3D.getRecoHit2D());
61  } else {
62  stereoSegment2D.push_back(recoHit3D.getRecoHit2D());
63  }
64  if ((axialSegment2D.size() > 6) and (stereoSegment2D.size() > 6)) break;
65  }
66 
67  if (not((axialSegment2D.size() > 6) and (stereoSegment2D.size() > 6))) continue;
68 
69  CDCTrajectory3D trajectory3D =
70  axialStereoFusionFitter.reconstructFuseTrajectories(axialSegment2D, stereoSegment2D, preliminaryTrajectory3D);
71  trajectory3D.setLocalOrigin({localOrigin, 0});
72 
73  // Copy only the covariance matrix, chi2 and ndf over for a conservative introduction for the moment.
74  UncertainHelix preliminaryUncertainHelix = preliminaryTrajectory3D.getLocalHelix();
75  UncertainHelix uncertainHelix = trajectory3D.getLocalHelix();
76 
77  preliminaryUncertainHelix.setHelixCovariance(uncertainHelix.helixCovariance());
78  preliminaryUncertainHelix.setChi2(uncertainHelix.chi2());
79  preliminaryUncertainHelix.setNDF(uncertainHelix.ndf());
80  preliminaryTrajectory3D.setLocalHelix(preliminaryUncertainHelix);
81 
82  track.setStartTrajectory3D(preliminaryTrajectory3D);
83  }
84 }
Belle2::TrackFindingCDC::UncertainHelix::setHelixCovariance
void setHelixCovariance(const HelixCovariance &helixCovariance)
Setter for the whole covariance matrix of the perigee parameters.
Definition: UncertainHelix.h:249
Belle2::TrackFindingCDC::CDCRecoHit3D
Class representing a three dimensional reconstructed hit.
Definition: CDCRecoHit3D.h:62
Belle2::TrackFindingCDC::CDCTrack
Class representing a sequence of three dimensional reconstructed hits.
Definition: CDCTrack.h:51
Belle2::TrackFindingCDC::UncertainHelix::chi2
double chi2() const
Getter for the chi square value of the helix fit.
Definition: UncertainHelix.h:273
Belle2::TrackFindingCDC::Vector2D
A two dimensional vector which is equipped with functions for correct handeling of orientation relat...
Definition: Vector2D.h:37
Belle2::TrackFindingCDC::CDCSZFitter::getFitter
static const CDCSZFitter & getFitter()
Getter for a standard sz line fitter instance.
Definition: CDCSZFitter.cc:38
Belle2::TrackFindingCDC::UncertainHelix::setChi2
void setChi2(const double chi2)
Setter for the chi square value of the helix fit.
Definition: UncertainHelix.h:279
Belle2::TrackFindingCDC::UncertainHelix::ndf
std::size_t ndf() const
Getter for the number of degrees of freediom used in the helix fit.
Definition: UncertainHelix.h:285
Belle2::TrackFindingCDC::CDCAxialStereoFusion::reconstructFuseTrajectories
void reconstructFuseTrajectories(const CDCSegmentPair &segmentPair)
Combine the two trajectories of the segments in the pair and assign the resulting three dimensional t...
Definition: CDCAxialStereoFusion.cc:43
Belle2::TrackFindingCDC::CDCTrajectory3D::setLocalOrigin
double setLocalOrigin(const Vector3D &localOrigin)
Setter for the origin of the local coordinate system.
Definition: CDCTrajectory3D.cc:369
Belle2::TrackFindingCDC::UncertainHelix::helixCovariance
const HelixCovariance & helixCovariance() const
Getter for the whole covariance matrix of the perigee parameters.
Definition: UncertainHelix.h:255
Belle2::TrackFindingCDC::CDCTrajectorySZ
Linear trajectory in sz space.
Definition: CDCTrajectorySZ.h:41
Belle2::TrackFindingCDC::CDCTrajectory2D
Particle trajectory as it is seen in xy projection represented as a circle.
Definition: CDCTrajectory2D.h:46
Belle2::TrackFindingCDC::CDCTrajectory2D::getLocalOrigin
const Vector2D & getLocalOrigin() const
Getter for the origin of the local coordinate system.
Definition: CDCTrajectory2D.h:508
Belle2::TrackFindingCDC::CDCSZFitter::fitWithStereoHits
CDCTrajectorySZ fitWithStereoHits(const CDCTrack &track) const
Returns the fitted sz trajectory of the track with the z-information of all stereo hits of the number...
Definition: CDCSZFitter.cc:127
Belle2::TrackFindingCDC::TrackSZFitter::apply
void apply(std::vector< CDCTrack > &tracks) final
Fit the tracks.
Definition: TrackSZFitter.cc:38
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::TrackFindingCDC::CDCSZFitter
Class implementing the z coordinate over travel distance line fit.
Definition: CDCSZFitter.h:37
Belle2::TrackFindingCDC::TrackSZFitter::getDescription
std::string getDescription() final
Short description of the findlet.
Definition: TrackSZFitter.cc:33
Belle2::TrackFindingCDC::CDCAxialStereoFusion
Utility class implementing the Kalmanesk combination of to two dimensional trajectories to one three ...
Definition: CDCAxialStereoFusion.h:40
Belle2::TrackFindingCDC::UncertainHelix::setNDF
void setNDF(std::size_t ndf)
Setter for the number of degrees of freediom used in the helix fit.
Definition: UncertainHelix.h:291
Belle2::TrackFindingCDC::CDCSegment2D
A reconstructed sequence of two dimensional hits in one super layer.
Definition: CDCSegment2D.h:40
Belle2::TrackFindingCDC::UncertainHelix
A general helix class including a covariance matrix.
Definition: UncertainHelix.h:44
Belle2::TrackFindingCDC::CDCTrajectory3D::getLocalHelix
const UncertainHelix & getLocalHelix() const
Getter for the helix in local coordinates.
Definition: CDCTrajectory3D.h:341
Belle2::TrackFindingCDC::CDCTrajectory3D
Particle full three dimensional trajectory.
Definition: CDCTrajectory3D.h:47