Computes the inverse of a square matrix.
More...
|
| mtx_inverse_dbl |
|
| mtx_inverse_cmplx |
|
Computes the inverse of a square matrix.
- Syntax
subroutine mtx_inverse(real(real64) a(:,:), optional integer(int32) iwork, optional real(real64) work(:), optional integer olwork, optional class(errors) err)
subroutine mtx_inverse(complex(real64) a(:,:), optional integer(int32) iwork, optional complex(real64) work(:), optional integer olwork, optional class(errors) err)
- Parameters
-
[in,out] | a | On input, the N-by-N matrix to invert. On output, the inverted matrix. |
[out] | iwork | An optional N-element integer workspace array. |
[out] | work | An optional input, that if provided, prevents any local memory allocation. If not provided, the memory required is allocated within. If provided, the length of the array must be at least olwork . |
[out] | olwork | An optional output used to determine workspace size. If supplied, the routine determines the optimal size for work , and returns without performing any actual calculations. |
[in,out] | err | An optional errors-based object that if provided can be used to retrieve information relating to any errors encountered during execution. If not provided, a default implementation of the errors class is used internally to provide error handling. Possible errors and warning messages that may be encountered are as follows.
- LA_ARRAY_SIZE_ERROR: Occurs if
a is not square. Will also occur if incorrectly sized workspace arrays are provided.
- LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.
- LA_SINGULAR_MATRIX_ERROR: Occurs if the input matrix is singular.
|
- Notes
- This routine utilizes the LAPACK routines DGETRF to perform an LU factorization of the matrix, and DGETRI to invert the LU factored matrix (ZGETRF and ZGETRI in the complex case).
- See Also
-
- Usage
- The following example illustrates the inversion of a 3-by-3 matrix.
program example
use iso_fortran_env, only : real64, int32
implicit none
real(real64) :: a(3,3), ai(3,3), c(3,3)
integer(int32) :: i
a = reshape([1.0d0, 4.0d0, 7.0d0, 2.0d0, 5.0d0, 8.0d0, 3.0d0, 6.0d0, &
0.0d0], [3, 3])
ai = a
print '(A)', "Inverse:"
do i = 1, size(ai, 1)
print *, ai(i,:)
end do
print '(A)', "A * A**-1:"
do i = 1, size(c, 1)
print *, c(i,:)
end do
end program
Performs sparse matrix multiplication C = A * B.
Computes the inverse of a square matrix.
Provides a set of common linear algebra routines.
The above program produces the following output. Inverse:
-1.7777777777777777 0.88888888888888884 -0.11111111111111110
1.5555555555555556 -0.77777777777777779 0.22222222222222221
-0.11111111111111119 0.22222222222222227 -0.11111111111111112
A * A**-1:
0.99999999999999989 5.5511151231257827E-017 -4.1633363423443370E-017
5.5511151231257827E-017 1.0000000000000000 -8.3266726846886741E-017
1.7763568394002505E-015 -8.8817841970012523E-016 1.0000000000000000
Definition at line 3051 of file linalg.f90.
The documentation for this interface was generated from the following file: