pf_pfasst_setup Subroutine

public subroutine pf_pfasst_setup(pf)

Uses

  • proc~~pf_pfasst_setup~~UsesGraph proc~pf_pfasst_setup pf_pfasst_setup module~pf_mod_utils pf_mod_utils proc~pf_pfasst_setup->module~pf_mod_utils module~pf_mod_dtype pf_mod_dtype module~pf_mod_utils->module~pf_mod_dtype module~pf_mod_timer pf_mod_timer module~pf_mod_utils->module~pf_mod_timer iso_c_binding iso_c_binding module~pf_mod_dtype->iso_c_binding module~pf_mod_timer->module~pf_mod_dtype

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)

Arguments

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

Main pfasst structure


Calls

proc~~pf_pfasst_setup~~CallsGraph proc~pf_pfasst_setup pf_pfasst_setup proc~pf_level_setup pf_level_setup proc~pf_pfasst_setup->proc~pf_level_setup proc~pf_time_interpolation_matrix pf_time_interpolation_matrix proc~pf_pfasst_setup->proc~pf_time_interpolation_matrix proc~pf_init_sdcmats pf_init_sdcmats proc~pf_level_setup->proc~pf_init_sdcmats proc~myluq myLUq proc~pf_init_sdcmats->proc~myluq proc~pf_quadrature pf_quadrature proc~pf_init_sdcmats->proc~pf_quadrature proc~sdc_qnodes sdc_qnodes proc~pf_quadrature->proc~sdc_qnodes proc~sdc_qmats sdc_qmats proc~pf_quadrature->proc~sdc_qmats proc~poly_legendre poly_legendre proc~sdc_qnodes->proc~poly_legendre proc~poly_diff poly_diff proc~sdc_qnodes->proc~poly_diff proc~poly_roots poly_roots proc~sdc_qnodes->proc~poly_roots proc~poly_int poly_int proc~sdc_qmats->proc~poly_int proc~not_proper not_proper proc~sdc_qmats->proc~not_proper interface~poly_eval poly_eval proc~sdc_qmats->interface~poly_eval proc~qsort qsort proc~poly_roots->proc~qsort proc~poly_eval_complex poly_eval_complex proc~poly_roots->proc~poly_eval_complex interface~poly_eval->interface~poly_eval interface~poly_eval->proc~poly_eval_complex

Contents

Source Code


Source Code

  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