pf_level_setup Subroutine

public subroutine pf_level_setup(pf, lev)

Uses

  • proc~~pf_level_setup~~UsesGraph proc~pf_level_setup pf_level_setup module~pf_mod_quadrature pf_mod_quadrature proc~pf_level_setup->module~pf_mod_quadrature module~pf_mod_dtype pf_mod_dtype module~pf_mod_quadrature->module~pf_mod_dtype iso_c_binding iso_c_binding module~pf_mod_quadrature->iso_c_binding module~pf_mod_dtype->iso_c_binding

Setup (allocate) PFASST level If the level is already setup, calling this again will allocate (or deallocate) tauQ appropriately. do some sanity checks (re)allocate tauQ (may to need create/destroy tauQ dynamically when doing AMR) skip the rest if we're already allocated allocate flat buffers for send, and recv allocate nodes, flags, and integration matrices make quadrature matrices Allocate and compute all the matrices initialize sweeper allocate solution and function arrays

Arguments

Type IntentOptional AttributesName
type(pf_pfasst_t), intent(in) :: pf

Main pfasst structure

class(pf_level_t), intent(inout), target:: lev

Level to set up


Calls

proc~~pf_level_setup~~CallsGraph proc~pf_level_setup pf_level_setup 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

Called by

proc~~pf_level_setup~~CalledByGraph proc~pf_level_setup pf_level_setup proc~pf_pfasst_setup pf_pfasst_setup proc~pf_pfasst_setup->proc~pf_level_setup

Contents

Source Code


Source Code

  subroutine pf_level_setup(pf, lev)
    use pf_mod_quadrature
    type(pf_pfasst_t), intent(in   )         :: pf   !!  Main pfasst structure
    class(pf_level_t), intent(inout), target :: lev  !!  Level to set up

    integer :: mpibuflen, nnodes, npieces, nnodes0
    integer :: i,ierr

    !> do some sanity checks
    if (lev%mpibuflen <= 0) stop "ERROR: Invalid mpibuflen, set before calling pfasst_setup."
    if (lev%nnodes <= 0) stop "ERROR: Invalid nnodes, should have been set in pfasst_create."
    if (lev%nsweeps <= 0) stop "ERROR: Invalid nsweeps, should have been set in pfasst_create."

    mpibuflen  = lev%mpibuflen
    nnodes = lev%nnodes

    lev%residual = -1.0_pfdp


    !> (re)allocate tauQ (may to need create/destroy tauQ dynamically  when doing AMR)
    if ((lev%index < pf%nlevels) .and. (.not. allocated(lev%tauQ))) then
       call lev%ulevel%factory%create_array(lev%tauQ, nnodes-1, lev%index,  lev%shape)
    else if ((lev%index >= pf%nlevels) .and. (allocated(lev%tauQ))) then
       deallocate(lev%tauQ)
    end if

    !> skip the rest if we're already allocated
    if (lev%allocated) return
    lev%allocated = .true.

    !> allocate flat buffers for send, and recv
    allocate(lev%send(mpibuflen),stat=ierr)
    if (ierr /= 0) stop "allocate error send buffer"
    allocate(lev%recv(mpibuflen),stat=ierr)
    if (ierr /= 0) stop "allocate error recieve buffer"

    !> allocate nodes, flags, and integration matrices
    allocate(lev%nodes(nnodes),stat=ierr)
    if (ierr /= 0) stop "allocate error nodes"
    allocate(lev%nflags(nnodes),stat=ierr)
    if (ierr /= 0) stop "allocate error nflags"
    lev%nflags=0
    !> make quadrature matrices
    if (btest(pf%qtype, 8)) then
       nnodes0=pf%levels(1)%nnodes
    else
       nnodes0=pf%levels(pf%nlevels)%nnodes
    end if
    !>  Allocate and compute all the matrices
    allocate(lev%sdcmats,stat=ierr)
    if (ierr /= 0) stop "allocate error sdcmats"


    call pf_init_sdcmats(lev%sdcmats,pf%qtype, nnodes,nnodes0,lev%nflags)

    lev%nodes = lev%sdcmats%qnodes

    !>  initialize sweeper
    lev%ulevel%sweeper%use_LUq=pf%use_LUq
    call lev%ulevel%sweeper%initialize(lev)

    
    if (pf%use_rk_stepper)  call lev%ulevel%stepper%initialize(lev)

    !> allocate solution and function arrays
    npieces = lev%ulevel%sweeper%npieces

    call lev%ulevel%factory%create_array(lev%Q, nnodes, lev%index,  lev%shape)
    call lev%ulevel%factory%create_array(lev%Fflt, nnodes*npieces, lev%index,  lev%shape)
    do i = 1, nnodes*npieces
       call lev%Fflt(i)%setval(0.0_pfdp, 0)
    end do
    lev%F(1:nnodes,1:npieces) => lev%Fflt
    call lev%ulevel%factory%create_array(lev%I, nnodes-1, lev%index,  lev%shape)
    call lev%ulevel%factory%create_array(lev%R, nnodes-1, lev%index,  lev%shape)

    !  Need space for old function values in imexR sweepers
    call lev%ulevel%factory%create_array(lev%pFflt, nnodes*npieces, lev%index, lev%shape)
    lev%pF(1:nnodes,1:npieces) => lev%pFflt
    if (lev%index < pf%nlevels) then
       call lev%ulevel%factory%create_array(lev%pQ, nnodes, lev%index,  lev%shape)
    end if
    call lev%ulevel%factory%create_single(lev%qend, lev%index,   lev%shape)
    call lev%ulevel%factory%create_single(lev%q0, lev%index,   lev%shape)
    
  end subroutine pf_level_setup