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