Setup both the PFASST object and the comm object loop over levels to set parameters Loop over levels setting interpolation and restriction matrices (in time)
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(pf_pfasst_t), | intent(inout), | target | :: | pf | Main pfasst structure |
subroutine pf_pfasst_setup(pf)
use pf_mod_utils
type(pf_pfasst_t), intent(inout), target :: pf !! Main pfasst structure
class(pf_level_t), pointer :: lev_fine, lev_coarse !! Pointers to level structures for brevity
integer :: l !! Level loop index
integer :: ierr !! error flag
!> loop over levels to set parameters
do l = 1, pf%nlevels
call pf_level_setup(pf, pf%levels(l))
end do
!> Loop over levels setting interpolation and restriction matrices (in time)
do l = pf%nlevels, 2, -1
lev_fine => pf%levels(l); lev_coarse => pf%levels(l-1)
allocate(lev_fine%tmat(lev_fine%nnodes,lev_coarse%nnodes),stat=ierr)
if (ierr /= 0) stop "allocate error tmat"
allocate(lev_fine%rmat(lev_coarse%nnodes,lev_fine%nnodes),stat=ierr)
if (ierr /= 0) stop "allocate error rmat"
! with the RK stepper, no need to interpolate and restrict in time
! we only copy the first node and last node betweem levels
if (pf%use_rk_stepper .eqv. .true.) then
lev_fine%tmat = 0.0_pfdp
lev_fine%rmat = 0.0_pfdp
lev_fine%tmat(1,1) = 1.0_pfdp
lev_fine%tmat(lev_fine%nnodes,lev_coarse%nnodes) = 1.0_pfdp
lev_fine%rmat(1,1) = 1.0_pfdp
lev_fine%rmat(lev_coarse%nnodes,lev_fine%nnodes) = 1.0_pfdp
else ! else compute the interpolation matrix
call pf_time_interpolation_matrix(lev_fine%nodes, lev_fine%nnodes, lev_coarse%nodes, lev_coarse%nnodes, lev_fine%tmat)
call pf_time_interpolation_matrix(lev_coarse%nodes, lev_coarse%nnodes, lev_fine%nodes, lev_fine%nnodes, lev_fine%rmat)
endif
end do
end subroutine pf_pfasst_setup