pf_apply_mat Subroutine

public 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 + 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

Arguments

Type IntentOptional AttributesName
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


Called by

proc~~pf_apply_mat~~CalledByGraph proc~pf_apply_mat pf_apply_mat proc~restrict_sdc restrict_sdc proc~restrict_sdc->proc~pf_apply_mat proc~interpolate_time_space interpolate_time_space proc~interpolate_time_space->proc~pf_apply_mat proc~pf_predictor_oc pf_predictor_oc proc~pf_predictor_oc->proc~interpolate_time_space proc~restrict_time_space_fas restrict_time_space_fas proc~pf_predictor_oc->proc~restrict_time_space_fas proc~pf_v_cycle pf_v_cycle proc~pf_v_cycle->proc~interpolate_time_space proc~pf_v_cycle->proc~restrict_time_space_fas proc~pf_predictor pf_predictor proc~pf_predictor->proc~interpolate_time_space proc~pf_predictor->proc~restrict_time_space_fas proc~restrict_time_space_fas->proc~restrict_sdc proc~pf_v_cycle_oc pf_v_cycle_oc proc~pf_v_cycle_oc->proc~interpolate_time_space proc~pf_v_cycle_oc->proc~restrict_time_space_fas proc~pf_block_run pf_block_run proc~pf_block_run->proc~pf_v_cycle proc~pf_block_run->proc~pf_predictor proc~pf_pfasst_block_oc pf_pfasst_block_oc proc~pf_pfasst_block_oc->proc~pf_predictor_oc proc~pf_pfasst_block_oc->proc~pf_v_cycle_oc proc~pf_pfasst_run pf_pfasst_run proc~pf_pfasst_run->proc~pf_block_run

Contents

Source Code


Source Code

  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