This interface provides methods to construct Circulant matrices.
Given a vector ,
the associated Circulant matrix is the following [n x n] matrix
integer, parameter :: n = 100
real(dp) :: c(n)
type(Circulant) :: A
call random_number(c)
A = Circulant(c)
Note
Only double precision is currently supported for this matrix type.
Construct a Circulant matrix from the rank-1 array c.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=dp), | intent(in) | :: | c(:) |
Generating vector. |
Corresponding Circulant matrix.
Convert a Circulant matrix to its dense representation.
B = dense(A)
A : Matrix of Circulant type.
It is an intent(in) argument.
B : Rank-2 array representation of the matrix .
Utility function to convert a Circulant matrix to a rank-2 array.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| type(Circulant), | intent(in) | :: | A |
Input diagonal matrix. |
Output dense rank-2 array.
This interface overloads the eig interface from stdlib_linalg to
compute the eigenvalues and eigenvectors of a real-valued matrix
whose type is Circulant.
call eig(A, lambda [, left] [, right])
A : real-valued matrix of Circulant.
It is an intent(in) argument.
lambda : Rank-1 real array returning the eigenvalues of A
in increasing order. It is an intent(out) argument.
left (optional) : complex rank-2 array of the same kind as
A returning the left eigenvectors of A.
It is an intent(out) argument.
right (optional) : complex rank-2 array of the same kind as
A returning the right eigenvectors of A.
It is an intent(out) argument.
Note
Eigenvalues of a circulant matrix can be efficiently computed using
the Fast Fourier Transform of the generating vector c. Likewise,
its eigenvectors are simply the corresponding Fourier modes.
Utility function to compute the eigenvalues and eigenvectors of a
Circulant matrix.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| type(Circulant), | intent(in) | :: | A |
Input matrix. |
||
| complex(kind=dp), | intent(out) | :: | lambda(:) |
Eigenvalues. |
||
| complex(kind=dp), | intent(out), | optional | :: | left(:,:) |
Eigenvectors. |
|
| complex(kind=dp), | intent(out), | optional | :: | right(:,:) |
Eigenvectors. |
This interface overloads the eigvals interface from stdlib_linalg
to compute the eigenvalues of a real-valued matrix whose
type is Circulant.
lambda = eigvals(A)
A : real-valued matrix of Circulant type.
It is an intent(in) argument.
lambda : Vector of eigenvalues in increasing order.
Note
Eigenvalues of a circulant matrix can be efficiently computed using
the Fast Fourier Transform of the generating vector c.
Utility function to compute the eigenvalues of a real Circulant
matrix.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| type(Circulant), | intent(in) | :: | A |
Input matrix. |
Eigenvalues.
This interface overloads the Fortran intrinsic matmul for a
Circulant matrix, both for matrix-vector and matrix-matrix products.
For a matrix-matrix product , only the matrix
has to be a Circulant matrix. Both and need to be
standard Fortran rank-2 arrays. All the underlying functions are
defined as pure.
y = matmul(A, x)
Compute the matrix-vector product for a Circulant matrix .
Both x and y are rank-1 arrays with the same kind as A.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| type(Circulant), | intent(in) | :: | A |
Input matrix. |
||
| real(kind=dp), | intent(in) | :: | x(:) |
Input vector. |
Output vector.
Compute the matrix-matrix product for a Circulant matrix A.
Both X and Y are rank-2 arrays with the same kind as A.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| type(Circulant), | intent(in) | :: | A |
Input matrix. |
||
| real(kind=dp), | intent(in) | :: | x(:,:) |
Input matrix. |
Output matrix.
Utility function to perform a scalar multiplication with a
Circulant matrix.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| type(Circulant), | intent(in) | :: | A | |||
| real(kind=dp), | intent(in) | :: | alpha |
Utility function to perform a scalar multiplication with a
Circulant matrix.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=dp), | intent(in) | :: | alpha | |||
| type(Circulant), | intent(in) | :: | A |
Utility function to return the shape of a Circulant matrix.
Utility function to get the shape of a Circulant matrix.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| type(Circulant), | intent(in) | :: | A |
Input matrix. |
Shape of the matrix.
Utility function to return the size of a Circulant matrix along
a given dimension.
Utility function to return the size of Circulant matrix along a
given dimension.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| type(Circulant), | intent(in) | :: | A |
Input matrix. |
||
| integer(kind=ilp), | intent(in), | optional | :: | dim |
Queried dimension. |
Size of the matrix along the dimension dim.
This interface overloads the solve interface from stdlib_linalg
for solving a linear system where is a
Circulant matrix. It also enables to solve a linear system with
multiple right-hand sides.
x = solve(A, b)
A : Matrix of Circulant type.
It is an intent(in) argument.
b : Rank-1 or rank-2 array defining the right-hand side(s).
It is an intent(in) argument.
x : Solution of the linear system.
It has the same type and shape as b.
Note
Linear systems characterized by a circulant matrix can be solved
efficiently in operations using the
Fast Fourier Transform algorithm available via fftpack.
Solve the linear system , where A is Circulant and
B is a rank-2 array. The solution matrix X has the same
dimension and kind as the right-hand side matrix B.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| type(Circulant), | intent(in) | :: | A |
Coefficient matrix. |
||
| real(kind=dp), | intent(in) | :: | B(:,:) |
Right-hand side vectors. |
Solution vectors.
Solve the linear system where is Circulant and
b a standard rank-1 array. The solution vector x has the same
dimension and kind as the right-hand side vector b.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| type(Circulant), | intent(in) | :: | A |
Coefficient matrix. |
||
| real(kind=dp), | intent(in) | :: | b(:) |
Right-hand side vector. |
Solution vector.
This interface overloads the svd interface from stdlib_linalg to
compute the the singular value decomposition of a Circulant matrix
.
call svd(A, s, u, vt)
A : Matrix of Circulant type.
It is an intent(in) argument.
s : Rank-1 array real array returning the singular values of
A. It is an intent(out) argument.
u (optional) : Rank-2 array of the same kind as A returning
the left singular vectors of A as columns.
Its size should be [n, n].
It is an intent(out) argument.
vt (optional) : Rank-2 array of the same kind as A returning
the right singular vectors of A as rows. Its
size should be [n, n].
It is an intent(out) argument.
Note
Singular values and singular vectors of a Circulant matrix can be
efficiently computed based on the Fast Fourier transform.
Compute the singular value decomposition of a Circulant matrix.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| type(Circulant), | intent(in) | :: | A |
Input matrix. |
||
| real(kind=dp), | intent(out) | :: | s(:) |
Singular values in descending order. |
||
| real(kind=dp), | intent(out), | optional | :: | u(:,:) |
Left singular vectors as columns. |
|
| real(kind=dp), | intent(out), | optional | :: | vt(:,:) |
Right singular vectors as rows. |
This interface overloads the svdvals interface from stdlib_linalg
to compute the singular values of a Circulant matrix .
s = svdvals(A)
A : Matrix of Circulant type.
It is an intent(in) argument.
s : Vector of singular values sorted in decreasing order.
Compute the singular values of a Circulant matrix.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| type(Circulant), | intent(in) | :: | A |
Input matrix. |
Singular values in descending order.
This interface overloads the Fortran intrinsic procedure to define
the transpose operation of a Circulant matrix.
B = transpose(A)
A : Matrix of Circulant.
It is an intent(in) argument.
B : Resulting transposed matrix. It is of the same type as A.
Base type to define a Circulant matrix of size [n x n] with elements
given by the vector
This interface provides methods to construct Circulant matrices.
Given a vector ,
the associated Circulant matrix is the following [n x n] matrix
| private pure, module function construct (c) | Construct a |