58 integer :: nbatch_to_receive
59 integer :: nbatch_to_send
62 integer :: nblock_to_receive
63 integer :: nblock_to_send
65 type(MPI_Request),
allocatable :: send_req(:)
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
92 this%task_from = task_from
93 this%task_to = task_to
99 this%n_comms = max(this%nbatch_to_send, this%nbatch_to_receive)
107 type(states_elec_t),
intent(in) :: st
108 integer,
intent(in) :: task_from
109 integer,
intent(out) :: nblock_to_receive
111 integer :: st_start, st_end, kpt_start, kpt_end, ib
115 nbatch_to_receive = 0
116 nblock_to_receive = 0
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
131 nbatch_to_receive = nblock_to_receive * (kpt_end-kpt_start+1)
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
144 type(states_elec_t),
intent(in) :: st
145 integer,
intent(in) :: task_to
146 integer,
intent(out) :: nblock_to_send
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)
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
170 n_comms = this%n_comms
178 nbatch_to_receive = this%nbatch_to_receive
183 integer pure function states_elec_all_to_all_communications_get_nsend(this) result(nbatch_to_send)
186 nbatch_to_send = this%nbatch_to_send
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
198 integer :: block_id, ib, ik
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))
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.)
211 call states_elec_parallel_allocate_batch(st, psib, np, ib, ik, packed=.
true.)
220 type(states_elec_t),
intent(in) :: st
221 integer,
intent(in) :: icom
222 integer,
intent(out) :: ib
223 integer,
intent(out) :: ik
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)
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.)
243 type(states_elec_t),
intent(in) :: st
244 integer,
intent(in) :: icom
245 integer,
intent(out) :: ib
246 integer,
intent(out) :: ik
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))
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.)
268 type(states_elec_t),
intent(in) :: st
271 call profiling_in(
"ALL_TO_ALL_COMM")
273 if (
allocated(this%send_req))
then
275 call st%st_kpt_mpi_grp%wait(this%nbatch_to_send, this%send_req)
277 safe_deallocate_a(this%send_req)
281 call profiling_out(
"ALL_TO_ALL_COMM")
289#include "states_elec_all_to_all_communications_inc.F90"
292#include "complex.F90"
293#include "states_elec_all_to_all_communications_inc.F90"
This module implements batches of mesh functions.
This module is intended to contain "only mathematical" functions and procedures.
This module defines the meshes, which are used in Octopus.
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
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.