CMS 3D CMS Logo

MultifitComputations.h
Go to the documentation of this file.
1 #ifndef DataFormats_CaloRecHit_interface_MultifitComputations_h
2 #define DataFormats_CaloRecHit_interface_MultifitComputations_h
3 
4 #include <cmath>
5 #include <limits>
6 #include <type_traits>
7 
8 #include <Eigen/Dense>
9 
10 namespace calo {
11  namespace multifit {
12 
13  template <int NROWS, int NCOLS>
14  using ColMajorMatrix = Eigen::Matrix<float, NROWS, NCOLS, Eigen::ColMajor>;
15 
16  template <int NROWS, int NCOLS>
17  using RowMajorMatrix = Eigen::Matrix<float, NROWS, NCOLS, Eigen::RowMajor>;
18 
19  template <int SIZE, typename T = float>
20  using ColumnVector = Eigen::Matrix<T, SIZE, 1>;
21 
22  template <int SIZE, typename T = float>
23  using RowVector = Eigen::Matrix<T, 1, SIZE>;
24 
25  // FIXME: provide specialization for Row Major layout
26  template <typename T, int Stride, int Order = Eigen::ColMajor>
27  struct MapSymM {
28  using type = T;
30 
31  static constexpr int total = Stride * (Stride + 1) / 2;
32  static constexpr int stride = Stride;
33  T* data;
34 
35  EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC MapSymM(T* data) : data{data} {}
36 
37  EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC T const& operator()(int const row, int const col) const {
38  auto const tmp = (Stride - col) * (Stride - col + 1) / 2;
39  auto const index = total - tmp + row - col;
40  return data[index];
41  }
42 
43  template <typename U = T>
44  EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC typename std::enable_if<std::is_same<base_type, U>::value, base_type>::type&
45  operator()(int const row, int const col) {
46  auto const tmp = (Stride - col) * (Stride - col + 1) / 2;
47  auto const index = total - tmp + row - col;
48  return data[index];
49  }
50  };
51 
52  // FIXME: either use/modify/improve eigen or make this more generic
53  // this is a map for a pulse matrix to building a 2d matrix for each channel
54  // and hide indexing
55  template <typename T>
56  struct MapMForPM {
57  using type = T;
59 
61  EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC MapMForPM(type* data) : data{data} {}
62 
63  EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC base_type operator()(int const row, int const col) const {
64  auto const index = 2 - col + row;
65  return index >= 0 ? data[index] : 0;
66  }
67  };
68 
69  // simple/trivial cholesky decomposition impl
70  template <typename MatrixType1, typename MatrixType2>
71  EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC void compute_decomposition_unrolled(MatrixType1& L, MatrixType2 const& M) {
72  auto const sqrtm_0_0 = std::sqrt(M(0, 0));
73  L(0, 0) = sqrtm_0_0;
74  using T = typename MatrixType1::base_type;
75 
76 #pragma unroll
77  for (int i = 1; i < MatrixType1::stride; i++) {
78  T sumsq{0};
79  for (int j = 0; j < i; j++) {
80  T sumsq2{0};
81  auto const m_i_j = M(i, j);
82  for (int k = 0; k < j; ++k)
83  sumsq2 += L(i, k) * L(j, k);
84 
85  auto const value_i_j = (m_i_j - sumsq2) / L(j, j);
86  L(i, j) = value_i_j;
87 
88  sumsq += value_i_j * value_i_j;
89  }
90 
91  auto const l_i_i = std::sqrt(M(i, i) - sumsq);
92  L(i, i) = l_i_i;
93  }
94  }
95 
96  template <typename MatrixType1, typename MatrixType2>
97  EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC void compute_decomposition(MatrixType1& L,
98  MatrixType2 const& M,
99  int const N) {
100  auto const sqrtm_0_0 = std::sqrt(M(0, 0));
101  L(0, 0) = sqrtm_0_0;
102  using T = typename MatrixType1::base_type;
103 
104  for (int i = 1; i < N; i++) {
105  T sumsq{0};
106  for (int j = 0; j < i; j++) {
107  T sumsq2{0};
108  auto const m_i_j = M(i, j);
109  for (int k = 0; k < j; ++k)
110  sumsq2 += L(i, k) * L(j, k);
111 
112  auto const value_i_j = (m_i_j - sumsq2) / L(j, j);
113  L(i, j) = value_i_j;
114 
115  sumsq += value_i_j * value_i_j;
116  }
117 
118  auto const l_i_i = std::sqrt(M(i, i) - sumsq);
119  L(i, i) = l_i_i;
120  }
121  }
122 
123  template <typename MatrixType1, typename MatrixType2, typename VectorType>
124  EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC void compute_decomposition_forwardsubst_with_offsets(
125  MatrixType1& L,
126  MatrixType2 const& M,
127  float b[MatrixType1::stride],
128  VectorType const& Atb,
129  int const N,
130  ColumnVector<MatrixType1::stride, int> const& pulseOffsets) {
131  auto const real_0 = pulseOffsets(0);
132  auto const sqrtm_0_0 = std::sqrt(M(real_0, real_0));
133  L(0, 0) = sqrtm_0_0;
134  using T = typename MatrixType1::base_type;
135  b[0] = Atb(real_0) / sqrtm_0_0;
136 
137  for (int i = 1; i < N; i++) {
138  auto const i_real = pulseOffsets(i);
139  T sumsq{0};
140  T total = 0;
141  auto const atb = Atb(i_real);
142  for (int j = 0; j < i; j++) {
143  auto const j_real = pulseOffsets(j);
144  T sumsq2{0};
145  auto const m_i_j = M(std::max(i_real, j_real), std::min(i_real, j_real));
146  for (int k = 0; k < j; ++k)
147  sumsq2 += L(i, k) * L(j, k);
148 
149  auto const value_i_j = (m_i_j - sumsq2) / L(j, j);
150  L(i, j) = value_i_j;
151 
152  sumsq += value_i_j * value_i_j;
153  total += value_i_j * b[j];
154  }
155 
156  auto const l_i_i = std::sqrt(M(i_real, i_real) - sumsq);
157  L(i, i) = l_i_i;
158  b[i] = (atb - total) / l_i_i;
159  }
160  }
161 
162  template <typename MatrixType1, typename MatrixType2, typename VectorType>
163  EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC void update_decomposition_forwardsubst_with_offsets(
164  MatrixType1& L,
165  MatrixType2 const& M,
166  float b[MatrixType1::stride],
167  VectorType const& Atb,
168  int const N,
169  ColumnVector<MatrixType1::stride, int> const& pulseOffsets) {
170  using T = typename MatrixType1::base_type;
171  auto const i = N - 1;
172  auto const i_real = pulseOffsets(i);
173  T sumsq{0};
174  T total = 0;
175  for (int j = 0; j < i; j++) {
176  auto const j_real = pulseOffsets(j);
177  T sumsq2{0};
178  auto const m_i_j = M(std::max(i_real, j_real), std::min(i_real, j_real));
179  for (int k = 0; k < j; ++k)
180  sumsq2 += L(i, k) * L(j, k);
181 
182  auto const value_i_j = (m_i_j - sumsq2) / L(j, j);
183  L(i, j) = value_i_j;
184  sumsq += value_i_j * value_i_j;
185 
186  total += value_i_j * b[j];
187  }
188 
189  auto const l_i_i = std::sqrt(M(i_real, i_real) - sumsq);
190  L(i, i) = l_i_i;
191  b[i] = (Atb(i_real) - total) / l_i_i;
192  }
193 
194  template <typename MatrixType1, typename MatrixType2, typename MatrixType3>
195  EIGEN_DEVICE_FUNC void solve_forward_subst_matrix(MatrixType1& A,
196  MatrixType2 const& pulseMatrixView,
197  MatrixType3 const& matrixL) {
198  // FIXME: this assumes pulses are on columns and samples on rows
199  constexpr auto NPULSES = MatrixType2::ColsAtCompileTime;
200  constexpr auto NSAMPLES = MatrixType2::RowsAtCompileTime;
201 
202 #pragma unroll
203  for (int icol = 0; icol < NPULSES; icol++) {
204  float reg_b[NSAMPLES];
205  float reg_L[NSAMPLES];
206 
207 // preload a column and load column 0 of cholesky
208 #pragma unroll
209  for (int i = 0; i < NSAMPLES; i++) {
210 #ifdef __CUDA_ARCH__
211  // load through the read-only cache
212  reg_b[i] = __ldg(&pulseMatrixView.coeffRef(i, icol));
213 #else
214  reg_b[i] = pulseMatrixView.coeffRef(i, icol);
215 #endif // __CUDA_ARCH__
216  reg_L[i] = matrixL(i, 0);
217  }
218 
219  // compute x0 and store it
220  auto x_prev = reg_b[0] / reg_L[0];
221  A(0, icol) = x_prev;
222 
223 // iterate
224 #pragma unroll
225  for (int iL = 1; iL < NSAMPLES; iL++) {
226 // update accum
227 #pragma unroll
228  for (int counter = iL; counter < NSAMPLES; counter++)
229  reg_b[counter] -= x_prev * reg_L[counter];
230 
231 // load the next column of cholesky
232 #pragma unroll
233  for (int counter = iL; counter < NSAMPLES; counter++)
234  reg_L[counter] = matrixL(counter, iL);
235 
236  // compute the next x for M(iL, icol)
237  x_prev = reg_b[iL] / reg_L[iL];
238 
239  // store the result value
240  A(iL, icol) = x_prev;
241  }
242  }
243  }
244 
245  template <typename MatrixType1, typename MatrixType2>
246  EIGEN_DEVICE_FUNC void solve_forward_subst_vector(float reg_b[MatrixType1::RowsAtCompileTime],
247  MatrixType1 inputAmplitudesView,
248  MatrixType2 matrixL) {
249  constexpr auto NSAMPLES = MatrixType1::RowsAtCompileTime;
250 
251  float reg_b_tmp[NSAMPLES];
252  float reg_L[NSAMPLES];
253 
254 // preload a column and load column 0 of cholesky
255 #pragma unroll
256  for (int i = 0; i < NSAMPLES; i++) {
257  reg_b_tmp[i] = inputAmplitudesView(i);
258  reg_L[i] = matrixL(i, 0);
259  }
260 
261  // compute x0 and store it
262  auto x_prev = reg_b_tmp[0] / reg_L[0];
263  reg_b[0] = x_prev;
264 
265 // iterate
266 #pragma unroll
267  for (int iL = 1; iL < NSAMPLES; iL++) {
268 // update accum
269 #pragma unroll
270  for (int counter = iL; counter < NSAMPLES; counter++)
271  reg_b_tmp[counter] -= x_prev * reg_L[counter];
272 
273 // load the next column of cholesky
274 #pragma unroll
275  for (int counter = iL; counter < NSAMPLES; counter++)
276  reg_L[counter] = matrixL(counter, iL);
277 
278  // compute the next x for M(iL, icol)
279  x_prev = reg_b_tmp[iL] / reg_L[iL];
280 
281  // store the result value
282  reg_b[iL] = x_prev;
283  }
284  }
285 
286  template <typename MatrixType1, typename MatrixType2, typename MatrixType3, typename MatrixType4>
287  EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC void calculateChiSq(MatrixType1 const& matrixL,
288  MatrixType2 const& pulseMatrixView,
289  MatrixType3 const& resultAmplitudesVector,
290  MatrixType4 const& inputAmplitudesView,
291  float& chi2) {
292  // FIXME: this assumes pulses are on columns and samples on rows
293  constexpr auto NPULSES = MatrixType2::ColsAtCompileTime;
294  constexpr auto NSAMPLES = MatrixType2::RowsAtCompileTime;
295 
296  // replace pulseMatrixView * resultAmplitudesVector - inputAmplitudesView
297  // NOTE:
298  float accum[NSAMPLES];
299  {
300  float results[NPULSES];
301 
302  // preload results and permute according to the pulse offsets /////////////// ??? this is not done in ECAL
303 #pragma unroll
304  for (int counter = 0; counter < NPULSES; counter++) {
305  results[counter] = resultAmplitudesVector[counter];
306  }
307 
308  // load accum
309 #pragma unroll
310  for (int counter = 0; counter < NSAMPLES; counter++)
311  accum[counter] = -inputAmplitudesView(counter);
312 
313  // iterate
314  for (int icol = 0; icol < NPULSES; icol++) {
315  float pm_col[NSAMPLES];
316 
317  // preload a column of pulse matrix
318 #pragma unroll
319  for (int counter = 0; counter < NSAMPLES; counter++)
320 #ifdef __CUDA_ARCH__
321  pm_col[counter] = __ldg(&pulseMatrixView.coeffRef(counter, icol));
322 #else
323  pm_col[counter] = pulseMatrixView.coeffRef(counter, icol);
324 #endif
325 
326  // accum
327 #pragma unroll
328  for (int counter = 0; counter < NSAMPLES; counter++)
329  accum[counter] += results[icol] * pm_col[counter];
330  }
331  }
332 
333  // compute chi2 and check that there is no rotation
334  // chi2 = matrixDecomposition
335  // .matrixL()
336  // . solve(mapAccum)
337  // .solve(pulseMatrixView * resultAmplitudesVector - inputAmplitudesView)
338  // .squaredNorm();
339 
340  {
341  float reg_L[NSAMPLES];
342  float accumSum = 0;
343 
344  // preload a column and load column 0 of cholesky
345 #pragma unroll
346  for (int i = 0; i < NSAMPLES; i++) {
347  reg_L[i] = matrixL(i, 0);
348  }
349 
350  // compute x0 and store it
351  auto x_prev = accum[0] / reg_L[0];
352  accumSum += x_prev * x_prev;
353 
354  // iterate
355 #pragma unroll
356  for (int iL = 1; iL < NSAMPLES; iL++) {
357  // update accum
358 #pragma unroll
359  for (int counter = iL; counter < NSAMPLES; counter++)
360  accum[counter] -= x_prev * reg_L[counter];
361 
362  // load the next column of cholesky
363 #pragma unroll
364  for (int counter = iL; counter < NSAMPLES; counter++)
365  reg_L[counter] = matrixL(counter, iL);
366 
367  // compute the next x for M(iL, icol)
368  x_prev = accum[iL] / reg_L[iL];
369 
370  // store the result value
371  accumSum += x_prev * x_prev;
372  }
373 
374  chi2 = accumSum;
375  }
376  }
377 
378  // TODO: add active bxs
379  template <typename MatrixType, typename VectorType>
380  EIGEN_DEVICE_FUNC void fnnls(MatrixType const& AtA,
381  VectorType const& Atb,
382  VectorType& solution,
383  int& npassive,
386  double eps, // convergence condition
387  const int maxIterations, // maximum number of iterations
388  const int relaxationPeriod, // every "relaxationPeriod" iterations
389  const int relaxationFactor) { // multiply "eps" by "relaxationFactor"
390  // constants
391  constexpr auto NPULSES = VectorType::RowsAtCompileTime;
392 
393  // to keep track of where to terminate if converged
394  Eigen::Index w_max_idx_prev = 0;
395  float w_max_prev = 0;
396  bool recompute = false;
397 
398  // used throughout
399  VectorType s;
400  float reg_b[NPULSES];
401  //float matrixLStorage[MapSymM<float, NPULSES>::total];
402  //MapSymM<float, NPULSES> matrixL{matrixLStorage};
403 
404  int iter = 0;
405  while (true) {
406  if (iter > 0 || npassive == 0) {
407  auto const nactive = NPULSES - npassive;
408  // exit if there are no more pulses to constrain
409  if (nactive == 0)
410  break;
411 
412  // compute the gradient
413  //w.tail(nactive) = Atb.tail(nactive) - (AtA * solution).tail(nactive);
414  Eigen::Index w_max_idx;
415  float w_max = -std::numeric_limits<float>::max();
416  for (int icol = npassive; icol < NPULSES; icol++) {
417  auto const icol_real = pulseOffsets(icol);
418  auto const atb = Atb(icol_real);
419  float sum = 0;
420 #pragma unroll
421  for (int counter = 0; counter < NPULSES; counter++)
422  sum += counter > icol_real ? AtA(counter, icol_real) * solution(counter)
423  : AtA(icol_real, counter) * solution(counter);
424 
425  auto const w = atb - sum;
426  if (w > w_max) {
427  w_max = w;
428  w_max_idx = icol - npassive;
429  }
430  }
431 
432  // check for convergence
433  if (w_max < eps || (w_max_idx == w_max_idx_prev && w_max == w_max_prev))
434  break;
435 
436  if (iter >= maxIterations)
437  break;
438 
439  w_max_prev = w_max;
440  w_max_idx_prev = w_max_idx;
441 
442  // move index to the right part of the vector
443  w_max_idx += npassive;
444 
445  Eigen::numext::swap(pulseOffsets.coeffRef(npassive), pulseOffsets.coeffRef(w_max_idx));
446  ++npassive;
447  }
448 
449  // inner loop
450  while (true) {
451  if (npassive == 0)
452  break;
453 
454  //s.head(npassive)
455  //auto const& matrixL =
456  // AtA.topLeftCorner(npassive, npassive)
457  // .llt().matrixL();
458  //.solve(Atb.head(npassive));
459  if (recompute || iter == 0)
460  compute_decomposition_forwardsubst_with_offsets(matrixL, AtA, reg_b, Atb, npassive, pulseOffsets);
461  else
462  update_decomposition_forwardsubst_with_offsets(matrixL, AtA, reg_b, Atb, npassive, pulseOffsets);
463 
464  // run backward substituion
465  s(npassive - 1) = reg_b[npassive - 1] / matrixL(npassive - 1, npassive - 1);
466  for (int i = npassive - 2; i >= 0; --i) {
467  float total = 0;
468  for (int j = i + 1; j < npassive; j++)
469  total += matrixL(j, i) * s(j);
470 
471  s(i) = (reg_b[i] - total) / matrixL(i, i);
472  }
473 
474  // done if solution values are all positive
475  bool hasNegative = false;
476  bool hasNans = false;
477  for (int counter = 0; counter < npassive; counter++) {
478  auto const s_ii = s(counter);
479  hasNegative |= s_ii <= 0;
480  hasNans |= std::isnan(s_ii);
481  }
482 
483  // FIXME: temporary solution. my cholesky impl is unstable yielding nans
484  // this check removes nans - do not accept solution unless all values
485  // are stable
486  if (hasNans)
487  break;
488  if (!hasNegative) {
489  for (int i = 0; i < npassive; i++) {
490  auto const i_real = pulseOffsets(i);
491  solution(i_real) = s(i);
492  }
493  //solution.head(npassive) = s.head(npassive);
494  recompute = false;
495  break;
496  }
497 
498  // there were negative values -> have to recompute the whole decomp
499  recompute = true;
500 
502  Eigen::Index alpha_idx = 0, alpha_idx_real = 0;
503  for (int i = 0; i < npassive; i++) {
504  if (s[i] <= 0.) {
505  auto const i_real = pulseOffsets(i);
506  auto const ratio = solution[i_real] / (solution[i_real] - s[i]);
507  if (ratio < alpha) {
508  alpha = ratio;
509  alpha_idx = i;
510  alpha_idx_real = i_real;
511  }
512  }
513  }
514 
515  // upadte solution
516  for (int i = 0; i < npassive; i++) {
517  auto const i_real = pulseOffsets(i);
518  solution(i_real) += alpha * (s(i) - solution(i_real));
519  }
520  //solution.head(npassive) += alpha *
521  // (s.head(npassive) - solution.head(npassive));
522  solution[alpha_idx_real] = 0;
523  --npassive;
524 
525  Eigen::numext::swap(pulseOffsets.coeffRef(npassive), pulseOffsets.coeffRef(alpha_idx));
526  }
527 
528  // as in cpu
529  ++iter;
530  if (iter % relaxationPeriod == 0)
531  eps *= relaxationFactor;
532  }
533  }
534 
535  } // namespace multifit
536 } // namespace calo
537 
538 #endif // DataFormats_CaloRecHit_interface_MultifitComputations_h
calo::multifit::compute_decomposition_forwardsubst_with_offsets
EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC void compute_decomposition_forwardsubst_with_offsets(MatrixType1 &L, MatrixType2 const &M, float b[MatrixType1::stride], VectorType const &Atb, int const N, ColumnVector< MatrixType1::stride, int > const &pulseOffsets)
Definition: MultifitComputations.h:124
calo::multifit::MapSymM::MapSymM
EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC MapSymM(T *data)
Definition: MultifitComputations.h:35
calo::multifit::solve_forward_subst_matrix
EIGEN_DEVICE_FUNC void solve_forward_subst_matrix(MatrixType1 &A, MatrixType2 const &pulseMatrixView, MatrixType3 const &matrixL)
Definition: MultifitComputations.h:195
calo::multifit::MapMForPM::base_type
typename std::remove_cv< type >::type base_type
Definition: MultifitComputations.h:58
counter
Definition: counter.py:1
dttmaxenums::L
Definition: DTTMax.h:29
mps_fire.i
i
Definition: mps_fire.py:428
calo::multifit::solve_forward_subst_vector
EIGEN_DEVICE_FUNC void solve_forward_subst_vector(float reg_b[MatrixType1::RowsAtCompileTime], MatrixType1 inputAmplitudesView, MatrixType2 matrixL)
Definition: MultifitComputations.h:246
calo::multifit::MapSymM::stride
static constexpr int stride
Definition: MultifitComputations.h:32
min
T min(T a, T b)
Definition: MathUtil.h:58
zMuMuMuonUserData.alpha
alpha
zGenParticlesMatch = cms.InputTag(""),
Definition: zMuMuMuonUserData.py:9
calo::multifit::RowMajorMatrix
Eigen::Matrix< float, NROWS, NCOLS, Eigen::RowMajor > RowMajorMatrix
Definition: MultifitComputations.h:17
VectorType
Vec4< T > VectorType
Definition: extBasic3DVector.h:165
cuy.col
col
Definition: cuy.py:1010
edm::swap
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:117
calo::multifit::update_decomposition_forwardsubst_with_offsets
EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC void update_decomposition_forwardsubst_with_offsets(MatrixType1 &L, MatrixType2 const &M, float b[MatrixType1::stride], VectorType const &Atb, int const N, ColumnVector< MatrixType1::stride, int > const &pulseOffsets)
Definition: MultifitComputations.h:163
calo::multifit::MapSymM::operator()
EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC std::enable_if< std::is_same< base_type, U >::value, base_type >::type & operator()(int const row, int const col)
Definition: MultifitComputations.h:45
bookConverter.results
results
Definition: bookConverter.py:144
hltPixelTracks_cff.chi2
chi2
Definition: hltPixelTracks_cff.py:25
calo::multifit::ColMajorMatrix
Eigen::Matrix< float, NROWS, NCOLS, Eigen::ColMajor > ColMajorMatrix
Definition: MultifitComputations.h:14
createJobs.tmp
tmp
align.sh
Definition: createJobs.py:716
calo::multifit::MapSymM::base_type
typename std::remove_const< type >::type base_type
Definition: MultifitComputations.h:29
calo::multifit::MapMForPM
Definition: MultifitComputations.h:56
alignCSCRings.s
s
Definition: alignCSCRings.py:92
calo::multifit::MapSymM::type
T type
Definition: MultifitComputations.h:28
w
const double w
Definition: UKUtility.cc:23
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
N
#define N
Definition: blowfish.cc:9
calo::multifit::MapSymM
Definition: MultifitComputations.h:27
CommonMethods.isnan
def isnan(num)
Definition: CommonMethods.py:98
dqmdumpme.k
k
Definition: dqmdumpme.py:60
b
double b
Definition: hdecay.h:118
calo::multifit::MapMForPM::MapMForPM
EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC MapMForPM(type *data)
Definition: MultifitComputations.h:61
particleFlowDisplacedVertex_cfi.ratio
ratio
Definition: particleFlowDisplacedVertex_cfi.py:93
calo::multifit::MapMForPM::type
T type
Definition: MultifitComputations.h:57
calo::multifit::compute_decomposition
EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC void compute_decomposition(MatrixType1 &L, MatrixType2 const &M, int const N)
Definition: MultifitComputations.h:97
HLT_FULL_cff.maxIterations
maxIterations
Definition: HLT_FULL_cff.py:13304
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
type
type
Definition: SiPixelVCal_PayloadInspector.cc:37
gainCalibHelper::gainCalibPI::type
type
Definition: SiPixelGainCalibHelper.h:39
calo::multifit::compute_decomposition_unrolled
EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC void compute_decomposition_unrolled(MatrixType1 &L, MatrixType2 const &M)
Definition: MultifitComputations.h:71
A
calo::multifit::MapSymM::total
static constexpr int total
Definition: MultifitComputations.h:31
counter
static std::atomic< unsigned int > counter
Definition: SharedResourceNames.cc:17
calo::multifit::MapMForPM::data
type * data
Definition: MultifitComputations.h:60
calo::multifit::calculateChiSq
EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC void calculateChiSq(MatrixType1 const &matrixL, MatrixType2 const &pulseMatrixView, MatrixType3 const &resultAmplitudesVector, MatrixType4 const &inputAmplitudesView, float &chi2)
Definition: MultifitComputations.h:287
calo::multifit::ColumnVector
Eigen::Matrix< T, SIZE, 1 > ColumnVector
Definition: MultifitComputations.h:20
calo::multifit::MapSymM::operator()
EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC T const & operator()(int const row, int const col) const
Definition: MultifitComputations.h:37
MaterialEffects_cfi.A
A
Definition: MaterialEffects_cfi.py:11
calo::multifit::MapMForPM::operator()
EIGEN_ALWAYS_INLINE EIGEN_DEVICE_FUNC base_type operator()(int const row, int const col) const
Definition: MultifitComputations.h:63
T
long double T
Definition: Basic3DVectorLD.h:48
relativeConstraints.value
value
Definition: relativeConstraints.py:53
calo::multifit::RowVector
Eigen::Matrix< T, 1, SIZE > RowVector
Definition: MultifitComputations.h:23
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46
dqmMemoryStats.total
total
Definition: dqmMemoryStats.py:152
calo::multifit::MapSymM::data
T * data
Definition: MultifitComputations.h:33
calo
Definition: Common.h:9
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
cms::cudacompat::__ldg
T __ldg(T const *x)
Definition: cudaCompat.h:77
calo::multifit::fnnls
EIGEN_DEVICE_FUNC void fnnls(MatrixType const &AtA, VectorType const &Atb, VectorType &solution, int &npassive, ColumnVector< VectorType::RowsAtCompileTime, int > &pulseOffsets, MapSymM< float, VectorType::RowsAtCompileTime > &matrixL, double eps, const int maxIterations, const int relaxationPeriod, const int relaxationFactor)
Definition: MultifitComputations.h:380