call_hooks Subroutine

public subroutine call_hooks(pf, level_ind, hook)

Uses

  • proc~~call_hooks~~UsesGraph proc~call_hooks call_hooks module~pf_mod_timer pf_mod_timer proc~call_hooks->module~pf_mod_timer module~pf_mod_dtype pf_mod_dtype module~pf_mod_timer->module~pf_mod_dtype iso_c_binding iso_c_binding module~pf_mod_dtype->iso_c_binding

Subroutine to call hooks associated with the hook and level

Arguments

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

main pfasst structure

integer, intent(in) :: level_ind

which pfasst level to call hook

integer, intent(in) :: hook

which hook to call


Calls

proc~~call_hooks~~CallsGraph proc~call_hooks call_hooks proc~start_timer start_timer proc~call_hooks->proc~start_timer proc~end_timer end_timer proc~call_hooks->proc~end_timer

Called by

proc~~call_hooks~~CalledByGraph proc~call_hooks call_hooks proc~pf_check_convergence_block pf_check_convergence_block proc~pf_check_convergence_block->proc~call_hooks proc~imexq_oc_sweep imexQ_oc_sweep proc~imexq_oc_sweep->proc~call_hooks proc~pf_pfasst_block_oc pf_pfasst_block_oc proc~pf_pfasst_block_oc->proc~call_hooks proc~pf_check_convergence_oc pf_check_convergence_oc proc~pf_pfasst_block_oc->proc~pf_check_convergence_oc proc~pf_predictor_oc pf_predictor_oc proc~pf_pfasst_block_oc->proc~pf_predictor_oc proc~pf_v_cycle_oc pf_v_cycle_oc proc~pf_pfasst_block_oc->proc~pf_v_cycle_oc proc~mkrk_step mkrk_step proc~mkrk_step->proc~call_hooks proc~restrict_time_space_fas restrict_time_space_fas proc~restrict_time_space_fas->proc~call_hooks proc~interpolate_qend interpolate_qend proc~interpolate_qend->proc~call_hooks proc~pf_check_convergence_oc->proc~call_hooks proc~interpolate_time_space interpolate_time_space proc~interpolate_time_space->proc~call_hooks proc~magpicard_sweep magpicard_sweep proc~magpicard_sweep->proc~call_hooks proc~interpolate_q0 interpolate_q0 proc~interpolate_q0->proc~call_hooks proc~rk_step rk_step proc~rk_step->proc~call_hooks proc~pf_block_run pf_block_run proc~pf_block_run->proc~call_hooks proc~pf_block_run->proc~pf_check_convergence_block proc~pf_predictor pf_predictor proc~pf_block_run->proc~pf_predictor proc~pf_v_cycle pf_v_cycle proc~pf_block_run->proc~pf_v_cycle proc~pf_predictor->proc~call_hooks proc~pf_predictor->proc~restrict_time_space_fas proc~pf_predictor->proc~interpolate_time_space proc~pf_predictor->proc~interpolate_q0 proc~pf_predictor_oc->proc~call_hooks proc~pf_predictor_oc->proc~restrict_time_space_fas proc~pf_predictor_oc->proc~interpolate_time_space proc~misdcq_oc_sweep misdcQ_oc_sweep proc~misdcq_oc_sweep->proc~call_hooks proc~verlet_sweep verlet_sweep proc~verlet_sweep->proc~call_hooks proc~misdcq_sweep misdcQ_sweep proc~misdcq_sweep->proc~call_hooks proc~imexq_sweep imexQ_sweep proc~imexq_sweep->proc~call_hooks proc~imk_sweep imk_sweep proc~imk_sweep->proc~mkrk_step proc~imk_sweep->proc~rk_step proc~pf_v_cycle->proc~restrict_time_space_fas proc~pf_v_cycle->proc~interpolate_time_space proc~pf_v_cycle->proc~interpolate_q0 proc~pf_pfasst_run pf_pfasst_run proc~pf_pfasst_run->proc~pf_block_run proc~pf_v_cycle_oc->proc~restrict_time_space_fas proc~pf_v_cycle_oc->proc~interpolate_time_space

Contents

Source Code


Source Code

  subroutine call_hooks(pf, level_ind, hook)
    use pf_mod_timer
    type(pf_pfasst_t), intent(inout), target :: pf         !! main pfasst structure
    integer,           intent(in)            :: level_ind  !! which pfasst level to call hook
    integer,           intent(in)            :: hook       !! which hook to call

    integer :: i  !!  hook loop index
    integer :: l  !!  level loop index

    call start_timer(pf, THOOKS)

    pf%state%hook = hook

    if (level_ind == -1) then  ! Do to all levels
       do l = 1, pf%nlevels
          do i = 1, pf%nhooks(l,hook)
             call pf%hooks(l,hook,i)%proc(pf, pf%levels(l), pf%state)
          end do
       end do
    else  ! Do to just level level_ind
       do i = 1, pf%nhooks(level_ind,hook)
          call pf%hooks(level_ind,hook,i)%proc(pf, pf%levels(level_ind), pf%state)
       end do
    end if

    call end_timer(pf, THOOKS)
  end subroutine call_hooks