sample_species_state Subroutine

public subroutine sample_species_state(sim, spec, n, x, v, barrier_normal_energy, vmin_normal, apply_barrier_energy_shift, temperature_k_override, drift_velocity_override)

1粒子種ぶんの位置・速度サンプルをまとめて生成する。

Arguments

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

ボックス境界・バッチ時間などのシミュレーション設定。

type(particle_species_spec), intent(in) :: spec

1粒子種の注入設定。

integer(kind=i32), intent(in) :: n

生成粒子数。

real(kind=dp), intent(out) :: x(:,:)
real(kind=dp), intent(out) :: v(:,:)
real(kind=dp), intent(in), optional :: barrier_normal_energy

reservoir_face 法線方向のエネルギー障壁 2 q Δφ / m [m^2/s^2]。

real(kind=dp), intent(in), optional :: vmin_normal

reservoir_face 法線速度の下限 [m/s]。

logical, intent(in), optional :: apply_barrier_energy_shift

reservoir_face 法線速度へ障壁エネルギー変換を適用するか。

real(kind=dp), intent(in), optional :: temperature_k_override
real(kind=dp), intent(in), optional :: drift_velocity_override(3)

Calls

proc~~sample_species_state~~CallsGraph proc~sample_species_state sample_species_state proc~lower_ascii lower_ascii proc~sample_species_state->proc~lower_ascii proc~sample_reservoir_face_particles sample_reservoir_face_particles proc~sample_species_state->proc~sample_reservoir_face_particles proc~sample_shifted_maxwell_velocities sample_shifted_maxwell_velocities proc~sample_species_state->proc~sample_shifted_maxwell_velocities proc~sample_uniform_positions sample_uniform_positions proc~sample_species_state->proc~sample_uniform_positions proc~species_temperature_k~2 species_temperature_k proc~sample_species_state->proc~species_temperature_k~2 proc~sample_reservoir_face_particles->proc~sample_shifted_maxwell_velocities

Called by

proc~~sample_species_state~~CalledByGraph proc~sample_species_state sample_species_state proc~init_particle_batch_from_config init_particle_batch_from_config proc~init_particle_batch_from_config->proc~sample_species_state proc~init_particles_from_config init_particles_from_config proc~init_particles_from_config->proc~sample_species_state

Source Code

  subroutine sample_species_state( &
    sim, spec, n, x, v, barrier_normal_energy, vmin_normal, apply_barrier_energy_shift, &
    temperature_k_override, drift_velocity_override &
    )
    type(sim_config), intent(in) :: sim
    type(particle_species_spec), intent(in) :: spec
    integer(i32), intent(in) :: n
    real(dp), intent(out) :: x(:, :)
    real(dp), intent(out) :: v(:, :)
    real(dp), intent(in), optional :: barrier_normal_energy
    real(dp), intent(in), optional :: vmin_normal
    logical, intent(in), optional :: apply_barrier_energy_shift
    real(dp), intent(in), optional :: temperature_k_override
    real(dp), intent(in), optional :: drift_velocity_override(3)
    logical :: apply_shift
    real(dp) :: temperature_k_local, drift_velocity_local(3)

    if (n <= 0_i32) return
    apply_shift = .true.
    if (present(apply_barrier_energy_shift)) apply_shift = apply_barrier_energy_shift
    temperature_k_local = species_temperature_k(spec)
    if (present(temperature_k_override)) temperature_k_local = temperature_k_override
    drift_velocity_local = spec%drift_velocity
    if (present(drift_velocity_override)) drift_velocity_local = drift_velocity_override
    select case (trim(lower_ascii(spec%source_mode)))
    case ('volume_seed')
      call sample_uniform_positions(spec%pos_low, spec%pos_high, x)
      call sample_shifted_maxwell_velocities( &
        drift_velocity_local, spec%m_particle, v, temperature_k=temperature_k_local &
        )
    case ('reservoir_face')
      if (present(barrier_normal_energy) .and. present(vmin_normal)) then
        call sample_reservoir_face_particles( &
          sim%box_min, sim%box_max, spec%inject_face, spec%pos_low, spec%pos_high, drift_velocity_local, &
          spec%m_particle, temperature_k_local, sim%batch_duration, x, v, &
          barrier_normal_energy=barrier_normal_energy, vmin_normal=vmin_normal, position_jitter_dt=sim%dt, &
          apply_barrier_energy_shift=apply_shift &
          )
      else if (present(barrier_normal_energy)) then
        call sample_reservoir_face_particles( &
          sim%box_min, sim%box_max, spec%inject_face, spec%pos_low, spec%pos_high, drift_velocity_local, &
          spec%m_particle, temperature_k_local, sim%batch_duration, x, v, &
          barrier_normal_energy=barrier_normal_energy, position_jitter_dt=sim%dt, &
          apply_barrier_energy_shift=apply_shift &
          )
      else if (present(vmin_normal)) then
        call sample_reservoir_face_particles( &
          sim%box_min, sim%box_max, spec%inject_face, spec%pos_low, spec%pos_high, drift_velocity_local, &
          spec%m_particle, temperature_k_local, sim%batch_duration, x, v, &
          vmin_normal=vmin_normal, position_jitter_dt=sim%dt, apply_barrier_energy_shift=apply_shift &
          )
      else
        call sample_reservoir_face_particles( &
          sim%box_min, sim%box_max, spec%inject_face, spec%pos_low, spec%pos_high, drift_velocity_local, &
          spec%m_particle, temperature_k_local, sim%batch_duration, x, v, position_jitter_dt=sim%dt, &
          apply_barrier_energy_shift=apply_shift &
          )
      end if
    case ('photo_raycast')
      error stop 'sample_species_state does not support photo_raycast. Use sample_photo_species_state.'
    case default
      error stop 'Unknown particles.species.source_mode.'
    end select
  end subroutine sample_species_state