32 type(mpi_grp_t) :: mpi_grp
35 logical :: double_sided_comms
53 class(ring_pattern_t),
intent(inout) :: this
54 type(mpi_grp_t),
intent(in) :: mpi_grp
55 logical,
intent(in) :: double_sided_comms
59 this%double_sided_comms = double_sided_comms
60 this%mpi_grp = mpi_grp
62 if (double_sided_comms)
then
63 this%nsteps = int((mpi_grp%size+2)/2)-1
65 this%nsteps = mpi_grp%size-1
69 write(
message(1),
'(a,i4,a)')
"Debug: The ring pattern will perform ", this%nsteps,
" steps."
77 integer pure function ring_pattern_get_nsteps(this)
78 class(ring_pattern_t),
intent(in) :: this
80 ring_pattern_get_nsteps = this%nsteps
102 integer pure function ring_pattern_get_rank_to(this, istep) result(rank_to)
103 class(ring_pattern_t),
intent(in) :: this
104 integer,
intent(in) :: istep
106 logical :: last_step, even_number_tasks, rank_to_first_half
108 rank_to = mod(this%mpi_grp%rank + istep, this%mpi_grp%size)
110 if (this%double_sided_comms)
then
111 last_step = istep == this%nsteps
112 even_number_tasks = mod(this%mpi_grp%size, 2) == 0
113 rank_to_first_half = rank_to < this%mpi_grp%size/2
115 if (last_step .and. even_number_tasks .and. rank_to_first_half)
then
122 integer pure function ring_pattern_get_rank_from(this, istep) result(rank_fr)
123 class(ring_pattern_t),
intent(in) :: this
124 integer,
intent(in) :: istep
126 logical :: last_step, even_number_tasks, rank_first_half
128 rank_fr = this%mpi_grp%rank - istep
129 if (rank_fr < 0) rank_fr = this%mpi_grp%size + rank_fr
130 rank_fr = mod(rank_fr, this%mpi_grp%size)
132 if (this%double_sided_comms)
then
133 last_step = istep == this%nsteps
134 even_number_tasks = mod(this%mpi_grp%size, 2) == 0
135 rank_first_half = this%mpi_grp%rank < this%mpi_grp%size/2
136 if (last_step .and. even_number_tasks .and. rank_first_half)
then
type(debug_t), save, public debug
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
This module is an helper to perform ring-pattern communications among all states.
integer pure function ring_pattern_get_nsteps(this)
Returns the total number of steps in the ring.
subroutine ring_pattern_start(this, mpi_grp, double_sided_comms)
Starts the ring pattern scheme.
integer pure function ring_pattern_get_rank_from(this, istep)
Returns the rank from which we receive the information.
integer pure function ring_pattern_get_rank_to(this, istep)
Returns the rank where to send the information.