Octopus
partition_transfer_oct_m Module Reference

Data Types

type  partition_transfer_t
 The partition transfer object ensures that during a mesh transfer points are associated to the correct ranks, if domain parallelization is used. More...
 

Functions/Subroutines

subroutine, public partition_transfer_init (this, np, global_index, mpi_grp_in, mpi_grp_out, part_out, nsend, nrec, order_in, order_out, inverse)
 initialize the partition transfer object More...
 
subroutine, public partition_transfer_end (this)
 
subroutine, public dpartition_transfer (this, f_in, f_out)
 
subroutine, public zpartition_transfer (this, f_in, f_out)
 

Function/Subroutine Documentation

◆ partition_transfer_init()

subroutine, public partition_transfer_oct_m::partition_transfer_init ( type(partition_transfer_t), intent(out)  this,
integer, intent(in)  np,
integer(int64), dimension(:), intent(in)  global_index,
type(mpi_grp_t), intent(in), target  mpi_grp_in,
type(mpi_grp_t), intent(in), target  mpi_grp_out,
integer, dimension(:), intent(in)  part_out,
integer, intent(out)  nsend,
integer, intent(out)  nrec,
integer(int64), dimension(:), intent(out), allocatable  order_in,
integer(int64), dimension(:), intent(out), allocatable  order_out,
logical, intent(in), optional  inverse 
)

initialize the partition transfer object

Warning
Input and output groups may have a different number of processes. In that case the transfer will only work if some further constraints are met (see sanity checks below). One of the cases where it should work is when one of the groups is a subgroup of a Cartesian topology created from the other group. This is the case when one of the groups is mpi_world, the other is the parallelization in domains group, and there are no slaves. If the optional argument inverse is set to .true., the direction of the transfer is inverse: part_out is assumed to specify the partition on the incoming mpi group. This is useful if one needs to get points defined on the output group.
Parameters
[in]npthe number of local points in the input partition
[in]global_indexthe global indices of the points of the input partition
[in]mpi_grp_inmpi group of in_mesh
[in]mpi_grp_outmpi group of out_mesh
[in]part_outpoint -> partition
[out]nsendTotal number of points to be sent
[out]nrecTotal number of points to be received

Definition at line 164 of file partition_transfer.F90.

◆ partition_transfer_end()

subroutine, public partition_transfer_oct_m::partition_transfer_end ( type(partition_transfer_t), intent(inout)  this)

Definition at line 389 of file partition_transfer.F90.

◆ dpartition_transfer()

subroutine, public partition_transfer_oct_m::dpartition_transfer ( type(partition_transfer_t), intent(in)  this,
real(real64), dimension(:), intent(in), contiguous  f_in,
real(real64), dimension(:), intent(out), contiguous  f_out 
)

Definition at line 469 of file partition_transfer.F90.

◆ zpartition_transfer()

subroutine, public partition_transfer_oct_m::zpartition_transfer ( type(partition_transfer_t), intent(in)  this,
complex(real64), dimension(:), intent(in), contiguous  f_in,
complex(real64), dimension(:), intent(out), contiguous  f_out 
)

Definition at line 556 of file partition_transfer.F90.