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