sfresvec Subroutine

public subroutine sfresvec(q, a, b, vec, eigv, nvec, rnorm, ndsr)

Uses

  • proc~~sfresvec~~UsesGraph proc~sfresvec sfresvec module~precision precision proc~sfresvec->module~precision iso_fortran_env iso_fortran_env module~precision->iso_fortran_env

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(out), dimension(:,:) :: q
real(kind=dp), intent(in), dimension(:,:) :: a
real(kind=dp), intent(in), dimension(:,:) :: b
real(kind=dp), intent(inout), dimension(:,:) :: vec
real(kind=dp), intent(in), dimension(:) :: eigv
integer, intent(in) :: nvec
real(kind=dp), intent(out), dimension(:) :: rnorm
integer, intent(in) :: ndsr

Calls

proc~~sfresvec~~CallsGraph proc~sfresvec sfresvec proc~oqp_dgemm_i64 oqp_dgemm_i64 proc~sfresvec->proc~oqp_dgemm_i64 dgemm dgemm proc~oqp_dgemm_i64->dgemm interface~show_message show_message proc~oqp_dgemm_i64->interface~show_message

Called by

proc~~sfresvec~~CalledByGraph proc~sfresvec sfresvec proc~tdhf_mrsf_energy tdhf_mrsf_energy proc~tdhf_mrsf_energy->proc~sfresvec proc~tdhf_sf_energy tdhf_sf_energy proc~tdhf_sf_energy->proc~sfresvec proc~tdhf_mrsf_energy_c tdhf_mrsf_energy_C proc~tdhf_mrsf_energy_c->proc~tdhf_mrsf_energy proc~tdhf_sf_energy_c tdhf_sf_energy_C proc~tdhf_sf_energy_c->proc~tdhf_sf_energy

Source Code

  subroutine sfresvec(q,a,b,vec,eigv,nvec,rnorm,ndsr)
    use precision, only: dp

    implicit none

    real(kind=dp), intent(out), dimension(:,:) :: q
    real(kind=dp), intent(in), dimension(:,:) :: a, b
    real(kind=dp), intent(inout), dimension(:,:) :: vec
    real(kind=dp), intent(in), dimension(:) :: eigv
    integer, intent(in) :: nvec
    real(kind=dp), intent(out), dimension(:) :: rnorm
    integer, intent(in) :: ndsr

    integer :: ist, xvec_dim

    xvec_dim = ubound(q, 1)

    call dgemm('n','n',xvec_dim,ndsr,nvec, &
               1.0_dp,b,xvec_dim, &
                      vec,nvec, &
               0.0_dp,q,xvec_dim)

    do ist = 1, ndsr
      vec(:,ist) = -vec(:,ist)*eigv(ist)
    end do

    call dgemm('n','n',xvec_dim,ndsr,nvec, &
               1.0_dp,a,xvec_dim, &
                      vec,nvec, &
               1.0_dp,q,xvec_dim)

    do ist = 1, ndsr
      rnorm(ist) = dot_product(q(:,ist),q(:,ist))
    end do

  end subroutine sfresvec