(git:d18deda)
Loading...
Searching...
No Matches
cp_blacs_env.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2025 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief methods related to the blacs parallel environment
10!> \par History
11!> 08.2002 created [fawzi]
12!> 02.2004 modified to associate a blacs_env with a given para_env
13!> \author Fawzi Mohamed
14! **************************************************************************************************
18 USE kinds, ONLY: dp
19 USE machine, ONLY: m_flush
20 USE mathlib, ONLY: gcd
23#include "../base/base_uses.f90"
24
25 IMPLICIT NONE
26 PRIVATE
27
28 LOGICAL, PRIVATE, PARAMETER :: debug_this_module = .true.
29 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'cp_blacs_env'
30
31 ! Blacs type of distribution
32 INTEGER, PARAMETER, PUBLIC :: blacs_grid_square = 1, &
33 blacs_grid_row = 2, &
35
36 PUBLIC :: cp_blacs_env_type
38
39! **************************************************************************************************
40!> \brief represent a blacs multidimensional parallel environment
41!> (for the mpi corrispective see cp_paratypes/mp_para_cart_type)
42!> \param ref_count the reference count, when it is zero this object gets
43!> deallocated
44!> \param my_pid process id of the actual processor
45!> \param n_pid number of process ids
46!> \param the para_env associated (and compatible) with this blacs_env
47!> \param blacs2mpi: maps mepos(1)-mepos(2) of blacs to its mpi rank
48!> \param mpi2blacs(i,rank): maps the mpi rank to the mepos(i)
49!> \par History
50!> 08.2002 created [fawzi]
51!> \author Fawzi Mohamed
52! **************************************************************************************************
54 INTEGER :: my_pid = -1, n_pid = -1, ref_count = -1
55 TYPE(mp_para_env_type), POINTER :: para_env => null()
56 INTEGER, DIMENSION(:, :), POINTER :: blacs2mpi => null()
57 INTEGER, DIMENSION(:, :), POINTER :: mpi2blacs => null()
58 LOGICAL :: repeatable = .false.
59 CONTAINS
60 PROCEDURE, PUBLIC, pass, non_overridable :: create => cp_blacs_env_create_low
61 PROCEDURE, PUBLIC, pass, non_overridable :: retain => cp_blacs_env_retain
62 PROCEDURE, PUBLIC, pass, non_overridable :: release => cp_blacs_env_release_low
63 PROCEDURE, PUBLIC, pass, non_overridable :: get => get_blacs_info
64 PROCEDURE, PUBLIC, pass, non_overridable :: write => cp_blacs_env_write
65 END TYPE cp_blacs_env_type
66
67!***
68CONTAINS
69
70! **************************************************************************************************
71!> \brief Return informations about the specified BLACS context.
72!> \param blacs_env ...
73!> \param my_process_row ...
74!> \param my_process_column ...
75!> \param my_process_number ...
76!> \param number_of_process_rows ...
77!> \param number_of_process_columns ...
78!> \param number_of_processes ...
79!> \param para_env ...
80!> \param blacs2mpi ...
81!> \param mpi2blacs ...
82!> \date 19.06.2001
83!> \par History
84!> MM.YYYY moved here from qs_blacs (Joost VandeVondele)
85!> \author Matthias Krack
86!> \version 1.0
87! **************************************************************************************************
88 SUBROUTINE get_blacs_info(blacs_env, my_process_row, my_process_column, &
89 my_process_number, number_of_process_rows, &
90 number_of_process_columns, number_of_processes, &
91 para_env, blacs2mpi, mpi2blacs)
92 CLASS(cp_blacs_env_type), INTENT(IN) :: blacs_env
93 INTEGER, INTENT(OUT), OPTIONAL :: my_process_row, my_process_column, my_process_number, &
94 number_of_process_rows, number_of_process_columns, number_of_processes
95 TYPE(mp_para_env_type), OPTIONAL, POINTER :: para_env
96 INTEGER, DIMENSION(:, :), OPTIONAL, POINTER :: blacs2mpi, mpi2blacs
97
98 IF (PRESENT(my_process_row)) my_process_row = blacs_env%mepos(1)
99 IF (PRESENT(my_process_column)) my_process_column = blacs_env%mepos(2)
100 IF (PRESENT(my_process_number)) my_process_number = blacs_env%my_pid
101 IF (PRESENT(number_of_process_rows)) number_of_process_rows = blacs_env%num_pe(1)
102 IF (PRESENT(number_of_process_columns)) number_of_process_columns = blacs_env%num_pe(2)
103 IF (PRESENT(number_of_processes)) number_of_processes = blacs_env%n_pid
104 IF (PRESENT(para_env)) para_env => blacs_env%para_env
105 IF (PRESENT(blacs2mpi)) blacs2mpi => blacs_env%blacs2mpi
106 IF (PRESENT(mpi2blacs)) mpi2blacs => blacs_env%mpi2blacs
107
108 END SUBROUTINE get_blacs_info
109
110! **************************************************************************************************
111!> \brief allocates and initializes a type that represent a blacs context
112!> \param blacs_env the type to initialize
113!> \param para_env the para_env for which a blacs env should be created
114!> \param blacs_grid_layout ...
115!> \param blacs_repeatable ...
116!> \param row_major ...
117!> \param grid_2d ...
118!> \par History
119!> 08.2002 created [fawzi]
120!> \author Fawzi Mohamed
121! **************************************************************************************************
122 SUBROUTINE cp_blacs_env_create(blacs_env, para_env, blacs_grid_layout, blacs_repeatable, row_major, grid_2d)
123 TYPE(cp_blacs_env_type), INTENT(OUT), POINTER :: blacs_env
124 TYPE(mp_para_env_type), INTENT(INOUT), TARGET :: para_env
125 INTEGER, INTENT(IN), OPTIONAL :: blacs_grid_layout
126 LOGICAL, INTENT(IN), OPTIONAL :: blacs_repeatable, row_major
127 INTEGER, DIMENSION(:), INTENT(IN), OPTIONAL :: grid_2d
128
129 ALLOCATE (blacs_env)
130 CALL blacs_env%create(para_env, blacs_grid_layout, blacs_repeatable, row_major, grid_2d)
131
132 END SUBROUTINE
133
134! **************************************************************************************************
135!> \brief allocates and initializes a type that represent a blacs context
136!> \param blacs_env the type to initialize
137!> \param para_env the para_env for which a blacs env should be created
138!> \param blacs_grid_layout ...
139!> \param blacs_repeatable ...
140!> \param row_major ...
141!> \param grid_2d ...
142!> \par History
143!> 08.2002 created [fawzi]
144!> \author Fawzi Mohamed
145! **************************************************************************************************
146 SUBROUTINE cp_blacs_env_create_low(blacs_env, para_env, blacs_grid_layout, blacs_repeatable, row_major, grid_2d)
147 CLASS(cp_blacs_env_type), INTENT(OUT) :: blacs_env
148 TYPE(mp_para_env_type), TARGET, INTENT(INOUT) :: para_env
149 INTEGER, INTENT(IN), OPTIONAL :: blacs_grid_layout
150 LOGICAL, INTENT(IN), OPTIONAL :: blacs_repeatable, row_major
151 INTEGER, DIMENSION(:), INTENT(IN), &
152 OPTIONAL :: grid_2d
153
154 INTEGER :: ipcol, iprow
155#if defined(__parallel)
156 INTEGER :: gcd_max, ipe, jpe, &
157 my_blacs_grid_layout, &
158 npcol, npe, nprow
159 LOGICAL :: my_blacs_repeatable, &
160 my_row_major
161#endif
162
163#ifdef __parallel
164 ! get the number of cpus for this blacs grid
165 nprow = 1
166 npcol = 1
167 npe = para_env%num_pe
168 ! get the layout of this grid
169
170 IF (PRESENT(grid_2d)) THEN
171 nprow = grid_2d(1)
172 npcol = grid_2d(2)
173 END IF
174
175 IF (nprow*npcol .NE. npe) THEN
176 ! hard code for the time being the grid layout
177 my_blacs_grid_layout = blacs_grid_square
178 IF (PRESENT(blacs_grid_layout)) my_blacs_grid_layout = blacs_grid_layout
179 ! XXXXXX
180 SELECT CASE (my_blacs_grid_layout)
181 CASE (blacs_grid_square)
182 ! make the grid as 'square' as possible, where square is defined as nprow and npcol
183 ! having the largest possible gcd
184 gcd_max = -1
185 DO ipe = 1, ceiling(sqrt(real(npe, dp)))
186 jpe = npe/ipe
187 IF (ipe*jpe .NE. npe) cycle
188 IF (gcd(ipe, jpe) >= gcd_max) THEN
189 nprow = ipe
190 npcol = jpe
191 gcd_max = gcd(ipe, jpe)
192 END IF
193 END DO
194 CASE (blacs_grid_row)
195 nprow = 1
196 npcol = npe
197 CASE (blacs_grid_col)
198 nprow = npe
199 npcol = 1
200 END SELECT
201 END IF
202
203 my_row_major = .true.
204 IF (PRESENT(row_major)) my_row_major = row_major
205 IF (my_row_major) THEN
206 CALL blacs_env%gridinit(para_env, "Row-major", nprow, npcol)
207 ELSE
208 CALL blacs_env%gridinit(para_env, "Col-major", nprow, npcol)
209 END IF
210
211 ! We set the components of blacs_env here such that we can still use INTENT(OUT) with gridinit
212 blacs_env%my_pid = para_env%mepos
213 blacs_env%n_pid = para_env%num_pe
214 blacs_env%ref_count = 1
215
216 my_blacs_repeatable = .false.
217 IF (PRESENT(blacs_repeatable)) my_blacs_repeatable = blacs_repeatable
218 blacs_env%repeatable = my_blacs_repeatable
219 IF (blacs_env%repeatable) CALL blacs_env%set(15, 1)
220
221#else
222 ! In serial mode, we just have to setup the object
223 CALL blacs_env%gridinit(para_env, "Row-major", 1, 1)
224
225 blacs_env%ref_count = 1
226 blacs_env%my_pid = 0
227 blacs_env%n_pid = 1
228 mark_used(blacs_grid_layout)
229 mark_used(blacs_repeatable)
230 mark_used(grid_2d)
231 mark_used(row_major)
232#endif
233
234 CALL para_env%retain()
235 blacs_env%para_env => para_env
236
237 ! generate the mappings blacs2mpi and mpi2blacs
238 ALLOCATE (blacs_env%blacs2mpi(0:blacs_env%num_pe(1) - 1, 0:blacs_env%num_pe(2) - 1))
239 blacs_env%blacs2mpi = 0
240 blacs_env%blacs2mpi(blacs_env%mepos(1), blacs_env%mepos(2)) = para_env%mepos
241 CALL para_env%sum(blacs_env%blacs2mpi)
242 ALLOCATE (blacs_env%mpi2blacs(2, 0:para_env%num_pe - 1))
243 blacs_env%mpi2blacs = -1
244 DO ipcol = 0, blacs_env%num_pe(2) - 1
245 DO iprow = 0, blacs_env%num_pe(1) - 1
246 blacs_env%mpi2blacs(1, blacs_env%blacs2mpi(iprow, ipcol)) = iprow
247 blacs_env%mpi2blacs(2, blacs_env%blacs2mpi(iprow, ipcol)) = ipcol
248 END DO
249 END DO
250 END SUBROUTINE cp_blacs_env_create_low
251
252! **************************************************************************************************
253!> \brief retains the given blacs env
254!> \param blacs_env the blacs env to retain
255!> \par History
256!> 08.2002 created [fawzi]
257!> \author Fawzi Mohamed
258! **************************************************************************************************
259 SUBROUTINE cp_blacs_env_retain(blacs_env)
260 CLASS(cp_blacs_env_type), INTENT(INOUT) :: blacs_env
261
262 cpassert(blacs_env%ref_count > 0)
263 blacs_env%ref_count = blacs_env%ref_count + 1
264 END SUBROUTINE cp_blacs_env_retain
265
266! **************************************************************************************************
267!> \brief releases the given blacs_env
268!> \param blacs_env the blacs env to release
269!> \par History
270!> 08.2002 created [fawzi]
271!> \author Fawzi Mohamed
272! **************************************************************************************************
273 SUBROUTINE cp_blacs_env_release(blacs_env)
274 TYPE(cp_blacs_env_type), POINTER :: blacs_env
275
276 IF (ASSOCIATED(blacs_env)) THEN
277 cpassert(blacs_env%ref_count > 0)
278 blacs_env%ref_count = blacs_env%ref_count - 1
279 IF (blacs_env%ref_count < 1) THEN
280 CALL blacs_env%release()
281 DEALLOCATE (blacs_env)
282 END IF
283 END IF
284 NULLIFY (blacs_env)
285 END SUBROUTINE cp_blacs_env_release
286
287! **************************************************************************************************
288!> \brief releases the given blacs_env
289!> \param blacs_env the blacs env to release
290!> \par History
291!> 08.2002 created [fawzi]
292!> \author Fawzi Mohamed
293! **************************************************************************************************
294 SUBROUTINE cp_blacs_env_release_low(blacs_env)
295 CLASS(cp_blacs_env_type), INTENT(INOUT) :: blacs_env
296
297 CALL blacs_env%gridexit()
298 CALL mp_para_env_release(blacs_env%para_env)
299 DEALLOCATE (blacs_env%mpi2blacs)
300 DEALLOCATE (blacs_env%blacs2mpi)
301
302 END SUBROUTINE cp_blacs_env_release_low
303
304! **************************************************************************************************
305!> \brief writes the description of the given blacs env
306!> \param blacs_env the blacs environment to write
307!> \param unit_nr the unit number where to write the description of the
308!> blacs environment
309!> \par History
310!> 08.2002 created [fawzi]
311!> \author Fawzi Mohamed
312! **************************************************************************************************
313 SUBROUTINE cp_blacs_env_write(blacs_env, unit_nr)
314 CLASS(cp_blacs_env_type), INTENT(IN) :: blacs_env
315 INTEGER, INTENT(in) :: unit_nr
316
317 WRITE (unit=unit_nr, fmt="(' group=',i10,', ref_count=',i10,',')") &
318 blacs_env%get_handle(), blacs_env%ref_count
319 WRITE (unit=unit_nr, fmt="(' mepos=(',i8,',',i8,'),')") &
320 blacs_env%mepos(1), blacs_env%mepos(2)
321 WRITE (unit=unit_nr, fmt="(' num_pe=(',i8,',',i8,'),')") &
322 blacs_env%num_pe(1), blacs_env%num_pe(2)
323 IF (ASSOCIATED(blacs_env%blacs2mpi)) THEN
324 WRITE (unit=unit_nr, fmt="(' blacs2mpi=')", advance="no")
325 CALL cp_2d_i_write(blacs_env%blacs2mpi, unit_nr=unit_nr)
326 ELSE
327 WRITE (unit=unit_nr, fmt="(' blacs2mpi=*null*')")
328 END IF
329 IF (ASSOCIATED(blacs_env%para_env)) THEN
330 WRITE (unit=unit_nr, fmt="(' para_env=<cp_para_env id=',i6,'>,')") &
331 blacs_env%para_env%get_handle()
332 ELSE
333 WRITE (unit=unit_nr, fmt="(' para_env=*null*')")
334 END IF
335 WRITE (unit=unit_nr, fmt="(' my_pid=',i10,', n_pid=',i10,' }')") &
336 blacs_env%my_pid, blacs_env%n_pid
337 CALL m_flush(unit_nr)
338 END SUBROUTINE cp_blacs_env_write
339
340END MODULE cp_blacs_env
various utilities that regard array of different kinds: output, allocation,... maybe it is not a good...
subroutine, public cp_2d_i_write(array, unit_nr, el_format)
writes an array to the given unit
methods related to the blacs parallel environment
integer, parameter, public blacs_grid_row
integer, parameter, public blacs_grid_col
integer, parameter, public blacs_grid_square
subroutine, public cp_blacs_env_release(blacs_env)
releases the given blacs_env
subroutine get_blacs_info(blacs_env, my_process_row, my_process_column, my_process_number, number_of_process_rows, number_of_process_columns, number_of_processes, para_env, blacs2mpi, mpi2blacs)
Return informations about the specified BLACS context.
subroutine, public cp_blacs_env_create(blacs_env, para_env, blacs_grid_layout, blacs_repeatable, row_major, grid_2d)
allocates and initializes a type that represent a blacs context
wrappers for the actual blacs calls. all functionality needed in the code should actually be provide ...
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
Machine interface based on Fortran 2003 and POSIX.
Definition machine.F:17
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
Definition machine.F:130
Collection of simple mathematical functions and subroutines.
Definition mathlib.F:15
elemental integer function, public gcd(a, b)
computes the greatest common divisor of two number
Definition mathlib.F:1280
Interface to the message passing library MPI.
subroutine, public mp_para_env_release(para_env)
releases the para object (to be called when you don't want anymore the shared copy of this object)
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
stores all the informations relevant to an mpi environment