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