(git:6a2e663)
cp_cfm_basic_linalg Module Reference

Basic linear algebra operations for complex full matrices. More...

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. More...
 
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 . More...
 
subroutine, public cp_cfm_scale_and_add (alpha, matrix_a, beta, matrix_b)
 Scale and add two BLACS matrices (a = alpha*a + beta*b). More...
 
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). More...
 
subroutine, public cp_cfm_lu_decompose (matrix_a, determinant)
 Computes LU decomposition of a given matrix. More...
 
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. More...
 
subroutine, public cp_cfm_column_scale (matrix_a, scaling)
 Scales columns of the full matrix by corresponding factors. More...
 
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'. More...
 
subroutine, public cp_cfm_lu_invert (matrix, info_out)
 Inverts a matrix using LU decomposition. The input matrix will be overwritten. More...
 
subroutine, public cp_cfm_cholesky_decompose (matrix, n, info_out)
 Used to replace a symmetric positive definite matrix M with its Cholesky decomposition U: M = U^T * U, with U upper triangular. More...
 
subroutine, public cp_cfm_cholesky_invert (matrix, n, info_out)
 Used to replace Cholesky decomposition by the inverse. More...
 
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)) . More...
 
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.) More...
 
subroutine, public cp_cfm_triangular_invert (matrix_a, uplo, info_out)
 Inverts a triangular matrix. More...
 
subroutine, public cp_cfm_transpose (matrix, trans, matrixt)
 Transposes a BLACS distributed complex matrix. More...
 
real(kind=dp) function, public cp_cfm_norm (matrix, mode)
 Norm of matrix using (p)zlange. More...
 
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. More...
 
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. More...
 

Detailed Description

Basic linear algebra operations for complex full matrices.

Note
  • not all functionality implemented
History
Nearly literal copy of Fawzi's routines
Author
Joost VandeVondele

Function/Subroutine Documentation

◆ cp_cfm_det()

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.

Parameters
matrix_a...
det_a...
Author
A. Sinyavskiy (andre.nosp@m.y.si.nosp@m.nyavs.nosp@m.kiy@.nosp@m.chem..nosp@m.uzh..nosp@m.ch)

Definition at line 75 of file cp_cfm_basic_linalg.F.

Here is the call graph for this function:

◆ cp_cfm_schur_product()

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 .

Parameters
matrix_athe first input matrix
matrix_bthe second input matrix
matrix_cmatrix to store the result

Definition at line 156 of file cp_cfm_basic_linalg.F.

Here is the caller graph for this function:

◆ cp_cfm_scale_and_add()

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).

Parameters
alpha...
matrix_a...
beta...
matrix_b...
Date
11.06.2001
Author
Matthias Krack
Version
1.0
Note
Use explicit loops to avoid temporary arrays, as a compiler reasonably assumes that arrays matrix_alocal_data and matrix_blocal_data may overlap (they are referenced by pointers). In general case (alpha*a + beta*b) explicit loops appears to be up to two times more efficient than equivalent LAPACK calls (zscale, zaxpy). This is because using LAPACK calls implies two passes through each array, so data need to be retrieved twice if arrays are large enough to not fit into the processor's cache.

Definition at line 243 of file cp_cfm_basic_linalg.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cp_cfm_scale_and_add_fm()

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).

Parameters
alpha...
matrix_a...
beta...
matrix_b...
Date
01.08.2014
Author
JGH
Version
1.0

Definition at line 354 of file cp_cfm_basic_linalg.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cp_cfm_lu_decompose()

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.

Parameters
matrix_afull matrix
determinantdeterminant
Date
11.06.2001
Author
Matthias Krack
Version
1.0
Note
The actual purpose right now is to efficiently compute the determinant of a given matrix. The original content of the matrix is destroyed.

Definition at line 461 of file cp_cfm_basic_linalg.F.

Here is the call graph for this function:

◆ cp_cfm_gemm()

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.

Parameters
transaform 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' ,
transbform of op2( matrix_b )
mnumber of rows of the matrix op1( matrix_a )
nnumber of columns of the matrix op2( matrix_b )
knumber of columns of the matrix op1( matrix_a ) as well as number of rows of the matrix op2( matrix_b )
alphascale factor
matrix_amatrix A
matrix_bmatrix B
betascale factor
matrix_cmatrix 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
Date
07.06.2001
Author
Matthias Krack
Version
1.0

Definition at line 564 of file cp_cfm_basic_linalg.F.

Here is the caller graph for this function:

◆ cp_cfm_column_scale()

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.

Parameters
matrix_amatrix to scale
scalingscale 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.
Author
Joost VandeVondele

Definition at line 648 of file cp_cfm_basic_linalg.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cp_cfm_solve()

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'.

Parameters
matrix_amatrix A (overwritten on exit)
general_a(input) matrix A_general, (output) matrix B
determinant(optional) determinant
Author
Florian Schiffmann

Definition at line 747 of file cp_cfm_basic_linalg.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cp_cfm_lu_invert()

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.

Parameters
matrixinput a general square non-singular matrix, outputs its inverse
info_outoptional, if present outputs the info from (p)zgetri
Author
Lianheng Tong

Definition at line 839 of file cp_cfm_basic_linalg.F.

Here is the caller graph for this function:

◆ cp_cfm_cholesky_decompose()

subroutine, public cp_cfm_basic_linalg::cp_cfm_cholesky_decompose ( type(cp_cfm_type), intent(in)  matrix,
integer, intent(in), optional  n,
integer, intent(out), optional  info_out 
)

Used to replace a symmetric positive definite matrix M with its Cholesky decomposition U: M = U^T * U, with U upper triangular.

Parameters
matrixthe matrix to replace with its Cholesky decomposition
nthe number of row (and columns) of the matrix & (defaults to the min(size(matrix)))
info_outif present, outputs info from (p)zpotrf
History
05.2002 created [JVdV] 12.2002 updated, added n optional parm [fawzi] 09.2021 removed CPASSERT(info == 0) since there is already check of info [Jan Wilhelm]
Author
Joost

Definition at line 925 of file cp_cfm_basic_linalg.F.

Here is the caller graph for this function:

◆ cp_cfm_cholesky_invert()

subroutine, public cp_cfm_basic_linalg::cp_cfm_cholesky_invert ( type(cp_cfm_type), intent(in)  matrix,
integer, intent(in), optional  n,
integer, intent(out), optional  info_out 
)

Used to replace Cholesky decomposition by the inverse.

Parameters
matrix: the matrix to invert (must be an upper triangular matrix), and is the output of Cholesky decomposition
n: size of the matrix to invert (defaults to the min(size(matrix)))
info_out: if present, outputs info of (p)zpotri
History
05.2002 created Lianheng Tong, based on cp_fm_cholesky_invert
Author
Lianheng Tong

Definition at line 981 of file cp_cfm_basic_linalg.F.

Here is the caller graph for this function:

◆ cp_cfm_trace()

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)) .

Parameters
matrix_aa complex matrix
matrix_banother complex matrix
tracevalue of the trace operator
History
  • 09.2017 created [Sergey Chulkov]
Author
Sergey Chulkov
Note
Based on the subroutine cp_fm_trace(). Note the transposition of matrix_a!

Definition at line 1036 of file cp_cfm_basic_linalg.F.

Here is the caller graph for this function:

◆ cp_cfm_triangular_multiply()

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.)

Parameters
triangular_matrixthe triangular matrix that multiplies the other
matrix_bthe matrix that gets multiplied and stores the result
sideon which side of matrix_b stays op(triangular_matrix) (defaults to 'L')
transa_tr...
invert_trif the triangular matrix should be inverted (defaults to false)
uplo_trif triangular_matrix is stored in the upper ('U') or lower ('L') triangle (defaults to 'U')
unit_diag_trif the diagonal elements of triangular_matrix should be assumed to be 1 (defaults to false)
n_rowsthe number of rows of the result (defaults to size(matrix_b,1))
n_colsthe number of columns of the result (defaults to size(matrix_b,2))
alpha...
History
08.2002 created [fawzi]
Author
Fawzi Mohamed
Note
needs an mpi env

Definition at line 1104 of file cp_cfm_basic_linalg.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cp_cfm_triangular_invert()

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.

Parameters
matrix_a...
uplo...
info_out...
Author
MI

Definition at line 1188 of file cp_cfm_basic_linalg.F.

Here is the caller graph for this function:

◆ cp_cfm_transpose()

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.

Parameters
matrixinput matrix
trans'T' for transpose, 'C' for Hermitian conjugate
matrixtoutput matrix
Author
Lianheng Tong

Definition at line 1238 of file cp_cfm_basic_linalg.F.

Here is the caller graph for this function:

◆ cp_cfm_norm()

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.

Parameters
matrixinput 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
Returns
the norm according to mode
Author
Lianheng Tong

Definition at line 1311 of file cp_cfm_basic_linalg.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cp_cfm_rot_rows()

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.

Parameters
matrix...
irow...
jrow...
cscosine of the rotation angle
snsinus of the rotation angle
Author
Ole Schuett

Definition at line 1382 of file cp_cfm_basic_linalg.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cp_cfm_rot_cols()

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.

Parameters
matrix...
icol...
jcol...
cscosine of the rotation angle
snsinus of the rotation angle
Author
Ole Schuett

Definition at line 1427 of file cp_cfm_basic_linalg.F.

Here is the call graph for this function:
Here is the caller graph for this function: