misdcQ_sweep Subroutine

public subroutine misdcQ_sweep(this, pf, level_index, t0, dt, nsweeps, flags)

Uses

  • proc~~misdcq_sweep~~UsesGraph proc~misdcq_sweep misdcQ_sweep module~pf_mod_hooks pf_mod_hooks proc~misdcq_sweep->module~pf_mod_hooks module~pf_mod_timer pf_mod_timer proc~misdcq_sweep->module~pf_mod_timer module~pf_mod_dtype pf_mod_dtype module~pf_mod_hooks->module~pf_mod_dtype module~pf_mod_timer->module~pf_mod_dtype iso_c_binding iso_c_binding module~pf_mod_dtype->iso_c_binding

Assign level pointer

Loop over sweeps

Arguments

Type IntentOptional AttributesName
class(pf_misdcQ_t), intent(inout) :: this
type(pf_pfasst_t), intent(inout), target:: pf

PFASST structure

integer, intent(in) :: level_index

which level to sweep on

real(kind=pfdp), intent(in) :: t0

Time at beginning of time step

real(kind=pfdp), intent(in) :: dt

time step size

integer, intent(in) :: nsweeps

number of sweeps to do

integer, intent(in), optional :: flags

Calls

proc~~misdcq_sweep~~CallsGraph proc~misdcq_sweep misdcQ_sweep proc~call_hooks call_hooks proc~misdcq_sweep->proc~call_hooks proc~pf_residual pf_residual proc~misdcq_sweep->proc~pf_residual proc~start_timer start_timer proc~misdcq_sweep->proc~start_timer proc~end_timer end_timer proc~misdcq_sweep->proc~end_timer proc~call_hooks->proc~start_timer proc~call_hooks->proc~end_timer proc~pf_residual->proc~start_timer proc~pf_residual->proc~end_timer

Contents

Source Code


Source Code

  subroutine misdcQ_sweep(this, pf, level_index, t0, dt, nsweeps,flags)
    use pf_mod_timer
    use pf_mod_hooks    
    class(pf_misdcQ_t),      intent(inout) :: this
    type(pf_pfasst_t),target,intent(inout) :: pf  !!  PFASST structure
    integer,          intent(in) :: level_index   !!  which level to sweep on
    real(pfdp),       intent(in) :: t0            !!  Time at beginning of time step
    real(pfdp),       intent(in) :: dt            !!  time step size
    integer,          intent(in) :: nsweeps       !!  number of sweeps to do
    integer, optional, intent(in   ) :: flags    
    !>  Local variables    
    class(pf_level_t), pointer:: lev
    integer                        :: m, n,k
    real(pfdp)                     :: t

    call start_timer(pf, TLEVEL+lev%index-1)

    lev => pf%levels(level_index)   !!  Assign level pointer

    do k = 1,nsweeps   !!  Loop over sweeps
       call call_hooks(pf, level_index, PF_PRE_SWEEP)    

       ! compute integrals and add fas correction
       
       do m = 1, lev%nnodes-1
          
          call lev%I(m)%setval(0.0_pfdp)
          call this%I3(m)%setval(0.0_pfdp)
          do n = 1, lev%nnodes
             call lev%I(m)%axpy(dt*this%QdiffE(m,n), lev%F(n,1))
             call lev%I(m)%axpy(dt*this%QdiffI(m,n), lev%F(n,2))
             call lev%I(m)%axpy(dt*lev%sdcmats%qmat(m,n),    lev%F(n,3))
             call this%I3(m)%axpy(dt*this%QtilI(m,n),     lev%F(n,3))
             !  Note we have to leave off the -dt*Qtil here and put it in after f2comp
          end do
          if (allocated(lev%tauQ)) then
             call lev%I(m)%axpy(1.0_pfdp, lev%tauQ(m))
          end if
       end do
       
       ! do the time-stepping
       if (k .eq. 1) then
          call lev%Q(1)%copy(lev%q0)
          call this%f_eval(lev%Q(1), t0, lev%index, lev%F(1,1),1)
          call this%f_eval(lev%Q(1), t0, lev%index, lev%F(1,2),2)
          call this%f_eval(lev%Q(1), t0, lev%index, lev%F(1,3),3)
       endif
       
       t = t0
       do m = 1, lev%nnodes-1
          t = t + dt*this%dtsdc(m)
          
          call this%rhs%setval(0.0_pfdp)
          do n = 1, m
             call this%rhs%axpy(dt*this%QtilE(m,n), lev%F(n,1))  
             call this%rhs%axpy(dt*this%QtilI(m,n), lev%F(n,2))  
          end do
          !  Add the tau term
          call this%rhs%axpy(1.0_pfdp, lev%I(m))
          !  Add the starting value
          call this%rhs%axpy(1.0_pfdp, lev%Q(1))
          
          call this%f_comp(lev%Q(m+1), t, dt*this%QtilI(m,m+1), this%rhs, lev%index, lev%F(m+1,2),2)
          
          !  Now we need to do the final subtraction for the f3 piece
          call this%rhs%copy(Lev%Q(m+1))       
          do n = 1, m
             call this%rhs%axpy(dt*this%QtilI(m,n), lev%F(n,3))  
          end do
          
          call this%rhs%axpy(-1.0_pfdp, this%I3(m))
          
          call this%f_comp(lev%Q(m+1), t, dt*this%QtilI(m,m+1), this%rhs, lev%index, lev%F(m+1,3),3)
          call this%f_eval(lev%Q(m+1), t, lev%index, lev%F(m+1,1),1)
          call this%f_eval(lev%Q(m+1), t, lev%index, lev%F(m+1,2),2)
       end do
       
       call pf_residual(pf, lev, dt)
       call lev%qend%copy(lev%Q(lev%nnodes))
       
       call call_hooks(pf, level_index, PF_POST_SWEEP)
    end do  !  End loop on sweeps
    call end_timer(pf, TLEVEL+lev%index-1)

  end subroutine misdcQ_sweep