spx_particle.f90 Source File


Contents

Source Code


Source Code

!> 粒子信息
!> @note 使用数组作为存储结构, SoA (Structure of Array) 结构效率更高
module spx_particle

    use spx_kinds, only: rk
    use spx_backup, only: backup_type
    implicit none

    private
    public :: particle_type

    !> 纯粒子域,用于高速存储与计算
    type particle_type

        real(rk) :: hsml  !! 光滑长度

        real(rk), allocatable :: loc(:, :)  !! 位置, m
        real(rk), allocatable :: vel(:, :)  !! 速度, m/s
        real(rk), allocatable :: acc(:, :)  !! 加速度, m/s^2

        real(rk), allocatable :: mass(:)    !! 质量, kg
        real(rk), allocatable :: rho(:)     !! 密度, kg/m^3
        real(rk), allocatable :: drho(:)    !! 密度变化率, kg/m^3/s
        real(rk), allocatable :: u(:)       !! 内能密度
        real(rk), allocatable :: du(:)      !! 内能密度变化率
        real(rk), allocatable :: p(:)       !! 压强, Pa
        real(rk), allocatable :: c(:)       !! 声速, m/s

        integer, allocatable :: itype(:)    !! 粒子类型

    contains

        procedure :: create_backup, restore_backup

    end type particle_type

contains

    !> 创建备份
    pure subroutine create_backup(particle, backup, istart, iend)
        class(particle_type), intent(in) :: particle  !! particle
        type(backup_type), intent(out) :: backup  !! backup
        integer, intent(in) :: istart, iend  !! start and end index

        allocate(backup%loc, source=particle%loc(:, istart:iend))
        allocate(backup%vel, source=particle%vel(:, istart:iend))
        allocate(backup%rho, source=particle%rho(istart:iend))
        allocate(backup%p, source=particle%p(istart:iend))

    end subroutine create_backup

    !> 从备份中恢复
    pure subroutine restore_backup(particle, backup, istart, iend)
        class(particle_type), intent(inout) :: particle  !! particle
        type(backup_type), intent(in) :: backup  !! backup
        integer, intent(in) :: istart, iend  !! start and end index

        particle%loc(:, istart:iend) = backup%loc
        particle%vel(:, istart:iend) = backup%vel
        particle%rho(istart:iend) = backup%rho
        particle%p(istart:iend) = backup%p

    end subroutine restore_backup

end module spx_particle