imexQ_sweep Subroutine

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

Uses

  • proc~~imexq_sweep~~UsesGraph proc~imexq_sweep imexQ_sweep module~pf_mod_hooks pf_mod_hooks proc~imexq_sweep->module~pf_mod_hooks module~pf_mod_timer pf_mod_timer proc~imexq_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

Perform nsweep SDC sweeps on level level_index and set qend appropriately. Assign level pointer

Loop over sweeps Loop over substeps Accumulate rhs Add the integral term

Add the starting value

Solve for the implicit piece Compute explicit function on new value

End substep loop

Arguments

Type IntentOptional AttributesName
class(pf_imexQ_t), intent(inout) :: this

Inputs

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~~imexq_sweep~~CallsGraph proc~imexq_sweep imexQ_sweep proc~call_hooks call_hooks proc~imexq_sweep->proc~call_hooks proc~pf_residual pf_residual proc~imexq_sweep->proc~pf_residual proc~start_timer start_timer proc~imexq_sweep->proc~start_timer proc~end_timer end_timer proc~imexq_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 imexQ_sweep(this, pf, level_index, t0, dt,nsweeps, flags)
    use pf_mod_timer
    use pf_mod_hooks

    !>  Inputs
    class(pf_imexQ_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(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    !!  points to current level

    integer     :: m, n,k   !!  Loop variables
    real(pfdp)  :: t        !!  Time at nodes
    
    lev => pf%levels(level_index)   !!  Assign level pointer

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

    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)
          if (this%explicit) then
             do n = 1, lev%nnodes
                call lev%I(m)%axpy(dt*this%QdiffE(m,n), lev%F(n,1))
             end do
          end if
          if (this%implicit) then
             do n = 1, lev%nnodes
                call lev%I(m)%axpy(dt*this%QdiffI(m,n), lev%F(n,2))
             end do
          end if
          if (allocated(lev%tauQ)) then
          call lev%I(m)%axpy(1.0_pfdp, lev%tauQ(m))
          end if
       end do
       !  Recompute the first function value if this is first sweep
       if (k .eq. 1) then
          call lev%Q(1)%copy(lev%q0)
          if (this%explicit) &
               call this%f_eval(lev%Q(1), t0, lev%index, lev%F(1,1),1)
          if (this%implicit) &
               call this%f_eval(lev%Q(1), t0, lev%index, lev%F(1,2),2)
       end if
       
       t = t0
       ! do the sub-stepping in sweep
       do m = 1, lev%nnodes-1  !!  Loop over substeps
          t = t + dt*this%dtsdc(m)

          !>  Accumulate rhs
          call this%rhs%setval(0.0_pfdp)
          do n = 1, m
             if (this%explicit) &
                  call this%rhs%axpy(dt*this%QtilE(m,n), lev%F(n,1))
             if (this%implicit) &
                  call this%rhs%axpy(dt*this%QtilI(m,n), lev%F(n,2))
          end do
          !>  Add the integral 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))

          !>  Solve for the implicit piece
          if (this%implicit) then
             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)
          else
             call lev%Q(m+1)%copy(this%rhs)
          end if
          !>  Compute explicit function on new value
          if (this%explicit) &
               call this%f_eval(lev%Q(m+1), t, lev%index, lev%F(m+1,1),1)

          
       end do  !!  End substep loop
       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 imexQ_sweep