pf_pfasst_run Subroutine

public subroutine pf_pfasst_run(pf, q0, dt, tend, nsteps, qend, flags)

This is the main interface to pfasst. It examines the parameters and decides which subroutine to call to execute the code correctly Set the number of time steps to do The user can either pass in the number of time steps or pass in the time step size and length of run Allocate stuff for holding results

Arguments

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

The complete PFASST structure

class(pf_encap_t), intent(in) :: q0

The initial condition

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

The time step for each processor

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

The final time of run

integer, intent(in), optional :: nsteps

The number of time steps

class(pf_encap_t), intent(inout), optional :: qend

The computed solution at tend

integer, intent(in), optional :: flags(:)

User defnined flags


Calls

proc~~pf_pfasst_run~~CallsGraph proc~pf_pfasst_run pf_pfasst_run proc~initialize_results initialize_results proc~pf_pfasst_run->proc~initialize_results proc~pf_block_run pf_block_run proc~pf_pfasst_run->proc~pf_block_run proc~pf_check_convergence_block pf_check_convergence_block proc~pf_block_run->proc~pf_check_convergence_block proc~pf_v_cycle pf_v_cycle proc~pf_block_run->proc~pf_v_cycle proc~pf_predictor pf_predictor proc~pf_block_run->proc~pf_predictor proc~end_timer end_timer proc~pf_block_run->proc~end_timer proc~pf_broadcast pf_broadcast proc~pf_block_run->proc~pf_broadcast proc~call_hooks call_hooks proc~pf_block_run->proc~call_hooks proc~start_timer start_timer proc~pf_block_run->proc~start_timer proc~pf_check_convergence_block->proc~call_hooks proc~pf_check_residual pf_check_residual proc~pf_check_convergence_block->proc~pf_check_residual proc~pf_send_status pf_send_status proc~pf_check_convergence_block->proc~pf_send_status proc~restrict_time_space_fas restrict_time_space_fas proc~pf_v_cycle->proc~restrict_time_space_fas proc~pf_post pf_post proc~pf_v_cycle->proc~pf_post proc~pf_residual pf_residual proc~pf_v_cycle->proc~pf_residual proc~interpolate_time_space interpolate_time_space proc~pf_v_cycle->proc~interpolate_time_space proc~pf_send pf_send proc~pf_v_cycle->proc~pf_send proc~interpolate_q0 interpolate_q0 proc~pf_v_cycle->proc~interpolate_q0 proc~pf_recv pf_recv proc~pf_v_cycle->proc~pf_recv proc~pf_predictor->proc~end_timer proc~pf_predictor->proc~call_hooks proc~pf_predictor->proc~start_timer proc~pf_predictor->proc~restrict_time_space_fas proc~pf_predictor->proc~pf_residual proc~pf_predictor->proc~interpolate_time_space proc~pf_predictor->proc~pf_send proc~pf_predictor->proc~interpolate_q0 proc~pf_predictor->proc~pf_recv end_timer end_timer proc~pf_broadcast->end_timer start_timer start_timer proc~pf_broadcast->start_timer proc~call_hooks->proc~end_timer proc~call_hooks->proc~start_timer proc~restrict_time_space_fas->proc~end_timer proc~restrict_time_space_fas->proc~call_hooks proc~restrict_time_space_fas->proc~start_timer proc~restrict_sdc restrict_sdc proc~restrict_time_space_fas->proc~restrict_sdc proc~pf_residual->proc~end_timer proc~pf_residual->proc~start_timer proc~interpolate_time_space->proc~end_timer proc~interpolate_time_space->proc~call_hooks proc~interpolate_time_space->proc~start_timer proc~pf_apply_mat pf_apply_mat proc~interpolate_time_space->proc~pf_apply_mat proc~pf_send->end_timer proc~pf_send->start_timer proc~interpolate_q0->proc~end_timer proc~interpolate_q0->proc~call_hooks proc~interpolate_q0->proc~start_timer proc~pf_recv->end_timer proc~pf_recv->start_timer proc~restrict_sdc->proc~pf_apply_mat

Contents

Source Code


Source Code

  subroutine pf_pfasst_run(pf, q0, dt, tend, nsteps, qend, flags)
    type(pf_pfasst_t), intent(inout), target   :: pf   !!  The complete PFASST structure
    class(pf_encap_t), intent(in   )           :: q0   !!  The initial condition
    real(pfdp),        intent(in   )           :: dt   !!  The time step for each processor
    real(pfdp),        intent(in   )           :: tend !!  The final time of run
    integer,           intent(in   ), optional :: nsteps  !!  The number of time steps
    class(pf_encap_t), intent(inout), optional :: qend    !!  The computed solution at tend
    integer,           intent(in   ), optional :: flags(:)!!  User defnined flags


    !  Local variables
    integer :: nproc  !!  Total number of processors
    integer :: nsteps_loc  !!  local number of time steps    
    real(pfdp) :: tend_loc !!  The final time of run

    ! make a local copy of nproc
    nproc = pf%comm%nproc

    !>  Set the number of time steps to do
    !!  The user can either pass in the number of time steps or
    !!  pass in the time step size and length of run
    if (present(nsteps)) then
      nsteps_loc = nsteps
      tend_loc=dble(nsteps_loc*dt)
    else
      nsteps_loc = ceiling(tend/dt)
      !  Do  sanity check on steps
      if (abs(real(nsteps_loc,pfdp)-tend/dt) > dt/100.0) then
        print *,'dt=',dt
        print *,'nsteps=',nsteps_loc
        print *,'tend=',tend
        stop "Invalid nsteps"
      end if
    end if
    pf%state%nsteps = nsteps_loc

    !>  Allocate stuff for holding results
    call initialize_results(pf%results,nsteps_loc, pf%niters, pf%comm%nproc, pf%nlevels,pf%rank)

    !  do sanity checks on Nproc
    if (mod(nsteps,nproc) > 0) stop "ERROR: nsteps must be multiple of nproc (pf_parallel.f90)."

    if (present(qend)) then
       call pf_block_run(pf, q0, dt, nsteps_loc,qend=qend,flags=flags)             
    else
       call pf_block_run(pf, q0, dt,  nsteps_loc,flags=flags)             
    end if

    if (pf%save_results) call pf%results%dump(pf%results)


    !  What we would like to do is check for
    !  1.  nlevels==1  and nprocs ==1 -> Serial SDC
    !      Predictor is either spreadQ or nothing
    !      Then we just call a loop on sweeps
    !      Communication is copy
    !  2.  nlevels > 1  and nprocs ==1 -> Serial MLSDC
    !      Predictor is needed to populate levels (or nothing)
    !      Then we just call a loop on MLSDC sweeps
    !      Communication is copy
    !  3.  nlevels == 1  and nprocs > 1 -> Pipelined SDC
    !      Predictor is just like PFASST, but on finest (only) level (or nothing)
    !  4.  nlevels > 1  and nprocs > 1 -> PFASST
  end subroutine pf_pfasst_run