core_update_state_impl Subroutine

public subroutine core_update_state_impl(plan, state, src_q)

ソース電荷から FMM state を更新する。

Arguments

Type IntentOptional Attributes Name
type(fmm_plan_type), intent(in) :: plan

構築済みの FMM 計画。

type(fmm_state_type), intent(inout) :: state

更新対象の state。

real(kind=dp), intent(in) :: src_q(:)

Calls

proc~~core_update_state_impl~~CallsGraph proc~core_update_state_impl core_update_state_impl proc~active_tree_level_bounds active_tree_level_bounds proc~core_update_state_impl->proc~active_tree_level_bounds proc~active_tree_level_node active_tree_level_node proc~core_update_state_impl->proc~active_tree_level_node proc~active_tree_max_depth active_tree_max_depth proc~core_update_state_impl->proc~active_tree_max_depth proc~active_tree_nnode active_tree_nnode proc~core_update_state_impl->proc~active_tree_nnode proc~initialize_fmm_state initialize_fmm_state proc~core_update_state_impl->proc~initialize_fmm_state

Source Code

  subroutine core_update_state_impl(plan, state, src_q)
    type(fmm_plan_type), intent(in) :: plan
    type(fmm_state_type), intent(inout) :: state
    real(dp), intent(in) :: src_q(:)
    integer(i32) :: n_target_nodes

    if (.not. plan%built) error stop 'FMM plan is not built.'
    if (size(src_q) /= plan%nsrc) error stop 'src_q size does not match plan.'

    state%ready = .false.
    call ensure_state_capacity(plan, state)
    n_target_nodes = merge(plan%target_nnode, plan%nnode, plan%target_tree_ready)
    !$omp parallel default(none) shared(plan, state, src_q, n_target_nodes)
    if (associated(state%multipole_active)) then
      !$omp workshare
      state%multipole_active = 0_i32
      !$omp end workshare
    end if
    if (associated(state%local_active)) then
      !$omp workshare
      state%local_active = 0_i32
      !$omp end workshare
    end if
    if (plan%nsrc > 0_i32) then
      !$omp workshare
      state%src_q = src_q
      !$omp end workshare
    end if
    if (plan%nsource_leaf <= 0_i32 .and. associated(state%multipole)) then
      !$omp workshare
      state%multipole = 0.0d0
      !$omp end workshare
    end if
    if (plan%m2l_pair_count <= 0_i32 .and. associated(state%local)) then
      !$omp workshare
      state%local = 0.0d0
      !$omp end workshare
    end if

    call p2m_leaf_moments(plan, state)
    call m2m_upward_pass(plan, state)
    if (n_target_nodes > 0_i32) then
      call m2l_accumulate(plan, state)
      call apply_periodic_target_correction(plan, state)
      call l2l_downward_pass(plan, state)
    end if
    !$omp end parallel

    state%ready = .true.
    state%update_count = state%update_count + 1_i32
  end subroutine core_update_state_impl