linalg 1.8.2
A linear algebra library that provides a user-friendly interface to several BLAS and LAPACK routines.
Loading...
Searching...
No Matches
linalg::diag_mtx_mult Interface Reference

Multiplies a diagonal matrix with another matrix or array. More...

Public Member Functions

 diag_mtx_mult_mtx
 
 diag_mtx_mult_mtx2
 
 diag_mtx_mult_mtx3
 
 diag_mtx_mult_mtx4
 
 diag_mtx_mult_mtx_cmplx
 
 diag_mtx_mult_mtx2_cmplx
 
 diag_mtx_mult_mtx_mix
 
 diag_mtx_mult_mtx2_mix
 
 diag_mtx_sparse_mult
 

Detailed Description

Multiplies a diagonal matrix with another matrix or array.

Syntax 1
Computes the matrix operation: C = alpha * A * op(B) + beta * C, or C = alpha * op(B) * A + beta * C.
subroutine diag_mtx_mult(logical lside, logical trans, real(real64) alpha, real(real64) a(:), real(real64) b(:,:), real(real64) beta, real(real64) c(:,:), optional class(errors) err)
subroutine diag_mtx_mult(logical lside, logical trans, real(real64) alpha, complex(real64) a(:), complex(real64) b(:,:), real(real64) beta, complex(real64) c(:,:), optional class(errors) err)
subroutine diag_mtx_mult(logical lside, logical trans, real(real64) alpha, complex(real64) a(:), real(real64) b(:,:), real(real64) beta, complex(real64) c(:,:), optional class(errors) err)
subroutine diag_mtx_mult(logical lside, logical trans, complex(real64) alpha, complex(real64) a(:), complex(real64) b(:,:), complex(real64) beta, complex(real64) c(:,:), optional class(errors) err)
subroutine diag_mtx_mult(logical lside, logical trans, complex(real64) alpha, real(real64) a(:), complex(real64) b(:,:), complex(real64) beta, complex(real64) c(:,:), optional class(errors) err)
Parameters
[in]lsideSet to true to apply matrix A from the left; else, set to false to apply matrix A from the left.
[in]transSet to true if op(B) = B^T; else, set to false for op(B) = B. In the complex case set to LA_TRANSPOSE if op(B) = B^T, set to LA_HERMITIAN_TRANSPOSE if op(B) = B^H, otherwise set to LA_NO_OPERATION if op(B) = B.
[in]alphaA scalar multiplier.
[in]aA K-element array containing the diagonal elements of A where K = MIN(M,P) if lside is true; else, if lside is false, K = MIN(N,P).
[in]bThe LDB-by-TDB matrix B where (LDB = leading dimension of B, and TDB = trailing dimension of B):
  • lside == true & trans == true: LDB = N, TDB = P
  • lside == true & trans == false: LDB = P, TDB = N
  • lside == false & trans == true: LDB = P, TDB = M
  • lside == false & trans == false: LDB = M, TDB = P
[in]betaA scalar multiplier.
[in,out]cOn input, the M-by-N matrix C. On output, the resulting M-by-N matrix.
[in,out]errAn 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 any of the input array sizes are incorrect.
Syntax 2
Computes the matrix operation: B = alpha * A * B, or B = alpha B * A.
subroutine diag_mtx_mult(logical lside, real(real64) alpha, real(real64) a(:), real(real64) b(:,:), optional class(errors) err)
subroutine diag_mtx_mult(logical lside, complex(real64) alpha, complex(real64) a(:), complex(real64) b(:,:), optional class(errors) err)
subroutine diag_mtx_mult(logical lside, complex(real64) alpha, real(real64) a(:), complex(real64) b(:,:), optional class(errors) err)
Parameters
[in]lsideSet to true to apply matrix A from the left; else, set to false to apply matrix A from the left.
[in]alphaA scalar multiplier.
[in]aA K-element array containing the diagonal elements of A where K = MIN(M,P) if lside is true; else, if lside is false, K = MIN(N,P).
[in]bOn input, the M-by-N matrix B. On output, the resulting M-by-N matrix.
[in,out]errAn 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 any of the input array sizes are incorrect.
Syntax 3
Computes the sparse-matrix operation: B = alpha * A * B or B = alpha * B * A where A is a diagonal matrix and B is a sparse matrix.
subroutine diag_mtx_mult(logical lside, real(real64) alpha, real(real64) a(:), class(csr_matrix) b, optional class(errors) err)
Parameters
[in]lsideSet to true to apply matrix A from the left; else, set to false to apply matrix A from the left.
[in]alphaA scalar multiplier.
[in]aA K-element array containing the diagonal elements of A where K = MIN(M,P) if lside is true; else, if lside is false, K = MIN(N,P).
[in]bOn input, the M-by-N matrix B. On output, the resulting M-by-N matrix.
[in,out]errAn 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 any of the input array sizes are incorrect.
Usage
The following example illustrates the use of the diagonal matrix multiplication routine to compute the \( S V^T \) component of a singular value decomposition.
program example
use iso_fortran_env, only : int32, real64
use linalg
implicit none
! Variables
real(real64) :: a(3,2), s(2), u(3,3), vt(2,2), ac(3,2)
integer(int32) :: i
! Initialize the 3-by-2 matrix A
! | 2 1 |
! A = |-3 1 |
! |-1 1 |
a = reshape([2.0d0, -3.0d0, -1.0d0, 1.0d0, 1.0d0, 1.0d0], [3, 2])
! Compute the singular value decomposition of A. Notice, V**T is returned
! instead of V. Also note, A is overwritten.
call svd(a, s, u, vt)
! Display the results
print '(A)', "U ="
do i = 1, size(u, 1)
print *, u(i,:)
end do
print '(A)', "S ="
print '(F9.5)', (s(i), i = 1, size(a, 2))
print '(A)', "V**T ="
do i = 1, size(vt, 1)
print *, vt(i,:)
end do
! Compute U * S * V**T
call diag_mtx_mult(.true., 1.0d0, s, vt) ! Compute: VT = S * V**T
ac = matmul(u(:,1:2), vt)
print '(A)', "U * S * V**T ="
do i = 1, size(ac, 1)
print *, ac(i,:)
end do
end program
Multiplies a diagonal matrix with another matrix or array.
Definition linalg.f90:553
Performs sparse matrix multiplication C = A * B.
Definition linalg.f90:5382
Computes the singular value decomposition of a matrix A. The SVD is defined as: , where is an M-by-M...
Definition linalg.f90:2182
Provides a set of common linear algebra routines.
Definition linalg.f90:145
The above program produces the following output.
U =
-0.47411577501825380 -0.81850539032073777 -0.32444284226152509
0.82566838523833064 -0.28535874325972488 -0.48666426339228758
0.30575472113569685 -0.49861740208412991 0.81110710565381272
S =
3.78845
1.62716
V**T =
-0.98483334211643059 0.17350299206578967
-0.17350299206578967 -0.98483334211643059
U * S * V**T =
1.9999999999999993 0.99999999999999956
-3.0000000000000000 1.0000000000000000
-1.0000000000000000 0.99999999999999967

Definition at line 553 of file linalg.f90.


The documentation for this interface was generated from the following file: