Octopus
ring_pattern.F90
Go to the documentation of this file.
1!! Copyright (C) 2023 N. Tancogne-Dejean
2!!
3!! This Source Code Form is subject to the terms of the Mozilla Public
4!! License, v. 2.0. If a copy of the MPL was not distributed with this
5!! file, You can obtain one at https://mozilla.org/MPL/2.0/.
6!!
7
8#include "global.h"
9
13!
15 use comm_oct_m
16 use debug_oct_m
17 use global_oct_m
19 use mpi_oct_m
21
22 implicit none
23
24 private
25
26 public :: &
28
30 private
31
32 type(mpi_grp_t) :: mpi_grp
33 integer :: nsteps
34 ! Each step might imply more than one communication step (e.g. multiple batches per rank).
35 logical :: double_sided_comms
36
37 contains
38 procedure :: start => ring_pattern_start
39 procedure :: get_nsteps => ring_pattern_get_nsteps
40 procedure :: get_rank_from => ring_pattern_get_rank_from
41 procedure :: get_rank_to => ring_pattern_get_rank_to
42 end type ring_pattern_t
43
44contains
45
52 subroutine ring_pattern_start(this, mpi_grp, 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
56
57 push_sub(ring_pattern_start)
58
59 this%double_sided_comms = double_sided_comms
60 this%mpi_grp = mpi_grp
61
62 if (this%mpi_grp%size <= 0) then
63 this%nsteps = 0
64 pop_sub(ring_pattern_start)
65 return
66 end if
67
68 if (double_sided_comms) then
69 this%nsteps = int((mpi_grp%size+2)/2)-1
70 else
71 this%nsteps = mpi_grp%size-1
72 end if
73
74 write(message(1), '(a,i6,a)') "Debug: The ring pattern will perform ", this%nsteps, " steps."
75 call messages_info(1, debug_only=.true.)
76
77 pop_sub(ring_pattern_start)
78 end subroutine ring_pattern_start
79
81 integer pure function ring_pattern_get_nsteps(this)
82 class(ring_pattern_t), intent(in) :: this
83
84 ring_pattern_get_nsteps = this%nsteps
85 end function ring_pattern_get_nsteps
86
106 integer pure function ring_pattern_get_rank_to(this, istep) result(rank_to)
107 class(ring_pattern_t), intent(in) :: this
108 integer, intent(in) :: istep
109
110 logical :: last_step, even_number_tasks, rank_to_first_half
111
112 rank_to = mod(this%mpi_grp%rank + istep, this%mpi_grp%size)
113
114 if (this%double_sided_comms) then
115 last_step = istep == this%nsteps
116 even_number_tasks = mod(this%mpi_grp%size, 2) == 0
117 rank_to_first_half = rank_to < this%mpi_grp%size/2
118
119 if (last_step .and. even_number_tasks .and. rank_to_first_half) then
120 rank_to = -1
121 end if
122 end if
123 end function ring_pattern_get_rank_to
124
126 integer pure function ring_pattern_get_rank_from(this, istep) result(rank_fr)
127 class(ring_pattern_t), intent(in) :: this
128 integer, intent(in) :: istep
129
130 logical :: last_step, even_number_tasks, rank_first_half
132 rank_fr = modulo(this%mpi_grp%rank - istep, this%mpi_grp%size)
134 if (this%double_sided_comms) then
135 last_step = istep == this%nsteps
136 even_number_tasks = mod(this%mpi_grp%size, 2) == 0
137 rank_first_half = this%mpi_grp%rank < this%mpi_grp%size/2
138 if (last_step .and. even_number_tasks .and. rank_first_half) then
139 rank_fr = -1
140 end if
141 end if
142 end function ring_pattern_get_rank_from
143
144end module ring_pattern_oct_m
146
147!! Local Variables:
148!! mode: f90
149!! coding: utf-8
150!! End:
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
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.
int true(void)