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)
77 multicomm_create_all_pairs, &
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 type(mpi_grp_t), public :: master_grp
129 type(mpi_grp_t), public :: base_grp
130
131 integer :: nthreads
132 integer :: node_type
133 logical :: have_slaves
134
135 type(MPI_Comm) :: full_comm = mpi_comm_undefined
136 integer :: full_comm_rank
137 type(MPI_Comm), public :: master_comm = mpi_comm_undefined
138 integer :: master_comm_rank
139 type(MPI_Comm), public :: slave_intercomm = mpi_comm_undefined
140
141 logical :: reorder_ranks
142 end type multicomm_t
143
146 private
147 type(mpi_grp_t) :: grp
148 integer :: rounds
149 integer, allocatable, public :: schedule(:, :)
150 end type multicomm_all_pairs_t
151
152contains
153
154 ! ---------------------------------------------------------
155 subroutine multicomm_all_pairs_copy(apout, apin)
156 type(multicomm_all_pairs_t), intent(inout) :: apout
157 type(multicomm_all_pairs_t), intent(in) :: apin
158
160
161 call mpi_grp_copy(apout%grp, apin%grp)
162 apout%rounds = apin%rounds
163 if (allocated(apin%schedule)) then
164 safe_allocate(apout%schedule(1:size(apin%schedule, 1), 1:size(apin%schedule, 2)))
165 apout%schedule = apin%schedule
166 else
167 if (allocated(apout%schedule)) then
168 safe_deallocate_a(apout%schedule)
169 end if
170 end if
171
173 end subroutine multicomm_all_pairs_copy
174
175 ! ---------------------------------------------------------
177 subroutine multicomm_init(mc, namespace, base_grp, mode_para, n_node, index_range, min_range)
178 type(multicomm_t), intent(out) :: mc
179 type(namespace_t), intent(in) :: namespace
180 type(mpi_grp_t), intent(in) :: base_grp
181 type(calc_mode_par_t), intent(in) :: mode_para
182 integer, intent(in) :: n_node
183 integer(int64), intent(inout) :: index_range(:)
184 integer, intent(in) :: min_range(:)
185
186 integer :: ii, num_slaves, slave_level, ipar
187 integer :: parse(1:p_strategy_max), default(1:p_strategy_max)
188 integer :: parallel_mask, default_mask
189
190 push_sub(multicomm_init)
192 mc%base_grp = base_grp
193
194 mc%n_node = n_node
195 parallel_mask = mode_para%parallel_mask()
196 default_mask = mode_para%default_parallel_mask()
197
198 call messages_print_with_emphasis(msg="Parallelization", namespace=namespace)
199
200 !%Variable ReorderRanks
201 !%Default no
202 !%Type logical
203 !%Section Execution::Parallelization
204 !%Description
205 !% This variable controls whether the ranks are reorganized to have a more
206 !% compact distribution with respect to domain parallelization which needs
207 !% to communicate most often. Depending on the system, this can improve
208 !% communication speeds.
209 !%End
210 call parse_variable(namespace, 'ReorderRanks', .false., mc%reorder_ranks)
211
212 call messages_obsolete_variable(namespace, 'ParallelizationStrategy')
213 call messages_obsolete_variable(namespace, 'ParallelizationGroupRanks')
215 do ipar = 1, p_strategy_max
216 default(ipar) = par_no
217 if (bitand(default_mask, ibset(0, ipar - 1)) /= 0) then
218 default(ipar) = par_auto
219 end if
220 end do
222 !%Variable ParDomains
223 !%Type integer
224 !%Default auto
225 !%Section Execution::Parallelization
226 !%Description
227 !% This variable controls the number of processors used for the
228 !% parallelization in domains.
229 !% The special value <tt>auto</tt>, the default, lets Octopus
230 !% decide how many processors will be assigned for this
231 !% strategy. To disable parallelization in domains, you can use
232 !% <tt>ParDomains = no</tt> (or set the number of processors to
233 !% 1).
234 !%
235 !% The total number of processors required is the multiplication
236 !% of the processors assigned to each parallelization strategy.
237 !%Option auto -1
238 !% The number of processors is assigned automatically.
239 !%Option no 0
240 !% This parallelization strategy is not used.
241 !%End
242 call parse_variable(namespace, 'ParDomains', default(p_strategy_domains), parse(p_strategy_domains))
244 !%Variable ParStates
245 !%Type integer
246 !%Section Execution::Parallelization
247 !%Description
248 !% This variable controls the number of processors used for the
249 !% parallelization in states. The special value <tt>auto</tt> lets
250 !% Octopus decide how many processors will be assigned for this
251 !% strategy. To disable parallelization in states, you can use
252 !% <tt>ParStates = no</tt> (or set the number of processors to 1).
253 !%
254 !% The default value depends on the <tt>CalculationMode</tt>. For
255 !% <tt>CalculationMode = td</tt> the default is <tt>auto</tt>, while
256 !% for for other modes the default is <tt>no</tt>.
257 !%
258 !% The total number of processors required is the multiplication
259 !% of the processors assigned to each parallelization strategy.
260 !%Option auto -1
261 !% The number of processors is assigned automatically.
262 !%Option no 0
263 !% This parallelization strategy is not used.
264 !%End
265 call parse_variable(namespace, 'ParStates', default(p_strategy_states), parse(p_strategy_states))
266
267 !%Variable ParKPoints
268 !%Type integer
269 !%Default auto
270 !%Section Execution::Parallelization
271 !%Description
272 !% This variable controls the number of processors used for the
273 !% parallelization in K-Points and/or spin.
274 !% The special value <tt>auto</tt> lets Octopus decide how many processors will be
275 !% assigned for this strategy. To disable parallelization in
276 !% KPoints, you can use <tt>ParKPoints = no</tt> (or set the
277 !% number of processors to 1).
278 !%
279 !% The total number of processors required is the multiplication
280 !% of the processors assigned to each parallelization strategy.
281 !%Option auto -1
282 !% The number of processors is assigned automatically.
283 !%Option no 0
284 !% This parallelization strategy is not used.
285 !%End
286 call parse_variable(namespace, 'ParKPoints', default(p_strategy_kpoints), parse(p_strategy_kpoints))
287
288 !%Variable ParOther
289 !%Type integer
290 !%Default auto
291 !%Section Execution::Parallelization
292 !%Description
293 !% This variable controls the number of processors used for the
294 !% 'other' parallelization mode, that is CalculatioMode
295 !% dependent. For <tt>CalculationMode = casida</tt>, it means
296 !% parallelization in electron-hole pairs.
297 !%
298 !% The special value <tt>auto</tt>,
299 !% the default, lets Octopus decide how many processors will be
300 !% assigned for this strategy. To disable parallelization in
301 !% Other, you can use <tt>ParOther = no</tt> (or set the
302 !% number of processors to 1).
303 !%
304 !% The total number of processors required is the multiplication
305 !% of the processors assigned to each parallelization strategy.
306 !%Option auto -1
307 !% The number of processors is assigned automatically.
308 !%Option no 0
309 !% This parallelization strategy is not used.
310 !%End
311 call parse_variable(namespace, 'ParOther', default(p_strategy_other), parse(p_strategy_other))
312
313 do ipar = 1, p_strategy_max
314 if (parse(ipar) == par_no) parse(ipar) = 1
315 end do
316
317 call strategy()
318
319 mc%have_slaves = .false.
320
321 if (mc%par_strategy /= p_strategy_serial) then
322 safe_allocate(mc%group_sizes(1:p_strategy_max))
323
324 mc%group_sizes = 1
325
326 do ipar = 1, p_strategy_max
327 if (multicomm_strategy_is_parallel(mc, ipar)) then
328 mc%group_sizes(ipar) = parse(ipar)
329 else if (parse(ipar) /= 1) then
330 call messages_write('Ignoring specification for ' // par_types(ipar))
331 call messages_new_line()
332 call messages_write('This parallelization strategy is not available.')
333 call messages_warning()
334 end if
335 end do
336
337 call assign_nodes()
338
339
340 !%Variable ParallelizationNumberSlaves
341 !%Type integer
342 !%Default 0
343 !%Section Execution::Parallelization
344 !%Description
345 !% Slaves are nodes used for task parallelization. The number of
346 !% such nodes is given by this variable multiplied by the number
347 !% of domains used in domain parallelization.
348 !%End
349 call parse_variable(namespace, 'ParallelizationNumberSlaves', 0, num_slaves)
350
351 ! the slaves must be defined at a certain parallelization level, for the moment this is state parallelization.
352 slave_level = p_strategy_states
353 mc%have_slaves = (num_slaves > 0)
354
355 if (mc%have_slaves) then
356 call messages_experimental('Task parallelization')
357 end if
358
359 ! clear parallel strategies that were available but will not be used
360 do ii = 1, p_strategy_max
361 if (mc%group_sizes(ii) == 1) mc%par_strategy = ibclr(mc%par_strategy, ii - 1)
362 end do
363
364 ! reset
365 call sanity_check()
366 end if
367
368 call group_comm_create()
369
370 call messages_print_with_emphasis(namespace=namespace)
371
372 pop_sub(multicomm_init)
373
374 contains
375
376 ! ---------------------------------------------------------
377 subroutine strategy()
378 integer :: jj, ipar
379
380 push_sub(multicomm_init.strategy)
381
382 if (base_grp%size > 1) then
383
384 mc%par_strategy = 0
385
386 do ipar = 1, p_strategy_max
387 if (parse(ipar) == par_auto .or. parse(ipar) > 1) then
388 mc%par_strategy = ibset(mc%par_strategy, ipar - 1)
389 end if
390 end do
391
392 if (mc%par_strategy /= bitand(mc%par_strategy, parallel_mask)) then
393 call messages_write('Parallelization strategies unavailable for this run mode are being discarded.')
394 call messages_warning()
395 end if
396
397 mc%par_strategy = bitand(mc%par_strategy, parallel_mask)
398
399 if (mc%par_strategy == p_strategy_serial) then
400 message(1) = "More than one node is available, but this run mode cannot run with the requested parallelization."
401 message(2) = "Please select a parallelization strategy compatible with"
402 jj = 2
403 do ii = 1, n_par_types
404 if (bitand(parallel_mask, 2**(ii - 1)) /= 0) then
405 jj = jj + 1
406 write(message(jj), '(2a)') " -> ", par_types(ii)
407 end if
408 end do
409 jj=jj+1
410 write(message(jj),'(a,i6)') "mc%par_strategy is : ",mc%par_strategy
411 call messages_fatal(jj, only_root_writes = .true.)
412 end if
413 else
414 mc%par_strategy = p_strategy_serial
415 end if
416
417 mc%nthreads = 1
418#if defined(HAVE_OPENMP)
419 !$omp parallel
420 !$omp master
421 mc%nthreads = omp_get_num_threads()
422 !$omp end master
423 !$omp end parallel
424#endif
425
426 if (mc%par_strategy == p_strategy_serial .and. mc%nthreads == 1) then
427 message(1) = "Info: Octopus will run in *serial*"
428 call messages_info(1, namespace=namespace)
429 else
430 write(message(1),'(a)') 'Info: Octopus will run in *parallel*'
431 write(message(2),'(a)') ''
432 write(message(3),'(a, i8)') ' Number of processes :', base_grp%size
433 write(message(4),'(a, i8)') ' Number of threads per process :', mc%nthreads
434 write(message(5),'(a)') ''
435 call messages_info(5, namespace=namespace)
436 end if
437
438 pop_sub(multicomm_init.strategy)
439 end subroutine strategy
440
441 ! ---------------------------------------------------------
442
443 subroutine assign_nodes()
444 integer :: ii, nn, kk, n_divisors, divisors(1:50)
445 integer(int64) :: n_group_max(1:p_strategy_max)
446
448
449 ! this is the maximum number of processors in each group
450 n_group_max(1:p_strategy_max) = max(index_range(1:p_strategy_max), 1_int64)
451 do kk = 1, p_strategy_max
452 if (.not. multicomm_strategy_is_parallel(mc, kk)) n_group_max(kk) = 1
453 end do
454
455 call messages_write('Debug info: Allowable group ranks:', new_line = .true.)
456 do kk = 1, p_strategy_max
457 call messages_write(par_types(kk), fmt = '2x,a12,":",1x')
458 call messages_write(n_group_max(kk), new_line = .true.)
459 end do
460 call messages_info(debug_only=.true.)
461
462 nn = mc%n_node
463
464 ! first loop, check the processors assigned by the user
465 do ipar = p_strategy_max, 1, -1
466
467 if (mc%group_sizes(ipar) == par_auto) cycle
468
469 if (mc%group_sizes(ipar) > n_group_max(ipar)) then
470 call messages_write('The number of processors specified for '//par_types(ipar)//'(')
471 call messages_write(mc%group_sizes(ipar))
472 call messages_write(')', new_line = .true.)
473 call messages_write('is larger than the degrees of freedom for that level (')
474 call messages_write(n_group_max(ipar))
475 call messages_write(').')
476 call messages_warning()
477 end if
478
479 if (mod(nn, mc%group_sizes(ipar)) /= 0) then
480 call messages_write('The number of processors specified for '//par_types(ipar)//'(')
481 call messages_write(mc%group_sizes(ipar))
482 call messages_write(')', new_line = .true.)
483 call messages_write('is not a divisor of the number of processors (')
484 call messages_write(mc%n_node)
485 call messages_write(').')
486 call messages_fatal()
487 end if
488
489 nn = nn/mc%group_sizes(ipar)
490
491 end do
492
493 ! second loop, now assign the rest automatically
494 do ipar = p_strategy_max, 1, -1
495
496 if (mc%group_sizes(ipar) /= par_auto) cycle
497
498 n_divisors = ubound(divisors, dim = 1)
499 call get_divisors(nn, n_divisors, divisors)
500
501 mc%group_sizes(ipar) = nn
502 do ii = 2, n_divisors
503 if (divisors(ii) > n_group_max(ipar)) then
504 mc%group_sizes(ipar) = divisors(ii - 1)
505 exit
506 end if
507 end do
508
509 nn = nn/mc%group_sizes(ipar)
510
511 end do
512
514 end subroutine assign_nodes
515
516
517 ! ---------------------------------------------------------
519 subroutine sanity_check()
520 real(real64) :: frac
521 integer :: ii, kk
522 integer(int64) :: jj, n_max
523 integer :: real_group_sizes(1:MAX_INDEX)
524
526
527 if (num_slaves > 0) then
528
529 if (mc%group_sizes(slave_level) < num_slaves + 1) then
530 message(1) = 'Too many nodes assigned to task parallelization.'
532 end if
533
534 write(message(1),'(a,i6)') 'Info: Number of slaves nodes :', &
535 num_slaves*product(mc%group_sizes(1:slave_level - 1))
536 call messages_info(1)
537
538 end if
539
540 ! print out some info
541 ii = 0
542 do kk = p_strategy_max, 1, -1
543 real_group_sizes(kk) = mc%group_sizes(kk)
544 if (.not. multicomm_strategy_is_parallel(mc, kk)) cycle
545 ii = ii + 1
546 if (kk == slave_level) real_group_sizes(kk) = real_group_sizes(kk) - num_slaves
547 write(message(ii),'(3a,i6,a,i12,a)') 'Info: Number of nodes in ', &
548 par_types(kk), ' group:', real_group_sizes(kk), ' (', index_range(kk), ')'
549 end do
550 call messages_info(ii)
551
552 ! do we have the correct number of processors
553 if (product(mc%group_sizes(1:p_strategy_max)) /= base_grp%size) then
554 write(message(1),'(a)') 'Inconsistent number of processors:'
555 write(message(2),'(a,i6)') ' MPI processes = ', base_grp%size
556 write(message(3),'(a,i6)') ' Required processes = ', product(mc%group_sizes(1:p_strategy_max))
557 message(4) = ''
558 message(5) = 'You probably have a problem in the ParDomains, ParStates, ParKPoints or ParOther.'
559 call messages_fatal(5, only_root_writes = .true.)
560 end if
561
562 if (any(real_group_sizes(1:p_strategy_max) > index_range(1:p_strategy_max))) then
563 message(1) = "Could not distribute nodes in parallel job. Most likely you are trying to"
564 message(2) = "use too many nodes for the job."
565 call messages_fatal(2, only_root_writes = .true.)
566 end if
567
568 if (any(index_range(1:p_strategy_max) / real_group_sizes(1:p_strategy_max) < min_range(1:p_strategy_max) .and. &
569 real_group_sizes(1:p_strategy_max) > 1)) then
570 message(1) = "I have fewer elements in a parallel group than recommended."
571 message(2) = "Maybe you should reduce the number of nodes."
572 call messages_warning(2)
573 end if
574
575 ! calculate fraction of idle time
576 frac = m_one
577 do ii = 1, p_strategy_max
578 n_max = ceiling(real(index_range(ii), real64) / real(real_group_sizes(ii)), real64)
579 jj = n_max*real_group_sizes(ii)
580 frac = frac*(m_one - real(jj - index_range(ii), real64) / real(jj, real64) )
581 end do
582
583 write(message(1), '(a,f5.2,a)') "Info: Octopus will waste at least ", &
584 (m_one - frac)*100.0_real64, "% of computer time."
585 if (frac < 0.8_real64) then
586 message(2) = "Usually a number of processors which is a multiple of small primes is best."
587 call messages_warning(2)
588 else
589 call messages_info(1)
590 end if
591
593 end subroutine sanity_check
594
595 ! ---------------------------------------------------------
596 subroutine group_comm_create()
597#if defined(HAVE_MPI)
598 logical :: dim_mask(MAX_INDEX)
599 integer :: i_strategy, irank
600 logical :: reorder, periodic_mask(MAX_INDEX)
601 integer :: coords(MAX_INDEX)
602 type(mpi_comm) :: new_comm
603 integer :: new_comm_size
604 character(len=6) :: node_type
605 type(mpi_grp_t) :: reorder_grp
606 type(mpi_group) :: base_group, reorder_group
607 integer :: ranks(base_grp%size)
608 integer :: ii, jj, kk, ll, nn
609 type(mpi_comm) :: reorder_comm
610#endif
611
613
614 mc%node_type = p_master
615
616 safe_allocate(mc%group_comm(1:p_strategy_max))
617 safe_allocate(mc%who_am_i(1:p_strategy_max))
618
619 mc%group_comm = mpi_comm_undefined
620 mc%who_am_i = 0
621
622#if defined(HAVE_MPI)
623 mc%full_comm = mpi_comm_null
624 mc%slave_intercomm = mpi_comm_null
625 if (mc%par_strategy /= p_strategy_serial) then
626 if (mc%reorder_ranks) then
627 ! first, reorder the ranks
628 ! this is done to get a column-major ordering of the ranks in the
629 ! Cartesian communicator, since they a ordered row-major otherwise
630 call mpi_comm_group(base_grp%comm, base_group)
631 ! now transpose the hypercube => get rank numbers in column-major order
632 nn = 1
633 do ii = 1, mc%group_sizes(1)
634 do jj = 1, mc%group_sizes(2)
635 do kk = 1, mc%group_sizes(3)
636 do ll = 1, mc%group_sizes(4)
637 ranks(nn) = (ll-1)*mc%group_sizes(3)*mc%group_sizes(2)*mc%group_sizes(1) &
638 + (kk-1)*mc%group_sizes(2)*mc%group_sizes(1) &
639 + (jj-1)*mc%group_sizes(1) + ii - 1
640 nn = nn + 1
641 end do
642 end do
643 end do
644 end do
645 call mpi_group_incl(base_group, base_grp%size, ranks, reorder_group)
646 ! now get the reordered communicator
647 call mpi_comm_create(base_grp%comm, reorder_group, reorder_comm)
648 call mpi_grp_init(reorder_grp, reorder_comm)
649 else
650 call mpi_grp_copy(reorder_grp, base_grp)
651 end if
652
653 ! Multilevel parallelization is organized in a hypercube. We
654 ! use an MPI Cartesian topology to generate the communicators
655 ! that correspond to each level.
656
657 ! create the topology
658 periodic_mask = .false.
659 reorder = .true.
660
661 ! The domain and states dimensions have to be periodic (2D torus)
662 ! in order to circulate matrix blocks.
665
666 ! We allow reordering of ranks.
667 call mpi_cart_create(reorder_grp%comm, p_strategy_max, mc%group_sizes, periodic_mask, reorder, mc%full_comm)
668
669 call mpi_comm_rank(mc%full_comm, mc%full_comm_rank)
670
671 ! get the coordinates of the current processor
672 call mpi_cart_coords(mc%full_comm, mc%full_comm_rank, p_strategy_max, coords)
673
674 ! find out what type of node this is
675 if (coords(slave_level) >= mc%group_sizes(slave_level) - num_slaves) then
676 mc%node_type = p_slave
677 end if
678
679 if (mc%node_type == p_master) then
680 mc%group_sizes(slave_level) = mc%group_sizes(slave_level) - num_slaves
681 else
682 mc%group_sizes(slave_level) = num_slaves
683 end if
685 call mpi_comm_split(mc%full_comm, mc%node_type, mc%full_comm_rank, new_comm)
686 assert(new_comm /= mpi_comm_null)
687 call mpi_comm_size(new_comm, new_comm_size)
688
689 reorder = .false.
690 if (product(mc%group_sizes(:)) /= new_comm_size) then
691 write(stderr,*) 'node ', base_grp%rank, ': mc%group_sizes = ', mc%group_sizes, ' new_comm_size = ', new_comm_size
692 call base_grp%barrier()
693 assert(product(mc%group_sizes(:)) == new_comm_size)
694 end if
695 call mpi_cart_create(new_comm, p_strategy_max, mc%group_sizes, periodic_mask, reorder, mc%master_comm)
696 assert(mc%master_comm /= mpi_comm_null)
697
698 call mpi_comm_free(new_comm)
699
700 call mpi_comm_rank(mc%master_comm, mc%master_comm_rank)
701
702 call mpi_grp_init(mc%master_grp, mc%master_comm)
703
704 ! The "lines" of the Cartesian grid.
705 ! Initialize all the communicators, even if they are not parallelized
706 do i_strategy = 1, p_strategy_max
707 dim_mask = .false.
708 dim_mask(i_strategy) = .true.
709 call mpi_cart_sub(mc%master_comm, dim_mask, mc%group_comm(i_strategy))
710 call mpi_comm_rank(mc%group_comm(i_strategy), mc%who_am_i(i_strategy))
711 end do
712
713 ! The domain-state "planes" of the grid (the ones with periodic dimensions).
714 dim_mask = .false.
715 dim_mask(p_strategy_domains) = .true.
716 dim_mask(p_strategy_states) = .true.
717 call mpi_cart_sub(mc%master_comm, dim_mask, mc%dom_st_comm)
718
719 ! The state-kpoints "planes" of the grid
720 dim_mask = .false.
721 dim_mask(p_strategy_states) = .true.
722 dim_mask(p_strategy_kpoints) = .true.
723 call mpi_cart_sub(mc%master_comm, dim_mask, mc%st_kpt_comm)
724
725 ! The domains-states-kpoints "cubes" of the grid
726 dim_mask = .false.
728 dim_mask(p_strategy_states) = .true.
729 dim_mask(p_strategy_kpoints) = .true.
730 call mpi_cart_sub(mc%master_comm, dim_mask, mc%dom_st_kpt_comm)
731
732 if (num_slaves > 0) call create_slave_intercommunicators()
733
734 call create_intranode_communicator(base_grp, mc%intranode_grp, mc%internode_grp)
735 else
736 ! we initialize these communicators so we can use them even in serial
737 mc%group_comm = base_grp%comm
738 mc%who_am_i = 0
739 mc%master_comm = base_grp%comm
740 mc%dom_st_comm = base_grp%comm
741 mc%st_kpt_comm = base_grp%comm
742 mc%dom_st_kpt_comm = base_grp%comm
743 call mpi_grp_copy(mc%master_grp, base_grp)
744 call mpi_grp_copy(mc%intranode_grp, base_grp)
745 call mpi_grp_copy(mc%internode_grp, base_grp)
746 end if
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, mc%base_grp%size - 1
760 if (mc%base_grp%rank == irank) then
761 write(message(1),'(5i10,5x,a)') mc%base_grp%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 mc%base_grp%barrier()
766 end do
767 end if
768#endif
769
771 end subroutine group_comm_create
772
773 ! -----------------------------------------------------
774#ifdef HAVE_MPI
775 subroutine create_slave_intercommunicators()
776 integer :: remote_leader
777 integer :: tag
778 integer :: coords(max_index)
779
780 push_sub(multicomm_init.create_slave_intercommunicators)
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)
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)
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)
798
799 pop_sub(multicomm_init.create_slave_intercommunicators)
800 end subroutine create_slave_intercommunicators
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))
821 end do
822 call mpi_comm_free(mc%dom_st_comm)
823 call mpi_comm_free(mc%st_kpt_comm)
824 call mpi_comm_free(mc%dom_st_kpt_comm)
825 call mpi_comm_free(mc%full_comm)
826 call mpi_comm_free(mc%master_comm)
827
828 if (multicomm_have_slaves(mc)) call mpi_comm_free(mc%slave_intercomm)
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
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
946 end subroutine multicomm_create_all_pairs
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
1038
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:
subroutine strategy()
Definition: multicomm.F90:473
subroutine group_comm_create()
Definition: multicomm.F90:685
subroutine sanity_check()
check if a balanced distribution of nodes will be used
Definition: multicomm.F90:608
subroutine assign_nodes()
Definition: multicomm.F90:532
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:158
real(real64), parameter, public m_one
Definition: global.F90:201
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:898
character(len=512), private msg
Definition: messages.F90:167
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:525
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1000
subroutine, public messages_new_line()
Definition: messages.F90:1089
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1040
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
This module contains some common usage patterns of MPI routines.
Definition: mpi_lib.F90:117
subroutine mpi_grp_copy(mpi_grp_out, mpi_grp_in)
MPI_THREAD_FUNNELED allows for calls to MPI from an OMP region if the thread is the team master.
Definition: mpi.F90:383
type(mpi_comm), parameter, public mpi_comm_undefined
used to indicate a communicator has not been initialized
Definition: mpi.F90:138
subroutine mpi_grp_init(grp, comm)
Initialize MPI group instance.
Definition: mpi.F90:341
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
logical pure function, public multicomm_strategy_is_parallel(mc, level)
Definition: multicomm.F90:728
integer, parameter, public p_slave
Definition: multicomm.F90:183
subroutine, public multicomm_divide_range(nobjs, nprocs, istart, ifinal, lsize, scalapack_compat)
This routine uses the one-factorization (or near-one-factorization of a complete graph to construct a...
Definition: multicomm.F90:748
subroutine, public multicomm_all_pairs_copy(apout, apin)
Definition: multicomm.F90:251
subroutine, public multicomm_end(mc)
Definition: multicomm.F90:706
integer, parameter, public par_no
Definition: multicomm.F90:187
logical pure function, public multicomm_have_slaves(this)
Definition: multicomm.F90:854
logical pure function, public multicomm_is_slave(this)
Definition: multicomm.F90:846
subroutine, public multicomm_init(mc, namespace, base_grp, mode_para, n_node, index_range, min_range)
create index and domain communicators
Definition: multicomm.F90:273
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:824
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:120
subroutine, public get_divisors(nn, n_divisors, divisors)
Definition: utils.F90:171
This is defined even when running serial.
Definition: mpi.F90:144
An all-pairs communication schedule for a given group.
Definition: multicomm.F90:240
Stores all communicators and groups.
Definition: multicomm.F90:208
int true(void)
subroutine parse()
Definition: xc.F90:487