ATLAS Offline Software
GaussianTrackDensity.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 #ifndef TRKVERTEXSEEDFINDERUTILIS_GAUSSIANTRACKDENSITY_H
6 #define TRKVERTEXSEEDFINDERUTILIS_GAUSSIANTRACKDENSITY_H
7 
10 
11 #include <map>
12 
14 namespace Trk
15 {
16 
17  class Track;
18  class GaussianTrackDensity;
19 
31  class GaussianTrackDensity final: public extends<AthAlgTool, IVertexTrackDensityEstimator>
32  {
33  public:
35  using base_class::base_class;
36 
41  virtual double
42  globalMaximum (const std::vector<const Track*>& vectorTrk) const override final;
43 
49  virtual double
50  globalMaximum (const std::vector<const Track*>& vectorTrk,
51  std::unique_ptr<ITrackDensity>& density) const override final;
52 
57  virtual double
58  globalMaximum (const std::vector<const TrackParameters*>& perigeeList) const override final;
59 
65  virtual double
66  globalMaximum (const std::vector<const TrackParameters*>& perigeeList,
67  std::unique_ptr<ITrackDensity>& density) const override final;
68 
69  virtual std::pair<double, double> globalMaximumWithWidth(
70  const std::vector<const TrackParameters*>& perigeeList)
71  const override final;
72 
73  private:
74  struct TrackEntry {
75  TrackEntry() = default;
76  TrackEntry(const TrackEntry&) = default;
77  TrackEntry(TrackEntry&&) = default;
78  TrackEntry& operator=(const TrackEntry&) = default;
80  ~TrackEntry() = default;
81 
82  TrackEntry(double c0, double c1, double c2, double zMin, double zMax)
83  : c_0(c0), c_1(c1), c_2(c2), lowerBound(zMin), upperBound(zMax) {}
84  explicit TrackEntry(double z)
85  : c_0(0), c_1(0), c_2(0), lowerBound(z), upperBound(z) {}
86  // Cached information for a single track
87  double c_0 = 0; // z-independent term in exponent
88  double c_1 = 1; // linear coefficient in exponent
89  double c_2 = 0; // quadratic coefficient in exponent
90  double lowerBound = 0;
91  double upperBound = 0;
92  };
93 
94  // helper to handle the evaluation of the parametrised track density
96  public:
97  // initialise with the z coordinate at which the density is to be evaluated
98  TrackDensityEval(double z_coordinate): m_z(z_coordinate){}
99  // add the contribution for one track to the density.
100  // will internally check if the z coordinate is within
101  // the bounds of the track
102  void addTrack (const TrackEntry & entry);
103  // retrieve the density and its derivatives
104  inline double density() const {return m_density;}
105  inline double firstDerivative() const {return m_firstDerivative;}
106  inline double secondDerivative() const {return m_secondDerivative;}
107  private:
108  double m_z;
109  double m_density{0};
110  double m_firstDerivative{0};
111  double m_secondDerivative{0};
112  };
113 
114  class TrackDensity final : public ITrackDensity {
115  public:
116  explicit TrackDensity(bool gaussStep) : m_gaussStep(gaussStep) {}
117  virtual ~TrackDensity() = default;
118 
123  virtual double trackDensity(double z) const override final;
124 
129  virtual void trackDensity(
130  double z, double& density, double& firstDerivative,
131  double& secondDerivative) const override final;
135  double globalMaximum() const;
140  std::pair<double, double> globalMaximumWithWidth() const;
147  void addTrack(const Perigee& itrk, const double d0SignificanceCut,
148  const double z0SignificanceCut);
149 
150  private:
151  inline void updateMaximum(double trialZ, double trialValue,
152  double secondDerivative, double& maxZ,
153  double& maxValue,
154  double& maxSecondDerivative) const {
155  if (trialValue > maxValue) {
156  maxZ = trialZ;
157  maxValue = trialValue;
158  maxSecondDerivative = secondDerivative;
159  }
160  }
161 
162  inline double stepSize(double y, double dy, double ddy) const {
163  return (m_gaussStep ? (y * dy) / (dy * dy - y * ddy) : -dy / ddy);
164  }
166  double m_maxRange = 0;
167  // Cache for track information
168  // functor to compare two Perigee values
169  struct pred_perigee {
170  inline bool operator()(const Perigee& left,
171  const Perigee& right) const {
172  return left.parameters()[Trk::z0] < right.parameters()[Trk::z0];
173  }
174  };
175  // functor to compare two TrackEntry values based on their upper limits (low
176  // to high)
178  inline bool operator()(const TrackEntry& left,
179  const TrackEntry& right) const {
180  return left.upperBound < right.upperBound;
181  }
182  };
183 
184  using trackMap =
186  SG::ArenaPoolSTLAllocator<std::pair<
188 
189  using lowerMap = std::map<
192  std::pair<const GaussianTrackDensity::TrackEntry, Perigee>>>;
193 
196  };
197 
203  double
204  globalMaximumImpl (const std::vector<const TrackParameters*>& perigeeList,
205  TrackDensity& density) const;
206 
212  std::pair<double,double>
213  globalMaximumWithWidthImpl (const std::vector<const TrackParameters*>& perigeeList,
214  TrackDensity& density) const;
215 
216 
222  void addTracks(const std::vector<const TrackParameters*>& perigeeList,
223  TrackDensity& density) const;
224 
225 
226 
227  // Cuts set by configurable properties
228 
229  // Maximum allowed d0 significance to use (in sigma)
230  Gaudi::Property<double> m_d0MaxSignificance{
231  this,
232  "MaxD0Significance",
233  3.5,
234  "Maximum radial impact parameter significance to use track"
235  };
236 
237  // Tracks within this many sigma(z) are added to weight; increasing cut
238  // trades CPU efficiency for improved smoothness in tails
239  Gaudi::Property<double> m_z0MaxSignificance{
240  this,
241  "MaxZ0Significance",
242  12.0,
243  "Maximum longitudinal impact parameter significance to include track in "
244  "weight"
245  };
246 
247  // Assumed shape of density function near peak; Gaussian (true) or parabolic
248  // (false)
249  Gaudi::Property<bool> m_gaussStep{
250  this,
251  "GaussianStep",
252  true,
253  "Peak search: True means assume Gaussian behavior, False means "
254  "Newton/parabolic"
255  };
256  };
257 }
258 #endif
Trk::y
@ y
Definition: ParamDefs.h:62
Trk::GaussianTrackDensity::TrackDensityEval::m_firstDerivative
double m_firstDerivative
Definition: GaussianTrackDensity.h:110
Trk::GaussianTrackDensity::globalMaximumWithWidth
virtual std::pair< double, double > globalMaximumWithWidth(const std::vector< const TrackParameters * > &perigeeList) const override final
Definition: GaussianTrackDensity.cxx:79
Trk::GaussianTrackDensity::TrackDensity::pred_perigee
Definition: GaussianTrackDensity.h:169
maxValue
#define maxValue(current, test)
Definition: CompoundLayerMaterialCreator.h:22
Trk::z
@ z
global position (cartesian)
Definition: ParamDefs.h:63
Trk::GaussianTrackDensity::TrackEntry::operator=
TrackEntry & operator=(TrackEntry &&)=default
Trk::GaussianTrackDensity::TrackDensityEval::TrackDensityEval
TrackDensityEval(double z_coordinate)
Definition: GaussianTrackDensity.h:98
Trk::GaussianTrackDensity::TrackDensity::pred_entry_by_max
Definition: GaussianTrackDensity.h:177
Trk::GaussianTrackDensity::globalMaximum
virtual double globalMaximum(const std::vector< const Track * > &vectorTrk) const override final
Find position of global maximum for density function.
Definition: GaussianTrackDensity.cxx:20
Trk::ParametersT
Dummy class used to allow special convertors to be called for surfaces owned by a detector element.
Definition: EMErrorDetail.h:25
Trk::GaussianTrackDensity::TrackDensity::trackDensity
virtual double trackDensity(double z) const override final
Evaluate the density function at the specified coordinate along the beamline.
Definition: GaussianTrackDensity.cxx:130
Trk::GaussianTrackDensity::TrackDensityEval::m_density
double m_density
Definition: GaussianTrackDensity.h:109
extractSporadic.c1
c1
Definition: extractSporadic.py:134
Trk::GaussianTrackDensity
Definition: GaussianTrackDensity.h:32
Trk::GaussianTrackDensity::TrackDensityEval::density
double density() const
Definition: GaussianTrackDensity.h:104
Trk::Perigee
ParametersT< 5, Charged, PerigeeSurface > Perigee
Definition: Tracking/TrkEvent/TrkParameters/TrkParameters/TrackParameters.h:29
Trk::z0
@ z0
Definition: ParamDefs.h:70
Trk::GaussianTrackDensity::TrackEntry::TrackEntry
TrackEntry(double z)
Definition: GaussianTrackDensity.h:84
Trk::GaussianTrackDensity::TrackDensity
Definition: GaussianTrackDensity.h:114
Trk::GaussianTrackDensity::TrackEntry::TrackEntry
TrackEntry(const TrackEntry &)=default
Trk::GaussianTrackDensity::TrackEntry::operator=
TrackEntry & operator=(const TrackEntry &)=default
Trk::GaussianTrackDensity::TrackEntry::TrackEntry
TrackEntry(double c0, double c1, double c2, double zMin, double zMax)
Definition: GaussianTrackDensity.h:82
Trk::GaussianTrackDensity::TrackDensity::globalMaximumWithWidth
std::pair< double, double > globalMaximumWithWidth() const
Return position of global maximum with Gaussian width for density function.
Definition: GaussianTrackDensity.cxx:160
hotSpotInTAG.c0
c0
Definition: hotSpotInTAG.py:192
ArenaPoolSTLAllocator.h
STL-style allocator wrapper for ArenaPoolAllocator.
Trk::GaussianTrackDensity::TrackDensity::trackMap
std::map< Perigee, GaussianTrackDensity::TrackEntry, pred_perigee, SG::ArenaPoolSTLAllocator< std::pair< const Perigee, GaussianTrackDensity::TrackEntry > >> trackMap
Definition: GaussianTrackDensity.h:187
Trk::GaussianTrackDensity::TrackEntry::c_1
double c_1
Definition: GaussianTrackDensity.h:88
Trk::GaussianTrackDensity::TrackDensity::pred_perigee::operator()
bool operator()(const Perigee &left, const Perigee &right) const
Definition: GaussianTrackDensity.h:170
Trk::GaussianTrackDensity::TrackEntry::lowerBound
double lowerBound
Definition: GaussianTrackDensity.h:90
AthAlgTool.h
Trk::GaussianTrackDensity::TrackDensity::lowerMap
std::map< GaussianTrackDensity::TrackEntry, Perigee, pred_entry_by_max, SG::ArenaPoolSTLAllocator< std::pair< const GaussianTrackDensity::TrackEntry, Perigee > >> lowerMap
Definition: GaussianTrackDensity.h:192
Trk::GaussianTrackDensity::TrackDensity::m_gaussStep
bool m_gaussStep
Definition: GaussianTrackDensity.h:165
Trk::GaussianTrackDensity::TrackDensityEval::firstDerivative
double firstDerivative() const
Definition: GaussianTrackDensity.h:105
Trk::GaussianTrackDensity::TrackDensity::m_trackMap
trackMap m_trackMap
Definition: GaussianTrackDensity.h:194
Trk::GaussianTrackDensity::TrackEntry::c_0
double c_0
Definition: GaussianTrackDensity.h:87
Trk::GaussianTrackDensity::TrackDensity::globalMaximum
double globalMaximum() const
Return position of global maximum for density function.
Definition: GaussianTrackDensity.cxx:210
Trk::GaussianTrackDensity::TrackDensity::TrackDensity
TrackDensity(bool gaussStep)
Definition: GaussianTrackDensity.h:116
Trk::GaussianTrackDensity::TrackDensityEval::addTrack
void addTrack(const TrackEntry &entry)
Definition: GaussianTrackDensity.cxx:257
GetAllXsec.entry
list entry
Definition: GetAllXsec.py:132
Trk::GaussianTrackDensity::TrackDensityEval::secondDerivative
double secondDerivative() const
Definition: GaussianTrackDensity.h:106
Trk
Ensure that the ATLAS eigen extensions are properly loaded.
Definition: FakeTrackBuilder.h:9
Trk::GaussianTrackDensity::addTracks
void addTracks(const std::vector< const TrackParameters * > &perigeeList, TrackDensity &density) const
Add a set of tracks to a density object.
Definition: GaussianTrackDensity.cxx:109
Trk::GaussianTrackDensity::TrackDensity::stepSize
double stepSize(double y, double dy, double ddy) const
Definition: GaussianTrackDensity.h:162
compileRPVLLRates.c2
c2
Definition: compileRPVLLRates.py:361
IVertexTrackDensityEstimator.h
Trk::GaussianTrackDensity::TrackEntry::TrackEntry
TrackEntry(TrackEntry &&)=default
makeTRTBarrelCans.dy
tuple dy
Definition: makeTRTBarrelCans.py:21
Trk::GaussianTrackDensity::TrackEntry::TrackEntry
TrackEntry()=default
Trk::GaussianTrackDensity::TrackDensity::m_lowerMap
lowerMap m_lowerMap
Definition: GaussianTrackDensity.h:195
Trk::GaussianTrackDensity::TrackDensityEval
Definition: GaussianTrackDensity.h:95
Trk::GaussianTrackDensity::m_d0MaxSignificance
Gaudi::Property< double > m_d0MaxSignificance
Definition: GaussianTrackDensity.h:230
SG::ArenaPoolSTLAllocator
STL-style allocator wrapper for ArenaPoolAllocator.
Definition: ArenaPoolSTLAllocator.h:106
Trk::GaussianTrackDensity::TrackEntry::upperBound
double upperBound
Definition: GaussianTrackDensity.h:91
Trk::GaussianTrackDensity::TrackEntry
Definition: GaussianTrackDensity.h:74
Trk::GaussianTrackDensity::TrackDensity::updateMaximum
void updateMaximum(double trialZ, double trialValue, double secondDerivative, double &maxZ, double &maxValue, double &maxSecondDerivative) const
Definition: GaussianTrackDensity.h:151
Track
Definition: TriggerChamberClusterOnTrackCreator.h:21
Trk::GaussianTrackDensity::TrackDensity::addTrack
void addTrack(const Perigee &itrk, const double d0SignificanceCut, const double z0SignificanceCut)
Add a track to the set being considered.
Definition: GaussianTrackDensity.cxx:222
Trk::GaussianTrackDensity::TrackDensity::~TrackDensity
virtual ~TrackDensity()=default
Trk::GaussianTrackDensity::TrackDensityEval::m_secondDerivative
double m_secondDerivative
Definition: GaussianTrackDensity.h:111
Trk::GaussianTrackDensity::TrackDensity::m_maxRange
double m_maxRange
Definition: GaussianTrackDensity.h:166
Trk::GaussianTrackDensity::globalMaximumImpl
double globalMaximumImpl(const std::vector< const TrackParameters * > &perigeeList, TrackDensity &density) const
Find position of global maximum for density function.
Definition: GaussianTrackDensity.cxx:90
Trk::GaussianTrackDensity::TrackEntry::~TrackEntry
~TrackEntry()=default
Trk::GaussianTrackDensity::m_gaussStep
Gaudi::Property< bool > m_gaussStep
Definition: GaussianTrackDensity.h:249
Trk::GaussianTrackDensity::TrackEntry::c_2
double c_2
Definition: GaussianTrackDensity.h:89
Trk::GaussianTrackDensity::m_z0MaxSignificance
Gaudi::Property< double > m_z0MaxSignificance
Definition: GaussianTrackDensity.h:239
Trk::GaussianTrackDensity::TrackDensity::pred_entry_by_max::operator()
bool operator()(const TrackEntry &left, const TrackEntry &right) const
Definition: GaussianTrackDensity.h:178
Trk::GaussianTrackDensity::TrackDensityEval::m_z
double m_z
Definition: GaussianTrackDensity.h:108
Trk::GaussianTrackDensity::globalMaximumWithWidthImpl
std::pair< double, double > globalMaximumWithWidthImpl(const std::vector< const TrackParameters * > &perigeeList, TrackDensity &density) const
Find position of global maximum with Gaussian width for density function.
Definition: GaussianTrackDensity.cxx:97