NORM
Computes the norm of an array.
Function Return Value
Norm of A. This is a scalar for the case where A is rank 1 or rank 2. For rank-3 arrays, the norms of each rank-2 array, in dimension 3, are computed. (Output)
Required Argument
A — An array of rank-1, rank-2, or rank-3, containing the values for which the norm is to be computed. It may be real, double, complex, or double complex. (Input)
Optional Arguments, Packaged Options
TYPE —Integer indicating the type of norm to be computed.
1 = l1
2 = l2 (default)
huge(1) = l∞
Use of the option number ?_reset_default_norm will switch the default from the l2 to the l1 or l∞ norms. (Input)
The option and derived type names are given in the following tables:
Option Names for NORM |
Option Value |
?_norm_for_lin_sol_svd |
1 |
?_reset_default_norm |
2 |
Name of Unallocated Option Array to Use for Setting Options |
Use |
Derived Type |
?_norm_options(:) |
Use when setting options for calls hereafter. |
?_options |
?_norm_options_once(:) |
Use when setting options for next call only. |
?_options |
For a description on how to use these options, see Matrix Optional Data Changes. See LIN_SOL_SVD in Chapter 1, “Linear Systems” for the specific options for these routines.
FORTRAN 90 Interface
NORM (A [, …])
Description
Computes the l2, l1, or l∞ norm. The l1 and l∞ norms are likely to be less expensive to compute than the l2 norm.
If the l2 norm is required, this function uses LIN_SOL_SVD (see Chapter 1, “Linear Systems”), to compute the largest singular value of A. For the other norms, Fortran 90 intrinsics are used.
Examples
Compute three norms of an array
use norm_int
real (kind(1e0)) A(5), n_1, n_2, n_inf
A = rand (A)
! I1
n_1 = norm(A, TYPE=1)
write (*,*) n_1
! I2
n_2 = norm(A)
write (*,*) n_2
! I infinity
n_inf = norm(A, TYPE=huge(1))
write (*,*) n_inf
end
Parallel Example (parallel_ex15.f90)
A “Polar Decomposition” of several matrices are computed. The box data type and the SVD() function are used. Orthogonality and small residuals are checked to verify that the results are correct.
use linear_operators
use mpi_setup_int
implicit none
! This is Parallel Example 15 using operators and
! functions for a polar decomposition.
integer, parameter :: n=33, nr=3
real(kind(1d0)) :: one=1d0, zero=0d0
real(kind(1d0)),dimension(n,n,nr) :: A, P, Q, &
S_D(n,nr), U_D, V_D
real(kind(1d0)) TEMP1(nr), TEMP2(nr)
! Setup for MPI:
mp_nprocs = mp_setup()
! Generate a random matrix.
if(mp_rank == 0) A = rand(A)
! Compute the singular value decomposition.
S_D = SVD(A, U=U_D, V=V_D)
! Compute the (left) orthogonal factor.
P = U_D .xt. V_D
! Compute the (right) self-adjoint factor.
Q = V_D .x. diag(S_D) .xt. V_D
! Check the results for orthogonality and
! small residuals.
TEMP1 = NORM(spread(EYE(n),3,nr) - (p .xt. p))
TEMP2 = NORM(A -(P .X. Q)) / NORM(A)
if (ALL(TEMP1 <= sqrt(epsilon(one))) .and. &
ALL(TEMP2 <= sqrt(epsilon(one)))) then
if(mp_rank == 0)&
write (*,*) 'Parallel Example 15 is correct.'
end if
! See to any error messages and exit MPI.
mp_nprocs = mp_setup('Final')
end