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 mesh_oct_m
37 use mpi_oct_m
43
44 implicit none
45
46 private
47
48 public :: &
50
52 private
53
54 integer :: task_from
55 integer :: task_to
56
57 integer :: nbatch_to_receive
58 integer :: nbatch_to_send
59 integer :: n_comms
60
61 integer :: nblock_to_receive
62 integer :: nblock_to_send
63 contains
68 procedure :: alloc_receive_batch => states_elec_all_to_all_communications_alloc_receive_batch
69 procedure :: get_send_indices => states_elec_all_to_all_communications_get_send_indices
70 procedure :: get_receive_indices => states_elec_all_to_all_communications_get_receive_indices
72
73contains
74
76 subroutine states_elec_all_to_all_communications_start(this, st, task_from, task_to)
77 class(states_elec_all_to_all_communications_t), intent(inout) :: this
78 type(states_elec_t), intent(in) :: st
79 integer, intent(in) :: task_from, task_to
80
82
83 this%task_from = task_from
84 this%task_to = task_to
85
86 this%nbatch_to_receive = states_elec_all_to_all_communications_eval_nreceive(st, task_from, this%nblock_to_receive)
87 this%nbatch_to_send = states_elec_all_to_all_communications_eval_nsend(st, task_to, this%nblock_to_send)
88
89 !Number of communications
90 this%n_comms = max(this%nbatch_to_send, this%nbatch_to_receive)
91
94
96 integer function states_elec_all_to_all_communications_eval_nreceive(st, task_from, nblock_to_receive) result(nbatch_to_receive)
97 type(states_elec_t), intent(in) :: st
98 integer, intent(in) :: task_from
99 integer, intent(out) :: nblock_to_receive
100
101 integer :: st_start, st_end, kpt_start, kpt_end, ib
102
104
105 nbatch_to_receive = 0
106 nblock_to_receive = 0
107
108 !What we receive for the wfn
109 if (task_from > -1) then
110 st_start = st%st_kpt_task(task_from, 1)
111 st_end = st%st_kpt_task(task_from, 2)
112 kpt_start = st%st_kpt_task(task_from, 3)
113 kpt_end = st%st_kpt_task(task_from, 4)
114
115 nblock_to_receive = 0
116 do ib = 1, st%group%nblocks
117 if (st%group%block_range(ib, 1) >= st_start .and. st%group%block_range(ib, 2) <= st_end) then
118 nblock_to_receive = nblock_to_receive + 1
119 end if
120 end do
121 nbatch_to_receive = nblock_to_receive * (kpt_end-kpt_start+1)
122 end if
123
124 if (debug%info) then
125 write(message(1), '(a,i2,a,i2,a,i2)') 'Task ', st%st_kpt_mpi_grp%rank, ' will receive ', &
126 nbatch_to_receive, ' batches from task ', task_from
127 call messages_info(1, all_nodes=.true.)
128 end if
129
132
134 integer function states_elec_all_to_all_communications_eval_nsend(st, task_to, nblock_to_send) result(nbatch_to_send)
135 type(states_elec_t), intent(in) :: st
136 integer, intent(in) :: task_to
137 integer, intent(out) :: nblock_to_send
138
140
141 nbatch_to_send = 0
142 nblock_to_send = 0
143
144 if (task_to > -1) then
145 nblock_to_send = (st%group%block_end-st%group%block_start+1)
146 nbatch_to_send = nblock_to_send*(st%d%kpt%end-st%d%kpt%start+1)
147 end if
149 if (debug%info) then
150 write(message(1), '(a,i2,a,i2,a,i2)') 'Task ', st%st_kpt_mpi_grp%rank, ' will send ', nbatch_to_send, &
151 ' batches to task ', task_to
152 call messages_info(1, all_nodes=.true.)
153 end if
159 integer pure function states_elec_all_to_all_communications_get_ncom(this) result(n_comms)
162 n_comms = this%n_comms
164
166 integer pure function states_elec_all_to_all_communications_get_nreceive(this) result(nbatch_to_receive)
167 class(states_elec_all_to_all_communications_t), intent(in) :: this
168
169 nbatch_to_receive = this%nbatch_to_receive
171
173 integer pure function states_elec_all_to_all_communications_get_nsend(this) result(nbatch_to_send)
174 class(states_elec_all_to_all_communications_t), intent(in) :: this
175
176 nbatch_to_send = this%nbatch_to_send
178
180 subroutine states_elec_all_to_all_communications_alloc_receive_batch(this, st, icom, np, psib)
181 class(states_elec_all_to_all_communications_t), intent(in) :: this
182 type(states_elec_t), intent(in) :: st
183 integer, intent(in) :: icom
184 integer, intent(in) :: np
185 class(wfs_elec_t), pointer, intent(out) :: psib
186
187 integer :: block_id, ib, ik
188
190
191 ! Given the icom, returns the id of the block of state communicated
192 block_id = mod(icom-1, this%nblock_to_receive)+1
193 ik = int((icom-block_id)/this%nblock_to_receive) + st%st_kpt_task(this%task_from, 3)
194 ib = block_id - 1 + st%group%iblock(st%st_kpt_task(this%task_from, 1))
195
196 if (debug%info) then
197 write(message(1), '(a,i2,a,i2,a,i2)') 'Task ', st%st_kpt_mpi_grp%rank, ' allocates memory for block ', &
198 ib, ' and k-point ', ik
199 call messages_info(1, all_nodes=.true.)
200 end if
201
202 call states_elec_parallel_allocate_batch(st, psib, np, ib, ik, packed=.true.)
203
206
208 subroutine states_elec_all_to_all_communications_get_send_indices(this, st, icom, ib, ik)
209 class(states_elec_all_to_all_communications_t), intent(in) :: this
210 type(states_elec_t), intent(in) :: st
211 integer, intent(in) :: icom
212 integer, intent(out) :: ib
213 integer, intent(out) :: ik
214
216
217 ! Given the icom, returns the id of the block of state communicated
218 ib = mod(icom-1, this%nblock_to_send) + 1
219 ik = int((icom-ib)/this%nblock_to_send) + st%d%kpt%start
220 ib = ib - 1 + st%group%iblock(st%st_start)
221
222 if (debug%info) then
223 write(message(1), '(a,i2,a,i2,a,i2)') 'Task ', st%st_kpt_mpi_grp%rank, ' will send the block ', &
224 ib, ' with k-point ', ik
225 call messages_info(1, all_nodes=.true.)
226 end if
230
232 subroutine states_elec_all_to_all_communications_get_receive_indices(this, st, icom, ib, ik)
233 class(states_elec_all_to_all_communications_t), intent(in) :: this
234 type(states_elec_t), intent(in) :: st
235 integer, intent(in) :: icom
236 integer, intent(out) :: ib
237 integer, intent(out) :: ik
238
240
241 ! Given the icom, returns the id of the block of state communicated
242 ib = mod(icom-1, this%nblock_to_receive)+1
243 ik = int((icom-ib)/this%nblock_to_receive) + st%st_kpt_task(this%task_from, 3)
244 ib = ib - 1 + st%group%iblock(st%st_kpt_task(this%task_from, 1))
245
246 if (debug%info) then
247 write(message(1), '(a,i2,a,i2,a,i2)') 'Task ', st%st_kpt_mpi_grp%rank, ' will receive the block ', &
248 ib, ' with k-point ', ik
249 call messages_info(1, all_nodes=.true.)
250 end if
251
254
256
257!! Local Variables:
258!! mode: f90
259!! coding: utf-8
260!! End:
This module implements batches of mesh functions.
Definition: batch.F90:133
type(debug_t), save, public debug
Definition: debug.F90:154
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
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.
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 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.
This module provides routines for communicating states when using states parallelization.
int true(void)