pf_residual Subroutine

public subroutine pf_residual(pf, lev, dt, flag)

Compute full residual at each node and measure it's size

Arguments

Type IntentOptional AttributesName
type(pf_pfasst_t), intent(inout) :: pf
class(pf_level_t), intent(inout) :: lev
real(kind=pfdp), intent(in) :: dt
integer, intent(in), optional :: flag

Calls

proc~~pf_residual~~CallsGraph proc~pf_residual pf_residual proc~start_timer start_timer proc~pf_residual->proc~start_timer proc~end_timer end_timer proc~pf_residual->proc~end_timer

Called by

proc~~pf_residual~~CalledByGraph proc~pf_residual pf_residual proc~pf_v_cycle pf_v_cycle proc~pf_v_cycle->proc~pf_residual proc~pf_predictor pf_predictor proc~pf_predictor->proc~pf_residual proc~mkrk_step mkrk_step proc~mkrk_step->proc~pf_residual proc~imk_actually_sweep imk_actually_sweep proc~imk_actually_sweep->proc~pf_residual proc~pf_v_cycle_oc pf_v_cycle_oc proc~pf_v_cycle_oc->proc~pf_residual proc~magpicard_sweep magpicard_sweep proc~magpicard_sweep->proc~pf_residual proc~rk_step rk_step proc~rk_step->proc~pf_residual proc~imexq_oc_sweep imexQ_oc_sweep proc~imexq_oc_sweep->proc~pf_residual proc~imexq_sweep imexQ_sweep proc~imexq_sweep->proc~pf_residual proc~pf_predictor_oc pf_predictor_oc proc~pf_predictor_oc->proc~pf_residual proc~misdcq_oc_sweep misdcQ_oc_sweep proc~misdcq_oc_sweep->proc~pf_residual proc~verlet_sweep verlet_sweep proc~verlet_sweep->proc~pf_residual proc~misdcq_sweep misdcQ_sweep proc~misdcq_sweep->proc~pf_residual proc~imk_sweep imk_sweep proc~imk_sweep->proc~mkrk_step proc~imk_sweep->proc~imk_actually_sweep proc~imk_sweep->proc~rk_step 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_v_cycle_oc proc~pf_pfasst_block_oc->proc~pf_predictor_oc proc~pf_pfasst_run pf_pfasst_run proc~pf_pfasst_run->proc~pf_block_run

Contents

Source Code


Source Code

  subroutine pf_residual(pf, lev, dt, flag)
    type(pf_pfasst_t),  intent(inout) :: pf
    class(pf_level_t),  intent(inout) :: lev
    real(pfdp),         intent(in)    :: dt
    integer, optional,  intent(in)    :: flag

    real(pfdp) :: res_norms(lev%nnodes-1)    !!  Holds norms of residual
    real(pfdp) :: sol_norms(lev%nnodes)      !!  Holds norms of solution ! for adjoint: need sol at t0 as well, not only t0+dt
    integer :: m
    
    call start_timer(pf, TRESIDUAL)

    call lev%ulevel%sweeper%residual(lev, dt, flag)

    ! compute max residual norm
    sol_norms(1) = lev%Q(1)%norm(flag) ! for adjoint
    do m = 1, lev%nnodes-1
       res_norms(m) = lev%R(m)%norm(flag)
       sol_norms(m+1) = lev%Q(m+1)%norm(flag) ! only the value at lev%nnodes is needed for forward integration, right? 
    end do

    !    lev%residual = res_norms(lev%nnodes-1)
    m = lev%nnodes  ! for usual forward integration
    if(present(flag)) then
      if(flag==2) m = 1
    end if
    lev%residual = maxval(res_norms)    
    if (sol_norms(m) > 0.0d0) then
       lev%residual_rel = lev%residual/sol_norms(m)
    else
       lev%residual_rel = 0.0d0
    end if

    call end_timer(pf, TRESIDUAL)

  end subroutine pf_residual