(git:34ef472)
|
basic linear algebra operations for full matrices More...
Functions/Subroutines | |
subroutine, public | cp_fm_det (matrix_a, det_a) |
Computes the determinant (with a correct sign even in parallel environment!) of a real square matrix. More... | |
subroutine, public | cp_fm_scale_and_add (alpha, matrix_a, beta, matrix_b) |
calc A <- alpha*A + beta*B optimized for alpha == 1.0 (just add beta*B) and beta == 0.0 (just scale A) More... | |
subroutine, public | cp_fm_geadd (alpha, trans, matrix_a, beta, matrix_b) |
interface to BLACS geadd: matrix_b = beta*matrix_b + alpha*opt(matrix_a) where opt(matrix_a) can be either: 'N': matrix_a 'T': matrix_a^T 'C': matrix_a^H (Hermitian conjugate) note that this is a level three routine, use cp_fm_scale_and_add if that is sufficient for your needs More... | |
subroutine, public | cp_fm_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) |
computes matrix_c = beta * matrix_c + alpha * ( matrix_a ** transa ) * ( matrix_b ** transb ) More... | |
subroutine, public | cp_fm_symm (side, uplo, m, n, alpha, matrix_a, matrix_b, beta, matrix_c) |
computes matrix_c = beta * matrix_c + alpha * matrix_a * matrix_b computes matrix_c = beta * matrix_c + alpha * matrix_b * matrix_a where matrix_a is symmetric More... | |
real(kind=dp) function, public | cp_fm_frobenius_norm (matrix_a) |
computes the Frobenius norm of matrix_a More... | |
subroutine, public | cp_fm_syrk (uplo, trans, k, alpha, matrix_a, ia, ja, beta, matrix_c) |
performs a rank-k update of a symmetric matrix_c matrix_c = beta * matrix_c + alpha * matrix_a * transpose ( matrix_a ) More... | |
subroutine, public | cp_fm_schur_product (matrix_a, matrix_b, matrix_c) |
computes the schur product of two matrices c_ij = a_ij * b_ij More... | |
subroutine, public | cp_fm_triangular_multiply (triangular_matrix, matrix_b, side, transpose_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 transpose_tr=.false. and invert_tr=.false.) triangular_matrix^T (if transpose_tr=.true. and invert_tr=.false.) triangular_matrix^(-1) (if transpose_tr=.false. and invert_tr=.true.) triangular_matrix^(-T) (if transpose_tr=.true. and invert_tr=.true.) More... | |
subroutine, public | cp_fm_scale (alpha, matrix_a) |
scales a matrix matrix_a = alpha * matrix_b More... | |
subroutine, public | cp_fm_transpose (matrix, matrixt) |
transposes a matrix matrixt = matrix ^ T More... | |
subroutine, public | cp_fm_upper_to_full (matrix, work) |
given an upper triangular matrix computes the corresponding full matrix More... | |
subroutine, public | cp_fm_column_scale (matrixa, scaling) |
scales column i of matrix a with scaling(i) More... | |
subroutine, public | cp_fm_row_scale (matrixa, scaling) |
scales row i of matrix a with scaling(i) More... | |
subroutine, public | cp_fm_invert (matrix_a, matrix_inverse, det_a, eps_svd, eigval) |
Inverts a cp_fm_type matrix, optionally returning the determinant of the input matrix. More... | |
subroutine, public | cp_fm_triangular_invert (matrix_a, uplo_tr) |
inverts a triangular matrix More... | |
subroutine, public | cp_fm_qr_factorization (matrix_a, matrix_r, nrow_fact, ncol_fact, first_row, first_col) |
performs a QR factorization of the input rectangular matrix A or of a submatrix of A the computed upper triangular matrix R is in output in the submatrix sub(A) of size NxN M and M give the dimension of the submatrix that has to be factorized (MxN) with M>N More... | |
subroutine, public | cp_fm_solve (matrix_a, general_a) |
computes the the solution to A*b=A_general using lu decomposition pay attention, both matrices are overwritten, a_general contais the result More... | |
subroutine, public | cp_complex_fm_gemm (transa, transb, m, n, k, alpha, A_re, A_im, B_re, B_im, beta, C_re, C_im, a_first_col, a_first_row, b_first_col, b_first_row, c_first_col, c_first_row) |
Convenience function. Computes the matrix multiplications needed for the multiplication of complex matrices. C = beta * C + alpha * ( A ** transa ) * ( B ** transb ) More... | |
real(kind=dp) function, public | cp_fm_norm (matrix, mode) |
norm of matrix using (p)dlange More... | |
subroutine, public | cp_fm_pdgeqpf (matrix, tau, nrow, ncol, first_row, first_col) |
compute a QR factorization with column pivoting of a M-by-N distributed matrix sub( A ) = A(IA:IA+M-1,JA:JA+N-1) More... | |
subroutine, public | cp_fm_pdorgqr (matrix, tau, nrow, first_row, first_col) |
generates an M-by-N real distributed matrix Q denoting A(IA:IA+M-1,JA:JA+N-1) with orthonormal columns, which is defined as the first N columns of a product of K elementary reflectors of order M More... | |
subroutine, public | cp_fm_rot_rows (matrix, irow, jrow, cs, sn) |
Applies a planar rotation defined by cs and sn to the i'th and j'th rows. More... | |
subroutine, public | cp_fm_rot_cols (matrix, icol, jcol, cs, sn) |
Applies a planar rotation defined by cs and sn to the i'th and j'th columnns. More... | |
subroutine, public | cp_fm_gram_schmidt_orthonorm (matrix_a, B, nrows, ncols, start_row, start_col, do_norm, do_print) |
Orthonormalizes selected rows and columns of a full matrix, matrix_a. More... | |
subroutine, public | cp_fm_potrf (fm_matrix, n) |
Cholesky decomposition. More... | |
subroutine, public | cp_fm_potri (fm_matrix, n) |
Invert trianguar matrix. More... | |
subroutine, public | cp_fm_cholesky_restore (fm_matrix, neig, fm_matrixb, fm_matrixout, op, pos, transa) |
... More... | |
basic linear algebra operations for full matrices
subroutine, public cp_fm_basic_linalg::cp_fm_det | ( | type(cp_fm_type), intent(in) | matrix_a, |
real(kind=dp), intent(out) | det_a | ||
) |
Computes the determinant (with a correct sign even in parallel environment!) of a real square matrix.
Definition at line 96 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_scale_and_add | ( | real(kind=dp), intent(in) | alpha, |
type(cp_fm_type), intent(in) | matrix_a, | ||
real(kind=dp), intent(in), optional | beta, | ||
type(cp_fm_type), intent(in), optional | matrix_b | ||
) |
calc A <- alpha*A + beta*B optimized for alpha == 1.0 (just add beta*B) and beta == 0.0 (just scale A)
alpha | ... |
matrix_a | ... |
beta | ... |
matrix_b | ... |
Definition at line 166 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_geadd | ( | real(kind=dp), intent(in) | alpha, |
character, intent(in) | trans, | ||
type(cp_fm_type), intent(in) | matrix_a, | ||
real(kind=dp), intent(in) | beta, | ||
type(cp_fm_type), intent(in) | matrix_b | ||
) |
interface to BLACS geadd: matrix_b = beta*matrix_b + alpha*opt(matrix_a) where opt(matrix_a) can be either: 'N': matrix_a 'T': matrix_a^T 'C': matrix_a^H (Hermitian conjugate) note that this is a level three routine, use cp_fm_scale_and_add if that is sufficient for your needs
alpha | : complex scalar |
trans | : 'N' normal, 'T' transposed |
matrix_a | : input matrix_a |
beta | : complex scalar |
matrix_b | : input matrix_b, upon out put the updated matrix_b |
Definition at line 273 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_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, | ||
real(kind=dp), intent(in) | alpha, | ||
type(cp_fm_type), intent(in) | matrix_a, | ||
type(cp_fm_type), intent(in) | matrix_b, | ||
real(kind=dp), intent(in) | beta, | ||
type(cp_fm_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 | ||
) |
computes matrix_c = beta * matrix_c + alpha * ( matrix_a ** transa ) * ( matrix_b ** transb )
transa | : 'N' -> normal 'T' -> transpose alpha,beta :: can be 0.0_dp and 1.0_dp |
transb | ... |
m | ... |
n | ... |
k | ... |
alpha | ... |
matrix_a | : m x k matrix ( ! for transa = 'N') |
matrix_b | : k x n matrix ( ! for transb = 'N') |
beta | ... |
matrix_c | : m x n matrix |
a_first_col | ... |
a_first_row | ... |
b_first_col | : the k x n matrix starts at col b_first_col of matrix_b (avoid usage) |
b_first_row | ... |
c_first_col | ... |
c_first_row | ... |
Definition at line 437 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_symm | ( | character(len=1), intent(in) | side, |
character(len=1), intent(in) | uplo, | ||
integer, intent(in) | m, | ||
integer, intent(in) | n, | ||
real(kind=dp), intent(in) | alpha, | ||
type(cp_fm_type), intent(in) | matrix_a, | ||
type(cp_fm_type), intent(in) | matrix_b, | ||
real(kind=dp), intent(in) | beta, | ||
type(cp_fm_type), intent(in) | matrix_c | ||
) |
computes matrix_c = beta * matrix_c + alpha * matrix_a * matrix_b computes matrix_c = beta * matrix_c + alpha * matrix_b * matrix_a where matrix_a is symmetric
side | : 'L' -> matrix_a is on the left 'R' -> matrix_a is on the right alpha,beta :: can be 0.0_dp and 1.0_dp |
uplo | ... |
m | ... |
n | ... |
alpha | ... |
matrix_a | : m x m matrix |
matrix_b | : m x n matrix |
beta | ... |
matrix_c | : m x n matrix |
Definition at line 575 of file cp_fm_basic_linalg.F.
real(kind=dp) function, public cp_fm_basic_linalg::cp_fm_frobenius_norm | ( | type(cp_fm_type), intent(in) | matrix_a | ) |
computes the Frobenius norm of matrix_a
computes the Frobenius norm of matrix_a
matrix_a | : m x n matrix |
Definition at line 628 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_syrk | ( | character(len=1), intent(in) | uplo, |
character(len=1), intent(in) | trans, | ||
integer, intent(in) | k, | ||
real(kind=dp), intent(in) | alpha, | ||
type(cp_fm_type), intent(in) | matrix_a, | ||
integer, intent(in) | ia, | ||
integer, intent(in) | ja, | ||
real(kind=dp), intent(in) | beta, | ||
type(cp_fm_type), intent(in) | matrix_c | ||
) |
performs a rank-k update of a symmetric matrix_c matrix_c = beta * matrix_c + alpha * matrix_a * transpose ( matrix_a )
uplo | : 'U' ('L') |
trans | : 'N' ('T') |
k | : number of cols to use in matrix_a ia,ja :: 1,1 (could be used for selecting subblock of a) |
alpha | ... |
matrix_a | ... |
ia | ... |
ja | ... |
beta | ... |
matrix_c | ... |
Definition at line 674 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_schur_product | ( | type(cp_fm_type), intent(in) | matrix_a, |
type(cp_fm_type), intent(in) | matrix_b, | ||
type(cp_fm_type), intent(in) | matrix_c | ||
) |
computes the schur product of two matrices c_ij = a_ij * b_ij
matrix_a | ... |
matrix_b | ... |
matrix_c | ... |
Definition at line 727 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_triangular_multiply | ( | type(cp_fm_type), intent(in) | triangular_matrix, |
type(cp_fm_type), intent(in) | matrix_b, | ||
character, intent(in), optional | side, | ||
logical, intent(in), optional | transpose_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, | ||
real(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 transpose_tr=.false. and invert_tr=.false.) triangular_matrix^T (if transpose_tr=.true. and invert_tr=.false.) triangular_matrix^(-1) (if transpose_tr=.false. and invert_tr=.true.) triangular_matrix^(-T) (if transpose_tr=.true. 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') |
transpose_tr | if the triangular matrix should be transposed (defaults to false) |
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 1589 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_scale | ( | real(kind=dp), intent(in) | alpha, |
type(cp_fm_type), intent(in) | matrix_a | ||
) |
scales a matrix matrix_a = alpha * matrix_b
alpha | ... |
matrix_a | ... |
Definition at line 1679 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_transpose | ( | type(cp_fm_type), intent(in) | matrix, |
type(cp_fm_type), intent(in) | matrixt | ||
) |
transposes a matrix matrixt = matrix ^ T
matrix | ... |
matrixt | ... |
Definition at line 1709 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_upper_to_full | ( | type(cp_fm_type), intent(in) | matrix, |
type(cp_fm_type), intent(in) | work | ||
) |
given an upper triangular matrix computes the corresponding full matrix
matrix | the upper triangular matrix as input, the full matrix as output |
work | a matrix of the same size as matrix |
Definition at line 1758 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_column_scale | ( | type(cp_fm_type), intent(in) | matrixa, |
real(kind=dp), dimension(:), intent(in) | scaling | ||
) |
scales column i of matrix a with scaling(i)
matrixa | ... |
scaling | : an array used for scaling the columns, SIZE(scaling) determines the number of columns to be scaled |
Definition at line 1881 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_row_scale | ( | type(cp_fm_type), intent(in) | matrixa, |
real(kind=dp), dimension(:), intent(in) | scaling | ||
) |
scales row i of matrix a with scaling(i)
matrixa | ... |
scaling | : an array used for scaling the columns, |
Definition at line 1944 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_invert | ( | type(cp_fm_type), intent(in) | matrix_a, |
type(cp_fm_type), intent(in) | matrix_inverse, | ||
real(kind=dp), intent(out), optional | det_a, | ||
real(kind=dp), intent(in), optional | eps_svd, | ||
real(kind=dp), dimension(:), intent(inout), optional, pointer | eigval | ||
) |
Inverts a cp_fm_type matrix, optionally returning the determinant of the input matrix.
matrix_a | the matrix to invert |
matrix_inverse | the inverse of matrix_a |
det_a | the determinant of matrix_a |
eps_svd | optional parameter to active SVD based inversion, singular values below eps_svd are screened |
eigval | optionally return matrix eigenvalues/singular values |
Definition at line 2014 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_triangular_invert | ( | type(cp_fm_type), intent(in) | matrix_a, |
character, intent(in), optional | uplo_tr | ||
) |
inverts a triangular matrix
matrix_a | ... |
uplo_tr | ... |
Definition at line 2192 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_qr_factorization | ( | type(cp_fm_type), intent(in) | matrix_a, |
type(cp_fm_type), intent(in) | matrix_r, | ||
integer, intent(in), optional | nrow_fact, | ||
integer, intent(in), optional | ncol_fact, | ||
integer, intent(in), optional | first_row, | ||
integer, intent(in), optional | first_col | ||
) |
performs a QR factorization of the input rectangular matrix A or of a submatrix of A the computed upper triangular matrix R is in output in the submatrix sub(A) of size NxN M and M give the dimension of the submatrix that has to be factorized (MxN) with M>N
matrix_a | ... |
matrix_r | ... |
nrow_fact | ... |
ncol_fact | ... |
first_row | ... |
first_col | ... |
Definition at line 2241 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_solve | ( | type(cp_fm_type), intent(in) | matrix_a, |
type(cp_fm_type), intent(in) | general_a | ||
) |
computes the the solution to A*b=A_general using lu decomposition pay attention, both matrices are overwritten, a_general contais the result
matrix_a | ... |
general_a | ... |
Definition at line 2323 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_complex_fm_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, | ||
real(kind=dp), intent(in) | alpha, | ||
type(cp_fm_type), intent(in) | A_re, | ||
type(cp_fm_type), intent(in) | A_im, | ||
type(cp_fm_type), intent(in) | B_re, | ||
type(cp_fm_type), intent(in) | B_im, | ||
real(kind=dp), intent(in) | beta, | ||
type(cp_fm_type), intent(in) | C_re, | ||
type(cp_fm_type), intent(in) | C_im, | ||
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 | ||
) |
Convenience function. Computes the matrix multiplications needed for the multiplication of complex matrices. C = beta * C + alpha * ( A ** transa ) * ( B ** transb )
transa | : 'N' -> normal 'T' -> transpose alpha,beta :: can be 0.0_dp and 1.0_dp |
transb | ... |
m | ... |
n | ... |
k | ... |
alpha | ... |
A_re | m x k matrix ( ! for transa = 'N'), real part |
A_im | m x k matrix ( ! for transa = 'N'), imaginary part |
B_re | k x n matrix ( ! for transa = 'N'), real part |
B_im | k x n matrix ( ! for transa = 'N'), imaginary part |
beta | ... |
C_re | m x n matrix, real part |
C_im | m x n matrix, imaginary part |
a_first_col | ... |
a_first_row | ... |
b_first_col | : the k x n matrix starts at col b_first_col of matrix_b (avoid usage) |
b_first_row | ... |
c_first_col | ... |
c_first_row | ... |
Definition at line 2392 of file cp_fm_basic_linalg.F.
real(kind=dp) function, public cp_fm_basic_linalg::cp_fm_norm | ( | type(cp_fm_type), intent(in) | matrix, |
character, intent(in) | mode | ||
) |
norm of matrix using (p)dlange
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 2579 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_pdgeqpf | ( | type(cp_fm_type), intent(in) | matrix, |
real(kind=dp), dimension(:), pointer | tau, | ||
integer, intent(in) | nrow, | ||
integer, intent(in) | ncol, | ||
integer, intent(in) | first_row, | ||
integer, intent(in) | first_col | ||
) |
compute a QR factorization with column pivoting of a M-by-N distributed matrix sub( A ) = A(IA:IA+M-1,JA:JA+N-1)
matrix | : input M-by-N distributed matrix sub( A ) which is to be factored |
tau | : scalar factors TAU of the elementary reflectors. TAU is tied to the distributed matrix A |
nrow | ... |
ncol | ... |
first_row | ... |
first_col | ... |
Definition at line 2731 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_pdorgqr | ( | type(cp_fm_type), intent(in) | matrix, |
real(kind=dp), dimension(:), pointer | tau, | ||
integer, intent(in) | nrow, | ||
integer, intent(in) | first_row, | ||
integer, intent(in) | first_col | ||
) |
generates an M-by-N real distributed matrix Q denoting A(IA:IA+M-1,JA:JA+N-1) with orthonormal columns, which is defined as the first N columns of a product of K elementary reflectors of order M
matrix | : On entry, the j-th column must contain the vector which defines the elementary reflector as returned from PDGEQRF On exit it contains the M-by-N distributed matrix Q |
tau | : contains the scalar factors TAU of elementary reflectors as returned by PDGEQRF |
nrow | ... |
first_row | ... |
first_col | ... |
Definition at line 2803 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_rot_rows | ( | type(cp_fm_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.
cs | cosine of the rotation angle |
sn | sinus of the rotation angle |
irow | ... |
jrow | ... |
Definition at line 2863 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_rot_cols | ( | type(cp_fm_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.
cs | cosine of the rotation angle |
sn | sinus of the rotation angle |
icol | ... |
jcol | ... |
Definition at line 2905 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_gram_schmidt_orthonorm | ( | type(cp_fm_type), intent(in) | matrix_a, |
real(kind=dp), dimension(:, :), intent(out) | B, | ||
integer, intent(in), optional | nrows, | ||
integer, intent(in), optional | ncols, | ||
integer, intent(in), optional | start_row, | ||
integer, intent(in), optional | start_col, | ||
logical, intent(in), optional | do_norm, | ||
logical, intent(in), optional | do_print | ||
) |
Orthonormalizes selected rows and columns of a full matrix, matrix_a.
matrix_a | ... |
B | ... |
nrows | number of rows of matrix_a, optional, defaults to size(matrix_a,1) |
ncols | number of columns of matrix_a, optional, defaults to size(matrix_a, 2) |
start_row | starting index of rows, optional, defaults to 1 |
start_col | starting index of columns, optional, defaults to 1 |
do_norm | ... |
do_print | ... |
Definition at line 2950 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_potrf | ( | type(cp_fm_type) | fm_matrix, |
integer, intent(in) | n | ||
) |
Cholesky decomposition.
fm_matrix | ... |
n | ... |
Definition at line 3084 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_potri | ( | type(cp_fm_type) | fm_matrix, |
integer, intent(in) | n | ||
) |
Invert trianguar matrix.
fm_matrix | the matrix to invert (must be an upper triangular matrix) |
n | size of the matrix to invert |
Definition at line 3121 of file cp_fm_basic_linalg.F.
subroutine, public cp_fm_basic_linalg::cp_fm_cholesky_restore | ( | type(cp_fm_type) | fm_matrix, |
integer, intent(in) | neig, | ||
type(cp_fm_type) | fm_matrixb, | ||
type(cp_fm_type) | fm_matrixout, | ||
character(len=*), intent(in) | op, | ||
character(len=*), intent(in) | pos, | ||
character(len=*), intent(in) | transa | ||
) |
...
fm_matrix | ... |
neig | ... |
fm_matrixb | ... |
fm_matrixout | ... |
op | ... |
pos | ... |
transa | ... |
Definition at line 3161 of file cp_fm_basic_linalg.F.