GNU Octave  9.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
oct-norm.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // Copyright (C) 2008-2024 The Octave Project Developers
4 //
5 // See the file COPYRIGHT.md in the top-level directory of this
6 // distribution or <https://octave.org/copyright/>.
7 //
8 // This file is part of Octave.
9 //
10 // Octave is free software: you can redistribute it and/or modify it
11 // under the terms of the GNU General Public License as published by
12 // the Free Software Foundation, either version 3 of the License, or
13 // (at your option) any later version.
14 //
15 // Octave is distributed in the hope that it will be useful, but
16 // WITHOUT ANY WARRANTY; without even the implied warranty of
17 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 // GNU General Public License for more details.
19 //
20 // You should have received a copy of the GNU General Public License
21 // along with Octave; see the file COPYING. If not, see
22 // <https://www.gnu.org/licenses/>.
23 //
24 ////////////////////////////////////////////////////////////////////////
25 
26 #if defined (HAVE_CONFIG_H)
27 # include "config.h"
28 #endif
29 
30 #include <cmath>
31 
32 #include <algorithm>
33 #include <limits>
34 #include <vector>
35 
36 #include "Array.h"
37 #include "CColVector.h"
38 #include "CMatrix.h"
39 #include "CRowVector.h"
40 #include "CSparse.h"
41 #include "MArray.h"
42 #include "dColVector.h"
43 #include "dDiagMatrix.h"
44 #include "dMatrix.h"
45 #include "dRowVector.h"
46 #include "dSparse.h"
47 #include "fCColVector.h"
48 #include "fCMatrix.h"
49 #include "fCRowVector.h"
50 #include "fColVector.h"
51 #include "fDiagMatrix.h"
52 #include "fMatrix.h"
53 #include "fRowVector.h"
54 #include "lo-error.h"
55 #include "lo-ieee.h"
56 #include "lo-mappers.h"
57 #include "mx-cm-s.h"
58 #include "mx-fcm-fs.h"
59 #include "mx-fs-fcm.h"
60 #include "mx-s-cm.h"
61 #include "oct-cmplx.h"
62 #include "oct-norm.h"
63 #include "quit.h"
64 #include "svd.h"
65 
67 
68 // Theory: norm accumulator is an object that has an accum method able
69 // to handle both real and complex element, and a cast operator
70 // returning the intermediate norm. Reference: Higham, N. "Estimating
71 // the Matrix p-Norm." Numer. Math. 62, 539-555, 1992.
72 
73 // norm accumulator for the p-norm
74 template <typename R>
75 class norm_accumulator_p
76 {
77 public:
78 
79  norm_accumulator_p (R pp) : m_p(pp), m_scl(0), m_sum(1) { }
80 
81  OCTAVE_DEFAULT_CONSTRUCT_COPY_MOVE_DELETE (norm_accumulator_p)
82 
83  template <typename U>
84  void accum (U val)
85  {
86  octave_quit ();
87  R t = std::abs (val);
88  if (m_scl == t) // we need this to handle Infs properly
89  m_sum += 1;
90  else if (m_scl < t)
91  {
92  m_sum *= std::pow (m_scl/t, m_p);
93  m_sum += 1;
94  m_scl = t;
95  }
96  else if (t != 0)
97  m_sum += std::pow (t/m_scl, m_p);
98  }
99 
100  operator R () { return m_scl * std::pow (m_sum, 1/m_p); }
101 
102 private:
103  R m_p, m_scl, m_sum;
104 };
105 
106 // norm accumulator for the minus p-pseudonorm
107 template <typename R>
108 class norm_accumulator_mp
109 {
110 public:
111 
112  norm_accumulator_mp (R pp) : m_p(pp), m_scl(0), m_sum(1) { }
113 
114  OCTAVE_DEFAULT_CONSTRUCT_COPY_MOVE_DELETE (norm_accumulator_mp)
115 
116  template <typename U>
117  void accum (U val)
118  {
119  octave_quit ();
120  R t = 1 / std::abs (val);
121  if (m_scl == t)
122  m_sum += 1;
123  else if (m_scl < t)
124  {
125  m_sum *= std::pow (m_scl/t, m_p);
126  m_sum += 1;
127  m_scl = t;
128  }
129  else if (t != 0)
130  m_sum += std::pow (t/m_scl, m_p);
131  }
132 
133  operator R () { return m_scl * std::pow (m_sum, -1/m_p); }
134 
135 private:
136  R m_p, m_scl, m_sum;
137 };
138 
139 // norm accumulator for the 2-norm (euclidean)
140 template <typename R>
141 class norm_accumulator_2
142 {
143 public:
144 
145  norm_accumulator_2 () : m_scl(0), m_sum(1) { }
146 
147  OCTAVE_DEFAULT_COPY_MOVE_DELETE (norm_accumulator_2)
148 
149  void accum (R val)
150  {
151  R t = std::abs (val);
152  if (m_scl == t)
153  m_sum += 1;
154  else if (m_scl < t)
155  {
156  m_sum *= pow2 (m_scl/t);
157  m_sum += 1;
158  m_scl = t;
159  }
160  else if (t != 0)
161  m_sum += pow2 (t/m_scl);
162  }
163 
164  void accum (std::complex<R> val)
165  {
166  accum (val.real ());
167  accum (val.imag ());
168  }
169 
170  operator R () { return m_scl * std::sqrt (m_sum); }
171 
172 private:
173  static inline R pow2 (R x) { return x*x; }
174 
175  //--------
176 
177  R m_scl, m_sum;
178 };
179 
180 // norm accumulator for the 1-norm (city metric)
181 template <typename R>
182 class norm_accumulator_1
183 {
184 public:
185 
186  norm_accumulator_1 () : m_sum (0) { }
187 
188  OCTAVE_DEFAULT_COPY_MOVE_DELETE (norm_accumulator_1)
189 
190  template <typename U>
191  void accum (U val)
192  {
193  m_sum += std::abs (val);
194  }
195 
196  operator R () { return m_sum; }
197 
198 private:
199  R m_sum;
200 };
201 
202 // norm accumulator for the inf-norm (max metric)
203 template <typename R>
204 class norm_accumulator_inf
205 {
206 public:
207 
208  norm_accumulator_inf () : m_max (0) { }
209 
210  OCTAVE_DEFAULT_COPY_MOVE_DELETE (norm_accumulator_inf)
211 
212  template <typename U>
213  void accum (U val)
214  {
215  if (math::isnan (val))
216  m_max = numeric_limits<R>::NaN ();
217  else
218  m_max = std::max (m_max, std::abs (val));
219  }
220 
221  operator R () { return m_max; }
222 
223 private:
224  R m_max;
225 };
226 
227 // norm accumulator for the -inf pseudonorm (min abs value)
228 template <typename R>
229 class norm_accumulator_minf
230 {
231 public:
232 
233  norm_accumulator_minf () : m_min (numeric_limits<R>::Inf ()) { }
234 
235  OCTAVE_DEFAULT_COPY_MOVE_DELETE (norm_accumulator_minf)
236 
237  template <typename U>
238  void accum (U val)
239  {
240  if (math::isnan (val))
241  m_min = numeric_limits<R>::NaN ();
242  else
243  m_min = std::min (m_min, std::abs (val));
244  }
245 
246  operator R () { return m_min; }
247 
248 private:
249  R m_min;
250 };
251 
252 // norm accumulator for the 0-pseudonorm (hamming distance)
253 template <typename R>
254 class norm_accumulator_0
255 {
256 public:
257 
258  norm_accumulator_0 () : m_num (0) { }
259 
260  OCTAVE_DEFAULT_COPY_MOVE_DELETE (norm_accumulator_0)
261 
262  template <typename U>
263  void accum (U val)
264  {
265  if (val != static_cast<U> (0)) ++m_num;
266  }
267 
268  operator R () { return m_num; }
269 
270 private:
271  unsigned int m_num;
272 };
273 
274 // OK, we're armed :) Now let's go for the fun
275 
276 template <typename T, typename R, typename ACC>
277 inline void
278 vector_norm (const Array<T>& v, R& res, ACC acc)
279 {
280  for (octave_idx_type i = 0; i < v.numel (); i++)
281  acc.accum (v(i));
282 
283  res = acc;
284 }
285 
286 // dense versions
287 template <typename T, typename R, typename ACC>
288 void
289 column_norms (const MArray<T>& m, MArray<R>& res, ACC acc)
290 {
291  res = MArray<R> (dim_vector (1, m.columns ()));
292  for (octave_idx_type j = 0; j < m.columns (); j++)
293  {
294  ACC accj = acc;
295  for (octave_idx_type i = 0; i < m.rows (); i++)
296  accj.accum (m(i, j));
297 
298  res.xelem (j) = accj;
299  }
300 }
301 
302 template <typename T, typename R, typename ACC>
303 void
304 row_norms (const MArray<T>& m, MArray<R>& res, ACC acc)
305 {
306  res = MArray<R> (dim_vector (m.rows (), 1));
307  std::vector<ACC> acci (m.rows (), acc);
308  for (octave_idx_type j = 0; j < m.columns (); j++)
309  {
310  for (octave_idx_type i = 0; i < m.rows (); i++)
311  acci[i].accum (m(i, j));
312  }
313 
314  for (octave_idx_type i = 0; i < m.rows (); i++)
315  res.xelem (i) = acci[i];
316 }
317 
318 // sparse versions
319 template <typename T, typename R, typename ACC>
320 void
321 column_norms (const MSparse<T>& m, MArray<R>& res, ACC acc)
322 {
323  res = MArray<R> (dim_vector (1, m.columns ()));
324  for (octave_idx_type j = 0; j < m.columns (); j++)
325  {
326  ACC accj = acc;
327  for (octave_idx_type k = m.cidx (j); k < m.cidx (j+1); k++)
328  accj.accum (m.data (k));
329 
330  res.xelem (j) = accj;
331  }
332 }
333 
334 template <typename T, typename R, typename ACC>
335 void
336 row_norms (const MSparse<T>& m, MArray<R>& res, ACC acc)
337 {
338  res = MArray<R> (dim_vector (m.rows (), 1));
339  std::vector<ACC> acci (m.rows (), acc);
340  for (octave_idx_type j = 0; j < m.columns (); j++)
341  {
342  for (octave_idx_type k = m.cidx (j); k < m.cidx (j+1); k++)
343  acci[m.ridx (k)].accum (m.data (k));
344  }
345 
346  for (octave_idx_type i = 0; i < m.rows (); i++)
347  res.xelem (i) = acci[i];
348 }
349 
350 // now the dispatchers
351 #define DEFINE_DISPATCHER(FCN_NAME, ARG_TYPE, RES_TYPE) \
352 template <typename T, typename R> \
353 RES_TYPE FCN_NAME (const ARG_TYPE& v, R p) \
354 { \
355  RES_TYPE res; \
356  if (p == 2) \
357  FCN_NAME (v, res, norm_accumulator_2<R> ()); \
358  else if (p == 1) \
359  FCN_NAME (v, res, norm_accumulator_1<R> ()); \
360  else if (lo_ieee_isinf (p)) \
361  { \
362  if (p > 0) \
363  FCN_NAME (v, res, norm_accumulator_inf<R> ()); \
364  else \
365  FCN_NAME (v, res, norm_accumulator_minf<R> ()); \
366  } \
367  else if (p == 0) \
368  FCN_NAME (v, res, norm_accumulator_0<R> ()); \
369  else if (p > 0) \
370  FCN_NAME (v, res, norm_accumulator_p<R> (p)); \
371  else \
372  FCN_NAME (v, res, norm_accumulator_mp<R> (p)); \
373  return res; \
374 }
375 
381 
382 // The approximate subproblem in Higham's method. Find lambda and mu such
383 // that norm ([lambda, mu], p) == 1 and norm (y*lambda + col*mu, p) is
384 // maximized.
385 // Real version. As in Higham's paper.
386 template <typename ColVectorT, typename R>
387 static void
388 higham_subp (const ColVectorT& y, const ColVectorT& col,
389  octave_idx_type nsamp, R p, R& lambda, R& mu)
390 {
391  R nrm = 0;
392  for (octave_idx_type i = 0; i < nsamp; i++)
393  {
394  octave_quit ();
395  R fi = i * static_cast<R> (M_PI) / nsamp;
396  R lambda1 = cos (fi);
397  R mu1 = sin (fi);
398  R lmnr = std::pow (std::pow (std::abs (lambda1), p) +
399  std::pow (std::abs (mu1), p), 1/p);
400  lambda1 /= lmnr; mu1 /= lmnr;
401  R nrm1 = vector_norm (lambda1 * y + mu1 * col, p);
402  if (nrm1 > nrm)
403  {
404  lambda = lambda1;
405  mu = mu1;
406  nrm = nrm1;
407  }
408  }
409 }
410 
411 // Complex version. Higham's paper does not deal with complex case, so we
412 // use a simple extension. First, guess the magnitudes as in real version,
413 // then try to rotate lambda to improve further.
414 template <typename ColVectorT, typename R>
415 static void
416 higham_subp (const ColVectorT& y, const ColVectorT& col,
417  octave_idx_type nsamp, R p,
418  std::complex<R>& lambda, std::complex<R>& mu)
419 {
420  typedef std::complex<R> CR;
421  R nrm = 0;
422  lambda = 1.0;
423  CR lamcu = lambda / std::abs (lambda);
424  // Probe magnitudes
425  for (octave_idx_type i = 0; i < nsamp; i++)
426  {
427  octave_quit ();
428  R fi = i * static_cast<R> (M_PI) / nsamp;
429  R lambda1 = cos (fi);
430  R mu1 = sin (fi);
431  R lmnr = std::pow (std::pow (std::abs (lambda1), p) +
432  std::pow (std::abs (mu1), p), 1/p);
433  lambda1 /= lmnr; mu1 /= lmnr;
434  R nrm1 = vector_norm (lambda1 * lamcu * y + mu1 * col, p);
435  if (nrm1 > nrm)
436  {
437  lambda = lambda1 * lamcu;
438  mu = mu1;
439  nrm = nrm1;
440  }
441  }
442  R lama = std::abs (lambda);
443  // Probe orientation
444  for (octave_idx_type i = 0; i < nsamp; i++)
445  {
446  octave_quit ();
447  R fi = i * static_cast<R> (M_PI) / nsamp;
448  lamcu = CR (cos (fi), sin (fi));
449  R nrm1 = vector_norm (lama * lamcu * y + mu * col, p);
450  if (nrm1 > nrm)
451  {
452  lambda = lama * lamcu;
453  nrm = nrm1;
454  }
455  }
456 }
457 
458 // the p-dual element (should work for both real and complex)
459 template <typename T, typename R>
460 inline T
461 elem_dual_p (T x, R p)
462 {
463  return math::signum (x) * std::pow (std::abs (x), p-1);
464 }
465 
466 // the VectorT is used for vectors, but actually it has to be
467 // a Matrix type to allow all the operations. For instance SparseMatrix
468 // does not support multiplication with column/row vectors.
469 // the dual vector
470 template <typename VectorT, typename R>
471 VectorT
472 dual_p (const VectorT& x, R p, R q)
473 {
474  VectorT res (x.dims ());
475  for (octave_idx_type i = 0; i < x.numel (); i++)
476  res.xelem (i) = elem_dual_p (x(i), p);
477  return res / vector_norm (res, q);
478 }
479 
480 // Higham's hybrid method
481 template <typename MatrixT, typename VectorT, typename R>
482 R
483 higham (const MatrixT& m, R p, R tol, int maxiter,
484  VectorT& x)
485 {
486  x.resize (m.columns (), 1);
487  // the OSE part
488  VectorT y(m.rows (), 1, 0), z(m.rows (), 1);
489  typedef typename VectorT::element_type RR;
490  RR lambda = 0;
491  RR mu = 1;
492  for (octave_idx_type k = 0; k < m.columns (); k++)
493  {
494  octave_quit ();
495  VectorT col (m.column (k));
496  if (k > 0)
497  higham_subp (y, col, 4*k, p, lambda, mu);
498  for (octave_idx_type i = 0; i < k; i++)
499  x(i) *= lambda;
500  x(k) = mu;
501  y = lambda * y + mu * col;
502  }
503 
504  // the PM part
505  x = x / vector_norm (x, p);
506  R q = p/(p-1);
507 
508  R gamma = 0, gamma1;
509  int iter = 0;
510  while (iter < maxiter)
511  {
512  octave_quit ();
513  y = m*x;
514  gamma1 = gamma;
515  gamma = vector_norm (y, p);
516  z = dual_p (y, p, q);
517  z = z.hermitian ();
518  z = z * m;
519 
520  if (iter > 0 && (vector_norm (z, q) <= gamma
521  || (gamma - gamma1) <= tol*gamma))
522  break;
523 
524  z = z.hermitian ();
525  x = dual_p (z, q, p);
526  iter++;
527  }
528 
529  return gamma;
530 }
531 
532 // derive column vector and SVD types
533 
534 static const char *p_less1_gripe = "xnorm: p must be >= 1";
535 
536 // Static constant to control the maximum number of iterations. 100 seems to
537 // be a good value. Eventually, we can provide a means to change this
538 // constant from Octave.
539 static int max_norm_iter = 100;
540 
541 // version with SVD for dense matrices
542 template <typename MatrixT, typename VectorT, typename R>
543 R
544 svd_matrix_norm (const MatrixT& m, R p, VectorT)
545 {
546  R res = 0;
547  if (p == 2)
548  {
549  math::svd<MatrixT> fact (m, math::svd<MatrixT>::Type::sigma_only);
550  res = fact.singular_values () (0, 0);
551  }
552  else if (p == 1)
553  res = xcolnorms (m, static_cast<R> (1)).max ();
554  else if (lo_ieee_isinf (p) && p > 1)
555  res = xrownorms (m, static_cast<R> (1)).max ();
556  else if (p > 1)
557  {
558  VectorT x;
559  const R sqrteps = std::sqrt (std::numeric_limits<R>::epsilon ());
560  res = higham (m, p, sqrteps, max_norm_iter, x);
561  }
562  else
563  (*current_liboctave_error_handler) ("%s", p_less1_gripe);
564 
565  return res;
566 }
567 
568 // SVD-free version for sparse matrices
569 template <typename MatrixT, typename VectorT, typename R>
570 R
571 matrix_norm (const MatrixT& m, R p, VectorT)
572 {
573  R res = 0;
574  if (p == 1)
575  res = xcolnorms (m, static_cast<R> (1)).max ();
576  else if (lo_ieee_isinf (p) && p > 1)
577  res = xrownorms (m, static_cast<R> (1)).max ();
578  else if (p > 1)
579  {
580  VectorT x;
581  const R sqrteps = std::sqrt (std::numeric_limits<R>::epsilon ());
582  res = higham (m, p, sqrteps, max_norm_iter, x);
583  }
584  else
585  (*current_liboctave_error_handler) ("%s", p_less1_gripe);
586 
587  return res;
588 }
589 
590 // and finally, here's what we've promised in the header file
591 
592 #define DEFINE_XNORM_FCNS(PREFIX, RTYPE) \
593 RTYPE xnorm (const PREFIX##ColumnVector& x, RTYPE p) \
594 { \
595  return vector_norm (x, p); \
596 } \
597 RTYPE xnorm (const PREFIX##RowVector& x, RTYPE p) \
598 { \
599  return vector_norm (x, p); \
600 } \
601 RTYPE xnorm (const PREFIX##Matrix& x, RTYPE p) \
602 { \
603  return svd_matrix_norm (x, p, PREFIX##Matrix ()); \
604 } \
605 RTYPE xfrobnorm (const PREFIX##Matrix& x) \
606 { \
607  return vector_norm (x, static_cast<RTYPE> (2)); \
608 }
609 
612 DEFINE_XNORM_FCNS(Float, float)
614 
615 // this is needed to avoid copying the sparse matrix for xfrobnorm
616 template <typename T, typename R>
617 inline void array_norm_2 (const T *v, octave_idx_type n, R& res)
618 {
619  norm_accumulator_2<R> acc;
620  for (octave_idx_type i = 0; i < n; i++)
621  acc.accum (v[i]);
622 
623  res = acc;
624 }
625 
626 #define DEFINE_XNORM_SPARSE_FCNS(PREFIX, RTYPE) \
627 RTYPE xnorm (const Sparse##PREFIX##Matrix& x, RTYPE p) \
628 { \
629  return matrix_norm (x, p, PREFIX##Matrix ()); \
630 } \
631 RTYPE xfrobnorm (const Sparse##PREFIX##Matrix& x) \
632 { \
633  RTYPE res; \
634  array_norm_2 (x.data (), x.nnz (), res); \
635  return res; \
636 }
637 
640 
641 #define DEFINE_COLROW_NORM_FCNS(PREFIX, RPREFIX, RTYPE) \
642 RPREFIX##RowVector \
643 xcolnorms (const PREFIX##Matrix& m, RTYPE p) \
644 { \
645  return column_norms (m, p); \
646 } \
647 RPREFIX##ColumnVector \
648 xrownorms (const PREFIX##Matrix& m, RTYPE p) \
649 { \
650  return row_norms (m, p); \
651 } \
652 
655 DEFINE_COLROW_NORM_FCNS(Float, Float, float)
657 
659 DEFINE_COLROW_NORM_FCNS(SparseComplex, , double)
660 
661 OCTAVE_END_NAMESPACE(octave)
#define Inf
Definition: Faddeeva.cc:260
#define NaN
Definition: Faddeeva.cc:261
charNDArray max(char d, const charNDArray &m)
Definition: chNDArray.cc:230
charNDArray min(char d, const charNDArray &m)
Definition: chNDArray.cc:207
T & xelem(octave_idx_type n)
Size of the specified dimension.
Definition: Array.h:524
octave_idx_type numel() const
Number of elements in the array.
Definition: Array.h:414
double max() const
Definition: dColVector.cc:261
Template for N-dimensional array classes with like-type math operators.
Definition: MArray.h:63
double max() const
Definition: dRowVector.cc:223
Definition: Sparse.h:49
Vector representing the dimensions (size) of an Array.
Definition: dim-vector.h:94
OCTAVE_BEGIN_NAMESPACE(octave) static octave_value daspk_fcn
function gamma(X)
Definition: gamma.f:3
OCTAVE_NORETURN liboctave_error_handler current_liboctave_error_handler
Definition: lo-error.c:41
#define lo_ieee_isinf(x)
Definition: lo-ieee.h:124
double signum(double x)
Definition: lo-mappers.h:229
bool isnan(bool)
Definition: lo-mappers.h:178
F77_RET_T const F77_DBLE * x
T octave_idx_type m
Definition: mx-inlines.cc:781
octave_idx_type n
Definition: mx-inlines.cc:761
std::complex< double > Complex
Definition: oct-cmplx.h:33
std::complex< float > FloatComplex
Definition: oct-cmplx.h:34
octave_int< T > pow(const octave_int< T > &a, const octave_int< T > &b)
VectorT dual_p(const VectorT &x, R p, R q)
Definition: oct-norm.cc:472
void array_norm_2(const T *v, octave_idx_type n, R &res)
Definition: oct-norm.cc:617
R matrix_norm(const MatrixT &m, R p, VectorT)
Definition: oct-norm.cc:571
#define DEFINE_COLROW_NORM_FCNS(PREFIX, RPREFIX, RTYPE)
Definition: oct-norm.cc:641
RowVector xcolnorms(const Matrix &m, double p)
Definition: oct-norm.cc:653
void column_norms(const MArray< T > &m, MArray< R > &res, ACC acc)
Definition: oct-norm.cc:289
T elem_dual_p(T x, R p)
Definition: oct-norm.cc:461
ColumnVector xrownorms(const Matrix &m, double p)
Definition: oct-norm.cc:653
#define DEFINE_XNORM_SPARSE_FCNS(PREFIX, RTYPE)
Definition: oct-norm.cc:626
void vector_norm(const Array< T > &v, R &res, ACC acc)
Definition: oct-norm.cc:278
void row_norms(const MArray< T > &m, MArray< R > &res, ACC acc)
Definition: oct-norm.cc:304
R svd_matrix_norm(const MatrixT &m, R p, VectorT)
Definition: oct-norm.cc:544
R higham(const MatrixT &m, R p, R tol, int maxiter, VectorT &x)
Definition: oct-norm.cc:483
#define DEFINE_XNORM_FCNS(PREFIX, RTYPE)
Definition: oct-norm.cc:592
#define DEFINE_DISPATCHER(FCN_NAME, ARG_TYPE, RES_TYPE)
Definition: oct-norm.cc:351