GNU Octave  6.2.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-2021 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 "quit.h"
63 #include "svd.h"
64 
65 // Theory: norm accumulator is an object that has an accum method able
66 // to handle both real and complex element, and a cast operator
67 // returning the intermediate norm. Reference: Higham, N. "Estimating
68 // the Matrix p-Norm." Numer. Math. 62, 539-555, 1992.
69 
70 // norm accumulator for the p-norm
71 template <typename R>
73 {
74  R p,scl,sum;
75 public:
76  norm_accumulator_p () { } // we need this one for Array
77  norm_accumulator_p (R pp) : p(pp), scl(0), sum(1) { }
78 
79  template <typename U>
80  void accum (U val)
81  {
82  octave_quit ();
83  R t = std::abs (val);
84  if (scl == t) // we need this to handle Infs properly
85  sum += 1;
86  else if (scl < t)
87  {
88  sum *= std::pow (scl/t, p);
89  sum += 1;
90  scl = t;
91  }
92  else if (t != 0)
93  sum += std::pow (t/scl, p);
94  }
95  operator R () { return scl * std::pow (sum, 1/p); }
96 };
97 
98 // norm accumulator for the minus p-pseudonorm
99 template <typename R>
101 {
102  R p,scl,sum;
103 public:
104  norm_accumulator_mp () { } // we need this one for Array
105  norm_accumulator_mp (R pp) : p(pp), scl(0), sum(1) { }
106 
107  template <typename U>
108  void accum (U val)
109  {
110  octave_quit ();
111  R t = 1 / std::abs (val);
112  if (scl == t)
113  sum += 1;
114  else if (scl < t)
115  {
116  sum *= std::pow (scl/t, p);
117  sum += 1;
118  scl = t;
119  }
120  else if (t != 0)
121  sum += std::pow (t/scl, p);
122  }
123  operator R () { return scl * std::pow (sum, -1/p); }
124 };
125 
126 // norm accumulator for the 2-norm (euclidean)
127 template <typename R>
129 {
130  R scl,sum;
131  static R pow2 (R x) { return x*x; }
132 public:
133  norm_accumulator_2 () : scl(0), sum(1) { }
134 
135  void accum (R val)
136  {
137  R t = std::abs (val);
138  if (scl == t)
139  sum += 1;
140  else if (scl < t)
141  {
142  sum *= pow2 (scl/t);
143  sum += 1;
144  scl = t;
145  }
146  else if (t != 0)
147  sum += pow2 (t/scl);
148  }
149 
150  void accum (std::complex<R> val)
151  {
152  accum (val.real ());
153  accum (val.imag ());
154  }
155 
156  operator R () { return scl * std::sqrt (sum); }
157 };
158 
159 // norm accumulator for the 1-norm (city metric)
160 template <typename R>
162 {
163  R sum;
164 public:
165  norm_accumulator_1 () : sum (0) { }
166  template <typename U>
167  void accum (U val)
168  {
169  sum += std::abs (val);
170  }
171  operator R () { return sum; }
172 };
173 
174 // norm accumulator for the inf-norm (max metric)
175 template <typename R>
177 {
178  R max;
179 public:
181  template <typename U>
182  void accum (U val)
183  {
184  if (octave::math::isnan (val))
186  else
187  max = std::max (max, std::abs (val));
188  }
189  operator R () { return max; }
190 };
191 
192 // norm accumulator for the -inf pseudonorm (min abs value)
193 template <typename R>
195 {
196  R min;
197 public:
198  norm_accumulator_minf () : min (octave::numeric_limits<R>::Inf ()) { }
199  template <typename U>
200  void accum (U val)
201  {
202  if (octave::math::isnan (val))
204  else
205  min = std::min (min, std::abs (val));
206  }
207  operator R () { return min; }
208 };
209 
210 // norm accumulator for the 0-pseudonorm (hamming distance)
211 template <typename R>
213 {
214  unsigned int num;
215 public:
216  norm_accumulator_0 () : num (0) { }
217  template <typename U>
218  void accum (U val)
219  {
220  if (val != static_cast<U> (0)) ++num;
221  }
222  operator R () { return num; }
223 };
224 
225 // OK, we're armed :) Now let's go for the fun
226 
227 template <typename T, typename R, typename ACC>
228 inline void vector_norm (const Array<T>& v, R& res, ACC acc)
229 {
230  for (octave_idx_type i = 0; i < v.numel (); i++)
231  acc.accum (v(i));
232 
233  res = acc;
234 }
235 
236 // dense versions
237 template <typename T, typename R, typename ACC>
238 void column_norms (const MArray<T>& m, MArray<R>& res, ACC acc)
239 {
240  res = MArray<R> (dim_vector (1, m.columns ()));
241  for (octave_idx_type j = 0; j < m.columns (); j++)
242  {
243  ACC accj = acc;
244  for (octave_idx_type i = 0; i < m.rows (); i++)
245  accj.accum (m(i, j));
246 
247  res.xelem (j) = accj;
248  }
249 }
250 
251 template <typename T, typename R, typename ACC>
252 void row_norms (const MArray<T>& m, MArray<R>& res, ACC acc)
253 {
254  res = MArray<R> (dim_vector (m.rows (), 1));
255  std::vector<ACC> acci (m.rows (), acc);
256  for (octave_idx_type j = 0; j < m.columns (); j++)
257  {
258  for (octave_idx_type i = 0; i < m.rows (); i++)
259  acci[i].accum (m(i, j));
260  }
261 
262  for (octave_idx_type i = 0; i < m.rows (); i++)
263  res.xelem (i) = acci[i];
264 }
265 
266 // sparse versions
267 template <typename T, typename R, typename ACC>
268 void column_norms (const MSparse<T>& m, MArray<R>& res, ACC acc)
269 {
270  res = MArray<R> (dim_vector (1, m.columns ()));
271  for (octave_idx_type j = 0; j < m.columns (); j++)
272  {
273  ACC accj = acc;
274  for (octave_idx_type k = m.cidx (j); k < m.cidx (j+1); k++)
275  accj.accum (m.data (k));
276 
277  res.xelem (j) = accj;
278  }
279 }
280 
281 template <typename T, typename R, typename ACC>
282 void row_norms (const MSparse<T>& m, MArray<R>& res, ACC acc)
283 {
284  res = MArray<R> (dim_vector (m.rows (), 1));
285  std::vector<ACC> acci (m.rows (), acc);
286  for (octave_idx_type j = 0; j < m.columns (); j++)
287  {
288  for (octave_idx_type k = m.cidx (j); k < m.cidx (j+1); k++)
289  acci[m.ridx (k)].accum (m.data (k));
290  }
291 
292  for (octave_idx_type i = 0; i < m.rows (); i++)
293  res.xelem (i) = acci[i];
294 }
295 
296 // now the dispatchers
297 #define DEFINE_DISPATCHER(FUNC_NAME, ARG_TYPE, RES_TYPE) \
298  template <typename T, typename R> \
299  RES_TYPE FUNC_NAME (const ARG_TYPE& v, R p) \
300  { \
301  RES_TYPE res; \
302  if (p == 2) \
303  FUNC_NAME (v, res, norm_accumulator_2<R> ()); \
304  else if (p == 1) \
305  FUNC_NAME (v, res, norm_accumulator_1<R> ()); \
306  else if (lo_ieee_isinf (p)) \
307  { \
308  if (p > 0) \
309  FUNC_NAME (v, res, norm_accumulator_inf<R> ()); \
310  else \
311  FUNC_NAME (v, res, norm_accumulator_minf<R> ()); \
312  } \
313  else if (p == 0) \
314  FUNC_NAME (v, res, norm_accumulator_0<R> ()); \
315  else if (p > 0) \
316  FUNC_NAME (v, res, norm_accumulator_p<R> (p)); \
317  else \
318  FUNC_NAME (v, res, norm_accumulator_mp<R> (p)); \
319  return res; \
320  }
321 
327 
328 // The approximate subproblem in Higham's method. Find lambda and mu such that
329 // norm ([lambda, mu], p) == 1 and norm (y*lambda + col*mu, p) is maximized.
330 // Real version. As in Higham's paper.
331 template <typename ColVectorT, typename R>
332 static void
333 higham_subp (const ColVectorT& y, const ColVectorT& col,
334  octave_idx_type nsamp, R p, R& lambda, R& mu)
335 {
336  R nrm = 0;
337  for (octave_idx_type i = 0; i < nsamp; i++)
338  {
339  octave_quit ();
340  R fi = i * static_cast<R> (M_PI) / nsamp;
341  R lambda1 = cos (fi);
342  R mu1 = sin (fi);
343  R lmnr = std::pow (std::pow (std::abs (lambda1), p) +
344  std::pow (std::abs (mu1), p), 1/p);
345  lambda1 /= lmnr; mu1 /= lmnr;
346  R nrm1 = vector_norm (lambda1 * y + mu1 * col, p);
347  if (nrm1 > nrm)
348  {
349  lambda = lambda1;
350  mu = mu1;
351  nrm = nrm1;
352  }
353  }
354 }
355 
356 // Complex version. Higham's paper does not deal with complex case, so we use
357 // a simple extension. First, guess the magnitudes as in real version, then
358 // try to rotate lambda to improve further.
359 template <typename ColVectorT, typename R>
360 static void
361 higham_subp (const ColVectorT& y, const ColVectorT& col,
362  octave_idx_type nsamp, R p,
363  std::complex<R>& lambda, std::complex<R>& mu)
364 {
365  typedef std::complex<R> CR;
366  R nrm = 0;
367  lambda = 1.0;
368  CR lamcu = lambda / std::abs (lambda);
369  // Probe magnitudes
370  for (octave_idx_type i = 0; i < nsamp; i++)
371  {
372  octave_quit ();
373  R fi = i * static_cast<R> (M_PI) / nsamp;
374  R lambda1 = cos (fi);
375  R mu1 = sin (fi);
376  R lmnr = std::pow (std::pow (std::abs (lambda1), p) +
377  std::pow (std::abs (mu1), p), 1/p);
378  lambda1 /= lmnr; mu1 /= lmnr;
379  R nrm1 = vector_norm (lambda1 * lamcu * y + mu1 * col, p);
380  if (nrm1 > nrm)
381  {
382  lambda = lambda1 * lamcu;
383  mu = mu1;
384  nrm = nrm1;
385  }
386  }
387  R lama = std::abs (lambda);
388  // Probe orientation
389  for (octave_idx_type i = 0; i < nsamp; i++)
390  {
391  octave_quit ();
392  R fi = i * static_cast<R> (M_PI) / nsamp;
393  lamcu = CR (cos (fi), sin (fi));
394  R nrm1 = vector_norm (lama * lamcu * y + mu * col, p);
395  if (nrm1 > nrm)
396  {
397  lambda = lama * lamcu;
398  nrm = nrm1;
399  }
400  }
401 }
402 
403 // the p-dual element (should work for both real and complex)
404 template <typename T, typename R>
405 inline T elem_dual_p (T x, R p)
406 {
407  return octave::math::signum (x) * std::pow (std::abs (x), p-1);
408 }
409 
410 // the VectorT is used for vectors, but actually it has to be
411 // a Matrix type to allow all the operations. For instance SparseMatrix
412 // does not support multiplication with column/row vectors.
413 // the dual vector
414 template <typename VectorT, typename R>
415 VectorT dual_p (const VectorT& x, R p, R q)
416 {
417  VectorT res (x.dims ());
418  for (octave_idx_type i = 0; i < x.numel (); i++)
419  res.xelem (i) = elem_dual_p (x(i), p);
420  return res / vector_norm (res, q);
421 }
422 
423 // Higham's hybrid method
424 template <typename MatrixT, typename VectorT, typename R>
425 R higham (const MatrixT& m, R p, R tol, int maxiter,
426  VectorT& x)
427 {
428  x.resize (m.columns (), 1);
429  // the OSE part
430  VectorT y(m.rows (), 1, 0), z(m.rows (), 1);
431  typedef typename VectorT::element_type RR;
432  RR lambda = 0;
433  RR mu = 1;
434  for (octave_idx_type k = 0; k < m.columns (); k++)
435  {
436  octave_quit ();
437  VectorT col (m.column (k));
438  if (k > 0)
439  higham_subp (y, col, 4*k, p, lambda, mu);
440  for (octave_idx_type i = 0; i < k; i++)
441  x(i) *= lambda;
442  x(k) = mu;
443  y = lambda * y + mu * col;
444  }
445 
446  // the PM part
447  x = x / vector_norm (x, p);
448  R q = p/(p-1);
449 
450  R gamma = 0, gamma1;
451  int iter = 0;
452  while (iter < maxiter)
453  {
454  octave_quit ();
455  y = m*x;
456  gamma1 = gamma;
457  gamma = vector_norm (y, p);
458  z = dual_p (y, p, q);
459  z = z.hermitian ();
460  z = z * m;
461 
462  if (iter > 0 && (vector_norm (z, q) <= gamma
463  || (gamma - gamma1) <= tol*gamma))
464  break;
465 
466  z = z.hermitian ();
467  x = dual_p (z, q, p);
468  iter++;
469  }
470 
471  return gamma;
472 }
473 
474 // derive column vector and SVD types
475 
476 static const char *p_less1_gripe = "xnorm: p must be >= 1";
477 
478 // Static constant to control the maximum number of iterations. 100 seems to
479 // be a good value. Eventually, we can provide a means to change this
480 // constant from Octave.
481 static int max_norm_iter = 100;
482 
483 // version with SVD for dense matrices
484 template <typename MatrixT, typename VectorT, typename R>
485 R svd_matrix_norm (const MatrixT& m, R p, VectorT)
486 {
487  R res = 0;
488  if (p == 2)
489  {
492  res = fact.singular_values () (0,0);
493  }
494  else if (p == 1)
495  res = xcolnorms (m, 1).max ();
496  else if (lo_ieee_isinf (p) && p > 1)
497  res = xrownorms (m, 1).max ();
498  else if (p > 1)
499  {
500  VectorT x;
501  const R sqrteps = std::sqrt (std::numeric_limits<R>::epsilon ());
502  res = higham (m, p, sqrteps, max_norm_iter, x);
503  }
504  else
506 
507  return res;
508 }
509 
510 // SVD-free version for sparse matrices
511 template <typename MatrixT, typename VectorT, typename R>
512 R matrix_norm (const MatrixT& m, R p, VectorT)
513 {
514  R res = 0;
515  if (p == 1)
516  res = xcolnorms (m, 1).max ();
517  else if (lo_ieee_isinf (p) && p > 1)
518  res = xrownorms (m, 1).max ();
519  else if (p > 1)
520  {
521  VectorT x;
522  const R sqrteps = std::sqrt (std::numeric_limits<R>::epsilon ());
523  res = higham (m, p, sqrteps, max_norm_iter, x);
524  }
525  else
527 
528  return res;
529 }
530 
531 // and finally, here's what we've promised in the header file
532 
533 #define DEFINE_XNORM_FUNCS(PREFIX, RTYPE) \
534  OCTAVE_API RTYPE xnorm (const PREFIX##ColumnVector& x, RTYPE p) \
535  { \
536  return vector_norm (x, p); \
537  } \
538  OCTAVE_API RTYPE xnorm (const PREFIX##RowVector& x, RTYPE p) \
539  { \
540  return vector_norm (x, p); \
541  } \
542  OCTAVE_API RTYPE xnorm (const PREFIX##Matrix& x, RTYPE p) \
543  { \
544  return svd_matrix_norm (x, p, PREFIX##Matrix ()); \
545  } \
546  OCTAVE_API RTYPE xfrobnorm (const PREFIX##Matrix& x) \
547  { \
548  return vector_norm (x, static_cast<RTYPE> (2)); \
549  }
550 
553 DEFINE_XNORM_FUNCS(Float, float)
555 
556 // this is needed to avoid copying the sparse matrix for xfrobnorm
557 template <typename T, typename R>
558 inline void array_norm_2 (const T *v, octave_idx_type n, R& res)
559 {
561  for (octave_idx_type i = 0; i < n; i++)
562  acc.accum (v[i]);
563 
564  res = acc;
565 }
566 
567 #define DEFINE_XNORM_SPARSE_FUNCS(PREFIX, RTYPE) \
568  OCTAVE_API RTYPE xnorm (const Sparse##PREFIX##Matrix& x, RTYPE p) \
569  { \
570  return matrix_norm (x, p, PREFIX##Matrix ()); \
571  } \
572  OCTAVE_API RTYPE xfrobnorm (const Sparse##PREFIX##Matrix& x) \
573  { \
574  RTYPE res; \
575  array_norm_2 (x.data (), x.nnz (), res); \
576  return res; \
577  }
578 
581 
582 #define DEFINE_COLROW_NORM_FUNCS(PREFIX, RPREFIX, RTYPE) \
583  extern OCTAVE_API RPREFIX##RowVector \
584  xcolnorms (const PREFIX##Matrix& m, RTYPE p) \
585  { \
586  return column_norms (m, p); \
587  } \
588  extern OCTAVE_API RPREFIX##ColumnVector \
589  xrownorms (const PREFIX##Matrix& m, RTYPE p) \
590  { \
591  return row_norms (m, p); \
592  } \
593 
596 DEFINE_COLROW_NORM_FUNCS(Float, Float, float)
598 
600 DEFINE_COLROW_NORM_FUNCS(SparseComplex, , double)
#define Inf
Definition: Faddeeva.cc:247
#define NaN
Definition: Faddeeva.cc:248
charNDArray max(char d, const charNDArray &m)
Definition: chNDArray.cc:230
charNDArray min(char d, const charNDArray &m)
Definition: chNDArray.cc:207
N Dimensional Array with copy-on-write semantics.
Definition: Array.h:128
T & xelem(octave_idx_type n)
Size of the specified dimension.
Definition: Array.h:469
octave_idx_type numel(void) const
Number of elements in the array.
Definition: Array.h:377
double max(void) const
Definition: dColVector.cc:261
Template for N-dimensional array classes with like-type math operators.
Definition: MArray.h:63
double max(void) const
Definition: dRowVector.cc:223
Definition: Sparse.h:49
Vector representing the dimensions (size) of an Array.
Definition: dim-vector.h:95
unsigned int num
Definition: oct-norm.cc:214
void accum(U val)
Definition: oct-norm.cc:218
void accum(U val)
Definition: oct-norm.cc:167
void accum(std::complex< R > val)
Definition: oct-norm.cc:150
static R pow2(R x)
Definition: oct-norm.cc:131
void accum(R val)
Definition: oct-norm.cc:135
void accum(U val)
Definition: oct-norm.cc:182
void accum(U val)
Definition: oct-norm.cc:200
void accum(U val)
Definition: oct-norm.cc:108
norm_accumulator_mp(R pp)
Definition: oct-norm.cc:105
norm_accumulator_p(R pp)
Definition: oct-norm.cc:77
void accum(U val)
Definition: oct-norm.cc:80
DM_T singular_values(void) const
Definition: svd.h:88
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:120
F77_RET_T const F77_DBLE * x
T octave_idx_type m
Definition: mx-inlines.cc:773
octave_idx_type n
Definition: mx-inlines.cc:753
double signum(double x)
Definition: lo-mappers.h:222
bool isnan(bool)
Definition: lo-mappers.h:178
static double fi[256]
Definition: randmtzig.cc:447
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:415
void array_norm_2(const T *v, octave_idx_type n, R &res)
Definition: oct-norm.cc:558
R matrix_norm(const MatrixT &m, R p, VectorT)
Definition: oct-norm.cc:512
#define DEFINE_XNORM_SPARSE_FUNCS(PREFIX, RTYPE)
Definition: oct-norm.cc:567
OCTAVE_API RowVector xcolnorms(const Matrix &m, double p)
Definition: oct-norm.cc:594
#define DEFINE_COLROW_NORM_FUNCS(PREFIX, RPREFIX, RTYPE)
Definition: oct-norm.cc:582
void column_norms(const MArray< T > &m, MArray< R > &res, ACC acc)
Definition: oct-norm.cc:238
T elem_dual_p(T x, R p)
Definition: oct-norm.cc:405
static const char * p_less1_gripe
Definition: oct-norm.cc:476
static void higham_subp(const ColVectorT &y, const ColVectorT &col, octave_idx_type nsamp, R p, R &lambda, R &mu)
Definition: oct-norm.cc:333
void vector_norm(const Array< T > &v, R &res, ACC acc)
Definition: oct-norm.cc:228
void row_norms(const MArray< T > &m, MArray< R > &res, ACC acc)
Definition: oct-norm.cc:252
R svd_matrix_norm(const MatrixT &m, R p, VectorT)
Definition: oct-norm.cc:485
#define DEFINE_DISPATCHER(FUNC_NAME, ARG_TYPE, RES_TYPE)
Definition: oct-norm.cc:297
R higham(const MatrixT &m, R p, R tol, int maxiter, VectorT &x)
Definition: oct-norm.cc:425
static int max_norm_iter
Definition: oct-norm.cc:481
OCTAVE_API ColumnVector xrownorms(const Matrix &m, double p)
Definition: oct-norm.cc:594
#define DEFINE_XNORM_FUNCS(PREFIX, RTYPE)
Definition: oct-norm.cc:533
static T abs(T x)
Definition: pr-output.cc:1678