Octopus
multicomm.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2015 M. Marques, A. Castro, A. Rubio, G. Bertsch, X. Andrade
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
51
52module multicomm_oct_m
53 use blacs_oct_m
55 use debug_oct_m
56 use global_oct_m
57 use, intrinsic :: iso_fortran_env
59 use mpi_oct_m
62#if defined(HAVE_OPENMP)
63 use omp_lib
64#endif
65 use parser_oct_m
67 use utils_oct_m
68
69 implicit none
70
71 private
72
73 public :: &
76#if defined(HAVE_MPI)
78#endif
87
88 integer, public, parameter :: &
89 P_MASTER = 1, &
90 p_slave = 2
91
92 integer, public, parameter :: &
93 PAR_AUTO = -1, &
94 par_no = 0
95
96 integer, parameter :: n_par_types = 4
97 character(len=11), parameter :: par_types(0:n_par_types) = &
98 (/ &
99 "serial ", &
100 "ParDomains", &
101 "ParStates ", &
102 "ParKPoints", &
103 "ParOther " &
104 /)
105
106 integer, parameter :: MAX_INDEX = 5
107
113 type multicomm_t
114 private
115 integer :: n_node
116
117 integer, public :: par_strategy
118
119 integer, allocatable :: group_sizes(:)
120 integer, allocatable, public :: who_am_i(:)
121 type(MPI_Comm), allocatable, public :: group_comm(:)
122 type(MPI_Comm), public :: dom_st_comm = mpi_comm_undefined
123 type(MPI_Comm), public :: st_kpt_comm = mpi_comm_undefined
124 type(MPI_Comm), public :: dom_st_kpt_comm = mpi_comm_undefined
125 type(mpi_grp_t), public :: intranode_grp
126 type(mpi_grp_t), public :: internode_grp
127
128 integer :: nthreads
129 integer :: node_type
130 logical :: have_slaves
131
132 type(MPI_Comm) :: full_comm = mpi_comm_undefined
133 integer :: full_comm_rank
134 type(MPI_Comm), public :: master_comm = mpi_comm_undefined
135 integer :: master_comm_rank
136 type(MPI_Comm), public :: slave_intercomm = mpi_comm_undefined
137
138 logical :: reorder_ranks
139 end type multicomm_t
140
143 private
144 type(mpi_grp_t) :: grp
145 integer :: rounds
146 integer, allocatable, public :: schedule(:, :)
147 end type multicomm_all_pairs_t
148
149contains
150
151 ! ---------------------------------------------------------
152 subroutine multicomm_all_pairs_copy(apout, apin)
153 type(multicomm_all_pairs_t), intent(inout) :: apout
154 type(multicomm_all_pairs_t), intent(in) :: apin
155
157
158 call mpi_grp_copy(apout%grp, apin%grp)
159 apout%rounds = apin%rounds
160 if (allocated(apin%schedule)) then
161 safe_allocate(apout%schedule(1:size(apin%schedule, 1), 1:size(apin%schedule, 2)))
162 apout%schedule = apin%schedule
163 else
164 if (allocated(apout%schedule)) then
165 safe_deallocate_a(apout%schedule)
166 end if
167 end if
168
170 end subroutine multicomm_all_pairs_copy
171
172 ! ---------------------------------------------------------
174 subroutine multicomm_init(mc, namespace, base_grp, mode_para, n_node, index_range, min_range)
175 type(multicomm_t), intent(out) :: mc
176 type(namespace_t), intent(in) :: namespace
177 type(mpi_grp_t), intent(in) :: base_grp
178 type(calc_mode_par_t), intent(in) :: mode_para
179 integer, intent(in) :: n_node
180 integer(int64), intent(inout) :: index_range(:)
181 integer, intent(in) :: min_range(:)
182
183 integer :: ii, num_slaves, slave_level, ipar
184 integer :: parse(1:p_strategy_max), default(1:p_strategy_max)
185 integer :: parallel_mask, default_mask
186
187 push_sub(multicomm_init)
188
189 mc%n_node = n_node
190 parallel_mask = mode_para%parallel_mask()
191 default_mask = mode_para%default_parallel_mask()
192
193 call messages_print_with_emphasis(msg="Parallelization", namespace=namespace)
194
195 !%Variable ReorderRanks
196 !%Default no
197 !%Type logical
198 !%Section Execution::Parallelization
199 !%Description
200 !% This variable controls whether the ranks are reorganized to have a more
201 !% compact distribution with respect to domain parallelization which needs
202 !% to communicate most often. Depending on the system, this can improve
203 !% communication speeds.
204 !%End
205 call parse_variable(namespace, 'ReorderRanks', .false., mc%reorder_ranks)
207 call messages_obsolete_variable(namespace, 'ParallelizationStrategy')
208 call messages_obsolete_variable(namespace, 'ParallelizationGroupRanks')
209
210 do ipar = 1, p_strategy_max
211 default(ipar) = par_no
212 if (bitand(default_mask, ibset(0, ipar - 1)) /= 0) then
213 default(ipar) = par_auto
214 end if
215 end do
217 !%Variable ParDomains
218 !%Type integer
219 !%Default auto
220 !%Section Execution::Parallelization
221 !%Description
222 !% This variable controls the number of processors used for the
223 !% parallelization in domains.
224 !% The special value <tt>auto</tt>, the default, lets Octopus
225 !% decide how many processors will be assigned for this
226 !% strategy. To disable parallelization in domains, you can use
227 !% <tt>ParDomains = no</tt> (or set the number of processors to
228 !% 1).
229 !%
230 !% The total number of processors required is the multiplication
231 !% of the processors assigned to each parallelization strategy.
232 !%Option auto -1
233 !% The number of processors is assigned automatically.
234 !%Option no 0
235 !% This parallelization strategy is not used.
236 !%End
237 call parse_variable(namespace, 'ParDomains', default(p_strategy_domains), parse(p_strategy_domains))
239 !%Variable ParStates
240 !%Type integer
241 !%Section Execution::Parallelization
242 !%Description
243 !% This variable controls the number of processors used for the
244 !% parallelization in states. The special value <tt>auto</tt> lets
245 !% Octopus decide how many processors will be assigned for this
246 !% strategy. To disable parallelization in states, you can use
247 !% <tt>ParStates = no</tt> (or set the number of processors to 1).
248 !%
249 !% The default value depends on the <tt>CalculationMode</tt>. For
250 !% <tt>CalculationMode = td</tt> the default is <tt>auto</tt>, while
251 !% for for other modes the default is <tt>no</tt>.
252 !%
253 !% The total number of processors required is the multiplication
254 !% of the processors assigned to each parallelization strategy.
255 !%Option auto -1
256 !% The number of processors is assigned automatically.
257 !%Option no 0
258 !% This parallelization strategy is not used.
259 !%End
260 call parse_variable(namespace, 'ParStates', default(p_strategy_states), parse(p_strategy_states))
261
262 !%Variable ParKPoints
263 !%Type integer
264 !%Default auto
265 !%Section Execution::Parallelization
266 !%Description
267 !% This variable controls the number of processors used for the
268 !% parallelization in K-Points and/or spin.
269 !% The special value <tt>auto</tt> lets Octopus decide how many processors will be
270 !% assigned for this strategy. To disable parallelization in
271 !% KPoints, you can use <tt>ParKPoints = no</tt> (or set the
272 !% number of processors to 1).
273 !%
274 !% The total number of processors required is the multiplication
275 !% of the processors assigned to each parallelization strategy.
276 !%Option auto -1
277 !% The number of processors is assigned automatically.
278 !%Option no 0
279 !% This parallelization strategy is not used.
280 !%End
281 call parse_variable(namespace, 'ParKPoints', default(p_strategy_kpoints), parse(p_strategy_kpoints))
282
283 !%Variable ParOther
284 !%Type integer
285 !%Default auto
286 !%Section Execution::Parallelization
287 !%Description
288 !% This variable controls the number of processors used for the
289 !% 'other' parallelization mode, that is CalculatioMode
290 !% dependent. For <tt>CalculationMode = casida</tt>, it means
291 !% parallelization in electron-hole pairs.
292 !%
293 !% The special value <tt>auto</tt>,
294 !% the default, lets Octopus decide how many processors will be
295 !% assigned for this strategy. To disable parallelization in
296 !% Other, you can use <tt>ParOther = no</tt> (or set the
297 !% number of processors to 1).
298 !%
299 !% The total number of processors required is the multiplication
300 !% of the processors assigned to each parallelization strategy.
301 !%Option auto -1
302 !% The number of processors is assigned automatically.
303 !%Option no 0
304 !% This parallelization strategy is not used.
305 !%End
306 call parse_variable(namespace, 'ParOther', default(p_strategy_other), parse(p_strategy_other))
307
308 do ipar = 1, p_strategy_max
309 if (parse(ipar) == par_no) parse(ipar) = 1
310 end do
311
312 call strategy()
313
314 mc%have_slaves = .false.
315
316 if (mc%par_strategy /= p_strategy_serial) then
317 safe_allocate(mc%group_sizes(1:p_strategy_max))
318
319 mc%group_sizes = 1
320
321 do ipar = 1, p_strategy_max
322 if (multicomm_strategy_is_parallel(mc, ipar)) then
323 mc%group_sizes(ipar) = parse(ipar)
324 else if (parse(ipar) /= 1) then
325 call messages_write('Ignoring specification for ' // par_types(ipar))
326 call messages_new_line()
327 call messages_write('This parallelization strategy is not available.')
328 call messages_warning()
329 end if
330 end do
331
332 call assign_nodes()
333
334
335 !%Variable ParallelizationNumberSlaves
336 !%Type integer
337 !%Default 0
338 !%Section Execution::Parallelization
339 !%Description
340 !% Slaves are nodes used for task parallelization. The number of
341 !% such nodes is given by this variable multiplied by the number
342 !% of domains used in domain parallelization.
343 !%End
344 call parse_variable(namespace, 'ParallelizationNumberSlaves', 0, num_slaves)
345
346 ! the slaves must be defined at a certain parallelization level, for the moment this is state parallelization.
347 slave_level = p_strategy_states
348 mc%have_slaves = (num_slaves > 0)
349
350 if (mc%have_slaves) then
351 call messages_experimental('Task parallelization')
352 end if
353
354 ! clear parallel strategies that were available but will not be used
355 do ii = 1, p_strategy_max
356 if (mc%group_sizes(ii) == 1) mc%par_strategy = ibclr(mc%par_strategy, ii - 1)
357 end do
358
359 ! reset
360 call sanity_check()
361 end if
362
363 call group_comm_create()
364
365 call messages_print_with_emphasis(namespace=namespace)
366
367 pop_sub(multicomm_init)
368
369 contains
370
371 ! ---------------------------------------------------------
372 subroutine strategy()
373 integer :: jj, ipar
374
375 push_sub(multicomm_init.strategy)
376
377 if (base_grp%size > 1) then
378
379 mc%par_strategy = 0
380
381 do ipar = 1, p_strategy_max
382 if (parse(ipar) == par_auto .or. parse(ipar) > 1) then
383 mc%par_strategy = ibset(mc%par_strategy, ipar - 1)
384 end if
385 end do
386
387 if (mc%par_strategy /= bitand(mc%par_strategy, parallel_mask)) then
388 call messages_write('Parallelization strategies unavailable for this run mode are being discarded.')
389 call messages_warning()
390 end if
391
392 mc%par_strategy = bitand(mc%par_strategy, parallel_mask)
393
394 if (mc%par_strategy == p_strategy_serial) then
395 message(1) = "More than one node is available, but this run mode cannot run with the requested parallelization."
396 message(2) = "Please select a parallelization strategy compatible with"
397 jj = 2
398 do ii = 1, n_par_types
399 if (bitand(parallel_mask, 2**(ii - 1)) /= 0) then
400 jj = jj + 1
401 write(message(jj), '(2a)') " -> ", par_types(ii)
402 end if
403 end do
404 jj=jj+1
405 write(message(jj),'(a,i6)') "mc%par_strategy is : ",mc%par_strategy
406 call messages_fatal(jj, only_root_writes = .true.)
407 end if
408 else
409 mc%par_strategy = p_strategy_serial
410 end if
411
412 mc%nthreads = 1
413#if defined(HAVE_OPENMP)
414 !$omp parallel
415 !$omp master
416 mc%nthreads = omp_get_num_threads()
417 !$omp end master
418 !$omp end parallel
419#endif
420
421 if (mc%par_strategy == p_strategy_serial .and. mc%nthreads == 1) then
422 message(1) = "Info: Octopus will run in *serial*"
423 call messages_info(1, namespace=namespace)
424 else
425 write(message(1),'(a)') 'Info: Octopus will run in *parallel*'
426 write(message(2),'(a)') ''
427 write(message(3),'(a, i8)') ' Number of processes :', base_grp%size
428 write(message(4),'(a, i8)') ' Number of threads per process :', mc%nthreads
429 write(message(5),'(a)') ''
430 call messages_info(5, namespace=namespace)
431 end if
432
433 pop_sub(multicomm_init.strategy)
434 end subroutine strategy
435
436 ! ---------------------------------------------------------
437
438 subroutine assign_nodes()
439 integer :: ii, nn, kk, n_divisors, divisors(1:50)
440 integer(int64) :: n_group_max(1:p_strategy_max)
441
443
444 ! this is the maximum number of processors in each group
445 n_group_max(1:p_strategy_max) = max(index_range(1:p_strategy_max), 1_int64)
446 do kk = 1, p_strategy_max
447 if (.not. multicomm_strategy_is_parallel(mc, kk)) n_group_max(kk) = 1
448 end do
449
450 call messages_write('Debug info: Allowable group ranks:', new_line = .true.)
451 do kk = 1, p_strategy_max
452 call messages_write(par_types(kk), fmt = '2x,a12,":",1x')
453 call messages_write(n_group_max(kk), new_line = .true.)
454 end do
455 call messages_info(debug_only=.true.)
456
457 nn = mc%n_node
458
459 ! first loop, check the processors assigned by the user
460 do ipar = p_strategy_max, 1, -1
461
462 if (mc%group_sizes(ipar) == par_auto) cycle
463
464 if (mc%group_sizes(ipar) > n_group_max(ipar)) then
465 call messages_write('The number of processors specified for '//par_types(ipar)//'(')
466 call messages_write(mc%group_sizes(ipar))
467 call messages_write(')', new_line = .true.)
468 call messages_write('is larger than the degrees of freedom for that level (')
469 call messages_write(n_group_max(ipar))
470 call messages_write(').')
471 call messages_warning()
472 end if
473
474 if (mod(nn, mc%group_sizes(ipar)) /= 0) then
475 call messages_write('The number of processors specified for '//par_types(ipar)//'(')
476 call messages_write(mc%group_sizes(ipar))
477 call messages_write(')', new_line = .true.)
478 call messages_write('is not a divisor of the number of processors (')
479 call messages_write(mc%n_node)
480 call messages_write(').')
481 call messages_fatal()
482 end if
483
484 nn = nn/mc%group_sizes(ipar)
485
486 end do
487
488 ! second loop, now assign the rest automatically
489 do ipar = p_strategy_max, 1, -1
490
491 if (mc%group_sizes(ipar) /= par_auto) cycle
492
493 n_divisors = ubound(divisors, dim = 1)
494 call get_divisors(nn, n_divisors, divisors)
495
496 mc%group_sizes(ipar) = nn
497 do ii = 2, n_divisors
498 if (divisors(ii) > n_group_max(ipar)) then
499 mc%group_sizes(ipar) = divisors(ii - 1)
500 exit
501 end if
502 end do
503
504 nn = nn/mc%group_sizes(ipar)
505
506 end do
507
509 end subroutine assign_nodes
510
511
512 ! ---------------------------------------------------------
514 subroutine sanity_check()
515 real(real64) :: frac
516 integer :: ii, kk
517 integer(int64) :: jj, n_max
518 integer :: real_group_sizes(1:MAX_INDEX)
519
521
522 if (num_slaves > 0) then
523
524 if (mc%group_sizes(slave_level) < num_slaves + 1) then
525 message(1) = 'Too many nodes assigned to task parallelization.'
526 call messages_fatal(1)
527 end if
528
529 write(message(1),'(a,i6)') 'Info: Number of slaves nodes :', &
530 num_slaves*product(mc%group_sizes(1:slave_level - 1))
531 call messages_info(1)
532
533 end if
534
535 ! print out some info
536 ii = 0
537 do kk = p_strategy_max, 1, -1
538 real_group_sizes(kk) = mc%group_sizes(kk)
539 if (.not. multicomm_strategy_is_parallel(mc, kk)) cycle
540 ii = ii + 1
541 if (kk == slave_level) real_group_sizes(kk) = real_group_sizes(kk) - num_slaves
542 write(message(ii),'(3a,i6,a,i12,a)') 'Info: Number of nodes in ', &
543 par_types(kk), ' group:', real_group_sizes(kk), ' (', index_range(kk), ')'
544 end do
545 call messages_info(ii)
546
547 ! do we have the correct number of processors
548 if (product(mc%group_sizes(1:p_strategy_max)) /= base_grp%size) then
549 write(message(1),'(a)') 'Inconsistent number of processors:'
550 write(message(2),'(a,i6)') ' MPI processes = ', base_grp%size
551 write(message(3),'(a,i6)') ' Required processes = ', product(mc%group_sizes(1:p_strategy_max))
552 message(4) = ''
553 message(5) = 'You probably have a problem in the ParDomains, ParStates, ParKPoints or ParOther.'
554 call messages_fatal(5, only_root_writes = .true.)
555 end if
556
557 if (any(real_group_sizes(1:p_strategy_max) > index_range(1:p_strategy_max))) then
558 message(1) = "Could not distribute nodes in parallel job. Most likely you are trying to"
559 message(2) = "use too many nodes for the job."
560 call messages_fatal(2, only_root_writes = .true.)
561 end if
562
563 if (any(index_range(1:p_strategy_max) / real_group_sizes(1:p_strategy_max) < min_range(1:p_strategy_max) .and. &
564 real_group_sizes(1:p_strategy_max) > 1)) then
565 message(1) = "I have fewer elements in a parallel group than recommended."
566 message(2) = "Maybe you should reduce the number of nodes."
567 call messages_warning(2)
568 end if
569
570 ! calculate fraction of idle time
571 frac = m_one
572 do ii = 1, p_strategy_max
573 n_max = ceiling(real(index_range(ii), real64) / real(real_group_sizes(ii)), real64)
574 jj = n_max*real_group_sizes(ii)
575 frac = frac*(m_one - real(jj - index_range(ii), real64) / real(jj, real64) )
576 end do
577
578 write(message(1), '(a,f5.2,a)') "Info: Octopus will waste at least ", &
579 (m_one - frac)*100.0_real64, "% of computer time."
580 if (frac < 0.8_real64) then
581 message(2) = "Usually a number of processors which is a multiple of small primes is best."
582 call messages_warning(2)
583 else
584 call messages_info(1)
585 end if
586
588 end subroutine sanity_check
589
590 ! ---------------------------------------------------------
591 subroutine group_comm_create()
592#if defined(HAVE_MPI)
593 logical :: dim_mask(MAX_INDEX)
594 integer :: i_strategy, irank
595 logical :: reorder, periodic_mask(MAX_INDEX)
596 integer :: coords(MAX_INDEX)
597 type(mpi_comm) :: new_comm
598 integer :: new_comm_size
599 character(len=6) :: node_type
600 type(mpi_grp_t) :: reorder_grp
601 type(mpi_group) :: base_group, reorder_group
602 integer :: ranks(base_grp%size)
603 integer :: ii, jj, kk, ll, nn
604 type(mpi_comm) :: reorder_comm
605#endif
606
608
609 mc%node_type = p_master
610
611 safe_allocate(mc%group_comm(1:p_strategy_max))
612 safe_allocate(mc%who_am_i(1:p_strategy_max))
613
614 mc%group_comm = mpi_comm_undefined
615 mc%who_am_i = 0
616
617#if defined(HAVE_MPI)
618 mc%full_comm = mpi_comm_null
619 mc%slave_intercomm = mpi_comm_null
620 if (mc%par_strategy /= p_strategy_serial) then
621 if (mc%reorder_ranks) then
622 ! first, reorder the ranks
623 ! this is done to get a column-major ordering of the ranks in the
624 ! Cartesian communicator, since they a ordered row-major otherwise
625 call mpi_comm_group(base_grp%comm, base_group, mpi_err)
626 if (mpi_err /= mpi_success) then
627 message(1) = "Error in getting MPI group!"
628 call messages_fatal(1)
629 end if
630 ! now transpose the hypercube => get rank numbers in column-major order
631 nn = 1
632 do ii = 1, mc%group_sizes(1)
633 do jj = 1, mc%group_sizes(2)
634 do kk = 1, mc%group_sizes(3)
635 do ll = 1, mc%group_sizes(4)
636 ranks(nn) = (ll-1)*mc%group_sizes(3)*mc%group_sizes(2)*mc%group_sizes(1) &
637 + (kk-1)*mc%group_sizes(2)*mc%group_sizes(1) &
638 + (jj-1)*mc%group_sizes(1) + ii - 1
639 nn = nn + 1
640 end do
641 end do
642 end do
643 end do
644 call mpi_group_incl(base_group, base_grp%size, ranks, reorder_group, mpi_err)
645 if (mpi_err /= mpi_success) then
646 message(1) = "Error in creating MPI group!"
647 call messages_fatal(1)
648 end if
649 ! now get the reordered communicator
650 call mpi_comm_create(base_grp%comm, reorder_group, reorder_comm, mpi_err)
651 if (mpi_err /= mpi_success) then
652 message(1) = "Error in creating reordered communicator!"
653 call messages_fatal(1)
654 end if
655 call mpi_grp_init(reorder_grp, reorder_comm)
656 else
657 call mpi_grp_copy(reorder_grp, base_grp)
658 end if
659
660 ! Multilevel parallelization is organized in a hypercube. We
661 ! use an MPI Cartesian topology to generate the communicators
662 ! that correspond to each level.
663
664 ! create the topology
665 periodic_mask = .false.
666 reorder = .true.
667
668 ! The domain and states dimensions have to be periodic (2D torus)
669 ! in order to circulate matrix blocks.
672
673 ! We allow reordering of ranks.
674 call mpi_cart_create(reorder_grp%comm, p_strategy_max, mc%group_sizes, periodic_mask, reorder, mc%full_comm, mpi_err)
675
676 call mpi_comm_rank(mc%full_comm, mc%full_comm_rank, mpi_err)
678 ! get the coordinates of the current processor
679 call mpi_cart_coords(mc%full_comm, mc%full_comm_rank, p_strategy_max, coords, mpi_err)
680
681 ! find out what type of node this is
682 if (coords(slave_level) >= mc%group_sizes(slave_level) - num_slaves) then
683 mc%node_type = p_slave
684 end if
685
686 if (mc%node_type == p_master) then
687 mc%group_sizes(slave_level) = mc%group_sizes(slave_level) - num_slaves
688 else
689 mc%group_sizes(slave_level) = num_slaves
690 end if
691
692 call mpi_comm_split(mc%full_comm, mc%node_type, mc%full_comm_rank, new_comm, mpi_err)
693 assert(new_comm /= mpi_comm_null)
694 call mpi_comm_size(new_comm, new_comm_size, mpi_err)
695
696 reorder = .false.
697 if (product(mc%group_sizes(:)) /= new_comm_size) then
698 write(stderr,*) 'node ', mpi_world%rank, ': mc%group_sizes = ', mc%group_sizes, ' new_comm_size = ', new_comm_size
699 call mpi_world%barrier()
700 assert(product(mc%group_sizes(:)) == new_comm_size)
701 end if
702 call mpi_cart_create(new_comm, p_strategy_max, mc%group_sizes, periodic_mask, reorder, mc%master_comm, mpi_err)
703 assert(mc%master_comm /= mpi_comm_null)
704
705 call mpi_comm_free(new_comm, mpi_err)
706
707 call mpi_comm_rank(mc%master_comm, mc%master_comm_rank, mpi_err)
708
709 ! The "lines" of the Cartesian grid.
710 ! Initialize all the communicators, even if they are not parallelized
711 do i_strategy = 1, p_strategy_max
712 dim_mask = .false.
713 dim_mask(i_strategy) = .true.
714 call mpi_cart_sub(mc%master_comm, dim_mask, mc%group_comm(i_strategy), mpi_err)
715 call mpi_comm_rank(mc%group_comm(i_strategy), mc%who_am_i(i_strategy), mpi_err)
716 end do
717
718 ! The domain-state "planes" of the grid (the ones with periodic dimensions).
719 dim_mask = .false.
720 dim_mask(p_strategy_domains) = .true.
721 dim_mask(p_strategy_states) = .true.
722 call mpi_cart_sub(mc%master_comm, dim_mask, mc%dom_st_comm, mpi_err)
723
724 ! The state-kpoints "planes" of the grid
725 dim_mask = .false.
726 dim_mask(p_strategy_states) = .true.
727 dim_mask(p_strategy_kpoints) = .true.
728 call mpi_cart_sub(mc%master_comm, dim_mask, mc%st_kpt_comm, mpi_err)
729
730 ! The domains-states-kpoints "cubes" of the grid
731 dim_mask = .false.
732 dim_mask(p_strategy_domains) = .true.
733 dim_mask(p_strategy_states) = .true.
734 dim_mask(p_strategy_kpoints) = .true.
735 call mpi_cart_sub(mc%master_comm, dim_mask, mc%dom_st_kpt_comm, mpi_err)
736
737 if (num_slaves > 0) call create_slave_intercommunicators()
738
739 call create_intranode_communicator(base_grp, mc%intranode_grp, mc%internode_grp)
740 else
741 ! we initialize these communicators so we can use them even in serial
742 mc%group_comm = base_grp%comm
743 mc%who_am_i = 0
744 mc%master_comm = base_grp%comm
745 mc%dom_st_comm = base_grp%comm
746 mc%st_kpt_comm = base_grp%comm
747 mc%dom_st_kpt_comm = base_grp%comm
748 call mpi_grp_copy(mc%intranode_grp, base_grp)
749 call mpi_grp_copy(mc%internode_grp, base_grp)
750 end if
751
752 ! This is temporary debugging information.
753 if (debug%info .and. mc%par_strategy /= p_strategy_serial) then
754 write(message(1),'(a)') 'Debug: MPI Task Assignment to MPI Groups'
755 write(message(2),'(5a10)') 'World', 'Domains', 'States', 'K-Points', 'Other'
756 call messages_info(1)
757
758 if (mc%node_type == p_slave) then
759 node_type = "slave"
760 else
761 node_type = "master"
762 end if
763 do irank = 0, mpi_world%size - 1
764 if (mpi_world%rank == irank) then
765 write(message(1),'(5i10,5x,a)') mpi_world%rank, mc%who_am_i(p_strategy_domains), mc%who_am_i(p_strategy_states), &
766 mc%who_am_i(p_strategy_kpoints), mc%who_am_i(p_strategy_other), trim(node_type)
767 call messages_info(1, all_nodes = .true.)
768 end if
769 call mpi_world%barrier()
770 end do
771 end if
772#endif
773
775 end subroutine group_comm_create
776
777 ! -----------------------------------------------------
778#ifdef HAVE_MPI
780 integer :: remote_leader
781 integer :: tag
782 integer :: coords(MAX_INDEX)
783
785
786 ! create the intercommunicators to communicate with slaves
787
788 ! get the coordinates of the current processor
789 call mpi_cart_coords(mc%full_comm, mc%full_comm_rank, p_strategy_max, coords, mpi_err)
790
791 !now get the rank of the remote_leader
792 if (mc%node_type == p_slave) then
793 coords(slave_level) = 0
794 else
795 coords(slave_level) = mc%group_sizes(slave_level)
796 end if
797 call mpi_cart_rank(mc%full_comm, coords, remote_leader, mpi_err)
798
799 ! now create the intercommunicator
800 tag = coords(p_strategy_domains)
801 call mpi_intercomm_create(mc%group_comm(slave_level), 0, base_grp%comm, remote_leader, tag, mc%slave_intercomm, mpi_err)
802
805#endif
806 end subroutine multicomm_init
807
808
809 ! ---------------------------------------------------------
810 subroutine multicomm_end(mc)
811 type(multicomm_t), intent(inout) :: mc
812
813#if defined(HAVE_MPI)
814 integer :: ii
815#endif
816
817 push_sub(multicomm_end)
818
819 if (mc%par_strategy /= p_strategy_serial) then
820#if defined(HAVE_MPI)
821 ! Delete communicators.
822 do ii = 1, p_strategy_max
823 ! initialized even if not parallelized
824 call mpi_comm_free(mc%group_comm(ii), mpi_err)
825 end do
826 call mpi_comm_free(mc%dom_st_comm, mpi_err)
827 call mpi_comm_free(mc%st_kpt_comm, mpi_err)
828 call mpi_comm_free(mc%dom_st_kpt_comm, mpi_err)
829 call mpi_comm_free(mc%full_comm, mpi_err)
830 call mpi_comm_free(mc%master_comm, mpi_err)
831
832 if (multicomm_have_slaves(mc)) call mpi_comm_free(mc%slave_intercomm, mpi_err)
833
834#endif
835 end if
836
837 safe_deallocate_a(mc%group_sizes)
838 safe_deallocate_a(mc%group_comm)
839 safe_deallocate_a(mc%who_am_i)
840
841 pop_sub(multicomm_end)
842 end subroutine multicomm_end
843
844
845 ! ---------------------------------------------------------
846 logical pure function multicomm_strategy_is_parallel(mc, level) result(rr)
847 type(multicomm_t), intent(in) :: mc
848 integer, intent(in) :: level
849
850 rr = bitand(mc%par_strategy, 2**(level - 1)) /= 0
851
853
854
855 ! ---------------------------------------------------------
862#if defined(HAVE_MPI)
863
864 subroutine multicomm_create_all_pairs(mpi_grp, ap)
865 type(mpi_grp_t), intent(in) :: mpi_grp
866 type(multicomm_all_pairs_t), intent(out) :: ap
867
868 integer :: grp_size, rounds, ir, in
869
870 push_sub(create_all_pairs)
871
872 ap%grp = mpi_grp
873 grp_size = mpi_grp%size
874
875 ! Number of rounds.
876 if (mod(grp_size, 2) == 0) then
877 rounds = grp_size - 1
878 else
879 rounds = grp_size
880 end if
881 ap%rounds = rounds
882
883 ! Calculate schedule.
884 safe_allocate(ap%schedule(0:grp_size - 1, 1:rounds))
885 do ir = 1, rounds
886 do in = 0, grp_size - 1
887 ap%schedule(in, ir) = get_partner(in + 1, ir) - 1
888 end do
889 end do
890
891 pop_sub(create_all_pairs)
892
893 contains
894
895 ! ---------------------------------------------------------
897 integer pure function get_partner(in, ir)
898 integer, intent(in) :: in, ir
899
900 ! No PUSH SUB, called too often.
901
902 if (mod(grp_size, 2) == 0) then
903 get_partner = get_partner_even(grp_size, in - 1, ir - 1) + 1
904 else
905 get_partner = get_partner_odd(grp_size, in - 1, ir - 1) + 1
906 end if
907
908 end function get_partner
909
910 ! ---------------------------------------------------------
911 integer pure function get_partner_even(grp_size, ii, rr) result(pp)
912 integer, intent(in) :: grp_size, ii, rr
913
914 integer :: mm
915
916 ! No PUSH SUB, called too often.
917
918 mm = grp_size / 2
919
920 if (ii == 0) then
921 pp = rr + 1
922 elseif (ii == rr + 1) then
923 pp = 0
924 else
925 ! I never know when to use which remainder function, but here
926 ! it has to be the modulo one. Do not change that!
927 pp = modulo(2 * rr - ii + 1, 2 * mm - 1) + 1
928 end if
929
930 end function get_partner_even
931
932 ! ---------------------------------------------------------
933 integer pure function get_partner_odd(grp_size, ii, rr) result(pp)
934 integer, intent(in) :: grp_size, ii, rr
935
936 integer :: mm
937
938 ! No PUSH SUB, called too often.
939
940 mm = (grp_size + 1) / 2
941
942 pp = get_partner_even(grp_size + 1, ii, rr)
943
944 if (pp == 2 * mm - 1) then
945 pp = ii
946 end if
947
948 end function get_partner_odd
949
951#endif
952
956 subroutine multicomm_divide_range(nobjs, nprocs, istart, ifinal, lsize, scalapack_compat)
957 integer, intent(in) :: nobjs
958 integer, intent(in) :: nprocs
959 integer, intent(out) :: istart(:)
960 integer, intent(out) :: ifinal(:)
961 integer, optional, intent(out) :: lsize(:)
962 logical, optional, intent(in) :: scalapack_compat
963
964 integer :: ii, jj, rank
965 logical :: scalapack_compat_
966 integer :: nbl, size
967
968 scalapack_compat_ = optional_default(scalapack_compat, .false.)
969#ifndef HAVE_SCALAPACK
970 scalapack_compat_ = .false.
971#endif
972
973 if (scalapack_compat_) then
974 nbl = nobjs/nprocs
975 if (mod(nobjs, nprocs) /= 0) nbl = nbl + 1
976
977 istart(1) = 1
978 do rank = 1, nprocs
979#ifdef HAVE_SCALAPACK
980 size = numroc(nobjs, nbl, rank - 1, 0, nprocs)
981#endif
982 if (size > 0) then
983 if (rank > 1) istart(rank) = ifinal(rank - 1) + 1
984 ifinal(rank) = istart(rank) + size - 1
985 else
986 istart(rank) = 1
987 ifinal(rank) = 0
988 end if
989 end do
990 else
991
992 if (nprocs <= nobjs) then
993
994 ! procs are assigned to groups by round robin
995 do rank = 0, nprocs - 1
996 jj = nobjs / nprocs
997 ii = nobjs - jj*nprocs
998 if (ii > 0 .and. rank < ii) then
999 jj = jj + 1
1000 istart(rank + 1) = rank*jj + 1
1001 ifinal(rank + 1) = istart(rank + 1) + jj - 1
1002 else
1003 ifinal(rank + 1) = nobjs - (nprocs - rank - 1)*jj
1004 istart(rank + 1) = ifinal(rank + 1) - jj + 1
1005 end if
1006 end do
1007
1008 else
1009 do ii = 1, nprocs
1010 if (ii <= nobjs) then
1011 istart(ii) = ii
1012 ifinal(ii) = ii
1013 else
1014 istart(ii) = 1
1015 ifinal(ii) = 0
1016 end if
1017 end do
1018 end if
1019 end if
1020
1021 if (present(lsize)) then
1022 lsize(1:nprocs) = ifinal(1:nprocs) - istart(1:nprocs) + 1
1023 assert(sum(lsize(1:nprocs)) == nobjs)
1024 end if
1025
1026 end subroutine multicomm_divide_range
1027
1028 ! ---------------------------------------------------------
1031 ! THREADSAFE
1032 subroutine multicomm_divide_range_omp(nobjs, ini, nobjs_loc)
1033 integer, intent(in) :: nobjs
1034 integer, intent(out) :: ini
1035 integer, intent(out) :: nobjs_loc
1036
1037 integer :: rank
1038#ifdef HAVE_OPENMP
1039 integer, allocatable :: istart(:), ifinal(:), lsize(:)
1040 integer :: nthreads
1041#endif
1043 ! no push_sub, threadsafe
1044 rank = 1
1045#ifdef HAVE_OPENMP
1046 nthreads = omp_get_num_threads()
1047 allocate(istart(1:nthreads))
1048 allocate(ifinal(1:nthreads))
1049 allocate(lsize(1:nthreads))
1050 call multicomm_divide_range(nobjs, nthreads, istart, ifinal, lsize)
1051 rank = 1 + omp_get_thread_num()
1052 ini = istart(rank)
1053 nobjs_loc = lsize(rank)
1054 deallocate(istart)
1055 deallocate(ifinal)
1056 deallocate(lsize)
1057#else
1058 ini = 1
1059 nobjs_loc = nobjs
1060#endif
1061
1062 end subroutine multicomm_divide_range_omp
1063
1064 ! ---------------------------------------------------------
1065
1066 logical pure function multicomm_is_slave(this) result(slave)
1067 type(multicomm_t), intent(in) :: this
1068
1069 slave = this%node_type == p_slave
1070 end function multicomm_is_slave
1071
1072 ! ---------------------------------------------------------
1073
1074 logical pure function multicomm_have_slaves(this) result(have_slaves)
1075 type(multicomm_t), intent(in) :: this
1076
1077 have_slaves = this%have_slaves
1078 end function multicomm_have_slaves
1079
1080end module multicomm_oct_m
1081
1082
1083!! Local Variables:
1084!! mode: f90
1085!! coding: utf-8
1086!! End:
integer pure function get_partner(in, ir)
Those are from the paper cited above.
Definition: multicomm.F90:984
subroutine strategy()
Definition: multicomm.F90:466
subroutine group_comm_create()
Definition: multicomm.F90:678
integer pure function get_partner_even(grp_size, ii, rr)
Definition: multicomm.F90:998
subroutine create_slave_intercommunicators()
Definition: multicomm.F90:866
subroutine sanity_check()
check if a balanced distribution of nodes will be used
Definition: multicomm.F90:601
subroutine assign_nodes()
Definition: multicomm.F90:525
integer pure function get_partner_odd(grp_size, ii, rr)
Definition: multicomm.F90:1020
This module contains interfaces for BLACS routines Interfaces are from http:
Definition: blacs.F90:27
This module handles the calculation mode.
integer, parameter, public p_strategy_max
integer, parameter, public p_strategy_kpoints
parallelization in k-points
integer, parameter, public p_strategy_other
something else like e-h pairs
integer, parameter, public p_strategy_domains
parallelization in domains
integer, parameter, public p_strategy_serial
single domain, all states, k-points on a single processor
integer, parameter, public p_strategy_states
parallelization in states
type(debug_t), save, public debug
Definition: debug.F90:156
real(real64), parameter, public m_one
Definition: global.F90:189
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:921
character(len=512), private msg
Definition: messages.F90:166
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:538
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1046
subroutine, public messages_new_line()
Definition: messages.F90:1135
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1086
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
This module contains some common usage patterns of MPI routines.
Definition: mpi_lib.F90:115
subroutine, public create_intranode_communicator(base_grp, intranode_grp, internode_grp)
Definition: mpi_lib.F90:209
subroutine mpi_grp_copy(mpi_grp_out, mpi_grp_in)
Definition: mpi.F90:403
type(mpi_comm), parameter, public mpi_comm_undefined
used to indicate a communicator has not been initialized
Definition: mpi.F90:136
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
subroutine mpi_grp_init(grp, comm)
Initialize MPI group instance.
Definition: mpi.F90:346
integer, public mpi_err
used to store return values of mpi calls
Definition: mpi.F90:269
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
logical pure function, public multicomm_strategy_is_parallel(mc, level)
Definition: multicomm.F90:933
integer, parameter, public p_slave
Definition: multicomm.F90:181
subroutine, public multicomm_divide_range(nobjs, nprocs, istart, ifinal, lsize, scalapack_compat)
Divide the range of numbers [1, nobjs] between nprocs processors.
Definition: multicomm.F90:1043
subroutine, public multicomm_all_pairs_copy(apout, apin)
Definition: multicomm.F90:246
subroutine, public multicomm_end(mc)
Definition: multicomm.F90:897
integer, parameter, public par_no
Definition: multicomm.F90:185
logical pure function, public multicomm_have_slaves(this)
Definition: multicomm.F90:1149
logical pure function, public multicomm_is_slave(this)
Definition: multicomm.F90:1141
subroutine, public multicomm_init(mc, namespace, base_grp, mode_para, n_node, index_range, min_range)
create index and domain communicators
Definition: multicomm.F90:268
subroutine, public multicomm_create_all_pairs(mpi_grp, ap)
This routine uses the one-factorization (or near-one-factorization of a complete graph to construct a...
Definition: multicomm.F90:951
subroutine, public multicomm_divide_range_omp(nobjs, ini, nobjs_loc)
Function to divide the range of numbers from 1 to nobjs between all available threads with OpenMP.
Definition: multicomm.F90:1119
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:118
subroutine, public get_divisors(nn, n_divisors, divisors)
Definition: utils.F90:168
This is defined even when running serial.
Definition: mpi.F90:142
An all-pairs communication schedule for a given group.
Definition: multicomm.F90:235
Stores all communicators and groups.
Definition: multicomm.F90:206
int true(void)
subroutine parse()
Definition: xc.F90:381