![]() |
(git:d18deda)
|
Basic linear algebra operations for complex full matrices. More...
Data Types | |
interface | cp_cfm_scale |
Functions/Subroutines | |
subroutine, public | cp_cfm_det (matrix_a, det_a) |
Computes the determinant (with a correct sign even in parallel environment!) of a complex square matrix. | |
subroutine, public | cp_cfm_schur_product (matrix_a, matrix_b, matrix_c) |
Computes the element-wise (Schur) product of two matrices: C = A \circ B . | |
subroutine, public | cp_cfm_scale_and_add (alpha, matrix_a, beta, matrix_b) |
Scale and add two BLACS matrices (a = alpha*a + beta*b). | |
subroutine, public | cp_cfm_scale_and_add_fm (alpha, matrix_a, beta, matrix_b) |
Scale and add two BLACS matrices (a = alpha*a + beta*b). where b is a real matrix (adapted from cp_cfm_scale_and_add). | |
subroutine, public | cp_cfm_lu_decompose (matrix_a, determinant) |
Computes LU decomposition of a given matrix. | |
subroutine, public | cp_cfm_gemm (transa, transb, m, n, k, alpha, matrix_a, matrix_b, beta, matrix_c, a_first_col, a_first_row, b_first_col, b_first_row, c_first_col, c_first_row) |
Performs one of the matrix-matrix operations: matrix_c = alpha * op1( matrix_a ) * op2( matrix_b ) + beta*matrix_c. | |
subroutine, public | cp_cfm_column_scale (matrix_a, scaling) |
Scales columns of the full matrix by corresponding factors. | |
subroutine, public | cp_cfm_solve (matrix_a, general_a, determinant) |
Solve the system of linear equations A*b=A_general using LU decomposition. Pay attention that both matrices are overwritten on exit and that the result is stored into the matrix 'general_a'. | |
subroutine, public | cp_cfm_lu_invert (matrix, info_out) |
Inverts a matrix using LU decomposition. The input matrix will be overwritten. | |
subroutine, public | cp_cfm_trace (matrix_a, matrix_b, trace) |
Returns the trace of matrix_a^T matrix_b, i.e sum_{i,j}(matrix_a(i,j)*matrix_b(i,j)) . | |
subroutine, public | cp_cfm_triangular_multiply (triangular_matrix, matrix_b, side, transa_tr, invert_tr, uplo_tr, unit_diag_tr, n_rows, n_cols, alpha) |
Multiplies in place by a triangular matrix: matrix_b = alpha op(triangular_matrix) matrix_b or (if side='R') matrix_b = alpha matrix_b op(triangular_matrix) op(triangular_matrix) is: triangular_matrix (if transa="N" and invert_tr=.false.) triangular_matrix^T (if transa="T" and invert_tr=.false.) triangular_matrix^H (if transa="C" and invert_tr=.false.) triangular_matrix^(-1) (if transa="N" and invert_tr=.true.) triangular_matrix^(-T) (if transa="T" and invert_tr=.true.) triangular_matrix^(-H) (if transa="C" and invert_tr=.true.) | |
subroutine, public | cp_cfm_triangular_invert (matrix_a, uplo, info_out) |
Inverts a triangular matrix. | |
subroutine, public | cp_cfm_transpose (matrix, trans, matrixt) |
Transposes a BLACS distributed complex matrix. | |
real(kind=dp) function, public | cp_cfm_norm (matrix, mode) |
Norm of matrix using (p)zlange. | |
subroutine, public | cp_cfm_rot_rows (matrix, irow, jrow, cs, sn) |
Applies a planar rotation defined by cs and sn to the i'th and j'th rows. | |
subroutine, public | cp_cfm_rot_cols (matrix, icol, jcol, cs, sn) |
Applies a planar rotation defined by cs and sn to the i'th and j'th columnns. | |
subroutine, public | cp_cfm_uplo_to_full (matrix, workspace, uplo) |
... | |
Basic linear algebra operations for complex full matrices.
subroutine, public cp_cfm_basic_linalg::cp_cfm_det | ( | type(cp_cfm_type), intent(in) | matrix_a, |
complex(kind=dp), intent(out) | det_a | ||
) |
Computes the determinant (with a correct sign even in parallel environment!) of a complex square matrix.
matrix_a | ... |
det_a | ... |
Definition at line 74 of file cp_cfm_basic_linalg.F.
subroutine, public cp_cfm_basic_linalg::cp_cfm_schur_product | ( | type(cp_cfm_type), intent(in) | matrix_a, |
type(cp_cfm_type), intent(in) | matrix_b, | ||
type(cp_cfm_type), intent(in) | matrix_c | ||
) |
Computes the element-wise (Schur) product of two matrices: C = A \circ B .
matrix_a | the first input matrix |
matrix_b | the second input matrix |
matrix_c | matrix to store the result |
Definition at line 152 of file cp_cfm_basic_linalg.F.
subroutine, public cp_cfm_basic_linalg::cp_cfm_scale_and_add | ( | complex(kind=dp), intent(in) | alpha, |
type(cp_cfm_type), intent(in) | matrix_a, | ||
complex(kind=dp), intent(in), optional | beta, | ||
type(cp_cfm_type), intent(in), optional | matrix_b | ||
) |
Scale and add two BLACS matrices (a = alpha*a + beta*b).
alpha | ... |
matrix_a | ... |
beta | ... |
matrix_b | ... |
Definition at line 239 of file cp_cfm_basic_linalg.F.
subroutine, public cp_cfm_basic_linalg::cp_cfm_scale_and_add_fm | ( | complex(kind=dp), intent(in) | alpha, |
type(cp_cfm_type), intent(in) | matrix_a, | ||
complex(kind=dp), intent(in) | beta, | ||
type(cp_fm_type), intent(in) | matrix_b | ||
) |
Scale and add two BLACS matrices (a = alpha*a + beta*b). where b is a real matrix (adapted from cp_cfm_scale_and_add).
alpha | ... |
matrix_a | ... |
beta | ... |
matrix_b | ... |
Definition at line 350 of file cp_cfm_basic_linalg.F.
subroutine, public cp_cfm_basic_linalg::cp_cfm_lu_decompose | ( | type(cp_cfm_type), intent(in) | matrix_a, |
complex(kind=dp), intent(out) | determinant | ||
) |
Computes LU decomposition of a given matrix.
matrix_a | full matrix |
determinant | determinant |
Definition at line 457 of file cp_cfm_basic_linalg.F.
subroutine, public cp_cfm_basic_linalg::cp_cfm_gemm | ( | character(len=1), intent(in) | transa, |
character(len=1), intent(in) | transb, | ||
integer, intent(in) | m, | ||
integer, intent(in) | n, | ||
integer, intent(in) | k, | ||
complex(kind=dp), intent(in) | alpha, | ||
type(cp_cfm_type), intent(in) | matrix_a, | ||
type(cp_cfm_type), intent(in) | matrix_b, | ||
complex(kind=dp), intent(in) | beta, | ||
type(cp_cfm_type), intent(in) | matrix_c, | ||
integer, intent(in), optional | a_first_col, | ||
integer, intent(in), optional | a_first_row, | ||
integer, intent(in), optional | b_first_col, | ||
integer, intent(in), optional | b_first_row, | ||
integer, intent(in), optional | c_first_col, | ||
integer, intent(in), optional | c_first_row | ||
) |
Performs one of the matrix-matrix operations: matrix_c = alpha * op1( matrix_a ) * op2( matrix_b ) + beta*matrix_c.
transa | form of op1( matrix_a ): op1( matrix_a ) = matrix_a, when transa == 'N' , op1( matrix_a ) = matrix_a^T, when transa == 'T' , op1( matrix_a ) = matrix_a^H, when transa == 'C' , |
transb | form of op2( matrix_b ) |
m | number of rows of the matrix op1( matrix_a ) |
n | number of columns of the matrix op2( matrix_b ) |
k | number of columns of the matrix op1( matrix_a ) as well as number of rows of the matrix op2( matrix_b ) |
alpha | scale factor |
matrix_a | matrix A |
matrix_b | matrix B |
beta | scale factor |
matrix_c | matrix C |
a_first_col | (optional) the first column of the matrix_a to multiply |
a_first_row | (optional) the first row of the matrix_a to multiply |
b_first_col | (optional) the first column of the matrix_b to multiply |
b_first_row | (optional) the first row of the matrix_b to multiply |
c_first_col | (optional) the first column of the matrix_c |
c_first_row | (optional) the first row of the matrix_c |
Definition at line 560 of file cp_cfm_basic_linalg.F.
subroutine, public cp_cfm_basic_linalg::cp_cfm_column_scale | ( | type(cp_cfm_type), intent(in) | matrix_a, |
complex(kind=dp), dimension(:), intent(in) | scaling | ||
) |
Scales columns of the full matrix by corresponding factors.
matrix_a | matrix to scale |
scaling | scale factors for every column. The actual number of scaled columns is limited by the number of scale factors given or by the actual number of columns whichever is smaller. |
Definition at line 632 of file cp_cfm_basic_linalg.F.
subroutine, public cp_cfm_basic_linalg::cp_cfm_solve | ( | type(cp_cfm_type), intent(in) | matrix_a, |
type(cp_cfm_type), intent(in) | general_a, | ||
complex(kind=dp), optional | determinant | ||
) |
Solve the system of linear equations A*b=A_general using LU decomposition. Pay attention that both matrices are overwritten on exit and that the result is stored into the matrix 'general_a'.
matrix_a | matrix A (overwritten on exit) |
general_a | (input) matrix A_general, (output) matrix B |
determinant | (optional) determinant |
Definition at line 731 of file cp_cfm_basic_linalg.F.
subroutine, public cp_cfm_basic_linalg::cp_cfm_lu_invert | ( | type(cp_cfm_type), intent(in) | matrix, |
integer, intent(out), optional | info_out | ||
) |
Inverts a matrix using LU decomposition. The input matrix will be overwritten.
matrix | input a general square non-singular matrix, outputs its inverse |
info_out | optional, if present outputs the info from (p)zgetri |
Definition at line 823 of file cp_cfm_basic_linalg.F.
subroutine, public cp_cfm_basic_linalg::cp_cfm_trace | ( | type(cp_cfm_type), intent(in) | matrix_a, |
type(cp_cfm_type), intent(in) | matrix_b, | ||
complex(kind=dp), intent(out) | trace | ||
) |
Returns the trace of matrix_a^T matrix_b, i.e sum_{i,j}(matrix_a(i,j)*matrix_b(i,j)) .
matrix_a | a complex matrix |
matrix_b | another complex matrix |
trace | value of the trace operator |
Definition at line 908 of file cp_cfm_basic_linalg.F.
subroutine, public cp_cfm_basic_linalg::cp_cfm_triangular_multiply | ( | type(cp_cfm_type), intent(in) | triangular_matrix, |
type(cp_cfm_type), intent(in) | matrix_b, | ||
character, intent(in), optional | side, | ||
character, intent(in), optional | transa_tr, | ||
logical, intent(in), optional | invert_tr, | ||
character, intent(in), optional | uplo_tr, | ||
logical, intent(in), optional | unit_diag_tr, | ||
integer, intent(in), optional | n_rows, | ||
integer, intent(in), optional | n_cols, | ||
complex(kind=dp), intent(in), optional | alpha | ||
) |
Multiplies in place by a triangular matrix: matrix_b = alpha op(triangular_matrix) matrix_b or (if side='R') matrix_b = alpha matrix_b op(triangular_matrix) op(triangular_matrix) is: triangular_matrix (if transa="N" and invert_tr=.false.) triangular_matrix^T (if transa="T" and invert_tr=.false.) triangular_matrix^H (if transa="C" and invert_tr=.false.) triangular_matrix^(-1) (if transa="N" and invert_tr=.true.) triangular_matrix^(-T) (if transa="T" and invert_tr=.true.) triangular_matrix^(-H) (if transa="C" and invert_tr=.true.)
triangular_matrix | the triangular matrix that multiplies the other |
matrix_b | the matrix that gets multiplied and stores the result |
side | on which side of matrix_b stays op(triangular_matrix) (defaults to 'L') |
transa_tr | ... |
invert_tr | if the triangular matrix should be inverted (defaults to false) |
uplo_tr | if triangular_matrix is stored in the upper ('U') or lower ('L') triangle (defaults to 'U') |
unit_diag_tr | if the diagonal elements of triangular_matrix should be assumed to be 1 (defaults to false) |
n_rows | the number of rows of the result (defaults to size(matrix_b,1)) |
n_cols | the number of columns of the result (defaults to size(matrix_b,2)) |
alpha | ... |
Definition at line 976 of file cp_cfm_basic_linalg.F.
subroutine, public cp_cfm_basic_linalg::cp_cfm_triangular_invert | ( | type(cp_cfm_type), intent(in) | matrix_a, |
character, intent(in), optional | uplo, | ||
integer, intent(out), optional | info_out | ||
) |
Inverts a triangular matrix.
matrix_a | ... |
uplo | ... |
info_out | ... |
Definition at line 1060 of file cp_cfm_basic_linalg.F.
subroutine, public cp_cfm_basic_linalg::cp_cfm_transpose | ( | type(cp_cfm_type), intent(in) | matrix, |
character, intent(in) | trans, | ||
type(cp_cfm_type), intent(in) | matrixt | ||
) |
Transposes a BLACS distributed complex matrix.
matrix | input matrix |
trans | 'T' for transpose, 'C' for Hermitian conjugate |
matrixt | output matrix |
Definition at line 1110 of file cp_cfm_basic_linalg.F.
real(kind=dp) function, public cp_cfm_basic_linalg::cp_cfm_norm | ( | type(cp_cfm_type), intent(in) | matrix, |
character, intent(in) | mode | ||
) |
Norm of matrix using (p)zlange.
matrix | input a general matrix |
mode | 'M' max abs element value, '1' or 'O' one norm, i.e. maximum column sum, 'I' infinity norm, i.e. maximum row sum, 'F' or 'E' Frobenius norm, i.e. sqrt of sum of all squares of elements |
Definition at line 1185 of file cp_cfm_basic_linalg.F.
subroutine, public cp_cfm_basic_linalg::cp_cfm_rot_rows | ( | type(cp_cfm_type), intent(in) | matrix, |
integer, intent(in) | irow, | ||
integer, intent(in) | jrow, | ||
real(dp), intent(in) | cs, | ||
real(dp), intent(in) | sn | ||
) |
Applies a planar rotation defined by cs and sn to the i'th and j'th rows.
matrix | ... |
irow | ... |
jrow | ... |
cs | cosine of the rotation angle |
sn | sinus of the rotation angle |
Definition at line 1256 of file cp_cfm_basic_linalg.F.
subroutine, public cp_cfm_basic_linalg::cp_cfm_rot_cols | ( | type(cp_cfm_type), intent(in) | matrix, |
integer, intent(in) | icol, | ||
integer, intent(in) | jcol, | ||
real(dp), intent(in) | cs, | ||
real(dp), intent(in) | sn | ||
) |
Applies a planar rotation defined by cs and sn to the i'th and j'th columnns.
matrix | ... |
icol | ... |
jcol | ... |
cs | cosine of the rotation angle |
sn | sinus of the rotation angle |
Definition at line 1301 of file cp_cfm_basic_linalg.F.
subroutine, public cp_cfm_basic_linalg::cp_cfm_uplo_to_full | ( | type(cp_cfm_type), intent(in) | matrix, |
type(cp_cfm_type), intent(in), optional | workspace, | ||
character, intent(in), optional | uplo | ||
) |
...
matrix | ... |
workspace | ... |
uplo | triangular format; defaults to 'U' |
Definition at line 1346 of file cp_cfm_basic_linalg.F.