Loading [MathJax]/jax/output/HTML-CSS/config.js
 (git:b77b4be)
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros
cp_fm_struct.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 represent the structure of a full matrix
10!> \par History
11!> 08.2002 created [fawzi]
12!> \author Fawzi Mohamed
13! **************************************************************************************************
21 USE machine, ONLY: m_cpuid_vlen,&
25#include "../base/base_uses.f90"
26
27 IMPLICIT NONE
28 PRIVATE
29
30 LOGICAL, PRIVATE, PARAMETER :: debug_this_module = .true.
31 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'cp_fm_struct'
32
33! the default blacs block sizes
34! consider using #ifdefs to give them the optimal values
35! these can be changed using scf_control
36! *** these are used by default
37 INTEGER, PRIVATE :: optimal_blacs_col_block_size = 64
38 INTEGER, PRIVATE :: optimal_blacs_row_block_size = 64
39 LOGICAL, PRIVATE :: force_block_size = .false.
40
47
48! **************************************************************************************************
49!> \brief keeps the information about the structure of a full matrix
50!> \param para_env the parallel environment of the matrices with this structure
51!> \param context the blacs context (parallel environment for scalapack),
52!> should be compatible with para_env
53!> \param descriptor the scalapack descriptor of the matrices, when using
54!> scalapack (ncol_block=descriptor(6), ncol_global=descriptor(4),
55!> nrow_block=descriptor(5), nrow_global=descriptor(3))
56!> \param ncol_block number of columns of a scalapack block
57!> \param nrow_block number of rows of a scalapack block
58!> \param nrow_global number of rows of the matrix
59!> \param ncol_global number of rows
60!> \param first_p_pos position of the first processor (for scalapack)
61!> \param row_indices real (global) indices of the rows (defined only for
62!> the local rows really used)
63!> \param col_indices real (global) indices of the cols (defined only for
64!> the local cols really used)
65!> \param nrow_locals nrow_locals(i) number of local rows of the matrix really
66!> used on the processors with context%mepos(1)==i
67!> \param ncol_locals ncol_locals(i) number of local rows of the matrix really
68!> used on the processors with context%mepos(2)==i
69!> \param ref_count reference count (see doc/ReferenceCounting.html)
70!> \param local_leading_dimension leading dimension of the data that is
71!> stored on this processor
72!>
73!> readonly attributes:
74!> \param nrow_local number of local rows really used on the actual processor
75!> \param ncol_local number of local cols really used on the actual processor
76!> \note
77!> use cp_fm_struct_get to extract information from this structure
78!> \par History
79!> 08.2002 created [fawzi]
80!> \author Fawzi Mohamed
81! **************************************************************************************************
83 TYPE(mp_para_env_type), POINTER :: para_env => null()
84 TYPE(cp_blacs_env_type), POINTER :: context => null()
85 INTEGER, DIMENSION(9) :: descriptor = -1
86 INTEGER :: nrow_block = -1, ncol_block = -1, nrow_global = -1, ncol_global = -1
87 INTEGER, DIMENSION(2) :: first_p_pos = -1
88 INTEGER, DIMENSION(:), POINTER :: row_indices => null(), col_indices => null(), &
89 nrow_locals => null(), ncol_locals => null()
90 INTEGER :: ref_count = -1, local_leading_dimension = -1
91 CONTAINS
92 PROCEDURE, pass(struct), non_overridable :: g2p_row => cp_fm_indxg2p_row
93 PROCEDURE, pass(struct), non_overridable :: g2p_col => cp_fm_indxg2p_col
94 PROCEDURE, pass(struct), non_overridable :: g2l_row => cp_fm_indxg2l_row
95 PROCEDURE, pass(struct), non_overridable :: g2l_col => cp_fm_indxg2l_col
96 PROCEDURE, pass(struct), non_overridable :: l2g_row => cp_fm_indxl2g_row
97 PROCEDURE, pass(struct), non_overridable :: l2g_col => cp_fm_indxl2g_col
98 END TYPE cp_fm_struct_type
99! **************************************************************************************************
101 TYPE(cp_fm_struct_type), POINTER :: struct => null()
102 END TYPE cp_fm_struct_p_type
103
104CONTAINS
105
106! **************************************************************************************************
107!> \brief allocates and initializes a full matrix structure
108!> \param fmstruct the pointer that will point to the new structure
109!> \param para_env the parallel environment
110!> \param context the blacs context of this matrix
111!> \param nrow_global the number of row of the full matrix
112!> \param ncol_global the number of columns of the full matrix
113!> \param nrow_block the number of rows of a block of the matrix,
114!> omit or set to -1 to use the built-in defaults
115!> \param ncol_block the number of columns of a block of the matrix,
116!> omit or set to -1 to use the built-in defaults
117!> \param descriptor the scalapack descriptor of the matrix (if not given
118!> a new one is allocated
119!> \param first_p_pos ...
120!> \param local_leading_dimension the leading dimension of the locally stored
121!> data block
122!> \param template_fmstruct a matrix structure where to take the default values
123!> \param square_blocks ...
124!> \param force_block ...
125!> \par History
126!> 08.2002 created [fawzi]
127!> \author Fawzi Mohamed
128! **************************************************************************************************
129 SUBROUTINE cp_fm_struct_create(fmstruct, para_env, context, nrow_global, &
130 ncol_global, nrow_block, ncol_block, descriptor, first_p_pos, &
131 local_leading_dimension, template_fmstruct, square_blocks, force_block)
132
133 TYPE(cp_fm_struct_type), POINTER :: fmstruct
134 TYPE(mp_para_env_type), TARGET, OPTIONAL :: para_env
135 INTEGER, INTENT(in), OPTIONAL :: nrow_global, ncol_global
136 INTEGER, INTENT(in), OPTIONAL :: nrow_block, ncol_block
137 INTEGER, INTENT(in), OPTIONAL :: local_leading_dimension
138 TYPE(cp_blacs_env_type), TARGET, OPTIONAL :: context
139 INTEGER, DIMENSION(9), INTENT(in), OPTIONAL :: descriptor
140 INTEGER, OPTIONAL, DIMENSION(2) :: first_p_pos
141 TYPE(cp_fm_struct_type), POINTER, OPTIONAL :: template_fmstruct
142 LOGICAL, OPTIONAL, INTENT(in) :: square_blocks
143 LOGICAL, OPTIONAL, INTENT(in) :: force_block
144
145 INTEGER :: i, nmax_block, vlen
146#if defined(__parallel)
147 INTEGER :: iunit, stat
148 INTEGER, EXTERNAL :: numroc
149 TYPE(cp_logger_type), POINTER :: logger
150#endif
151
152 LOGICAL :: my_square_blocks, my_force_block
153
154 ALLOCATE (fmstruct)
155
156 IF (.NOT. PRESENT(template_fmstruct)) THEN
157 cpassert(PRESENT(context))
158 cpassert(PRESENT(nrow_global))
159 cpassert(PRESENT(ncol_global))
160 fmstruct%local_leading_dimension = 1
161 fmstruct%nrow_block = 0 ! populate default later
162 fmstruct%ncol_block = 0 ! populate default later
163 ELSE
164 fmstruct%context => template_fmstruct%context
165 fmstruct%para_env => template_fmstruct%para_env
166 fmstruct%descriptor = template_fmstruct%descriptor
167 fmstruct%nrow_block = template_fmstruct%nrow_block
168 fmstruct%nrow_global = template_fmstruct%nrow_global
169 fmstruct%ncol_block = template_fmstruct%ncol_block
170 fmstruct%ncol_global = template_fmstruct%ncol_global
171 fmstruct%first_p_pos = template_fmstruct%first_p_pos
172 fmstruct%local_leading_dimension = &
173 template_fmstruct%local_leading_dimension
174 END IF
175
176 ! allow to request default block size (zero or negative value)
177 IF (PRESENT(nrow_block)) fmstruct%nrow_block = nrow_block
178 IF (PRESENT(ncol_block)) fmstruct%ncol_block = ncol_block
179 IF (0 >= fmstruct%nrow_block) THEN
180 fmstruct%nrow_block = optimal_blacs_row_block_size
181 END IF
182 IF (0 >= fmstruct%ncol_block) THEN
183 fmstruct%ncol_block = optimal_blacs_col_block_size
184 END IF
185 cpassert(0 < fmstruct%nrow_block .AND. 0 < fmstruct%ncol_block)
186
187 IF (PRESENT(context)) THEN
188 fmstruct%context => context
189 fmstruct%para_env => context%para_env
190 END IF
191 IF (PRESENT(para_env)) fmstruct%para_env => para_env
192 CALL fmstruct%context%retain()
193 CALL fmstruct%para_env%retain()
194
195 IF (PRESENT(nrow_global)) THEN
196 fmstruct%nrow_global = nrow_global
197 fmstruct%local_leading_dimension = 1
198 END IF
199 IF (PRESENT(ncol_global)) THEN
200 fmstruct%ncol_global = ncol_global
201 END IF
202
203 my_force_block = force_block_size
204 IF (PRESENT(force_block)) my_force_block = force_block
205 IF (.NOT. my_force_block) THEN
206 vlen = m_cpuid_vlen()
207 nmax_block = (fmstruct%nrow_global + fmstruct%context%num_pe(1) - 1)/ &
208 (fmstruct%context%num_pe(1))
209 IF (1 < vlen) THEN ! flooring not ceiling (OOB)
210 fmstruct%nrow_block = fmstruct%nrow_block/vlen*vlen
211 nmax_block = nmax_block/vlen*vlen
212 END IF
213 fmstruct%nrow_block = max(min(fmstruct%nrow_block, nmax_block), 1)
214
215 nmax_block = (fmstruct%ncol_global + fmstruct%context%num_pe(2) - 1)/ &
216 (fmstruct%context%num_pe(2))
217 IF (1 < vlen) THEN ! flooring not ceiling (OOB)
218 fmstruct%ncol_block = fmstruct%ncol_block/vlen*vlen
219 nmax_block = nmax_block/vlen*vlen
220 END IF
221 fmstruct%ncol_block = max(min(fmstruct%ncol_block, nmax_block), 1)
222 END IF
223
224 ! square matrix -> square blocks (otherwise, e.g., PDPOTRF fails)
225 my_square_blocks = fmstruct%nrow_global == fmstruct%ncol_global
226 ! however, requesting non-square blocks takes precedence
227 IF (PRESENT(square_blocks)) my_square_blocks = square_blocks
228 IF (my_square_blocks) THEN
229 fmstruct%nrow_block = min(fmstruct%nrow_block, fmstruct%ncol_block)
230 fmstruct%ncol_block = fmstruct%nrow_block
231 END IF
232
233 ALLOCATE (fmstruct%nrow_locals(0:(fmstruct%context%num_pe(1) - 1)), &
234 fmstruct%ncol_locals(0:(fmstruct%context%num_pe(2) - 1)))
235 IF (.NOT. PRESENT(template_fmstruct)) &
236 fmstruct%first_p_pos = (/0, 0/)
237 IF (PRESENT(first_p_pos)) fmstruct%first_p_pos = first_p_pos
238
239 fmstruct%nrow_locals = 0
240 fmstruct%ncol_locals = 0
241#if defined(__parallel)
242 fmstruct%nrow_locals(fmstruct%context%mepos(1)) = &
243 numroc(fmstruct%nrow_global, fmstruct%nrow_block, &
244 fmstruct%context%mepos(1), fmstruct%first_p_pos(1), &
245 fmstruct%context%num_pe(1))
246 fmstruct%ncol_locals(fmstruct%context%mepos(2)) = &
247 numroc(fmstruct%ncol_global, fmstruct%ncol_block, &
248 fmstruct%context%mepos(2), fmstruct%first_p_pos(2), &
249 fmstruct%context%num_pe(2))
250 CALL fmstruct%para_env%sum(fmstruct%nrow_locals)
251 CALL fmstruct%para_env%sum(fmstruct%ncol_locals)
252 fmstruct%nrow_locals(:) = fmstruct%nrow_locals(:)/fmstruct%context%num_pe(2)
253 fmstruct%ncol_locals(:) = fmstruct%ncol_locals(:)/fmstruct%context%num_pe(1)
254
255 IF (sum(fmstruct%ncol_locals) .NE. fmstruct%ncol_global .OR. &
256 sum(fmstruct%nrow_locals) .NE. fmstruct%nrow_global) THEN
257 ! try to collect some output if this is going to happen again
258 ! this seems to trigger on blanc, but should really never happen
259 logger => cp_get_default_logger()
260 iunit = cp_logger_get_default_unit_nr(logger, local=.true.)
261 WRITE (iunit, *) "mepos", fmstruct%context%mepos(1:2), "numpe", fmstruct%context%num_pe(1:2)
262 WRITE (iunit, *) "ncol_global", fmstruct%ncol_global
263 WRITE (iunit, *) "nrow_global", fmstruct%nrow_global
264 WRITE (iunit, *) "ncol_locals", fmstruct%ncol_locals
265 WRITE (iunit, *) "nrow_locals", fmstruct%nrow_locals
266 CALL m_flush(iunit)
267 END IF
268
269 IF (sum(fmstruct%ncol_locals) .NE. fmstruct%ncol_global) &
270 cpabort("sum of local cols not equal global cols")
271 IF (sum(fmstruct%nrow_locals) .NE. fmstruct%nrow_global) &
272 cpabort("sum of local row not equal global rows")
273#else
274 ! block = full matrix
275 fmstruct%nrow_block = fmstruct%nrow_global
276 fmstruct%ncol_block = fmstruct%ncol_global
277 fmstruct%nrow_locals(fmstruct%context%mepos(1)) = fmstruct%nrow_global
278 fmstruct%ncol_locals(fmstruct%context%mepos(2)) = fmstruct%ncol_global
279#endif
280
281 fmstruct%local_leading_dimension = max(fmstruct%local_leading_dimension, &
282 fmstruct%nrow_locals(fmstruct%context%mepos(1)))
283 IF (PRESENT(local_leading_dimension)) THEN
284 IF (max(1, fmstruct%nrow_locals(fmstruct%context%mepos(1))) > local_leading_dimension) &
285 CALL cp_abort(__location__, "local_leading_dimension too small ("// &
286 cp_to_string(local_leading_dimension)//"<"// &
287 cp_to_string(fmstruct%local_leading_dimension)//")")
288 fmstruct%local_leading_dimension = local_leading_dimension
289 END IF
290
291 NULLIFY (fmstruct%row_indices, fmstruct%col_indices)
292
293 ! the max should go away
294 ALLOCATE (fmstruct%row_indices(max(fmstruct%nrow_locals(fmstruct%context%mepos(1)), 1)))
295 DO i = 1, SIZE(fmstruct%row_indices)
296#ifdef __parallel
297 fmstruct%row_indices(i) = fmstruct%l2g_row(i, fmstruct%context%mepos(1))
298#else
299 fmstruct%row_indices(i) = i
300#endif
301 END DO
302 ALLOCATE (fmstruct%col_indices(max(fmstruct%ncol_locals(fmstruct%context%mepos(2)), 1)))
303 DO i = 1, SIZE(fmstruct%col_indices)
304#ifdef __parallel
305 fmstruct%col_indices(i) = fmstruct%l2g_col(i, fmstruct%context%mepos(2))
306#else
307 fmstruct%col_indices(i) = i
308#endif
309 END DO
310
311 fmstruct%ref_count = 1
312
313 IF (PRESENT(descriptor)) THEN
314 fmstruct%descriptor = descriptor
315 ELSE
316 fmstruct%descriptor = 0
317#if defined(__parallel)
318 ! local leading dimension needs to be at least 1
319 CALL descinit(fmstruct%descriptor, fmstruct%nrow_global, &
320 fmstruct%ncol_global, fmstruct%nrow_block, &
321 fmstruct%ncol_block, fmstruct%first_p_pos(1), &
322 fmstruct%first_p_pos(2), fmstruct%context, &
323 fmstruct%local_leading_dimension, stat)
324 cpassert(stat == 0)
325#endif
326 END IF
327 END SUBROUTINE cp_fm_struct_create
328
329! **************************************************************************************************
330!> \brief retains a full matrix structure
331!> \param fmstruct the structure to retain
332!> \par History
333!> 08.2002 created [fawzi]
334!> \author Fawzi Mohamed
335! **************************************************************************************************
336 SUBROUTINE cp_fm_struct_retain(fmstruct)
337 TYPE(cp_fm_struct_type), INTENT(INOUT) :: fmstruct
338
339 cpassert(fmstruct%ref_count > 0)
340 fmstruct%ref_count = fmstruct%ref_count + 1
341 END SUBROUTINE cp_fm_struct_retain
342
343! **************************************************************************************************
344!> \brief releases a full matrix structure
345!> \param fmstruct the structure to release
346!> \par History
347!> 08.2002 created [fawzi]
348!> \author Fawzi Mohamed
349! **************************************************************************************************
350 SUBROUTINE cp_fm_struct_release(fmstruct)
351 TYPE(cp_fm_struct_type), POINTER :: fmstruct
352
353 IF (ASSOCIATED(fmstruct)) THEN
354 cpassert(fmstruct%ref_count > 0)
355 fmstruct%ref_count = fmstruct%ref_count - 1
356 IF (fmstruct%ref_count < 1) THEN
357 CALL cp_blacs_env_release(fmstruct%context)
358 CALL mp_para_env_release(fmstruct%para_env)
359 IF (ASSOCIATED(fmstruct%row_indices)) THEN
360 DEALLOCATE (fmstruct%row_indices)
361 END IF
362 IF (ASSOCIATED(fmstruct%col_indices)) THEN
363 DEALLOCATE (fmstruct%col_indices)
364 END IF
365 IF (ASSOCIATED(fmstruct%nrow_locals)) THEN
366 DEALLOCATE (fmstruct%nrow_locals)
367 END IF
368 IF (ASSOCIATED(fmstruct%ncol_locals)) THEN
369 DEALLOCATE (fmstruct%ncol_locals)
370 END IF
371 DEALLOCATE (fmstruct)
372 END IF
373 END IF
374 NULLIFY (fmstruct)
375 END SUBROUTINE cp_fm_struct_release
376
377! **************************************************************************************************
378!> \brief returns true if the two matrix structures are equivalent, false
379!> otherwise.
380!> \param fmstruct1 one of the full matrix structures to compare
381!> \param fmstruct2 the second of the full matrix structures to compare
382!> \return ...
383!> \par History
384!> 08.2002 created [fawzi]
385!> \author Fawzi Mohamed
386! **************************************************************************************************
387 FUNCTION cp_fm_struct_equivalent(fmstruct1, fmstruct2) RESULT(res)
388 TYPE(cp_fm_struct_type), POINTER :: fmstruct1, fmstruct2
389 LOGICAL :: res
390
391 INTEGER :: i
392
393 cpassert(ASSOCIATED(fmstruct1))
394 cpassert(ASSOCIATED(fmstruct2))
395 cpassert(fmstruct1%ref_count > 0)
396 cpassert(fmstruct2%ref_count > 0)
397 IF (ASSOCIATED(fmstruct1, fmstruct2)) THEN
398 res = .true.
399 ELSE
400 res = (fmstruct1%context == fmstruct2%context) .AND. &
401 (fmstruct1%nrow_global == fmstruct2%nrow_global) .AND. &
402 (fmstruct1%ncol_global == fmstruct2%ncol_global) .AND. &
403 (fmstruct1%nrow_block == fmstruct2%nrow_block) .AND. &
404 (fmstruct1%ncol_block == fmstruct2%ncol_block) .AND. &
405 (fmstruct1%local_leading_dimension == &
406 fmstruct2%local_leading_dimension)
407 DO i = 1, 9
408 res = res .AND. (fmstruct1%descriptor(i) == fmstruct1%descriptor(i))
409 END DO
410 END IF
411 END FUNCTION cp_fm_struct_equivalent
412
413! **************************************************************************************************
414!> \brief returns the values of various attributes of the matrix structure
415!> \param fmstruct the structure you want info about
416!> \param para_env ...
417!> \param context ...
418!> \param descriptor ...
419!> \param ncol_block ...
420!> \param nrow_block ...
421!> \param nrow_global ...
422!> \param ncol_global ...
423!> \param first_p_pos ...
424!> \param row_indices ...
425!> \param col_indices ...
426!> \param nrow_local ...
427!> \param ncol_local ...
428!> \param nrow_locals ...
429!> \param ncol_locals ...
430!> \param local_leading_dimension ...
431!> \par History
432!> 08.2002 created [fawzi]
433!> \author Fawzi Mohamed
434! **************************************************************************************************
435 SUBROUTINE cp_fm_struct_get(fmstruct, para_env, context, &
436 descriptor, ncol_block, nrow_block, nrow_global, &
437 ncol_global, first_p_pos, row_indices, &
438 col_indices, nrow_local, ncol_local, nrow_locals, ncol_locals, &
439 local_leading_dimension)
440 TYPE(cp_fm_struct_type), INTENT(IN) :: fmstruct
441 TYPE(mp_para_env_type), OPTIONAL, POINTER :: para_env
442 TYPE(cp_blacs_env_type), OPTIONAL, POINTER :: context
443 INTEGER, DIMENSION(9), INTENT(OUT), OPTIONAL :: descriptor
444 INTEGER, INTENT(out), OPTIONAL :: ncol_block, nrow_block, nrow_global, &
445 ncol_global
446 INTEGER, DIMENSION(2), INTENT(out), OPTIONAL :: first_p_pos
447 INTEGER, DIMENSION(:), OPTIONAL, POINTER :: row_indices, col_indices
448 INTEGER, INTENT(out), OPTIONAL :: nrow_local, ncol_local
449 INTEGER, DIMENSION(:), OPTIONAL, POINTER :: nrow_locals, ncol_locals
450 INTEGER, INTENT(out), OPTIONAL :: local_leading_dimension
451
452 IF (PRESENT(para_env)) para_env => fmstruct%para_env
453 IF (PRESENT(context)) context => fmstruct%context
454 IF (PRESENT(descriptor)) descriptor = fmstruct%descriptor
455 IF (PRESENT(ncol_block)) ncol_block = fmstruct%ncol_block
456 IF (PRESENT(nrow_block)) nrow_block = fmstruct%nrow_block
457 IF (PRESENT(nrow_global)) nrow_global = fmstruct%nrow_global
458 IF (PRESENT(ncol_global)) ncol_global = fmstruct%ncol_global
459 IF (PRESENT(first_p_pos)) first_p_pos = fmstruct%first_p_pos
460 IF (PRESENT(nrow_locals)) nrow_locals => fmstruct%nrow_locals
461 IF (PRESENT(ncol_locals)) ncol_locals => fmstruct%ncol_locals
462 IF (PRESENT(local_leading_dimension)) local_leading_dimension = &
463 fmstruct%local_leading_dimension
464
465 IF (PRESENT(nrow_local)) nrow_local = fmstruct%nrow_locals(fmstruct%context%mepos(1))
466 IF (PRESENT(ncol_local)) ncol_local = fmstruct%ncol_locals(fmstruct%context%mepos(2))
467
468 IF (PRESENT(row_indices)) row_indices => fmstruct%row_indices
469 IF (PRESENT(col_indices)) col_indices => fmstruct%col_indices
470 END SUBROUTINE cp_fm_struct_get
471
472! **************************************************************************************************
473!> \brief Write nicely formatted info about the FM struct to the given I/O unit
474!> \param fmstruct a cp_fm_struct_type instance
475!> \param io_unit the I/O unit to use for writing
476! **************************************************************************************************
477 SUBROUTINE cp_fm_struct_write_info(fmstruct, io_unit)
478 TYPE(cp_fm_struct_type), INTENT(IN) :: fmstruct
479 INTEGER, INTENT(IN) :: io_unit
480
481 INTEGER, PARAMETER :: oblock_size = 8
482
483 CHARACTER(len=30) :: fm
484 INTEGER :: oblock
485
486 WRITE (fm, "(A,I2,A)") "(A,I5,A,I5,A,", oblock_size, "I6)"
487
488 WRITE (io_unit, '(A,I12)') "CP_FM_STRUCT | No. of matrix columns: ", fmstruct%ncol_global
489 WRITE (io_unit, '(A,I12)') "CP_FM_STRUCT | No. of matrix rows: ", fmstruct%nrow_global
490 WRITE (io_unit, '(A,I12)') "CP_FM_STRUCT | No. of block columns: ", fmstruct%ncol_block
491 WRITE (io_unit, '(A,I12)') "CP_FM_STRUCT | No. of block rows: ", fmstruct%nrow_block
492
493 WRITE (io_unit, '(A)') "CP_FM_STRUCT | Number of local columns: "
494 DO oblock = 0, (SIZE(fmstruct%ncol_locals) - 1)/oblock_size
495 WRITE (io_unit, fm) "CP_FM_STRUCT | CPUs ", &
496 oblock*oblock_size, "..", (oblock + 1)*oblock_size - 1, ": ", &
497 fmstruct%ncol_locals(oblock*oblock_size:min(SIZE(fmstruct%ncol_locals), (oblock + 1)*oblock_size) - 1)
498 END DO
499
500 WRITE (io_unit, '(A)') "CP_FM_STRUCT | Number of local rows: "
501 DO oblock = 0, (SIZE(fmstruct%nrow_locals) - 1)/oblock_size
502 WRITE (io_unit, fm) "CP_FM_STRUCT | CPUs ", &
503 oblock*oblock_size, "..", (oblock + 1)*oblock_size - 1, ": ", &
504 fmstruct%nrow_locals(oblock*oblock_size:min(SIZE(fmstruct%nrow_locals), (oblock + 1)*oblock_size) - 1)
505 END DO
506 END SUBROUTINE cp_fm_struct_write_info
507
508! **************************************************************************************************
509!> \brief creates a struct with twice the number of blocks on each core.
510!> If matrix A has to be multiplied with B anc C, a
511!> significant speedup of pdgemm can be acchieved by joining the matrices
512!> in a new one with this structure (see arnoldi in rt_matrix_exp)
513!> \param fmstruct the struct to create
514!> \param struct struct of either A or B
515!> \param context ...
516!> \param col in which direction the matrix should be enlarged
517!> \param row in which direction the matrix should be enlarged
518!> \par History
519!> 06.2009 created [fschiff]
520!> \author Florian Schiffmann
521! **************************************************************************************************
522 SUBROUTINE cp_fm_struct_double(fmstruct, struct, context, col, row)
523 TYPE(cp_fm_struct_type), POINTER :: fmstruct
524 TYPE(cp_fm_struct_type), INTENT(INOUT) :: struct
525 TYPE(cp_blacs_env_type), INTENT(INOUT), TARGET :: context
526 LOGICAL, INTENT(in) :: col, row
527
528 INTEGER :: n_doubled_items_in_partially_filled_block, ncol_block, ncol_global, newdim_col, &
529 newdim_row, nfilled_blocks, nfilled_blocks_remain, nprocs_col, nprocs_row, nrow_block, &
530 nrow_global
531 TYPE(mp_para_env_type), POINTER :: para_env
532
533 CALL cp_fm_struct_get(struct, nrow_global=nrow_global, &
534 ncol_global=ncol_global, nrow_block=nrow_block, &
535 ncol_block=ncol_block)
536 newdim_row = nrow_global
537 newdim_col = ncol_global
538 nprocs_row = context%num_pe(1)
539 nprocs_col = context%num_pe(2)
540 para_env => struct%para_env
541
542 IF (col) THEN
543 IF (ncol_global == 0) THEN
544 newdim_col = 0
545 ELSE
546 ! ncol_block nfilled_blocks_remain * ncol_block
547 ! |<--->| |<--->|
548 ! |-----|-----|-----|-----|---|
549 ! | 0 | 1 | 2 | 0 | 1 | <- context%mepos(2)
550 ! |-----|-----|-----|-----|---|
551 ! |<--- nfilled_blocks -->|<-> -- items (columns) in partially filled blocks
552 ! | * ncol_block |
553 n_doubled_items_in_partially_filled_block = 2*mod(ncol_global, ncol_block)
554 nfilled_blocks = ncol_global/ncol_block
555 nfilled_blocks_remain = mod(nfilled_blocks, nprocs_col)
556 newdim_col = 2*(nfilled_blocks/nprocs_col)
557 IF (n_doubled_items_in_partially_filled_block > ncol_block) THEN
558 ! doubled number of columns in a partially filled block does not fit into a single block.
559 ! Due to cyclic distribution of ScaLAPACK blocks, an extra block for each core needs to be added
560 ! |-----|-----|-----|----| |-----|-----|-----|-----|-----|-----|-----|-----|-----|---|
561 ! | 0 | 1 | 2 | 0 | --> | 0 | 1 | 2 | 0 | 1 | 2 | 0 | 1 | 2 | 0|
562 ! |-----|-----|-----|----| |-----|-----|-----|-----|-----|-----|-----|-----|-----|---|
563 ! a a a b a1 a1 a1 a2 a2 a2 b1 empty empty b2
564 newdim_col = newdim_col + 1
565
566 ! the number of columns which does not fit into the added extra block
567 n_doubled_items_in_partially_filled_block = n_doubled_items_in_partially_filled_block - ncol_block
568 ELSE IF (nfilled_blocks_remain > 0) THEN
569 ! |-----|-----|-----|-----|--| |-----|-----|-----|-----|-----|-----|-----|-----|-----|-----|
570 ! | 0 | 1 | 2 | 0 | 1| -> | 0 | 1 | 2 | 0 | 1 | 2 | 0 | 1 | 2 | 0 |
571 ! |-----|-----|-----|-----|--| |-----|-----|-----|-----|-----|-----|-----|-----|-----|-----|
572 ! a a a b b a1 a1 a1 a2 a2 a2 b1 b1 b2 empty b2
573 newdim_col = newdim_col + 1
574 n_doubled_items_in_partially_filled_block = 0
575 END IF
576
577 newdim_col = (newdim_col*nprocs_col + nfilled_blocks_remain)*ncol_block + n_doubled_items_in_partially_filled_block
578 END IF
579 END IF
580
581 IF (row) THEN
582 IF (nrow_global == 0) THEN
583 newdim_row = 0
584 ELSE
585 n_doubled_items_in_partially_filled_block = 2*mod(nrow_global, nrow_block)
586 nfilled_blocks = nrow_global/nrow_block
587 nfilled_blocks_remain = mod(nfilled_blocks, nprocs_row)
588 newdim_row = 2*(nfilled_blocks/nprocs_row)
589 IF (n_doubled_items_in_partially_filled_block > nrow_block) THEN
590 newdim_row = newdim_row + 1
591 n_doubled_items_in_partially_filled_block = n_doubled_items_in_partially_filled_block - nrow_block
592 ELSE IF (nfilled_blocks_remain > 0) THEN
593 newdim_row = newdim_row + 1
594 n_doubled_items_in_partially_filled_block = 0
595 END IF
596
597 newdim_row = (newdim_row*nprocs_row + nfilled_blocks_remain)*nrow_block + n_doubled_items_in_partially_filled_block
598 END IF
599 END IF
600
601 ! square_blocks=.FALSE. ensures that matrix blocks of the doubled matrix will have
602 ! nrow_block x ncol_block shape even in case of a square doubled matrix
603 CALL cp_fm_struct_create(fmstruct=fmstruct, para_env=para_env, &
604 context=context, &
605 nrow_global=newdim_row, &
606 ncol_global=newdim_col, &
607 ncol_block=ncol_block, &
608 nrow_block=nrow_block, &
609 square_blocks=.false.)
610
611 END SUBROUTINE cp_fm_struct_double
612! **************************************************************************************************
613!> \brief allows to modify the default settings for matrix creation
614!> \param nrow_block ...
615!> \param ncol_block ...
616!> \param force_block ...
617! **************************************************************************************************
618 SUBROUTINE cp_fm_struct_config(nrow_block, ncol_block, force_block)
619 INTEGER, INTENT(IN), OPTIONAL :: nrow_block, ncol_block
620 LOGICAL, INTENT(IN), OPTIONAL :: force_block
621
622 INTEGER :: vlen
623
624 vlen = m_cpuid_vlen()
625 IF (PRESENT(ncol_block)) THEN
626 IF (0 < ncol_block) THEN
627 optimal_blacs_col_block_size = (ncol_block + vlen - 1)/vlen*vlen
628 END IF
629 END IF
630 IF (PRESENT(nrow_block)) THEN
631 IF (0 < nrow_block) THEN
632 optimal_blacs_row_block_size = (nrow_block + vlen - 1)/vlen*vlen
633 END IF
634 END IF
635 IF (PRESENT(force_block)) force_block_size = force_block
636
637 END SUBROUTINE cp_fm_struct_config
638
639! **************************************************************************************************
640!> \brief ...
641!> \return ...
642! **************************************************************************************************
643 FUNCTION cp_fm_struct_get_nrow_block() RESULT(res)
644 INTEGER :: res
645
646 res = optimal_blacs_row_block_size
647 END FUNCTION cp_fm_struct_get_nrow_block
648
649! **************************************************************************************************
650!> \brief ...
651!> \return ...
652! **************************************************************************************************
653 FUNCTION cp_fm_struct_get_ncol_block() RESULT(res)
654 INTEGER :: res
655
656 res = optimal_blacs_col_block_size
657 END FUNCTION cp_fm_struct_get_ncol_block
658
659! **************************************************************************************************
660!> \brief wrapper to scalapack function INDXG2P that computes the row process
661!> coordinate which possesses the entry of a distributed matrix specified
662!> by a global index INDXGLOB.
663!> \param struct ...
664!> \param INDXGLOB ...
665!> \return ...
666!> \author Mauro Del Ben [MDB] - 12.2012, modified by F. Stein
667! **************************************************************************************************
668 FUNCTION cp_fm_indxg2p_row(struct, INDXGLOB) RESULT(G2P)
669 CLASS(cp_fm_struct_type), INTENT(IN) :: struct
670 INTEGER, INTENT(IN) :: indxglob
671 INTEGER :: g2p
672
673#if defined(__parallel)
674 INTEGER :: number_of_process_rows
675 INTEGER, EXTERNAL :: indxg2p
676#endif
677
678#if defined(__parallel)
679
680 CALL struct%context%get(number_of_process_rows=number_of_process_rows)
681
682 g2p = indxg2p(indxglob, struct%nrow_block, 0, struct%first_p_pos(1), number_of_process_rows)
683
684#else
685 mark_used(struct)
686 mark_used(indxglob)
687
688 g2p = 0
689
690#endif
691
692 END FUNCTION cp_fm_indxg2p_row
693
694! **************************************************************************************************
695!> \brief wrapper to scalapack function INDXG2P that computes the col process
696!> coordinate which possesses the entry of a distributed matrix specified
697!> by a global index INDXGLOB.
698!> \param struct ...
699!> \param INDXGLOB ...
700!> \return ...
701!> \author Mauro Del Ben [MDB] - 12.2012, modified by F. Stein
702! **************************************************************************************************
703 FUNCTION cp_fm_indxg2p_col(struct, INDXGLOB) RESULT(G2P)
704 CLASS(cp_fm_struct_type), INTENT(IN) :: struct
705 INTEGER, INTENT(IN) :: indxglob
706 INTEGER :: g2p
707
708#if defined(__parallel)
709 INTEGER :: number_of_process_columns
710 INTEGER, EXTERNAL :: indxg2p
711#endif
712
713#if defined(__parallel)
714
715 CALL struct%context%get(number_of_process_columns=number_of_process_columns)
716
717 g2p = indxg2p(indxglob, struct%ncol_block, 0, struct%first_p_pos(2), number_of_process_columns)
718
719#else
720 mark_used(struct)
721 mark_used(indxglob)
722
723 g2p = 0
724
725#endif
726
727 END FUNCTION cp_fm_indxg2p_col
728
729! **************************************************************************************************
730!> \brief wrapper to scalapack function INDXG2L that computes the local index
731!> of a distributed matrix entry pointed to by the global index INDXGLOB.
732!>
733!> Arguments
734!> =========
735!>
736!> INDXGLOB (global input) INTEGER
737!> The global index of the distributed matrix entry.
738!>
739!> NB (global input) INTEGER
740!> Block size, size of the blocks the distributed matrix is
741!> split into.
742!>
743!> IPROC (local dummy) INTEGER
744!> Dummy argument in this case in order to unify the calling
745!> sequence of the tool-routines.
746!>
747!> ISRCPROC (local dummy) INTEGER
748!> Dummy argument in this case in order to unify the calling
749!> sequence of the tool-routines.
750!>
751!> NPROCS (global input) INTEGER
752!> The total number processes over which the distributed
753!> matrix is distributed.
754!>
755!> \param struct ...
756!> \param INDXGLOB ...
757!> \return ...
758!> \author Mauro Del Ben [MDB] - 12.2012
759! **************************************************************************************************
760 FUNCTION cp_fm_indxg2l_row(struct, INDXGLOB) RESULT(G2L)
761 CLASS(cp_fm_struct_type), INTENT(IN) :: struct
762 INTEGER, INTENT(IN) :: indxglob
763 INTEGER :: g2l
764
765#if defined(__parallel)
766 INTEGER :: number_of_process_rows
767 INTEGER, EXTERNAL :: indxg2l
768#endif
769
770#if defined(__parallel)
771
772 CALL struct%context%get(number_of_process_rows=number_of_process_rows)
773
774 g2l = indxg2l(indxglob, struct%nrow_block, 0, struct%first_p_pos(1), number_of_process_rows)
775
776#else
777 mark_used(struct)
778
779 g2l = indxglob
780
781#endif
782
783 END FUNCTION cp_fm_indxg2l_row
784
785! **************************************************************************************************
786!> \brief wrapper to scalapack function INDXG2L that computes the local index
787!> of a distributed matrix entry pointed to by the global index INDXGLOB.
788!>
789!> Arguments
790!> =========
791!>
792!> INDXGLOB (global input) INTEGER
793!> The global index of the distributed matrix entry.
794!>
795!> NB (global input) INTEGER
796!> Block size, size of the blocks the distributed matrix is
797!> split into.
798!>
799!> IPROC (local dummy) INTEGER
800!> Dummy argument in this case in order to unify the calling
801!> sequence of the tool-routines.
802!>
803!> ISRCPROC (local dummy) INTEGER
804!> Dummy argument in this case in order to unify the calling
805!> sequence of the tool-routines.
806!>
807!> NPROCS (global input) INTEGER
808!> The total number processes over which the distributed
809!> matrix is distributed.
810!>
811!> \param struct ...
812!> \param INDXGLOB ...
813!> \return ...
814!> \author Mauro Del Ben [MDB] - 12.2012
815! **************************************************************************************************
816 FUNCTION cp_fm_indxg2l_col(struct, INDXGLOB) RESULT(G2L)
817 CLASS(cp_fm_struct_type), INTENT(IN) :: struct
818 INTEGER, INTENT(IN) :: indxglob
819 INTEGER :: g2l
820
821#if defined(__parallel)
822 INTEGER :: number_of_process_columns
823 INTEGER, EXTERNAL :: indxg2l
824#endif
825
826#if defined(__parallel)
827
828 CALL struct%context%get(number_of_process_columns=number_of_process_columns)
829
830 g2l = indxg2l(indxglob, struct%ncol_block, 0, struct%first_p_pos(2), number_of_process_columns)
831
832#else
833 mark_used(struct)
834
835 g2l = indxglob
836
837#endif
838
839 END FUNCTION cp_fm_indxg2l_col
840
841! **************************************************************************************************
842!> \brief wrapper to scalapack function INDXL2G that computes the global index
843!> of a distributed matrix entry pointed to by the local index INDXLOC
844!> of the process indicated by IPROC.
845!>
846!> Arguments
847!> =========
848!>
849!> INDXLOC (global input) INTEGER
850!> The local index of the distributed matrix entry.
851!>
852!> NB (global input) INTEGER
853!> Block size, size of the blocks the distributed matrix is
854!> split into.
855!>
856!> IPROC (local input) INTEGER
857!> The coordinate of the process whose local array row or
858!> column is to be determined.
859!>
860!> ISRCPROC (global input) INTEGER
861!> The coordinate of the process that possesses the first
862!> row/column of the distributed matrix.
863!>
864!> NPROCS (global input) INTEGER
865!> The total number processes over which the distributed
866!> matrix is distributed.
867!>
868!> \param struct ...
869!> \param INDXLOC ...
870!> \param IPROC ...
871!> \return ...
872!> \author Mauro Del Ben [MDB] - 12.2012
873! **************************************************************************************************
874 FUNCTION cp_fm_indxl2g_row(struct, INDXLOC, IPROC) RESULT(L2G)
875 CLASS(cp_fm_struct_type), INTENT(IN) :: struct
876 INTEGER, INTENT(IN) :: indxloc, iproc
877 INTEGER :: l2g
878
879#if defined(__parallel)
880 INTEGER :: number_of_process_rows
881 INTEGER, EXTERNAL :: indxl2g
882
883 CALL struct%context%get(number_of_process_rows=number_of_process_rows)
884
885 l2g = indxl2g(indxloc, struct%nrow_block, iproc, struct%first_p_pos(1), number_of_process_rows)
886
887#else
888 mark_used(struct)
889 mark_used(indxloc)
890 mark_used(iproc)
891
892 l2g = indxloc
893
894#endif
895
896 END FUNCTION cp_fm_indxl2g_row
897
898! **************************************************************************************************
899!> \brief wrapper to scalapack function INDXL2G that computes the global index
900!> of a distributed matrix entry pointed to by the local index INDXLOC
901!> of the process indicated by IPROC.
902!>
903!> Arguments
904!> =========
905!>
906!> INDXLOC (global input) INTEGER
907!> The local index of the distributed matrix entry.
908!>
909!> NB (global input) INTEGER
910!> Block size, size of the blocks the distributed matrix is
911!> split into.
912!>
913!> IPROC (local input) INTEGER
914!> The coordinate of the process whose local array row or
915!> column is to be determined.
916!>
917!> ISRCPROC (global input) INTEGER
918!> The coordinate of the process that possesses the first
919!> row/column of the distributed matrix.
920!>
921!> NPROCS (global input) INTEGER
922!> The total number processes over which the distributed
923!> matrix is distributed.
924!>
925!> \param struct ...
926!> \param INDXLOC ...
927!> \param IPROC ...
928!> \return ...
929!> \author Mauro Del Ben [MDB] - 12.2012
930! **************************************************************************************************
931 FUNCTION cp_fm_indxl2g_col(struct, INDXLOC, IPROC) RESULT(L2G)
932 CLASS(cp_fm_struct_type), INTENT(IN) :: struct
933 INTEGER, INTENT(IN) :: indxloc, iproc
934 INTEGER :: l2g
935
936#if defined(__parallel)
937 INTEGER :: number_of_process_columns
938 INTEGER, EXTERNAL :: indxl2g
939
940 CALL struct%context%get(number_of_process_columns=number_of_process_columns)
941
942 l2g = indxl2g(indxloc, struct%ncol_block, iproc, struct%first_p_pos(2), number_of_process_columns)
943
944#else
945 mark_used(struct)
946 mark_used(indxloc)
947 mark_used(iproc)
948
949 l2g = indxloc
950
951#endif
952
953 END FUNCTION cp_fm_indxl2g_col
954
955END MODULE cp_fm_struct
methods related to the blacs parallel environment
subroutine, public cp_blacs_env_release(blacs_env)
releases the given blacs_env
represent the structure of a full matrix
subroutine, public cp_fm_struct_create(fmstruct, para_env, context, nrow_global, ncol_global, nrow_block, ncol_block, descriptor, first_p_pos, local_leading_dimension, template_fmstruct, square_blocks, force_block)
allocates and initializes a full matrix structure
integer function, public cp_fm_struct_get_nrow_block()
...
integer function, public cp_fm_struct_get_ncol_block()
...
subroutine, public cp_fm_struct_config(nrow_block, ncol_block, force_block)
allows to modify the default settings for matrix creation
subroutine, public cp_fm_struct_get(fmstruct, para_env, context, descriptor, ncol_block, nrow_block, nrow_global, ncol_global, first_p_pos, row_indices, col_indices, nrow_local, ncol_local, nrow_locals, ncol_locals, local_leading_dimension)
returns the values of various attributes of the matrix structure
subroutine, public cp_fm_struct_double(fmstruct, struct, context, col, row)
creates a struct with twice the number of blocks on each core. If matrix A has to be multiplied with ...
logical function, public cp_fm_struct_equivalent(fmstruct1, fmstruct2)
returns true if the two matrix structures are equivalent, false otherwise.
subroutine, public cp_fm_struct_retain(fmstruct)
retains a full matrix structure
subroutine, public cp_fm_struct_write_info(fmstruct, io_unit)
Write nicely formatted info about the FM struct to the given I/O unit.
subroutine, public cp_fm_struct_release(fmstruct)
releases a full matrix structure
various routines to log and control the output. The idea is that decisions about where to log should ...
recursive integer function, public cp_logger_get_default_unit_nr(logger, local, skip_not_ionode)
asks the default unit number of the given logger. try to use cp_logger_get_unit_nr
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
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
pure integer function, public m_cpuid_vlen(cpuid, typesize)
Determine vector-length for a given CPUID.
Definition machine.F:276
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...
keeps the information about the structure of a full matrix
type of a logger, at the moment it contains just a print level starting at which level it should be l...
stores all the informations relevant to an mpi environment