TMatVec_offdiag Subroutine

private subroutine TMatVec_offdiag(eel, x, y)

Uses

  • proc~~tmatvec_offdiag~~UsesGraph proc~tmatvec_offdiag TMatVec_offdiag module~mod_memory mod_memory proc~tmatvec_offdiag->module~mod_memory module~mod_io mod_io module~mod_memory->module~mod_io module~mod_constants mod_constants module~mod_memory->module~mod_constants iso_c_binding iso_c_binding module~mod_memory->iso_c_binding module~mod_io->module~mod_constants module~mod_constants->iso_c_binding

Perform matrix vector multiplication y = [TMat-diag(TMat)]*x, where TMat is polarization matrix (precomputed and stored in memory) and x and y are column vectors

Arguments

Type IntentOptional Attributes Name
type(ommp_electrostatics_type), intent(in) :: eel

The electostatic data structure

real(kind=rp), intent(in), dimension(3*eel%pol_atoms) :: x

Input vector

real(kind=rp), intent(out), dimension(3*eel%pol_atoms) :: y

Output vector


Calls

proc~~tmatvec_offdiag~~CallsGraph proc~tmatvec_offdiag TMatVec_offdiag dgemm dgemm proc~tmatvec_offdiag->dgemm

Called by

proc~~tmatvec_offdiag~~CalledByGraph proc~tmatvec_offdiag TMatVec_offdiag proc~tmatvec_incore TMatVec_incore proc~tmatvec_incore->proc~tmatvec_offdiag

Contents

Source Code


Source Code

    subroutine TMatVec_offdiag(eel, x, y)
        !! Perform matrix vector multiplication y = [TMat-diag(TMat)]*x,
        !! where TMat is polarization matrix (precomputed and stored in memory)
        !! and x and y are column vectors 
        use mod_memory, only: mallocate, mfree
        
        implicit none
        
        type(ommp_electrostatics_type), intent(in) :: eel
        !! The electostatic data structure 
        real(rp), dimension(3*eel%pol_atoms), intent(in) :: x
        !! Input vector
        real(rp), dimension(3*eel%pol_atoms), intent(out) :: y
        !! Output vector
        
        integer(ip) :: i, n

        n = 3*eel%pol_atoms
       
        ! Compute the matrix vector product
        call dgemm('N', 'N', n, 1, n, 1.0_rp, eel%tmat, n, x, n, 0.0_rp, y, n)
        ! Subtract the product of diagonal 
        !$omp parallel do default(shared) private(i) 
        do i = 1, n
            y(i) = y(i) - eel%tmat(i,i) * x(i)
        end do
    
    end subroutine TMatVec_offdiag