apply_box_boundary Subroutine

public subroutine apply_box_boundary(cfg, x, v, alive, escaped_boundary)

1ステップの更新候補位置にボックス境界条件を適用し、生存/流出状態と位置速度を更新する。

Arguments

Type IntentOptional Attributes Name
type(sim_config), intent(in) :: cfg

境界条件とボックス範囲を含むシミュレーション設定。

real(kind=dp), intent(inout) :: x(3)
real(kind=dp), intent(inout) :: v(3)
logical, intent(inout) :: alive

粒子生存フラグ(流出時は .false. へ更新)。

logical, intent(out) :: escaped_boundary

この呼び出しで境界流出が発生したか。


Called by

proc~~apply_box_boundary~~CalledByGraph proc~apply_box_boundary apply_box_boundary proc~sample_photo_raycast_particles sample_photo_raycast_particles proc~sample_photo_raycast_particles->proc~apply_box_boundary proc~sample_photo_species_state sample_photo_species_state proc~sample_photo_species_state->proc~sample_photo_raycast_particles proc~init_particle_batch_from_config init_particle_batch_from_config proc~init_particle_batch_from_config->proc~sample_photo_species_state

Source Code

  subroutine apply_box_boundary(cfg, x, v, alive, escaped_boundary)
    type(sim_config), intent(in) :: cfg
    real(dp), intent(inout) :: x(3)
    real(dp), intent(inout) :: v(3)
    logical, intent(inout) :: alive
    logical, intent(out) :: escaped_boundary

    integer(i32) :: axis
    integer(i32) :: bc
    real(dp) :: span, eps

    escaped_boundary = .false.
    if (.not. cfg%use_box) return
    if (.not. alive) return

    eps = 1.0d-12
    do axis = 1, 3
      span = cfg%box_max(axis) - cfg%box_min(axis)
      if (span <= 0.0_dp) then
        alive = .false.
        escaped_boundary = .true.
        return
      end if

      do while (x(axis) < cfg%box_min(axis))
        bc = cfg%bc_low(axis)
        call apply_one_side_boundary( &
          axis, bc, cfg%box_min(axis), cfg%box_max(axis), span, eps, x, v, alive, escaped_boundary &
          )
        if ((.not. alive) .or. (.not. escaped_boundary .and. x(axis) >= cfg%box_min(axis))) exit
      end do
      if (.not. alive) return

      do while (x(axis) > cfg%box_max(axis))
        bc = cfg%bc_high(axis)
        call apply_one_side_boundary( &
          axis, bc, cfg%box_min(axis), cfg%box_max(axis), span, eps, x, v, alive, escaped_boundary &
          )
        if ((.not. alive) .or. (.not. escaped_boundary .and. x(axis) <= cfg%box_max(axis))) exit
      end do
      if (.not. alive) return
    end do
  end subroutine apply_box_boundary