CMS 3D CMS Logo

KFUpdator.cc
Go to the documentation of this file.
9 
10 // test of joseph form
11 #ifdef KU_JF_TEST
12 
13 #include <atomic>
14 #include <iostream>
15 namespace {
16  struct Stat {
17  Stat() : tot(0), nopd(0), jnopd(0), fnopd(0), inopd(0), ijnopd(0), ifnopd(0), dtot(0) {
18  for (int i = 0; i < 5; ++i)
19  idmm[i] = 0;
20  }
21 
22  std::atomic<long long> tot;
23  std::atomic<long long> nopd;
24  std::atomic<long long> jnopd;
25  std::atomic<long long> fnopd;
26  std::atomic<long long> inopd;
27  std::atomic<long long> ijnopd;
28  std::atomic<long long> ifnopd;
29 
30  std::atomic<unsigned long long> dtot;
31  std::atomic<unsigned long long> idmm[5];
32  double mm = 0;
33  double dmm[5] = {0};
34  bool ok = true;
35  ~Stat() {
36  double n = 1.e-3 / double(dtot);
37  std::cout << "KF " << tot * 1.e-9 << "G " << mm << " " << dmm[0] << '/' << dmm[1] << '/' << dmm[2] << '/'
38  << dmm[3] << '/' << dmm[4] << '\n'
39  << dtot << ' ' << idmm[0] << '/' << idmm[1] << '/' << idmm[2] << '/' << idmm[3] << '/' << idmm[4]
40  << '\n'
41  << std::sqrt(idmm[0] * n) << '/' << std::sqrt(idmm[1] * n) << '/' << std::sqrt(idmm[2] * n) << '/'
42  << std::sqrt(idmm[3] * n) << '/' << std::sqrt(idmm[4] * n) << " " << nopd << " " << jnopd << " "
43  << fnopd << " " << inopd << " " << ijnopd << " " << ifnopd << std::endl;
44  }
45  };
46 
47  Stat stat;
48 
49  bool isNopd(AlgebraicSymMatrix55 const& m) {
50  return m(0, 0) < 0 || m(1, 1) < 0 || m(2, 2) < 0 || m(3, 3) < 0 || m(4, 4) < 0;
51  }
52 } // namespace
53 #endif
54 
55 namespace {
56 
57  template <unsigned int D>
58  TrajectoryStateOnSurface lupdate(const TrajectoryStateOnSurface& tsos, const TrackingRecHit& aRecHit) {
59  typedef typename AlgebraicROOTObject<D, 5>::Matrix MatD5;
60  typedef typename AlgebraicROOTObject<5, D>::Matrix Mat5D;
61  typedef typename AlgebraicROOTObject<D, D>::SymMatrix SMatDD;
62  typedef typename AlgebraicROOTObject<D>::Vector VecD;
63  using ROOT::Math::SMatrixNoInit;
64  double pzSign = tsos.localParameters().pzSign();
65 
66  auto&& x = tsos.localParameters().vector();
67  auto&& C = tsos.localError().matrix();
68 
69  // projection matrix (assume element of "H" to be just 0 or 1)
71 
72  // Measurement matrix
73  VecD r, rMeas;
74  SMatDD V(SMatrixNoInit{}), VMeas(SMatrixNoInit{});
75 
76  KfComponentsHolder holder;
77  holder.template setup<D>(&r, &V, &pf, &rMeas, &VMeas, x, C);
78  aRecHit.getKfComponents(holder);
79 
80  r -= rMeas;
81 
82  // and covariance matrix of residuals
83  SMatDD R = V + VMeas;
84  bool ok = invertPosDefMatrix(R);
85 
86  // Compute Kalman gain matrix
88  Mat5D K = C * pf.project(R);
89  pf.projectAndSubtractFrom(M, K);
90 
91  // Compute local filtered state vector
92  AlgebraicVector5 fsv = x + K * r;
93 
94  // Compute covariance matrix of local filtered state vector
95 #define KU_JosephForm
96 #ifdef KU_JosephForm
97  AlgebraicSymMatrix55 fse = ROOT::Math::Similarity(M, C) + ROOT::Math::Similarity(K, V);
98 #else
100  ROOT::Math::AssignSym::Evaluate(fse, M * C);
101 #endif
102 
103 // test of joseph form
104 #ifdef KU_JF_TEST
105 
107  ROOT::Math::AssignSym::Evaluate(fse2, M * C);
108 
109  // std::cout << "Joseph Form \n" << fse << std::endl;
110  // std::cout << "Fast Form \n" << fse2 << std::endl;
111 
112  stat.tot++;
113  auto n1 = isNopd(fse);
114  auto n2 = isNopd(fse2);
115  if (n1 && n2)
116  stat.nopd++;
117  if (n1)
118  stat.jnopd++;
119  if (n2)
120  stat.fnopd++;
121 
122  AlgebraicSymMatrix55 dd = fse2 - fse;
123  auto dda = dd.Array();
124  auto fsa = fse.Array();
125  double ddd[15];
126  for (int i = 0; i < 15; ++i)
127  ddd[i] = std::abs(dda[i]) / std::abs(fsa[i]);
128  auto mm = *std::max_element(ddd, ddd + 15);
129  stat.mm = std::max(stat.mm, mm);
130  if (stat.ok && !(n1 || n2)) {
131  stat.dtot++;
132  for (int i = 0; i < 5; ++i) {
133  auto dmm = std::sqrt(fse2(i, i)) - std::sqrt(fse(i, i));
134  stat.dmm[i] = std::max(stat.dmm[i], std::abs(dmm));
135  stat.idmm[i] += (unsigned long long)(1000. * dmm * dmm);
137  stat.ok = false;
138  }
139  }
140 
141  AlgebraicSymMatrix55 ifse = fse;
142  invertPosDefMatrix(ifse);
143  AlgebraicSymMatrix55 ifse2 = fse2;
144  invertPosDefMatrix(ifse2);
145  n1 = isNopd(ifse);
146  n2 = isNopd(ifse2);
147  if (n1 && n2)
148  stat.inopd++;
149  if (n1)
150  stat.ijnopd++;
151  if (n2)
152  stat.ifnopd++;
153 #endif
154 
155  /*
156  // expanded similariy
157  AlgebraicSymMatrix55 fse;
158  ROOT::Math::AssignSym::Evaluate(fse, (M* C) * ( ROOT::Math::Transpose(M)));
159  AlgebraicSymMatrix55 tmp;
160  ROOT::Math::AssignSym::Evaluate(tmp, (K*V) * (ROOT::Math::Transpose(K)));
161  fse += tmp;
162  */
163 
164  if (ok) {
167  tsos.surface(),
168  &(tsos.globalParameters().magneticField()),
169  tsos.surfaceSide());
170  } else {
171  edm::LogError("KFUpdator") << " could not invert martix:\n" << (V + VMeas);
172  return TrajectoryStateOnSurface();
173  }
174  }
175 } // namespace
176 
178  switch (aRecHit.dimension()) {
179  case 1:
180  return lupdate<1>(tsos, aRecHit);
181  case 2:
182  return lupdate<2>(tsos, aRecHit);
183  case 3:
184  return lupdate<3>(tsos, aRecHit);
185  case 4:
186  return lupdate<4>(tsos, aRecHit);
187  case 5:
188  return lupdate<5>(tsos, aRecHit);
189  }
190  throw cms::Exception("Rec hit of invalid dimension (not 1,2,3,4,5)")
191  << "The value was " << aRecHit.dimension() << ", type is " << typeid(aRecHit).name() << "\n";
192 }
KfComponentsHolder.h
cms::cuda::V
cudaStream_t T uint32_t const T *__restrict__ const uint32_t *__restrict__ uint32_t int cudaStream_t V
Definition: HistoContainer.h:99
KFUpdator::update
TrajectoryStateOnSurface update(const TrajectoryStateOnSurface &, const TrackingRecHit &) const override
Definition: KFUpdator.cc:177
mps_fire.i
i
Definition: mps_fire.py:428
MessageLogger.h
dqmiodumpmetadata.n
n
Definition: dqmiodumpmetadata.py:28
ProjectMatrix.h
MeasurementExtractor.h
gather_cfg.cout
cout
Definition: gather_cfg.py:144
AlgebraicROOTObject::SymMatrix
ROOT::Math::SMatrix< double, D1, D1, ROOT::Math::MatRepSym< double, D1 > > SymMatrix
Definition: AlgebraicROOTObjects.h:68
TransientTrackingRecHit.h
LocalTrajectoryError::matrix
const AlgebraicSymMatrix55 & matrix() const
Definition: LocalTrajectoryError.h:60
convertSQLiteXML.ok
bool ok
Definition: convertSQLiteXML.py:98
LocalTrajectoryParameters
Definition: LocalTrajectoryParameters.h:25
AlgebraicMatrixID
ROOT::Math::SMatrixIdentity AlgebraicMatrixID
Definition: AlgebraicROOTObjects.h:72
Plane.h
TrajectoryStateOnSurface
Definition: TrajectoryStateOnSurface.h:16
AlgebraicMatrix55
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepStd< double, 5, 5 > > AlgebraicMatrix55
Definition: AlgebraicROOTObjects.h:55
hgcalPlots.stat
stat
Definition: hgcalPlots.py:1119
visualization-live-secondInstance_cfg.m
m
Definition: visualization-live-secondInstance_cfg.py:72
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
LocalTrajectoryParameters::vector
AlgebraicVector5 vector() const
Definition: LocalTrajectoryParameters.h:120
createTree.dd
string dd
Definition: createTree.py:154
TrajectoryStateOnSurface::localParameters
const LocalTrajectoryParameters & localParameters() const
Definition: TrajectoryStateOnSurface.h:73
LocalTrajectoryError
Definition: LocalTrajectoryError.h:20
LocalTrajectoryParameters::pzSign
float pzSign() const
Sign of the z-component of the momentum in the local frame.
Definition: LocalTrajectoryParameters.h:146
KFUpdator.h
AlgebraicROOTObject::Matrix
ROOT::Math::SMatrix< double, D1, D2, ROOT::Math::MatRepStd< double, D1, D2 > > Matrix
Definition: AlgebraicROOTObjects.h:69
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
invertPosDefMatrix
bool invertPosDefMatrix(ROOT::Math::SMatrix< T, N, N, ROOT::Math::MatRepSym< T, N > > &m)
Definition: invertPosDefMatrix.h:10
edm::LogError
Log< level::Error, false > LogError
Definition: MessageLogger.h:123
AlgebraicVector5
ROOT::Math::SVector< double, 5 > AlgebraicVector5
Definition: AlgebraicROOTObjects.h:14
KfComponentsHolder
Definition: KfComponentsHolder.h:13
TrackingRecHit
Definition: TrackingRecHit.h:21
alignCSCRings.r
r
Definition: alignCSCRings.py:93
gen::C
C
Definition: PomwigHadronizer.cc:78
TrackingRecHit::getKfComponents
virtual void getKfComponents(KfComponentsHolder &holder) const
Definition: TrackingRecHit.cc:23
TrackingRecHit::dimension
virtual int dimension() const =0
packedPFCandidateRefMixer_cfi.pf
pf
Definition: packedPFCandidateRefMixer_cfi.py:4
Exception
Definition: hltDiff.cc:246
genVertex_cff.x
x
Definition: genVertex_cff.py:12
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
GlobalTrajectoryParameters::magneticField
const MagneticField & magneticField() const
Definition: GlobalTrajectoryParameters.h:106
TrajectoryStateOnSurface::surface
const SurfaceType & surface() const
Definition: TrajectoryStateOnSurface.h:78
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
TrajectoryStateOnSurface::globalParameters
const GlobalTrajectoryParameters & globalParameters() const
Definition: TrajectoryStateOnSurface.h:64
TrajectoryStateOnSurface::localError
const LocalTrajectoryError & localError() const
Definition: TrajectoryStateOnSurface.h:77
AlgebraicSymMatrix55
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
Definition: AlgebraicROOTObjects.h:23
TrajectoryStateOnSurface::surfaceSide
SurfaceSide surfaceSide() const
Position relative to material, defined relative to momentum vector.
Definition: TrajectoryStateOnSurface.h:89
dttmaxenums::R
Definition: DTTMax.h:29
ProjectMatrix
Definition: ProjectMatrix.h:8
AlgebraicROOTObject::Vector
ROOT::Math::SVector< double, D1 > Vector
Definition: AlgebraicROOTObjects.h:67
invertPosDefMatrix.h