(git:374b731)
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-2024 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#ifndef __SCALAPACK
165 CALL cp_abort(__location__, &
166 "to USE the blacs environment "// &
167 "you need the blacs/scalapack library : recompile with -D__SCALAPACK (and link scalapack and blacs) ")
168#endif
169#endif
170
171#ifdef __SCALAPACK
172 ! get the number of cpus for this blacs grid
173 nprow = 1
174 npcol = 1
175 npe = para_env%num_pe
176 ! get the layout of this grid
177
178 IF (PRESENT(grid_2d)) THEN
179 nprow = grid_2d(1)
180 npcol = grid_2d(2)
181 END IF
182
183 IF (nprow*npcol .NE. npe) THEN
184 ! hard code for the time being the grid layout
185 my_blacs_grid_layout = blacs_grid_square
186 IF (PRESENT(blacs_grid_layout)) my_blacs_grid_layout = blacs_grid_layout
187 ! XXXXXX
188 SELECT CASE (my_blacs_grid_layout)
189 CASE (blacs_grid_square)
190 ! make the grid as 'square' as possible, where square is defined as nprow and npcol
191 ! having the largest possible gcd
192 gcd_max = -1
193 DO ipe = 1, ceiling(sqrt(real(npe, dp)))
194 jpe = npe/ipe
195 IF (ipe*jpe .NE. npe) cycle
196 IF (gcd(ipe, jpe) >= gcd_max) THEN
197 nprow = ipe
198 npcol = jpe
199 gcd_max = gcd(ipe, jpe)
200 END IF
201 END DO
202 CASE (blacs_grid_row)
203 nprow = 1
204 npcol = npe
205 CASE (blacs_grid_col)
206 nprow = npe
207 npcol = 1
208 END SELECT
209 END IF
210
211 my_row_major = .true.
212 IF (PRESENT(row_major)) my_row_major = row_major
213 IF (my_row_major) THEN
214 CALL blacs_env%gridinit(para_env, "Row-major", nprow, npcol)
215 ELSE
216 CALL blacs_env%gridinit(para_env, "Col-major", nprow, npcol)
217 END IF
218
219 ! We set the components of blacs_env here such that we can still use INTENT(OUT) with gridinit
220 blacs_env%my_pid = para_env%mepos
221 blacs_env%n_pid = para_env%num_pe
222 blacs_env%ref_count = 1
223
224 my_blacs_repeatable = .false.
225 IF (PRESENT(blacs_repeatable)) my_blacs_repeatable = blacs_repeatable
226 blacs_env%repeatable = my_blacs_repeatable
227 IF (blacs_env%repeatable) CALL blacs_env%set(15, 1)
228
229#else
230 ! In serial mode, we just have to setup the object
231 CALL blacs_env%gridinit(para_env, "Row-major", 1, 1)
232
233 blacs_env%ref_count = 1
234 blacs_env%my_pid = 0
235 blacs_env%n_pid = 1
236 mark_used(blacs_grid_layout)
237 mark_used(blacs_repeatable)
238 mark_used(grid_2d)
239 mark_used(row_major)
240#endif
241
242 CALL para_env%retain()
243 blacs_env%para_env => para_env
244
245 ! generate the mappings blacs2mpi and mpi2blacs
246 ALLOCATE (blacs_env%blacs2mpi(0:blacs_env%num_pe(1) - 1, 0:blacs_env%num_pe(2) - 1))
247 blacs_env%blacs2mpi = 0
248 blacs_env%blacs2mpi(blacs_env%mepos(1), blacs_env%mepos(2)) = para_env%mepos
249 CALL para_env%sum(blacs_env%blacs2mpi)
250 ALLOCATE (blacs_env%mpi2blacs(2, 0:para_env%num_pe - 1))
251 blacs_env%mpi2blacs = -1
252 DO ipcol = 0, blacs_env%num_pe(2) - 1
253 DO iprow = 0, blacs_env%num_pe(1) - 1
254 blacs_env%mpi2blacs(1, blacs_env%blacs2mpi(iprow, ipcol)) = iprow
255 blacs_env%mpi2blacs(2, blacs_env%blacs2mpi(iprow, ipcol)) = ipcol
256 END DO
257 END DO
258 END SUBROUTINE cp_blacs_env_create_low
259
260! **************************************************************************************************
261!> \brief retains the given blacs env
262!> \param blacs_env the blacs env to retain
263!> \par History
264!> 08.2002 created [fawzi]
265!> \author Fawzi Mohamed
266! **************************************************************************************************
267 SUBROUTINE cp_blacs_env_retain(blacs_env)
268 CLASS(cp_blacs_env_type), INTENT(INOUT) :: blacs_env
269
270 cpassert(blacs_env%ref_count > 0)
271 blacs_env%ref_count = blacs_env%ref_count + 1
272 END SUBROUTINE cp_blacs_env_retain
273
274! **************************************************************************************************
275!> \brief releases the given blacs_env
276!> \param blacs_env the blacs env to release
277!> \par History
278!> 08.2002 created [fawzi]
279!> \author Fawzi Mohamed
280! **************************************************************************************************
281 SUBROUTINE cp_blacs_env_release(blacs_env)
282 TYPE(cp_blacs_env_type), POINTER :: blacs_env
283
284 IF (ASSOCIATED(blacs_env)) THEN
285 cpassert(blacs_env%ref_count > 0)
286 blacs_env%ref_count = blacs_env%ref_count - 1
287 IF (blacs_env%ref_count < 1) THEN
288 CALL blacs_env%release()
289 DEALLOCATE (blacs_env)
290 END IF
291 END IF
292 NULLIFY (blacs_env)
293 END SUBROUTINE cp_blacs_env_release
294
295! **************************************************************************************************
296!> \brief releases the given blacs_env
297!> \param blacs_env the blacs env to release
298!> \par History
299!> 08.2002 created [fawzi]
300!> \author Fawzi Mohamed
301! **************************************************************************************************
302 SUBROUTINE cp_blacs_env_release_low(blacs_env)
303 CLASS(cp_blacs_env_type), INTENT(INOUT) :: blacs_env
304
305 CALL blacs_env%gridexit()
306 CALL mp_para_env_release(blacs_env%para_env)
307 DEALLOCATE (blacs_env%mpi2blacs)
308 DEALLOCATE (blacs_env%blacs2mpi)
309
310 END SUBROUTINE cp_blacs_env_release_low
311
312! **************************************************************************************************
313!> \brief writes the description of the given blacs env
314!> \param blacs_env the blacs environment to write
315!> \param unit_nr the unit number where to write the description of the
316!> blacs environment
317!> \par History
318!> 08.2002 created [fawzi]
319!> \author Fawzi Mohamed
320! **************************************************************************************************
321 SUBROUTINE cp_blacs_env_write(blacs_env, unit_nr)
322 CLASS(cp_blacs_env_type), INTENT(IN) :: blacs_env
323 INTEGER, INTENT(in) :: unit_nr
324
325 WRITE (unit=unit_nr, fmt="(' group=',i10,', ref_count=',i10,',')") &
326 blacs_env%get_handle(), blacs_env%ref_count
327 WRITE (unit=unit_nr, fmt="(' mepos=(',i8,',',i8,'),')") &
328 blacs_env%mepos(1), blacs_env%mepos(2)
329 WRITE (unit=unit_nr, fmt="(' num_pe=(',i8,',',i8,'),')") &
330 blacs_env%num_pe(1), blacs_env%num_pe(2)
331 IF (ASSOCIATED(blacs_env%blacs2mpi)) THEN
332 WRITE (unit=unit_nr, fmt="(' blacs2mpi=')", advance="no")
333 CALL cp_2d_i_write(blacs_env%blacs2mpi, unit_nr=unit_nr)
334 ELSE
335 WRITE (unit=unit_nr, fmt="(' blacs2mpi=*null*')")
336 END IF
337 IF (ASSOCIATED(blacs_env%para_env)) THEN
338 WRITE (unit=unit_nr, fmt="(' para_env=<cp_para_env id=',i6,'>,')") &
339 blacs_env%para_env%get_handle()
340 ELSE
341 WRITE (unit=unit_nr, fmt="(' para_env=*null*')")
342 END IF
343 WRITE (unit=unit_nr, fmt="(' my_pid=',i10,', n_pid=',i10,' }')") &
344 blacs_env%my_pid, blacs_env%n_pid
345 CALL m_flush(unit_nr)
346 END SUBROUTINE cp_blacs_env_write
347
348END 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:106
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:1291
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