particle_set_SetMPIDataTypeParams_BH Subroutine

private impure subroutine particle_set_SetMPIDataTypeParams_BH(this, types, lengths, displacement)

Uses

  • proc~~particle_set_setmpidatatypeparams_bh~~UsesGraph proc~particle_set_setmpidatatypeparams_bh particle_set_SetMPIDataTypeParams_BH mpi_f08 mpi_f08 proc~particle_set_setmpidatatypeparams_bh->mpi_f08

Sets up parameters used when creating the MPI derived type of particles with BH type.

Arguments

Type IntentOptional Attributes Name
class(particle_set), intent(inout) :: this

Set of Lagrangian objects

type(MPI_Datatype), intent(out), allocatable :: types(:)

Array of types

integer, intent(out), allocatable :: lengths(:)

Array of lengths

integer(kind=MPI_ADDRESS_KIND), intent(out), allocatable :: displacement(:)

Array of displacements


Calls

proc~~particle_set_setmpidatatypeparams_bh~~CallsGraph proc~particle_set_setmpidatatypeparams_bh particle_set_SetMPIDataTypeParams_BH mpi_get_address mpi_get_address proc~particle_set_setmpidatatypeparams_bh->mpi_get_address

Called by

proc~~particle_set_setmpidatatypeparams_bh~~CalledByGraph proc~particle_set_setmpidatatypeparams_bh particle_set_SetMPIDataTypeParams_BH proc~particle_set_setmpidatatypeparams particle_set%particle_set_SetMPIDataTypeParams proc~particle_set_setmpidatatypeparams->proc~particle_set_setmpidatatypeparams_bh

Source Code

    impure subroutine particle_set_SetMPIDataTypeParams_BH(this,types,lengths,displacement)
      !> Sets up parameters used when creating the MPI derived type of
      ! particles with BH type.
      use mpi_f08
      implicit none
      class(particle_set), intent(inout) :: this                               !! Set of Lagrangian objects
      type(MPI_Datatype),  intent(out),   &
                             allocatable :: types(:)                           !! Array of types
      integer,             intent(out),   &
                             allocatable :: lengths(:)                         !! Array of lengths
      integer(kind=MPI_ADDRESS_KIND),    &
                           intent(out),   &
                             allocatable :: displacement(:)                    !! Array of displacements
      ! Work variable
      integer :: N
      type(particle_BH_obj) :: sample
      integer :: ierr

      associate(mpi=>this%parallel)
        ! Create the MPI structure
        N=27
        allocate(types(N),displacement(N),lengths(N))
        types( 1) = mpi%INT8    ; lengths( 1) = 1                    ! id     --  integer8
        types( 2) = mpi%REAL_WP ; lengths( 2) = size(sample%p)       ! p      --  real*3
        types( 3) = mpi%INTEGER ; lengths( 3) = size(sample%c)       ! c      --  integer*3
        types( 4) = mpi%INTEGER ; lengths( 4) = 1                    ! s      --  integer
        types( 5) = mpi%REAL_WP ; lengths( 5) = 1                    ! d      --  real
        types( 6) = mpi%REAL_WP ; lengths( 6) = 1                    ! rho    --  real
        types( 7) = mpi%REAL_WP ; lengths( 7) = size(sample%v)       ! v      --  real*3
        types( 8) = mpi%REAL_WP ; lengths( 8) = size(sample%w)       ! w      --  real*3
        types( 9) = mpi%REAL_WP ; lengths( 9) = size(sample%Fh)      ! Fh     --  real*3
        types(10) = mpi%REAL_WP ; lengths(10) = size(sample%Th)      ! Th     --  real*3
        types(11) = mpi%REAL_WP ; lengths(11) = size(sample%Fc)      ! Fc     --  real*3
        types(12) = mpi%REAL_WP ; lengths(12) = size(sample%Tc)      ! Tc     --  real*3
        types(13) = mpi%REAL_WP ; lengths(13) = size(sample%pold)    ! pold   --  real*3
        types(14) = mpi%REAL_WP ; lengths(14) = size(sample%vold)    ! vold   --  real*3
        types(15) = mpi%REAL_WP ; lengths(15) = size(sample%wold)    ! wold   --  real*3
        types(16) = mpi%REAL_WP ; lengths(16) = size(sample%Fhold)   ! Fhold  --  real*3
        types(17) = mpi%REAL_WP ; lengths(17) = size(sample%Thold)   ! Thold  --  real*3
        types(18) = mpi%REAL_WP ; lengths(18) = size(sample%Fcold)   ! Fcold  --  real*3
        types(19) = mpi%REAL_WP ; lengths(19) = size(sample%Tcold)   ! Tcold  --  real*3
        types(20) = mpi%REAL_WP ; lengths(20) = size(sample%Fb)      ! Fb     --  real*3
        types(21) = mpi%REAL_WP ; lengths(21) = size(sample%vs)      ! vs     --  real*3
        types(22) = mpi%REAL_WP ; lengths(22) = size(sample%gx)      ! gx     --  real*3
        types(23) = mpi%REAL_WP ; lengths(23) = size(sample%gy)      ! gy     --  real*3
        types(24) = mpi%REAL_WP ; lengths(24) = size(sample%gz)      ! gz     --  real*3
        types(25) = mpi%REAL_WP ; lengths(25) = size(sample%Fix)     ! Fix    --  real*3
        types(26) = mpi%REAL_WP ; lengths(26) = size(sample%Fiy)     ! Fiy    --  real*3
        types(27) = mpi%REAL_WP ; lengths(27) = size(sample%Fiz)     ! Fiz    --  real*3
        ! Count the displacement for this structure
        call MPI_GET_ADDRESS(sample%id,      displacement( 1),ierr)
        call MPI_GET_ADDRESS(sample%p(1),    displacement( 2),ierr)
        call MPI_GET_ADDRESS(sample%c(1),    displacement( 3),ierr)
        call MPI_GET_ADDRESS(sample%s,       displacement( 4),ierr)
        call MPI_GET_ADDRESS(sample%d,       displacement( 5),ierr)
        call MPI_GET_ADDRESS(sample%rho,     displacement( 6),ierr)
        call MPI_GET_ADDRESS(sample%v(1),    displacement( 7),ierr)
        call MPI_GET_ADDRESS(sample%w(1),    displacement( 8),ierr)
        call MPI_GET_ADDRESS(sample%Fh(1),   displacement( 9),ierr)
        call MPI_GET_ADDRESS(sample%Th(1),   displacement(10),ierr)
        call MPI_GET_ADDRESS(sample%Fc(1),   displacement(11),ierr)
        call MPI_GET_ADDRESS(sample%Tc(1),   displacement(12),ierr)
        call MPI_GET_ADDRESS(sample%pold(1), displacement(13),ierr)
        call MPI_GET_ADDRESS(sample%vold(1), displacement(14),ierr)
        call MPI_GET_ADDRESS(sample%wold(1), displacement(15),ierr)
        call MPI_GET_ADDRESS(sample%Fhold(1),displacement(16),ierr)
        call MPI_GET_ADDRESS(sample%Thold(1),displacement(17),ierr)
        call MPI_GET_ADDRESS(sample%Fcold(1),displacement(18),ierr)
        call MPI_GET_ADDRESS(sample%Tcold(1),displacement(19),ierr)
        call MPI_GET_ADDRESS(sample%Fb(1),   displacement(20),ierr)
        call MPI_GET_ADDRESS(sample%vs(1),   displacement(21),ierr)
        call MPI_GET_ADDRESS(sample%gx(1),   displacement(22),ierr)
        call MPI_GET_ADDRESS(sample%gy(1),   displacement(23),ierr)
        call MPI_GET_ADDRESS(sample%gz(1),   displacement(24),ierr)
        call MPI_GET_ADDRESS(sample%Fix(1),  displacement(25),ierr)
        call MPI_GET_ADDRESS(sample%Fiy(1),  displacement(26),ierr)
        call MPI_GET_ADDRESS(sample%Fiz(1),  displacement(27),ierr)
        displacement = displacement - displacement( 1)
      end associate

      return
    end subroutine particle_set_SetMPIDataTypeParams_BH