mathlib Module


Uses

  • module~~mathlib~~UsesGraph module~mathlib mathlib module~oqp_linalg oqp_linalg module~mathlib->module~oqp_linalg module~precision precision module~mathlib->module~precision module~blas_wrap blas_wrap module~oqp_linalg->module~blas_wrap module~lapack_wrap lapack_wrap module~oqp_linalg->module~lapack_wrap iso_fortran_env iso_fortran_env module~precision->iso_fortran_env module~blas_wrap->module~precision module~mathlib_types mathlib_types module~blas_wrap->module~mathlib_types module~messages messages module~blas_wrap->module~messages module~lapack_wrap->module~precision module~lapack_wrap->module~mathlib_types module~lapack_wrap->module~messages module~messages->module~precision comm_IOFILE comm_IOFILE module~messages->comm_IOFILE comm_PAR comm_PAR module~messages->comm_PAR module~io_constants io_constants module~messages->module~io_constants

Used by

  • module~~mathlib~~UsedByGraph module~mathlib mathlib module~grd1 grd1 module~grd1->module~mathlib proc~build_pfon_density build_pfon_density proc~build_pfon_density->module~mathlib proc~dftder dftder proc~dftder->module~mathlib proc~eijden eijden proc~eijden->module~mathlib proc~electric_moments electric_moments proc~electric_moments->module~mathlib proc~form_rohf_fock form_rohf_fock proc~form_rohf_fock->module~mathlib proc~get_ab_initio_density get_ab_initio_density proc~get_ab_initio_density->module~mathlib proc~get_ab_initio_orbital get_ab_initio_orbital proc~get_ab_initio_orbital->module~mathlib proc~get_spin_square get_spin_square proc~get_spin_square->module~mathlib proc~get_td_transition_dipole get_td_transition_dipole proc~get_td_transition_dipole->module~mathlib proc~get_transition_dipole get_transition_dipole proc~get_transition_dipole->module~mathlib proc~guess_hcore guess_hcore proc~guess_hcore->module~mathlib proc~guess_huckel guess_huckel proc~guess_huckel->module~mathlib proc~huckel_guess huckel_guess proc~huckel_guess->module~mathlib proc~int2_td_data_t_parallel_stop int2_td_data_t%int2_td_data_t_parallel_stop proc~int2_td_data_t_parallel_stop->module~mathlib proc~mo_to_ao mo_to_ao proc~mo_to_ao->module~mathlib proc~oqp_resp_charges oqp_resp_charges proc~oqp_resp_charges->module~mathlib proc~oqp_tdhf_z_vector oqp_tdhf_z_vector proc~oqp_tdhf_z_vector->module~mathlib proc~run_population_analysis run_population_analysis proc~run_population_analysis->module~mathlib proc~scf_driver scf_driver proc~scf_driver->module~mathlib proc~sfdmat sfdmat proc~sfdmat->module~mathlib proc~symmetrize_matrices symmetrize_matrices proc~symmetrize_matrices->module~mathlib proc~tddft_fxc tddft_fxc proc~tddft_fxc->module~mathlib proc~tddft_gxc tddft_gxc proc~tddft_gxc->module~mathlib proc~tdhf_1e_grad tdhf_1e_grad proc~tdhf_1e_grad->module~mathlib proc~tdhf_1e_grad->module~grd1 proc~tdhf_energy tdhf_energy proc~tdhf_energy->module~mathlib proc~tdhf_gradient tdhf_gradient proc~tdhf_gradient->module~mathlib proc~tdhf_gradient->module~grd1 proc~tdhf_mrsf_energy tdhf_mrsf_energy proc~tdhf_mrsf_energy->module~mathlib proc~tdhf_sf_energy tdhf_sf_energy proc~tdhf_sf_energy->module~mathlib proc~tdhf_unrelaxed_density tdhf_unrelaxed_density proc~tdhf_unrelaxed_density->module~mathlib proc~utddft_fxc utddft_fxc proc~utddft_fxc->module~mathlib proc~hf_gradient hf_gradient proc~hf_gradient->module~grd1

Interfaces

public interface pack_matrix

  • public subroutine PACK_F90(A, AP, UPLO)

    @brief Fortran-90 routine for packing symmetric matrix to 1D array

    @date 6 October 2021 - Initial release - @author Igor S. Gerasimov

    @param[in] A - matrix for packing (N x N) @param[out] AP - packed matrix ( N*(N+1)/2 ) @param[in,optional] UPLO - format of packed matrix, U for upper and L for lower

    Arguments

    Type IntentOptional Attributes Name
    real(kind=dp), intent(in) :: A(:,:)
    real(kind=dp), intent(out) :: AP(:)
    character(len=1), intent(in), optional :: UPLO
  • private subroutine PACK_F77(A, N, AP, UPLO) bind(C, name="MTX_PACK")

    @brief Fortran-77 routine for packing symmetric matrix to 1D array

    @date 6 October 2021 - Initial release - @author Igor S. Gerasimov

    @param[in] A - matrix for packing (N x N) @param[in] N - shape of matrix A @param[out] AP - packed matrix ( N*(N+1)/2 ) @param[in] UPLO - format of packed matrix, U for upper and L for lower

    Arguments

    Type IntentOptional Attributes Name
    real(kind=dp), intent(in) :: A(N,*)
    integer, intent(in) :: N
    real(kind=dp), intent(out) :: AP(*)
    character(len=1), intent(in) :: UPLO

public interface unpack_matrix

  • public subroutine UNPACK_F90(AP, A, UPLO)

    @brief Fortran-90 routine for unpacking 1D array to symmetric matrix

    @details LAPACK returns only upper or lower filling of matrix

    @date 6 October 2021 - Initial release - @author Igor S. Gerasimov

    @param[in] AP - packed matrix (N x N) @param[out] A - matrix for unpacking ( N*(N+1)/2 ) @param[in,optional] UPLO - format of packed matrix, U for upper and L for lower

    Arguments

    Type IntentOptional Attributes Name
    real(kind=dp), intent(in) :: AP(*)
    real(kind=dp), intent(out) :: A(:,:)
    character(len=1), intent(in), optional :: UPLO
  • private subroutine UNPACK_F77(AP, A, N, UPLO) bind(C, name="MTX_UNPACK")

    @brief Fortran-77 routine for unpacking 1D array to symmetric matrix

    @details LAPACK returns only upper or lower filling of matrix

    @date 6 October 2021 - Initial release - @author Igor S. Gerasimov

    @param[in] AP - packed matrix @param[out] A - matrix for unpacking (N x N) @param[in] N - shape of matrix A @param[in] UPLO - format of packed matrix, U for upper and L for lower

    Arguments

    Type IntentOptional Attributes Name
    real(kind=dp), intent(in) :: AP(*)
    real(kind=dp), intent(out) :: A(N,*)
    integer, intent(in) :: N
    character(len=1), intent(in) :: UPLO

Functions

public function traceprod_sym_packed(a, b, n) result(res)

@brief Compute the trace of the product of two symmetric matrices in packed format @detail The trace is actually an inner product of matrices, assuming they are vectors @param[in] a first matrix @param[in] b second matrix @param[in] n dimension of matrices A and B

Arguments

Type IntentOptional Attributes Name
real(kind=dp) :: a(*)
real(kind=dp) :: b(*)
integer :: n

Return Value real(kind=dp)


Subroutines

public subroutine orb_to_dens(d, v, x, m, n, ldv)

@brief Compute density matrix from a set of orbitals and respective occupation numbers @detail Compute the transformation: D = V * X * V^T @param[out] d density matrix @param[in] v matrix of orbitals @param[in] x vector of occupation numbers @param[in] m number of columns in V @param[in] n dimension of D @param[in] ldv leading dimension of V

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(out) :: d(*)
real(kind=dp), intent(in) :: v(ldv,*)
real(kind=dp), intent(in) :: x(*)
integer, intent(in) :: m
integer, intent(in) :: n
integer, intent(in) :: ldv

public subroutine solve_linear_equations(a, b, n, nrhs, lda, ierr)

@brief Compute the solution to a real system of linear equations A * X = B @detai Wrapper for DSYSV from Lapack @param[in,out] A general matrix, destroyed on exit. @param[in,out] B RHS on entry. The solution on exit. @param[in] n size of the problem. @param[in] nrhs number of RHS vectors @param[in] lda leading dimension of matrix A @param[out] ierr Error flag, ierr=0 if no errors. Read DSYEV manual for details

Arguments

Type IntentOptional Attributes Name
real(kind=dp) :: a(*)
real(kind=dp) :: b(*)
integer, intent(in) :: n
integer, intent(in) :: nrhs
integer, intent(in) :: lda
integer, intent(inout) :: ierr

public subroutine symmetrize_matrix(a, n)

@brief Compute A = A + A^T of a square matrix @param[inout] a square NxN matrix @param[in] n matrix dimension

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(inout) :: a(n,*)
integer, intent(in) :: n

public subroutine antisymmetrize_matrix(a, n)

@brief Compute A = A - A^T of a square matrix @param[inout] a square NxN matrix @param[in] n matrix dimension

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(inout) :: a(n,*)
integer, intent(in) :: n

public subroutine triangular_to_full(a, n, uplo)

@brief Fill the upper/lower triangle of the symmetric matrix @param[inout] a square NxN matrix in triangular form @param[in] n matrix dimension @param[in] uplo U if A is upper triangular, L if lower triangular

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(inout) :: a(n,*)
integer, intent(in) :: n
character(len=1), intent(in) :: uplo

public subroutine orthogonal_transform_sym(n, m, a, u, ldu, b)

@brief Compute orthogonal transformation of a symmetric marix A in packed format: B = U^T * A * U @param[in] a Matrix to transform @param[in] u Orthogonal matrix U(ldu,m) @param[in] n dimension of matrix A @param[in] m dimension of matrix B @param[in] ldu leading dimension of matrix U @param[out] b Result @param[inout] wrk Scratch space @author Vladimir Mironov

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n
integer, intent(in) :: m
real(kind=8), intent(in) :: a(*)
real(kind=8), intent(in) :: u(n,*)
integer, intent(in) :: ldu
real(kind=8), intent(out) :: b(*)

public subroutine orthogonal_transform(trans, ld, u, a, b, wrk)

@brief Compute orthogonal transformation of a square marix @param[in] trans If trans='n' compute B = U^T * A * U If trans='t' compute B = U * A * U^T @param[in] ld Dimension of matrices @param[in] u Square orthogonal matrix @param[inout] a Matrix to transform, optionally output matrix @param[out] b Result, can be absent for in-place transform of matrix A @param[inout] wrk Scratch space, optional @author Vladimir Mironov

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: trans
integer :: ld
real(kind=dp), intent(in) :: u(*)
real(kind=dp), intent(in) :: a(*)
real(kind=dp), intent(out), optional :: b(*)
real(kind=dp), intent(inout), optional, target :: wrk(*)

public subroutine orthogonal_transform2(trans, m, n, u, ldu, a, lda, b, ldb, wrk)

@brief Compute orthogonal transformation of a square marix @param[in] trans If trans='n' compute B = U^T * A * U If trans='t' compute B = U * A * U^T @param[in] ld Dimension of matrices @param[in] u Square orthogonal matrix @param[in] a Matrix to transform @param[out] b Result @param[inout] wrk Scratch space @author Vladimir Mironov

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: trans
integer :: m
integer :: n
real(kind=dp), intent(in) :: u(*)
integer :: ldu
real(kind=dp), intent(in) :: a(*)
integer :: lda
real(kind=dp), intent(out) :: b(*)
integer :: ldb
real(kind=dp), intent(inout) :: wrk(*)

public subroutine matrix_invsqrt(s, q, nbf, qrnk, tol)

@brief Compute matrix inverse square root using SVD and removing linear dependency @detail This subroutine is used to obtain set of canonical orbitals by diagonalization of the basis set overlap matrix Q = S^{-1/2}, Q^T * S * Q = I @param[in] s Overlap matrix, symmetric packed format @param[out] q Matrix inverse square root, square matrix @param[in] nbf Dimeension of matrices S and Q, basis set size @param[out] qrnk Rank of matrix Q @param[in] tol optional, tolerance to remove linear dependency, default = 1.0e-8

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(in) :: s(*)
real(kind=dp), intent(out) :: q(nbf,*)
integer, intent(in) :: nbf
integer, intent(out), optional :: qrnk
real(kind=dp), intent(in), optional :: tol

public subroutine PACK_F90(A, AP, UPLO)

@brief Fortran-90 routine for packing symmetric matrix to 1D array

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(in) :: A(:,:)
real(kind=dp), intent(out) :: AP(:)
character(len=1), intent(in), optional :: UPLO

public subroutine UNPACK_F90(AP, A, UPLO)

@brief Fortran-90 routine for unpacking 1D array to symmetric matrix

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(in) :: AP(*)
real(kind=dp), intent(out) :: A(:,:)
character(len=1), intent(in), optional :: UPLO