Octopus
states_elec_all_to_all_communications.F90
Go to the documentation of this file.
1!! Copyright (C) 2023 N. Tancogne-Dejean
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
21
32 use batch_oct_m
33 use debug_oct_m
34 use global_oct_m
35 use math_oct_m
36 use mesh_oct_m
38 use mpi_oct_m
44
45 implicit none
46
47 private
48
49 public :: &
51
53 private
54
55 integer :: task_from
56 integer :: task_to
57
58 integer :: nbatch_to_receive
59 integer :: nbatch_to_send
60 integer :: n_comms
61
62 integer :: nblock_to_receive
63 integer :: nblock_to_send
64
65 type(MPI_Request), allocatable :: send_req(:)
66 contains
71 procedure :: alloc_receive_batch => states_elec_all_to_all_communications_alloc_receive_batch
72 procedure :: get_send_indices => states_elec_all_to_all_communications_get_send_indices
73 procedure :: get_receive_indices => states_elec_all_to_all_communications_get_receive_indices
74 procedure :: dpost_all_mpi_isend => dstates_elec_all_to_all_communications_post_all_mpi_isend
75 procedure :: zpost_all_mpi_isend => zstates_elec_all_to_all_communications_post_all_mpi_isend
76 procedure :: dmpi_recv_batch => dstates_elec_all_to_all_communications_mpi_recv_batch
77 procedure :: zmpi_recv_batch => zstates_elec_all_to_all_communications_mpi_recv_batch
78 procedure :: wait_all_isend => states_elec_all_to_all_communications_wait_all_isend
80
81contains
82
83 !------------------------------------------------------------
85 subroutine states_elec_all_to_all_communications_start(this, st, task_from, task_to)
86 class(states_elec_all_to_all_communications_t), intent(inout) :: this
87 type(states_elec_t), intent(in) :: st
88 integer, intent(in) :: task_from, task_to
89
91
92 this%task_from = task_from
93 this%task_to = task_to
94
95 this%nbatch_to_receive = states_elec_all_to_all_communications_eval_nreceive(st, task_from, this%nblock_to_receive)
96 this%nbatch_to_send = states_elec_all_to_all_communications_eval_nsend(st, task_to, this%nblock_to_send)
97
98 !Number of communications
99 this%n_comms = max(this%nbatch_to_send, this%nbatch_to_receive)
100
103
104 !------------------------------------------------------------
106 integer function states_elec_all_to_all_communications_eval_nreceive(st, task_from, nblock_to_receive) result(nbatch_to_receive)
107 type(states_elec_t), intent(in) :: st
108 integer, intent(in) :: task_from
109 integer, intent(out) :: nblock_to_receive
110
111 integer :: st_start, st_end, kpt_start, kpt_end, ib
112
114
115 nbatch_to_receive = 0
116 nblock_to_receive = 0
117
118 !What we receive for the wfn
119 if (task_from > -1) then
120 st_start = st%st_kpt_task(task_from, 1)
121 st_end = st%st_kpt_task(task_from, 2)
122 kpt_start = st%st_kpt_task(task_from, 3)
123 kpt_end = st%st_kpt_task(task_from, 4)
125 nblock_to_receive = 0
126 do ib = 1, st%group%nblocks
127 if (st%group%block_range(ib, 1) >= st_start .and. st%group%block_range(ib, 2) <= st_end) then
128 nblock_to_receive = nblock_to_receive + 1
129 end if
130 end do
131 nbatch_to_receive = nblock_to_receive * (kpt_end-kpt_start+1)
132 end if
133
134 write(message(1), '(a,i2,a,i2,a,i2)') 'Debug: Task ', st%st_kpt_mpi_grp%rank, ' will receive ', &
135 nbatch_to_receive, ' batches from task ', task_from
136 call messages_info(1, all_nodes=.true., debug_only=.true.)
137
140
141 !------------------------------------------------------------
143 integer function states_elec_all_to_all_communications_eval_nsend(st, task_to, nblock_to_send) result(nbatch_to_send)
144 type(states_elec_t), intent(in) :: st
145 integer, intent(in) :: task_to
146 integer, intent(out) :: nblock_to_send
147
150 nbatch_to_send = 0
151 nblock_to_send = 0
153 if (task_to > -1) then
154 nblock_to_send = (st%group%block_end-st%group%block_start+1)
155 nbatch_to_send = nblock_to_send*(st%d%kpt%end-st%d%kpt%start+1)
156 end if
157
158 write(message(1), '(a,i2,a,i2,a,i2)') 'Debug: Task ', st%st_kpt_mpi_grp%rank, ' will send ', nbatch_to_send, &
159 ' batches to task ', task_to
160 call messages_info(1, all_nodes=.true., debug_only=.true.)
165 !------------------------------------------------------------
167 integer pure function states_elec_all_to_all_communications_get_ncom(this) result(n_comms)
170 n_comms = this%n_comms
172
173 !------------------------------------------------------------
175 integer pure function states_elec_all_to_all_communications_get_nreceive(this) result(nbatch_to_receive)
176 class(states_elec_all_to_all_communications_t), intent(in) :: this
177
178 nbatch_to_receive = this%nbatch_to_receive
180
181 !------------------------------------------------------------
183 integer pure function states_elec_all_to_all_communications_get_nsend(this) result(nbatch_to_send)
184 class(states_elec_all_to_all_communications_t), intent(in) :: this
185
186 nbatch_to_send = this%nbatch_to_send
188
189 !------------------------------------------------------------
191 subroutine states_elec_all_to_all_communications_alloc_receive_batch(this, st, icom, np, psib)
192 class(states_elec_all_to_all_communications_t), intent(in) :: this
193 type(states_elec_t), intent(in) :: st
194 integer, intent(in) :: icom
195 integer, intent(in) :: np
196 type(wfs_elec_t), intent(out) :: psib
197
198 integer :: block_id, ib, ik
201
202 ! Given the icom, returns the id of the block of state communicated
203 block_id = mod(icom-1, this%nblock_to_receive)+1
204 ik = int((icom-block_id)/this%nblock_to_receive) + st%st_kpt_task(this%task_from, 3)
205 ib = block_id - 1 + st%group%iblock(st%st_kpt_task(this%task_from, 1))
206
207 write(message(1), '(a,i2,a,i2,a,i2)') 'Debug: Task ', st%st_kpt_mpi_grp%rank, ' allocates memory for block ', &
208 ib, ' and k-point ', ik
209 call messages_info(1, all_nodes=.true., debug_only=.true.)
210
211 call states_elec_parallel_allocate_batch(st, psib, np, ib, ik, packed=.true.)
212
215
216 !------------------------------------------------------------
218 subroutine states_elec_all_to_all_communications_get_send_indices(this, st, icom, ib, ik)
219 class(states_elec_all_to_all_communications_t), intent(in) :: this
220 type(states_elec_t), intent(in) :: st
221 integer, intent(in) :: icom
222 integer, intent(out) :: ib
223 integer, intent(out) :: ik
224
226
227 ! Given the icom, returns the id of the block of state communicated
228 ib = mod(icom-1, this%nblock_to_send) + 1
229 ik = int((icom-ib)/this%nblock_to_send) + st%d%kpt%start
230 ib = ib - 1 + st%group%iblock(st%st_start)
231
232 write(message(1), '(a,i2,a,i2,a,i2)') 'Debug: Task ', st%st_kpt_mpi_grp%rank, ' will send the block ', &
233 ib, ' with k-point ', ik
234 call messages_info(1, all_nodes=.true., debug_only=.true.)
235
238
239 !------------------------------------------------------------
241 subroutine states_elec_all_to_all_communications_get_receive_indices(this, st, icom, ib, ik)
242 class(states_elec_all_to_all_communications_t), intent(in) :: this
243 type(states_elec_t), intent(in) :: st
244 integer, intent(in) :: icom
245 integer, intent(out) :: ib
246 integer, intent(out) :: ik
247
249
250 ! Given the icom, returns the id of the block of state communicated
251 ib = mod(icom-1, this%nblock_to_receive)+1
252 ik = int((icom-ib)/this%nblock_to_receive) + st%st_kpt_task(this%task_from, 3)
253 ib = ib - 1 + st%group%iblock(st%st_kpt_task(this%task_from, 1))
254
255 if (debug%info) then
256 write(message(1), '(a,i2,a,i2,a,i2)') 'Task ', st%st_kpt_mpi_grp%rank, ' will receive the block ', &
257 ib, ' with k-point ', ik
258 call messages_info(1, all_nodes=.true.)
259 end if
263
264 !------------------------------------------------------------
267 class(states_elec_all_to_all_communications_t), intent(inout) :: this
268 type(states_elec_t), intent(in) :: st
269
271 call profiling_in("ALL_TO_ALL_COMM")
272
273 if (allocated(this%send_req)) then
274
275 call st%st_kpt_mpi_grp%wait(this%nbatch_to_send, this%send_req)
277 safe_deallocate_a(this%send_req)
278
279 end if
280
281 call profiling_out("ALL_TO_ALL_COMM")
285
286
287#include "undef.F90"
288#include "real.F90"
289#include "states_elec_all_to_all_communications_inc.F90"
290
291#include "undef.F90"
292#include "complex.F90"
293#include "states_elec_all_to_all_communications_inc.F90"
294#include "undef.F90"
295
297
298!! Local Variables:
299!! mode: f90
300!! coding: utf-8
301!! End:
This module implements batches of mesh functions.
Definition: batch.F90:133
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:624
This module provides routines for communicating all batches in a ring-pattern scheme.
subroutine states_elec_all_to_all_communications_get_receive_indices(this, st, icom, ib, ik)
Given the icom step, returns the block and k-point indices to be received.
integer function states_elec_all_to_all_communications_eval_nsend(st, task_to, nblock_to_send)
How many batches we will send from task_send.
subroutine states_elec_all_to_all_communications_wait_all_isend(this, st)
Do a MPI waitall for the isend requests.
subroutine zstates_elec_all_to_all_communications_post_all_mpi_isend(this, st, np, node_to)
Post all isend commands for all batches of a given task.
subroutine dstates_elec_all_to_all_communications_mpi_recv_batch(this, st, np, node_fr, icom, psib_receiv)
Allocate a batch and perform the MPI_Recv. On exist, the batch contains the received information.
integer pure function states_elec_all_to_all_communications_get_nreceive(this)
Returns the number of receiv calls.
subroutine states_elec_all_to_all_communications_alloc_receive_batch(this, st, icom, np, psib)
Given the icom step, allocate the receiv buffer (wfs_elec_t)
integer pure function states_elec_all_to_all_communications_get_ncom(this)
Returns the number of communications.
subroutine dstates_elec_all_to_all_communications_post_all_mpi_isend(this, st, np, node_to)
Post all isend commands for all batches of a given task.
subroutine states_elec_all_to_all_communications_get_send_indices(this, st, icom, ib, ik)
Given the icom step, returns the block and k-point indices to be sent.
integer function states_elec_all_to_all_communications_eval_nreceive(st, task_from, nblock_to_receive)
How many batches we will receive from task_from.
subroutine states_elec_all_to_all_communications_start(this, st, task_from, task_to)
Given a task to send to, and a task to receive from, initializes a states_elec_all_to_all_communicati...
integer pure function states_elec_all_to_all_communications_get_nsend(this)
Returns the number send calls.
subroutine zstates_elec_all_to_all_communications_mpi_recv_batch(this, st, np, node_fr, icom, psib_receiv)
Allocate a batch and perform the MPI_Recv. On exist, the batch contains the received information.
This module provides routines for communicating states when using states parallelization.
int true(void)