Apply a matrix (tmat or rmat) to src and add to dst. Mathematically this is dst= dst + amatsrc Where dst and src are vectors, mat is a matrix, and a is a scalar If the optional variable "zero" is provided and is true, then we compute dst= amatsrc
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(pf_encap_t), | intent(inout) | :: | dst(:) | destination vector |
||
real(kind=pfdp), | intent(in) | :: | a | scalar |
||
real(kind=pfdp), | intent(in) | :: | mat(:,:) | matrix |
||
class(pf_encap_t), | intent(in) | :: | src(:) | src vector |
||
logical, | intent(in), | optional | :: | zero | If false, don't zero out the the dst variable before computing |
|
integer, | intent(in), | optional | :: | flags | Used for choosing which variable to operate on Local variables |
subroutine pf_apply_mat(dst, a, mat, src, zero, flags)
!! Apply a matrix (tmat or rmat) to src and add to dst.
!! Mathematically this is
!! dst= dst + a*mat*src
!! Where dst and src are vectors, mat is a matrix, and a is a scalar
!! If the optional variable "zero" is provided and is true, then we compute
!! dst= a*mat*src
class(pf_encap_t), intent(inout) :: dst(:) !! destination vector
real(pfdp), intent(in) :: a !! scalar
real(pfdp), intent(in) :: mat(:, :) !! matrix
class(pf_encap_t), intent(in) :: src(:) !! src vector
logical, intent(in), optional :: zero !! If false, don't zero out the the dst variable before computing
integer, intent(in), optional :: flags !! Used for choosing which variable to operate on
!! Local variables
logical :: lzero !! local version of input parameter zero
integer :: which !! local version of flags
integer :: n, m !! size of mat
integer :: i, j !! loop variables
lzero = .true.; if (present(zero)) lzero = zero
which = 1; if(present(flags)) which = flags
n = size(mat, dim=1)
m = size(mat, dim=2)
do i = 1, n
if (lzero) call dst(i)%setval(0.0_pfdp, flags)
do j = 1, m
if (a*mat(i, j) /= 0.0_pfdp) call dst(i)%axpy(a * mat(i, j), src(j), flags)
end do
end do
end subroutine pf_apply_mat