TABLE OF CONTENTS
- ABINIT/distrb2
- ABINIT/distrb2_hf
- ABINIT/initmpi_atom
- ABINIT/initmpi_band
- ABINIT/initmpi_grid
- ABINIT/initmpi_img
- ABINIT/initmpi_pert
- ABINIT/initmpi_seq
- ABINIT/initmpi_world
- ABINIT/m_mpinfo
- ABINIT/proc_distrb_band
- ABINIT/proc_distrb_cycle
- ABINIT/proc_distrb_cycle_bands
- ABINIT/proc_distrb_kptband
- ABINIT/proc_distrb_nband
- m_mpinfo/clnmpi_atom
- m_mpinfo/clnmpi_grid
- m_mpinfo/clnmpi_img
- m_mpinfo/clnmpi_pert
- m_mpinfo/copy_mpi_enreg
- m_mpinfo/destroy_mpi_enreg
- m_mpinfo/init_mpi_enreg
- m_mpinfo/iwrite_fftdatar
- m_mpinfo/mpi_distrib_is_ok
- m_mpinfo/nullify_mpi_enreg
- m_mpinfo/pre_gather
- m_mpinfo/pre_scatter
- m_mpinfo/ptabs_fourdp
- m_mpinfo/ptabs_fourwf
- m_mpinfo/set_mpi_enreg_fft
- m_mpinfo/unset_mpi_enreg_fft
ABINIT/distrb2 [ Functions ]
NAME
distrb2
FUNCTION
Creates the tabs of repartition of processors for sharing the jobs on k-points, spins and bands.
INPUTS
mband = maximum number of bands nband(nkpt*nsppol) = number of bands per k point, for each spin nkpt = number of k-points nproc= number of processors available for this distribution nsppol = 1 for unpolarized, 2 for polarized
SIDE EFFECTS
mpi_enreg%proc_distrb(nkpt,mband,nsppol)=number of the processor that will treat each band in each k point. mpi_enreg%nproc_spkpt is set
NOTES
For the time being, the band parallelisation works only when the number of bands is identical for spin up and spin down at the same k point. The problem is the most clearly seen in the kpgio routine, where a different parallel repartition of k points for spin up and spin down would conflict with the present computation of k+G sphere, independent of the spin.
SOURCE
2296 subroutine distrb2(mband,mband_mem_out,nband,nkpt,nproc,nsppol,mpi_enreg) 2297 2298 !Arguments ------------------------------------ 2299 integer,intent(in) :: mband,nkpt,nproc,nsppol 2300 integer,intent(in) :: nband(nkpt*nsppol) 2301 integer,intent(out) :: mband_mem_out 2302 type(MPI_type),intent(inout) :: mpi_enreg 2303 2304 !Local variables------------------------------- 2305 integer :: maxproc_bandpool 2306 integer :: nproc_band 2307 integer :: inb1,ind,ind0,nband_k,proc_max,proc_min 2308 integer :: nband_k_sp2,minb_per_proc 2309 integer :: iiband,iikpt,iisppol,ikpt_this_proc,nb_per_proc,nproc_spkpt,temp_unit 2310 integer :: kpt_distrb(nkpt) 2311 logical,save :: first=.true.,has_file 2312 character(len=500) :: msg 2313 2314 !****************************************************************** 2315 2316 nproc_spkpt=mpi_enreg%nproc_spkpt 2317 if (mpi_enreg%paral_pert==1) nproc_spkpt=nproc 2318 ! if (mpi_enreg%paral_pert==1) nproc_spkpt=nproc/mpi_enreg%nproc_pert 2319 2320 mband_mem_out = 0 2321 2322 !Initialization of proc_distrb 2323 mpi_enreg%proc_distrb = nproc+1 2324 do iisppol=1,nsppol 2325 do iikpt=1,nkpt 2326 do iiband=1,nband(iikpt+(iisppol-1)*nkpt) 2327 mpi_enreg%proc_distrb(iikpt,iiband,iisppol)=nproc_spkpt-1 2328 end do 2329 end do 2330 end do 2331 !That s all for an empty communication space 2332 if (nproc==0) return 2333 2334 !Some checks 2335 !if (mpi_enreg%paralbd==0 .and. any(dtset%optdriver == [RUNL_GSTATE, RUNL_RESPFN])) then 2336 if (mpi_enreg%paralbd==0) then 2337 ! Check if nkpt and nproc_spkpt match 2338 if(nproc_spkpt>nkpt*nsppol) then 2339 ! Too many proc. with respect to nkpt 2340 write(msg,'(a,i0,a,i0,a,i0,4a)')& 2341 'nproc_spkpt= ',nproc_spkpt,' >= nkpt= ',nkpt,'* nsppol= ',nsppol,ch10,& 2342 'The number of processors is larger than nkpt*nsppol. This is a WASTE.',ch10,& 2343 ' Ignore this warning if this is not a GS run' 2344 ABI_WARNING(msg) 2345 else if (mod(nkpt*nsppol,nproc_spkpt) /= 0) then 2346 ! nkpt not a multiple of nproc_spkpt 2347 write(msg,'(a,i0,a,i0,5a)')& 2348 'nkpt*nsppol (', nkpt*nsppol, ') is not a multiple of nproc_spkpt (',nproc_spkpt, ')', ch10,& 2349 'The k-point parallelisation is INEFFICIENT. ',ch10,& 2350 'Ignore this warning if this is not a GS run.' 2351 ABI_WARNING(msg) 2352 end if 2353 end if 2354 2355 !Inquire whether there exist a file containing the processor distribution 2356 if (first) then 2357 ! Case first time: test file to do 2358 ! Open the file containing the k-point distribution 2359 first=.false.; has_file = file_exists("kpt_distrb") 2360 end if 2361 2362 !--------------------------------------------------------------------------- 2363 !Initialize the processor distribution, either from a file, or from an algorithm 2364 if (has_file) then 2365 if (open_file('kpt_distrb',msg,newunit=temp_unit,form='formatted',status='old') /= 0) then 2366 ABI_ERROR(msg) 2367 end if 2368 rewind(unit=temp_unit) 2369 if (mpi_enreg%paralbd == 1) then 2370 ! -> read bands distribution 2371 read(temp_unit,*) mpi_enreg%proc_distrb 2372 else 2373 read(temp_unit,*) kpt_distrb 2374 end if 2375 close(temp_unit) 2376 proc_max=0 2377 proc_min=nproc_spkpt 2378 ! -> determine the range of proc. requested 2379 if (mpi_enreg%paralbd == 1) then 2380 do iisppol=1,nsppol 2381 do iikpt=1,nkpt 2382 nband_k = nband(iikpt+(iisppol-1)*nkpt) 2383 proc_max=maxval(mpi_enreg%proc_distrb(iikpt,1:nband_k,iisppol)) 2384 proc_min=minval(mpi_enreg%proc_distrb(iikpt,1:nband_k,iisppol)) 2385 end do 2386 end do 2387 else 2388 proc_max=maxval(kpt_distrb(1:nkpt)) 2389 proc_min=minval(kpt_distrb(1:nkpt)) 2390 ! -> fill the tab proc_distrb with kpt_distrb 2391 do iisppol=1,nsppol 2392 do iikpt=1,nkpt 2393 nband_k = nband(iikpt+(iisppol-1)*nkpt) 2394 do iiband=1,nband_k 2395 mpi_enreg%proc_distrb(iikpt,iiband,iisppol)=kpt_distrb(iikpt) 2396 end do 2397 end do 2398 end do 2399 end if ! mpi_enreg%paralbd 2400 2401 if(proc_max>(nproc_spkpt-1)) then 2402 ! Too much proc. requested 2403 write(msg, '(a,a,a,i0,a,a,a)' )& 2404 'The number of processors mentioned in the kpt_distrb file',ch10,& 2405 'must be lower or equal to the actual number of processors =',nproc_spkpt-1,ch10,& 2406 'Action: change the kpt_distrb file, or increase the',' number of processors.' 2407 ABI_ERROR(msg) 2408 end if 2409 2410 if(proc_max/=(nproc_spkpt-1)) then 2411 ! Too few proc. used 2412 write(msg, '(a,i0,a,a,a,i0,a,a,a)' )& 2413 'Only ',proc_max+1,' processors are used (from kpt_distrb file),',ch10,& 2414 'when',nproc_spkpt,' processors are available.',ch10,& 2415 'Action: adjust number of processors and kpt_distrb file.' 2416 ABI_ERROR(msg) 2417 end if 2418 2419 if(proc_min<0) then 2420 write(msg, '(a,a,a)' )& 2421 'The number of processors must be bigger than 0 in kpt_distrb file.',ch10,& 2422 'Action: modify kpt_distrb file.' 2423 ABI_ERROR(msg) 2424 end if 2425 2426 !--------------------------------------------------------------------------- 2427 else 2428 ! 'kpt_distrb' file does not exist 2429 2430 if (mpi_enreg%paralbd==1) then 2431 2432 ! No possible band parallelization 2433 if (nproc<(nkpt*nsppol)) then 2434 2435 ! Does not allow a processor to treat different spins 2436 ! NB: for odd nproc this will happen anyway for the middle proc - will this not unbalance things? 2437 ind0=0 2438 inb1=(nkpt*nsppol)/nproc 2439 if (mod((nkpt*nsppol),nproc)/=0) inb1=inb1+1 2440 do iikpt=1,nkpt 2441 nband_k=nband(iikpt) 2442 nband_k_sp2=nband(iikpt+nkpt*(nsppol-1)) 2443 !ind=ind0/inb1 2444 ind=mod(ind0,nproc) 2445 do iiband=1,nband_k 2446 mpi_enreg%proc_distrb(iikpt,iiband,1)=ind 2447 if (nsppol==2 .and. iiband <= nband_k_sp2) then 2448 mpi_enreg%proc_distrb(iikpt,iiband,2)=ind 2449 end if 2450 !if (nsppol==2) mpi_enreg%proc_distrb(iikpt,iiband,2)=nproc-ind-1 2451 end do 2452 ind0=ind0+1 2453 end do 2454 2455 ! MT130831 : OLD CODING 2456 ! do iisppol=1,nsppol;do iikpt=1,nkpt 2457 ! ind=(iikpt+(iisppol-1)*nkpt-1)/inb1 2458 ! nband_k=nband(iikpt+(iisppol-1)*nkpt) 2459 ! do iiband=1,nband_k 2460 ! mpi_enreg%proc_distrb(iikpt,iiband,iisppol)=ind 2461 ! end do;end do;end do 2462 ! MT130831 : END OF OLD CODING 2463 2464 ! Possible band parallelization 2465 else 2466 ! Does not allow a processor to treat different spins 2467 ind0=0 2468 maxproc_bandpool=floor(nproc*one/(nkpt*nsppol)) 2469 do iikpt=1,nkpt 2470 nband_k=nband(iikpt) 2471 nband_k_sp2=nband(iikpt+nkpt*(nsppol-1)) 2472 minb_per_proc=floor(nband_k*one/maxproc_bandpool) 2473 if (mod(nband_k,maxproc_bandpool)/=0) minb_per_proc=minb_per_proc+1 2474 do nb_per_proc = minb_per_proc, nband_k 2475 if (mod(nband_k,nb_per_proc)==0) exit 2476 end do 2477 nproc_band = nband_k / nb_per_proc 2478 2479 mband_mem_out = max(mband_mem_out,nb_per_proc) 2480 do iiband=1,nband_k 2481 ind=mod((iiband-1)/nb_per_proc+ind0, nproc) 2482 !ind = mod( mod((iiband-1),nproc_band) + ind0, nproc ) 2483 mpi_enreg%proc_distrb(iikpt,iiband,1)=ind 2484 !TODO : could end up with 0 bands on certain procs with this configuration and nband(k) /= constant 2485 if (nsppol==2 .and. iiband <= nband_k_sp2) then 2486 mpi_enreg%proc_distrb(iikpt,iiband,2)=ind+nkpt*nproc_band 2487 end if 2488 !if (nsppol==2) mpi_enreg%proc_distrb(iikpt,iiband,2)=nproc-ind-1 2489 end do 2490 ind0=ind+1 ! take last proc associated to the bands at iikpt, and increment 2491 !ind0 = ind0 + nproc_band 2492 end do 2493 2494 ! MT130831 : OLD CODING 2495 ! ind0=0;maxproc_bandpool=nproc/(nkpt*nsppol) 2496 ! do iisppol=1,nsppol;do iikpt=1,nkpt 2497 ! nband_k=nband(iikpt+(iisppol-1)*nkpt) 2498 ! nb_per_proc=nband_k/maxproc_bandpool;if (mod(nband_k,maxproc_bandpool)/=0) nb_per_proc=nb_per_proc+1 2499 ! do iiband=1,nband_k 2500 ! ind=(iiband-1)/nb_per_proc+ind0 2501 ! mpi_enreg%proc_distrb(iikpt,iiband,iisppol)=ind 2502 ! end do 2503 ! ind0=ind+1 2504 ! end do;end do 2505 ! MT130831 : END OF OLD CODING 2506 2507 end if 2508 2509 ! XG060807 : OLD CODING 2510 ! ind=0 2511 ! do iisppol=1,nsppol;do iikpt=1,nkpt 2512 ! nband_k=nband(iikpt+(iisppol-1)*nkpt) 2513 ! do iiband=1,nband_k 2514 ! mpi_enreg%proc_distrb(iikpt,iiband,iisppol)=ind/nb_per_proc 2515 ! ind=ind+1 2516 ! end do;end do;end do 2517 ! XG060807 : END OF OLD CODING 2518 2519 elseif (mpi_enreg%paralbd==0) then 2520 2521 ! Does not allow a processor to treat different spins 2522 ind0=0 2523 nb_per_proc=(nsppol*nkpt)/nproc_spkpt 2524 if (mod((nsppol*nkpt),nproc_spkpt)/=0) nb_per_proc=nb_per_proc+1 2525 do iikpt=1,nkpt 2526 nband_k=nband(iikpt) 2527 ind=ind0/nb_per_proc 2528 do iiband=1,nband_k 2529 mpi_enreg%proc_distrb(iikpt,iiband,1)=ind 2530 if (nsppol==2) then 2531 !TODO: why are these bands ordered in the opposite direction wrt spin 1?? 2532 mpi_enreg%proc_distrb(iikpt,iiband,2)=nproc_spkpt-ind-1 2533 end if 2534 end do 2535 ind0=ind0+1 2536 end do 2537 2538 ! XG060807 : OLD CODING 2539 ! ind=0 2540 ! do iisppol=1,nsppol;do iikpt=1,nkpt 2541 ! nband_k = nband(iikpt+(iisppol-1)*nkpt) 2542 ! do iiband=1,nband_k 2543 ! Distribute k-points homogeneously 2544 ! proc_distrb(iikpt,iiband,iisppol)=mod(iikpt-1,nproc_spkpt) 2545 ! mpi_enreg%proc_distrb(iikpt,iiband,iisppol)=ind/nb_per_proc 2546 ! end do 2547 ! ind=ind + 1 2548 ! end do;end do 2549 ! XG060807 : END OF OLD CODING 2550 2551 end if ! mpi_enreg%paralbd 2552 2553 end if ! has_file 2554 2555 ! local indices on cpus, for each sppol kpt band 2556 mpi_enreg%my_kpttab(:)=0 2557 mpi_enreg%my_isppoltab(:)=0 2558 do iisppol=1,nsppol 2559 ikpt_this_proc=0 2560 do iikpt=1,nkpt 2561 nband_k=nband(iikpt+(iisppol-1)*nkpt) 2562 if(proc_distrb_cycle(mpi_enreg%proc_distrb,iikpt,1,nband_k,iisppol,mpi_enreg%me_kpt)) cycle 2563 ikpt_this_proc=ikpt_this_proc+1 2564 mpi_enreg%my_kpttab(iikpt)=ikpt_this_proc 2565 mpi_enreg%my_isppoltab(iisppol)=1 2566 end do 2567 end do 2568 2569 if (mband_mem_out == 0) mband_mem_out = mband 2570 2571 end subroutine distrb2
ABINIT/distrb2_hf [ Functions ]
NAME
distrb2_hf
FUNCTION
Create the tabs of repartition of processors for sharing the jobs on occupied states (labeled by k-points, bands and spin indices) for Hartree-Fock calculation.
INPUTS
nbandhf = maximum number of occupied bands nkpthf = number of k-points in full BZ nproc= number of processors available for this distribution nsppol = 1 for unpolarized, 2 for polarized
SIDE EFFECTS
mpi_enreg%proc_distrb(nkpthf,nbandhf,nsppol)=number of the processor that will treat each band in each k point. mpi_enreg%nproc_spkpt is set
NOTES
For the time being, the band parallelisation works only when the number of bands is identical for spin up and spin down at the same k point. The problem is the most clearly seen in the kpgio routine, where a different parallel repartition of k points for spin up and spin down would conflict with the present computation of k+G sphere, independent of the spin.
SOURCE
2603 subroutine distrb2_hf(nbandhf,nkpthf, nproc, nsppol, mpi_enreg) 2604 2605 !Arguments ------------------------------------ 2606 integer,intent(in) :: nbandhf,nkpthf,nproc,nsppol 2607 type(MPI_type),intent(inout) :: mpi_enreg 2608 2609 !Local variables------------------------------- 2610 integer :: ind,iiband,iikpt,iistep,nproc_hf 2611 character(len=500) :: msg 2612 2613 !****************************************************************** 2614 2615 nproc_hf=mpi_enreg%nproc_hf 2616 2617 !* Initialize distrb_hf (the array exists necessarily) 2618 do iiband=1,nbandhf 2619 do iikpt=1,nkpthf 2620 mpi_enreg%distrb_hf(iikpt,iiband,1)=nproc_hf-1 2621 end do 2622 end do 2623 2624 !* End of the routine for an empty communication space 2625 if (nproc==0) return 2626 2627 !*** Testing section *** 2628 2629 if (nsppol==2) then 2630 !* Check that the distribution over (spin,k point) allow to consider spin up and spin dn independently. 2631 if (mpi_enreg%nproc_spkpt/=1.and.mod(mpi_enreg%nproc_spkpt,2)/=0) then 2632 ABI_ERROR( 'The variable nproc_spkpt is not even but nsppol= 2') 2633 !* In this case, one processor will carry both spin. (this will create pbms for the calculation) 2634 end if 2635 !* Check that the number of band is the same for each spin, at each k-point. (It should be) 2636 !* do iikpt=1,nkpthf 2637 !* if (nband(iikpt)/=nband(iikpt+nkpthf)) then 2638 !* msg = ' WARNING - the number of bands is different for spin up or spin down. ' 2639 !* ABI_ERROR(msg) 2640 !* end if 2641 !* end do 2642 !* If one of this test is not good for one proc then other procs fall in deadlock, according to distrb2. 2643 !* What does it mean ??? 2644 end if 2645 2646 2647 !* Check if nkpthf and nproc_hf match 2648 if (nproc_hf>nkpthf*nbandhf) then 2649 !* There are too many processors with respect to nkpthf*nbandhf 2650 write(msg, '(a,a,i4,a,i4,a,i4,a,a)' ) ch10,& 2651 & 'nproc_hf=',nproc_hf,' >= nkpthf=',nkpthf,'* nbandhf=',nbandhf,ch10,& 2652 & 'The number of processors is larger than nkpthf*nbandhf. This is a waste.' 2653 ABI_WARNING(msg) 2654 2655 else if(mod(nkpthf*nbandhf,nproc_hf)/=0) then 2656 !* nkpthf*nbandhf is not a multiple of nproc_hf 2657 write(msg, '(2a,i5,a,i5,3a)' ) ch10,& 2658 & 'nkpthf*nbandhf (', nkpthf*nbandhf, ') is not a multiple of nproc_hf (',nproc_hf, ')', ch10,& 2659 & 'The parallelisation may not be efficient.' 2660 ABI_WARNING(msg) 2661 end if 2662 2663 !*** End of testing section *** 2664 2665 !*** Initialize the processor distribution from a simple algorithm *** 2666 2667 if (nproc_hf<nkpthf) then 2668 !* In this case, a parallelization over kpts only. 2669 iistep=nkpthf/nproc_hf 2670 if (mod(nkpthf,nproc_hf) /=0) iistep=iistep+1 2671 ind=0 2672 do iikpt=1,nkpthf 2673 !*** Only the first "nbandhf" bands are considered (they are assumed to be the only occupied ones) 2674 do iiband=1,nbandhf 2675 mpi_enreg%distrb_hf(iikpt,iiband,1)=ind/iistep 2676 end do 2677 ind=ind+1 2678 end do 2679 2680 else 2681 !* In this case, a parallelization over all the occupied states is possible. 2682 if (nproc_hf < nbandhf*nkpthf) then 2683 iistep=(nbandhf*nkpthf)/nproc_hf; 2684 if (mod((nbandhf*nkpthf),nproc_hf) /=0) iistep=iistep+1 2685 else 2686 iistep=1 2687 end if 2688 ind=0 2689 do iikpt=1,nkpthf 2690 !*** Only the first "nbandhf" bands are considered (they are assumed to be the only occupied ones) 2691 do iiband=1,nbandhf 2692 mpi_enreg%distrb_hf(iikpt,iiband,1)=ind/iistep 2693 ind=ind+1 2694 end do 2695 end do 2696 end if 2697 2698 !*** Initialization of processor distribution from a file (simple copy from distrb2, not yet implemented) *** 2699 2700 ! !* Inquire whether there exist a file containing the processor distribution 2701 ! if (first) then 2702 ! ! Case first time : test file to do 2703 ! ! Open the file containing the (k-points,bands) distribution 2704 ! open(unit=temp_unit,file='kpt_distrb_hf',form='formatted',status='old',iostat=ios) 2705 ! if(ios==0) then 2706 ! ! 'kpt_distrb_hf' file exists 2707 ! file_exist=1 2708 ! close(temp_unit) 2709 ! else 2710 ! file_exist=0 2711 ! end if 2712 ! first=.false. 2713 ! end if 2714 ! 2715 ! !* Initialize the processor distribution, either from a file, or from an algorithm 2716 ! if (file_exist == 1) then 2717 ! !* Read (k-points,bands) distribution out of the file 2718 ! open(unit=temp_unit,file='kpt_distrb_hf',form='formatted',status='old',iostat=ios) 2719 ! rewind(unit=temp_unit) 2720 ! read(temp_unit,*) mpi_enreg%distrb_hf 2721 ! close(temp_unit) 2722 ! !* Determine the range of processors requested 2723 ! proc_max=0 2724 ! proc_min=nproc_hf 2725 ! do iikpt=1,nkpthf 2726 ! mband_occ_k = mband_occ(iikpt+(iisppol-1)*nkpthf) 2727 ! proc_max=maxval(mpi_enreg%distrb_hf(iikpt,1:mband_occ_k,1)) 2728 ! proc_min=minval(mpi_enreg%distrb_hf(iikpt,1:mband_occ_k,1)) 2729 ! end do 2730 ! 2731 ! if(proc_max>(nproc_hf-1)) then 2732 ! !* Too much proc. requested 2733 ! write(msg, '(a,a,a,i4,a,a,a)' )& 2734 ! & ' The number of processors mentioned in the kpt_distrb file',ch10,& 2735 ! & ' must be lower or equal to the actual number of processors =',& 2736 ! & nproc_hf-1,ch10,& 2737 ! & ' Action: change the kpt_distrb file, or increase the',& 2738 ! & ' number of processors.' 2739 ! ABI_ERROR(msg) 2740 ! end if 2741 ! 2742 ! if(proc_max/=(nproc_hf-1)) then 2743 ! !* Too few proc. used 2744 ! write(msg, '(a,i4,a,a,a,i4,a,a,a)' )& 2745 ! & ' Only ',proc_max+1,' processors are used (from kpt_distrb file),',ch10,& 2746 ! & ' when',nproc_hf,' processors are available.',ch10,& 2747 ! & ' Action: adjust number of processors and kpt_distrb file.' 2748 ! ABI_ERROR(msg) 2749 ! end if 2750 ! 2751 ! if(proc_min<0) then 2752 ! write(msg, '(a,a,a)' )& 2753 ! & ' The number of processors must be bigger than 0 in kpt_distrb file.',ch10,& 2754 ! & ' Action: modify kpt_distrb file.' 2755 ! ABI_ERROR(msg) 2756 ! end if 2757 ! else 2758 ! !* The file does not exist... 2759 ! end if ! file_exist 2760 2761 !DEBUG 2762 !write(std_out,*)' distrb2_hf: exit ' 2763 !ENDDEBUG 2764 2765 end subroutine distrb2_hf
ABINIT/initmpi_atom [ Functions ]
NAME
initmpi_atom
FUNCTION
Initializes the mpi information for parallelism over atoms (PAW).
INPUTS
dtset <type(dataset_type)>=all input variables for this dataset mpi_enreg= information about MPI parallelization
OUTPUT
mpi_enreg= information about MPI parallelization comm_atom =communicator over atoms nproc_atom =size of the communicator over atoms my_natom =number of atoms treated by current proc my_atmtab(mpi_enreg%natom)=indexes of the atoms treated by current processor
SOURCE
1032 subroutine initmpi_atom(dtset,mpi_enreg) 1033 1034 !Arguments ------------------------------------ 1035 !scalars 1036 type(dataset_type),intent(in) :: dtset 1037 type(MPI_type),intent(inout) :: mpi_enreg 1038 1039 !Local variables------------------------------- 1040 !scalars 1041 logical :: my_atmtab_allocated,paral_atom 1042 character(len=500) :: msg 1043 integer :: iatom 1044 1045 ! *********************************************************************** 1046 1047 DBG_ENTER("COLL") 1048 1049 mpi_enreg%nproc_atom=1 1050 mpi_enreg%comm_atom=xmpi_comm_self 1051 mpi_enreg%my_natom=dtset%natom 1052 if (associated(mpi_enreg%my_atmtab))then 1053 ABI_FREE(mpi_enreg%my_atmtab) 1054 end if 1055 nullify(mpi_enreg%my_atmtab) 1056 1057 if (xmpi_paral==0) then 1058 mpi_enreg%nproc_atom=0 1059 ABI_MALLOC(mpi_enreg%my_atmtab,(0)) 1060 return 1061 end if 1062 1063 !Check compatibility 1064 if (dtset%paral_atom>0) then 1065 msg='' 1066 if (dtset%usepaw==0) msg= 'Parallelisation over atoms not compatible with usepaw=0 !' 1067 if (dtset%usedmft==1) msg=' Parallelisation over atoms not compatible with usedmft=1 !' 1068 if (dtset%usewvl==1) msg= 'Parallelisation over atoms not compatible with usewvl=1 !' 1069 if (dtset%prtden>1.and.dtset%paral_kgb<=0) & 1070 & msg= 'Parallelisation over atoms not compatible with prtden>1 (PAW AE densities) !' 1071 if (dtset%optdriver/=RUNL_GSTATE.and.dtset%optdriver/=RUNL_RESPFN) & 1072 & msg=' Parallelisation over atoms only compatible with GS or RF !' 1073 if (dtset%macro_uj/=0)msg=' Parallelisation over atoms not compatible with macro_uj!=0 !' 1074 if (msg/='') then 1075 ABI_ERROR(msg) 1076 end if 1077 end if 1078 1079 if (mpi_enreg%comm_atom==xmpi_comm_null) then 1080 mpi_enreg%nproc_atom=0;mpi_enreg%my_natom=0 1081 ABI_MALLOC(mpi_enreg%my_atmtab,(0)) 1082 return 1083 end if 1084 1085 if (dtset%paral_atom>0) then 1086 1087 ! Build correct atom communicator 1088 if (dtset%optdriver==RUNL_GSTATE.and.dtset%paral_kgb==1) then 1089 mpi_enreg%comm_atom=mpi_enreg%comm_kptband 1090 else 1091 mpi_enreg%comm_atom=mpi_enreg%comm_cell 1092 end if 1093 1094 ! Get number of processors sharing the atomic data distribution 1095 mpi_enreg%nproc_atom=xmpi_comm_size(mpi_enreg%comm_atom) 1096 1097 ! Get local number of atoms 1098 call get_my_natom(mpi_enreg%comm_atom,mpi_enreg%my_natom,dtset%natom) 1099 paral_atom=(mpi_enreg%my_natom/=dtset%natom) 1100 1101 ! Build atom table 1102 if (mpi_enreg%my_natom>0.and.paral_atom) then 1103 my_atmtab_allocated=.false. 1104 call get_my_atmtab(mpi_enreg%comm_atom,mpi_enreg%my_atmtab,my_atmtab_allocated, & 1105 & paral_atom,dtset%natom) 1106 else if (.not.paral_atom) then 1107 ABI_MALLOC(mpi_enreg%my_atmtab,(dtset%natom)) 1108 mpi_enreg%my_atmtab(1:dtset%natom)=(/(iatom, iatom=1,dtset%natom)/) 1109 else if (mpi_enreg%my_natom==0) then 1110 ABI_MALLOC(mpi_enreg%my_atmtab,(0)) 1111 end if 1112 1113 end if 1114 1115 DBG_EXIT("COLL") 1116 1117 end subroutine initmpi_atom
ABINIT/initmpi_band [ Functions ]
NAME
initmpi_band
FUNCTION
Initializes the mpi information for band parallelism (paralbd=1).
INPUTS
mpi_enreg= information about MPI parallelization nband(nkpt*nsppol)= number of bands per k point, for each spin nkpt= number of k-points nsppol= 1 for unpolarized, 2 for polarized
OUTPUT
mpi_enreg=information about MPI parallelization mpi_enreg%comm_band=communicator of BAND set
SOURCE
2047 subroutine initmpi_band(mkmem,mpi_enreg,nband,nkpt,nsppol) 2048 2049 !Arguments ------------------------------------ 2050 !scalars 2051 integer,intent(in) :: mkmem 2052 integer,intent(in) :: nkpt,nsppol 2053 integer,intent(in) :: nband(nkpt*nsppol) 2054 type(MPI_type),intent(inout) :: mpi_enreg 2055 2056 !Local variables------------------------------- 2057 !scalars 2058 integer :: ii,ikpt,iproc_min,iproc_max,irank,isppol 2059 integer :: me,nband_k,nproc,nb_per_proc,nrank,nstates,spacecomm 2060 integer :: maxproc_bandpool, mband 2061 character(len=500) :: msg 2062 !arrays 2063 integer,allocatable :: ranks(:) 2064 2065 ! *********************************************************************** 2066 2067 ! reinstate default just to be sure - can be switched inside a previous part of the same dtset! 2068 mpi_enreg%comm_band=xmpi_comm_self 2069 mpi_enreg%nproc_band=1 2070 2071 mband = maxval(nband) 2072 2073 ABI_UNUSED(mkmem) 2074 2075 ! Comm_kpt is supposed to treat spins, k-points and bands 2076 !MJV: I think we need to make a proper subcomm here, not treat bands inside the same comm... 2077 spacecomm=mpi_enreg%comm_kpt 2078 nproc=mpi_enreg%nproc_spkpt 2079 2080 ! make sure we have saturated kpt parallelization 2081 if (mpi_enreg%paralbd==1 .and. xmpi_paral==1 .and. nproc >= 2*nkpt*nsppol) then 2082 2083 ! number of procs per kpt/spin, on which we can distribute bands 2084 maxproc_bandpool=floor(nproc*one/(nkpt*nsppol)) 2085 me=mpi_enreg%me_kpt 2086 2087 !! total number of states/bands, over all k and spin 2088 nstates=sum(nband(1:nkpt*nsppol)) 2089 ! number of bands per proc in the band pool 2090 ! nb_per_proc=nstates/maxproc_bandpool 2091 2092 do nb_per_proc = mband / maxproc_bandpool, mband 2093 if (mod(mband,nb_per_proc)==0) exit 2094 end do 2095 2096 nrank=0 2097 2098 ! NB: do this for all procs even if mkmem == 0, otherwise the subcomm call below fails 2099 if (nb_per_proc<mband) then 2100 do isppol=1,nsppol 2101 do ikpt=1,nkpt 2102 ii=ikpt+(isppol-1)*nkpt 2103 nband_k=nband(ii) 2104 if (nb_per_proc<nband_k) then 2105 iproc_min=minval(mpi_enreg%proc_distrb(ikpt,:,isppol)) 2106 iproc_max=maxval(mpi_enreg%proc_distrb(ikpt,:,isppol)) 2107 if ((me>=iproc_min).and.(me<=iproc_max)) then 2108 nrank=iproc_max-iproc_min+1 2109 if (.not.allocated(ranks)) then 2110 ABI_MALLOC(ranks,(nrank)) 2111 if (nrank>0) ranks=(/((iproc_min+irank-1),irank=1,nrank)/) 2112 else if (nrank/=size(ranks)) then 2113 ! TODO MJV: still can not lift this restriction... 2114 ABI_BUG('Number of bands per proc should be the same for all k-points!') 2115 end if 2116 end if 2117 end if 2118 end do 2119 end do 2120 2121 if (.not.allocated(ranks)) then 2122 nrank = 0 2123 ABI_MALLOC(ranks,(0)) 2124 end if 2125 2126 ! ABI_CHECK(nrank*nkpt==nproc, ' band and k-point distribution should be rectangular: make sure nproc=nkpt*integer') 2127 2128 if (nrank*nkpt*nsppol < nproc) then 2129 write(unit=msg,fmt='(a,i6,2a,i6,a,i6,4a)') & 2130 'The number of processors nproc = ', nproc, ch10,& 2131 ' is not equal to nrank (=',nrank,& 2132 ') times nkpt*nsppol (',nkpt*nsppol,& 2133 ' , which may change with perturbation) !',ch10,& 2134 ' This is inefficient (load unbalancing). Adjust nband to have a divisor <= nproc/nkpt/nsppol',ch10 2135 ABI_WARNING(msg) 2136 end if 2137 ! NB: everyone in spacecomm has to call subcomm, even if it is a trivial call with self_comm for the subcomm 2138 mpi_enreg%comm_band=xmpi_subcomm(spacecomm,nrank,ranks, my_rank_in_group=mpi_enreg%me_band) 2139 mpi_enreg%nproc_band=nrank 2140 ! mpi_enreg%me_band=mod(me, nrank) 2141 2142 write(msg,'(4(a,i0))') 'P Present parallel dimensions: nkpt= ',nkpt,' nsppol ',nsppol,& 2143 ' nband per processor= ', nb_per_proc, ' npband= ',nrank 2144 call wrtout(std_out,msg) 2145 2146 ABI_FREE(ranks) 2147 end if 2148 2149 end if 2150 2151 end subroutine initmpi_band
ABINIT/initmpi_grid [ Functions ]
NAME
initmpi_grid
FUNCTION
Initializes the MPI information for the grid: * 2D if parallization KPT/FFT (paral_kgb == 0 & MPI) * 3D if parallization KPT/FFT/BAND (paral_kgb == 1 & MPI) * 2D in case of an Hartree-Fock calculation
INPUTS
OUTPUT
SOURCE
1178 subroutine initmpi_grid(mpi_enreg) 1179 1180 !Arguments ------------------------------------ 1181 type(MPI_type),intent(inout) :: mpi_enreg 1182 1183 !Local variables------------------------------- 1184 !scalars 1185 integer :: nproc,nproc_eff,spacecomm 1186 character(len=500) :: msg 1187 #if defined HAVE_MPI 1188 integer :: commcart_4d,dimcart,ierr,me_cart_4d 1189 integer :: commcart_2d,me_cart_2d 1190 logical :: reorder 1191 !arrays 1192 integer,allocatable :: coords(:),sizecart(:) 1193 logical,allocatable :: periode(:), keepdim(:) 1194 #endif 1195 1196 ! ********************************************************************* 1197 1198 DBG_ENTER("COLL") 1199 1200 !Select the correct "world" communicator" 1201 nproc=mpi_enreg%nproc_cell 1202 if(mpi_enreg%paral_pert==1) nproc=mpi_enreg%nproc_cell 1203 spacecomm=mpi_enreg%comm_cell 1204 1205 !Fake values for null communicator 1206 if (nproc==0) then 1207 mpi_enreg%nproc_fft = 0 1208 mpi_enreg%nproc_band = 0 1209 mpi_enreg%nproc_hf = 0 1210 mpi_enreg%nproc_spkpt = 0 1211 mpi_enreg%nproc_spinor = 0 1212 mpi_enreg%comm_fft = xmpi_comm_null 1213 mpi_enreg%comm_band = xmpi_comm_null 1214 mpi_enreg%comm_hf = xmpi_comm_null 1215 mpi_enreg%comm_kpt = xmpi_comm_null 1216 mpi_enreg%comm_kptband = xmpi_comm_null 1217 mpi_enreg%comm_spinor = xmpi_comm_null 1218 mpi_enreg%comm_bandspinor = xmpi_comm_null 1219 mpi_enreg%comm_spinorfft = xmpi_comm_null 1220 mpi_enreg%comm_bandfft = xmpi_comm_null 1221 mpi_enreg%comm_bandspinorfft = xmpi_comm_null 1222 mpi_enreg%bandpp = 1 1223 return 1224 end if 1225 1226 #if defined HAVE_MPI 1227 if (mpi_enreg%paral_hf==0) then 1228 ! either the option Fock exchange is not active or there is no parallelization on Fock exchange calculation. 1229 1230 if (mpi_enreg%nproc_spinor>1) mpi_enreg%paral_spinor=1 1231 1232 !Effective number of processors used for the grid 1233 nproc_eff=mpi_enreg%nproc_fft*mpi_enreg%nproc_band *mpi_enreg%nproc_spkpt*mpi_enreg%nproc_spinor 1234 if(nproc_eff/=nproc) then 1235 write(msg,'(4a,5(a,i0,a))') & 1236 ' The number of band*FFT*spin*kpt*spinor processors, npband*npfft*np_spkpt*npspinor should be',ch10,& 1237 ' equal to the total number of processors, nproc.',ch10,& 1238 ' However, npband =',mpi_enreg%nproc_band, ch10, & 1239 ' npfft =',mpi_enreg%nproc_fft, ch10, & 1240 ' np_spkpt =',mpi_enreg%nproc_spkpt, ch10, & 1241 ' npspinor =',mpi_enreg%nproc_spinor, ch10, & 1242 ' nproc =',nproc,ch10 1243 ABI_WARNING(msg) 1244 end if 1245 1246 ! Nothing to do if only 1 proc 1247 if (nproc_eff==1) return 1248 1249 ! Initialize the communicator for Hartree-Fock to xmpi_comm_self 1250 mpi_enreg%me_hf =0 1251 mpi_enreg%comm_hf=xmpi_comm_self 1252 1253 if(mpi_enreg%paral_kgb==0) then 1254 mpi_enreg%me_fft =0 1255 mpi_enreg%me_band=0 1256 mpi_enreg%me_kpt =mpi_enreg%me_cell 1257 mpi_enreg%me_spinor=0 1258 mpi_enreg%comm_fft=xmpi_comm_self 1259 mpi_enreg%comm_band=xmpi_comm_self 1260 mpi_enreg%comm_kpt=mpi_enreg%comm_cell 1261 mpi_enreg%comm_spinor=xmpi_comm_self 1262 mpi_enreg%comm_bandspinor=xmpi_comm_self 1263 mpi_enreg%comm_kptband=mpi_enreg%comm_cell 1264 mpi_enreg%comm_spinorfft=xmpi_comm_self 1265 mpi_enreg%comm_bandfft=xmpi_comm_self 1266 mpi_enreg%comm_bandspinorfft=xmpi_comm_self 1267 else 1268 ! CREATE THE 4D GRID 1269 ! ================================================== 1270 1271 ! Create the global cartesian 4D- communicator 1272 ! valgrind claims this is not deallocated in test v5/72 1273 ! Can someone knowledgable check? 1274 dimcart=4 1275 ABI_MALLOC(sizecart,(dimcart)) 1276 ABI_MALLOC(periode,(dimcart)) 1277 ! MT 2012-june: change the order of the indexes; not sure this is efficient 1278 ! (not efficient on TGCC-Curie). 1279 sizecart(1)=mpi_enreg%nproc_spkpt ! mpi_enreg%nproc_spkpt 1280 sizecart(2)=mpi_enreg%nproc_band ! mpi_enreg%nproc_band 1281 sizecart(3)=mpi_enreg%nproc_spinor ! mpi_enreg%nproc_spinor 1282 sizecart(4)=mpi_enreg%nproc_fft ! mpi_enreg%nproc_fft 1283 periode(:)=.false.;reorder=.false. 1284 call MPI_CART_CREATE(spacecomm,dimcart,sizecart,periode,reorder,commcart_4d,ierr) 1285 ABI_FREE(periode) 1286 ABI_FREE(sizecart) 1287 1288 ! Find the index and coordinates of the current processor 1289 call MPI_COMM_RANK(commcart_4d, me_cart_4d, ierr) 1290 ABI_MALLOC(coords,(dimcart)) 1291 call MPI_CART_COORDS(commcart_4d, me_cart_4d,dimcart,coords,ierr) 1292 mpi_enreg%me_kpt =coords(1) 1293 mpi_enreg%me_band=coords(2) 1294 mpi_enreg%me_spinor=coords(3) 1295 mpi_enreg%me_fft =coords(4) 1296 ABI_FREE(coords) 1297 1298 ABI_MALLOC(keepdim,(dimcart)) 1299 1300 ! Create the communicator for fft distribution 1301 keepdim(1)=.false. 1302 keepdim(2)=.false. 1303 keepdim(3)=.false. 1304 keepdim(4)=.true. 1305 call MPI_CART_SUB(commcart_4d, keepdim, mpi_enreg%comm_fft,ierr) 1306 1307 ! Create the communicator for band distribution 1308 keepdim(1)=.false. 1309 keepdim(2)=.true. 1310 keepdim(3)=.false. 1311 keepdim(4)=.false. 1312 call MPI_CART_SUB(commcart_4d, keepdim, mpi_enreg%comm_band,ierr) 1313 1314 ! Create the communicator for kpt distribution 1315 keepdim(1)=.true. 1316 keepdim(2)=.false. 1317 keepdim(3)=.false. 1318 keepdim(4)=.false. 1319 call MPI_CART_SUB(commcart_4d, keepdim, mpi_enreg%comm_kpt,ierr) 1320 1321 ! Create the communicator for spinor distribution 1322 keepdim(1)=.false. 1323 keepdim(2)=.false. 1324 keepdim(3)=.true. 1325 keepdim(4)=.false. 1326 call MPI_CART_SUB(commcart_4d, keepdim, mpi_enreg%comm_spinor,ierr) 1327 1328 ! Create the communicator for band-spinor distribution 1329 keepdim(1)=.false. 1330 keepdim(2)=.true. 1331 keepdim(3)=.true. 1332 keepdim(4)=.false. 1333 call MPI_CART_SUB(commcart_4d, keepdim, mpi_enreg%comm_bandspinor,ierr) 1334 if (ierr /= MPI_SUCCESS ) then 1335 call xmpi_abort(mpi_enreg%comm_world,ierr) 1336 end if 1337 1338 ! Create the communicator for kpt-band distribution 1339 keepdim(1)=.true. 1340 keepdim(2)=.true. 1341 keepdim(3)=.false. 1342 keepdim(4)=.false. 1343 call MPI_CART_SUB(commcart_4d, keepdim, mpi_enreg%comm_kptband,ierr) 1344 1345 ! Create the communicator for fft-spinor distribution 1346 keepdim(1)=.false. 1347 keepdim(2)=.false. 1348 keepdim(3)=.true. 1349 keepdim(4)=.true. 1350 call MPI_CART_SUB(commcart_4d, keepdim, mpi_enreg%comm_spinorfft,ierr) 1351 1352 ! Create the communicator for fft-band distribution 1353 keepdim(1)=.false. 1354 keepdim(2)=.true. 1355 keepdim(3)=.false. 1356 keepdim(4)=.true. 1357 call MPI_CART_SUB(commcart_4d, keepdim, mpi_enreg%comm_bandfft,ierr) 1358 1359 ! Create the communicator for fft-band-spinor distribution 1360 keepdim(1)=.false. 1361 keepdim(2)=.true. 1362 keepdim(3)=.true. 1363 keepdim(4)=.true. 1364 call MPI_CART_SUB(commcart_4d, keepdim, mpi_enreg%comm_bandspinorfft,ierr) 1365 1366 ABI_FREE(keepdim) 1367 call xmpi_comm_free(commcart_4d) 1368 end if 1369 1370 !Write some data 1371 !write(msg,'(a,4i5)') 'npfft, npband, npspinor and np_spkpt: ',& 1372 !mpi_enreg%nproc_fft,mpi_enreg%nproc_band, mpi_enreg%nproc_spinor,mpi_enreg%nproc_spkpt 1373 !call wrtout(std_out,msg,'COLL') 1374 !write(msg,'(a,4i5)') 'me_fft, me_band, me_spinor , me_kpt: ',& 1375 !mpi_enreg%me_fft,mpi_enreg%me_band,mpi_enreg%me_spinor, mpi_enreg%me_kpt 1376 !call wrtout(std_out,msg,'COLL') 1377 1378 else ! paral_hf==1 1379 !* Option Hartree-Fock is active and more than 1 processor is dedicated to the parallelization over occupied states. 1380 1381 !* Initialize the values of fft, band and spinor communicators, as in the case paral_kgb==0. 1382 mpi_enreg%me_fft =0 1383 mpi_enreg%me_band=0 1384 mpi_enreg%me_spinor=0 1385 mpi_enreg%comm_fft=xmpi_comm_self 1386 mpi_enreg%comm_band=xmpi_comm_self 1387 mpi_enreg%comm_spinor=xmpi_comm_self 1388 mpi_enreg%comm_bandspinor=xmpi_comm_self 1389 mpi_enreg%comm_kptband=mpi_enreg%comm_cell 1390 mpi_enreg%comm_spinorfft=xmpi_comm_self 1391 mpi_enreg%comm_bandfft=xmpi_comm_self 1392 mpi_enreg%comm_bandspinorfft=xmpi_comm_self 1393 1394 !* Create the global cartesian 2D- communicator 1395 dimcart=2 1396 ABI_MALLOC(sizecart,(dimcart)) 1397 ABI_MALLOC(periode,(dimcart)) 1398 sizecart(1)=mpi_enreg%nproc_spkpt ! mpi_enreg%nproc_spkpt 1399 sizecart(2)=mpi_enreg%nproc_hf ! mpi_enreg%nproc_hf 1400 periode(:)=.false.;reorder=.false. 1401 call MPI_CART_CREATE(spacecomm,dimcart,sizecart,periode,reorder,commcart_2d,ierr) 1402 ABI_FREE(periode) 1403 ABI_FREE(sizecart) 1404 1405 !* Find the index and coordinates of the current processor 1406 call MPI_COMM_RANK(commcart_2d, me_cart_2d, ierr) 1407 ABI_MALLOC(coords,(dimcart)) 1408 call MPI_CART_COORDS(commcart_2d, me_cart_2d,dimcart,coords,ierr) 1409 mpi_enreg%me_kpt =coords(1) 1410 mpi_enreg%me_hf=coords(2) 1411 ABI_FREE(coords) 1412 1413 ABI_MALLOC(keepdim,(dimcart)) 1414 1415 !* Create the communicator for kpt distribution 1416 keepdim(1)=.true. 1417 keepdim(2)=.false. 1418 call MPI_CART_SUB(commcart_2d, keepdim, mpi_enreg%comm_kpt,ierr) 1419 1420 !* Create the communicator for hf distribution 1421 keepdim(1)=.false. 1422 keepdim(2)=.true. 1423 call MPI_CART_SUB(commcart_2d, keepdim, mpi_enreg%comm_hf,ierr) 1424 1425 ABI_FREE(keepdim) 1426 call xmpi_comm_free(commcart_2d) 1427 1428 !* Write some data 1429 write(msg,'(a,2(1x,i0))') 'nphf and np_spkpt: ',mpi_enreg%nproc_hf, mpi_enreg%nproc_spkpt 1430 call wrtout(std_out,msg) 1431 write(msg,'(a,2(1x,i0))') 'me_hf, me_kpt: ',mpi_enreg%me_hf, mpi_enreg%me_kpt 1432 call wrtout(std_out,msg) 1433 end if 1434 #endif 1435 1436 DBG_EXIT("COLL") 1437 1438 end subroutine initmpi_grid
ABINIT/initmpi_img [ Functions ]
NAME
initmpi_img
FUNCTION
Initializes the mpi information for parallelism over images of the cell (npimage>1).
INPUTS
dtset <type(dataset_type)>=all input variables in this dataset mpi_enreg= information about MPI parallelization option= see below
OUTPUT
mpi_enreg%my_nimage= number of images of the cell treated by current proc ===== if option==1 or option==-1 mpi_enreg%my_imgtab= indexes of images of the cell treated by current proc ===== if option==2 or option==3 or option==-1 mpi_enreg%comm_cell=Communicator over all processors treating the same image mpi_enreg%nproc_cell=size of comm_cell mpi_enreg%me_cell=my rank in comm_cell ===== if option==3 or option==-1 mpi_enreg%comm_img=Communicator over all images mpi_enreg%nproc_img=size of comm_img mpi_enreg%me_img=my rank in comm_img mpi_enreg%distrb_img(:)=index of processor treating each image (in comm_img communicator)
SOURCE
1544 subroutine initmpi_img(dtset,mpi_enreg,option) 1545 1546 !Arguments ------------------------------------ 1547 integer,intent(in) :: option 1548 type(dataset_type),intent(in) :: dtset 1549 type(MPI_type),intent(inout) :: mpi_enreg 1550 1551 !Local variables------------------------------- 1552 integer :: imod,irank,iprocmax,iprocmin,jrank 1553 integer :: ndynimage_eff,nimage_eff,nproc_per_image,nrank 1554 logical,parameter :: debug=.false. 1555 integer,allocatable :: ranks(:) 1556 character(len=500) :: msg 1557 1558 !integer :: group_cell,ierr 1559 1560 ! *********************************************************************** 1561 1562 DBG_ENTER("COLL") 1563 if (option/=0) then 1564 mpi_enreg%comm_img=xmpi_comm_self 1565 mpi_enreg%comm_cell=mpi_enreg%comm_world 1566 end if 1567 1568 if (xmpi_paral==1.and.dtset%npimage>1.and.dtset%npimage<=mpi_enreg%nproc.and. dtset%optdriver==RUNL_GSTATE) then 1569 1570 ! Activate flag for parallelization over images 1571 mpi_enreg%paral_img=1 1572 1573 ndynimage_eff=dtset%ndynimage;if (dtset%ntimimage<=1) ndynimage_eff=0 1574 1575 ! Print several warnings 1576 if (option==0) then 1577 nimage_eff=max(ndynimage_eff,dtset%nimage-ndynimage_eff) 1578 if (dtset%npimage>nimage_eff) then 1579 write(unit=msg,fmt='(3a,i4,a,i4,4a)') & 1580 'The number of processors used for the parallelization',ch10,& 1581 ' over images (npimage=',dtset%npimage,& 1582 ') is greater than the number of dynamic (or static) images (',nimage_eff,') !',ch10,& 1583 ' This is inefficient.',ch10 1584 ABI_WARNING(msg) 1585 end if 1586 if (dtset%npimage>mpi_enreg%nproc) then 1587 write(unit=msg,fmt='(3a,i6,a,i4,4a)') & 1588 'The number of processors used for the parallelization',ch10,& 1589 ' over images (nproc=',mpi_enreg%nproc,& 1590 ') is smaller than npimage in input file (',dtset%npimage,& 1591 ')!',ch10,' This is unconsistent.',ch10 1592 ABI_ERROR(msg) 1593 end if 1594 if (mod(nimage_eff,dtset%npimage)/=0) then 1595 write(unit=msg,fmt='(3a,i4,a,i4,4a)') & 1596 'The number of processors used for the parallelization',ch10,& 1597 ' over images (npimage=',dtset%npimage,& 1598 ') does not divide the number of dynamic images (',nimage_eff,& 1599 ') !',ch10,' This is inefficient (charge unbalancing).',ch10 1600 ABI_WARNING(msg) 1601 end if 1602 end if 1603 1604 ! # of images treated by current proc 1605 nproc_per_image=mpi_enreg%nproc/dtset%npimage 1606 iprocmax=nproc_per_image*dtset%npimage-1 1607 if (mpi_enreg%me<=iprocmax) then 1608 mpi_enreg%my_nimage=(ndynimage_eff/dtset%npimage)+((dtset%nimage-ndynimage_eff)/dtset%npimage) 1609 imod=mod(ndynimage_eff,dtset%npimage)-1 1610 if (mpi_enreg%me/nproc_per_image<=imod) mpi_enreg%my_nimage=mpi_enreg%my_nimage+1 1611 imod=mod((dtset%nimage-ndynimage_eff),dtset%npimage)-1 1612 if (mpi_enreg%me/nproc_per_image<=imod) mpi_enreg%my_nimage=mpi_enreg%my_nimage+1 1613 else 1614 mpi_enreg%my_nimage=0 1615 end if 1616 if (option==1.or.option==-1) then 1617 ! Indexes of images treated by current proc 1618 if (mpi_enreg%me<=iprocmax) then 1619 ABI_MALLOC(mpi_enreg%my_imgtab,(mpi_enreg%my_nimage)) 1620 nrank=0 1621 imod=mpi_enreg%me/nproc_per_image+1;imod=mod(imod,dtset%npimage) 1622 ! Dynamic images 1623 irank=0 1624 do jrank=1,dtset%nimage 1625 if (dtset%dynimage(jrank)/=0.and.dtset%ntimimage>1) then 1626 irank=irank+1 1627 if (mod(irank,dtset%npimage)==imod) then 1628 nrank=nrank+1 1629 mpi_enreg%my_imgtab(nrank)=jrank 1630 end if 1631 end if 1632 end do 1633 ! Static images 1634 irank=0 1635 do jrank=1,dtset%nimage 1636 if (dtset%dynimage(jrank)==0.or.dtset%ntimimage<=1) then 1637 irank=irank+1 1638 if (mod(irank,dtset%npimage)==imod) then 1639 nrank=nrank+1 1640 mpi_enreg%my_imgtab(nrank)=jrank 1641 end if 1642 end if 1643 end do 1644 if (nrank/=mpi_enreg%my_nimage) then 1645 ABI_BUG('Error on nrank !') 1646 end if 1647 ! Sort images by increasing index (this step is MANDATORY !!) 1648 ABI_MALLOC(ranks,(nrank)) 1649 call sort_int(nrank,mpi_enreg%my_imgtab,ranks) 1650 ABI_FREE(ranks) 1651 else 1652 ABI_MALLOC(mpi_enreg%my_imgtab,(0)) 1653 end if 1654 end if 1655 if (option==2.or.option==3.or.option==-1) then 1656 ! Communicator over one image 1657 if (mpi_enreg%me<=iprocmax) then 1658 ABI_MALLOC(ranks,(nproc_per_image)) 1659 iprocmin=(mpi_enreg%me/nproc_per_image)*nproc_per_image 1660 ranks=(/((iprocmin+irank-1),irank=1,nproc_per_image)/) 1661 mpi_enreg%comm_cell=xmpi_subcomm(mpi_enreg%comm_world,nproc_per_image,ranks) 1662 ABI_FREE(ranks) 1663 mpi_enreg%me_cell=xmpi_comm_rank(mpi_enreg%comm_cell) 1664 mpi_enreg%nproc_cell=nproc_per_image 1665 if (mpi_enreg%me_cell==0.and.mod(mpi_enreg%me,nproc_per_image)/=0) then 1666 ABI_BUG('Error on me_cell !') 1667 end if 1668 else 1669 mpi_enreg%comm_img=xmpi_comm_null 1670 mpi_enreg%nproc_cell=0 1671 mpi_enreg%me_cell=-1 1672 end if 1673 end if 1674 if (option==3.or.option==-1) then 1675 ! Communicator over all images 1676 if (mpi_enreg%me<=iprocmax) then 1677 ABI_MALLOC(ranks,(dtset%npimage)) 1678 iprocmin=mod(mpi_enreg%me,nproc_per_image) 1679 ranks=(/((iprocmin+(irank-1)*nproc_per_image),irank=1,dtset%npimage)/) 1680 mpi_enreg%comm_img=xmpi_subcomm(mpi_enreg%comm_world,dtset%npimage,ranks) 1681 ABI_FREE(ranks) 1682 mpi_enreg%me_img=xmpi_comm_rank(mpi_enreg%comm_img) 1683 mpi_enreg%nproc_img=dtset%npimage 1684 if (iprocmin==0.and.mpi_enreg%me_img==0.and.mpi_enreg%me/=0) then 1685 ABI_BUG('Error on me_img!') 1686 end if 1687 ABI_MALLOC(mpi_enreg%distrb_img,(dtset%nimage)) 1688 ! Dynamic images 1689 nrank=0 1690 do irank=1,dtset%nimage 1691 if (dtset%dynimage(irank)/=0.and.dtset%ntimimage>1) then 1692 nrank=nrank+1 1693 mpi_enreg%distrb_img(irank)=mod(nrank,dtset%npimage)-1 1694 if (mpi_enreg%distrb_img(irank)==-1) mpi_enreg%distrb_img(irank)=dtset%npimage-1 1695 end if 1696 end do 1697 ! Static images 1698 nrank=0 1699 do irank=1,dtset%nimage 1700 if (dtset%dynimage(irank)==0.or.dtset%ntimimage<=1) then 1701 nrank=nrank+1 1702 mpi_enreg%distrb_img(irank)=mod(nrank,dtset%npimage)-1 1703 if (mpi_enreg%distrb_img(irank)==-1) mpi_enreg%distrb_img(irank)=dtset%npimage-1 1704 end if 1705 end do 1706 else 1707 mpi_enreg%comm_img=xmpi_comm_null 1708 mpi_enreg%nproc_img=0 1709 mpi_enreg%me_img=-1 1710 ABI_MALLOC(mpi_enreg%distrb_img,(0)) 1711 end if 1712 end if 1713 1714 ! if (debug) then 1715 ! write(200+mpi_enreg%me,*) "==================================" 1716 ! write(200+mpi_enreg%me,*) "DEBUGGING STATMENTS IN INITMPI_IMG" 1717 ! write(200+mpi_enreg%me,*) "==================================" 1718 ! write(200+mpi_enreg%me,*) "option =",option 1719 ! write(200+mpi_enreg%me,*) "MPI_UNDEFINED =",MPI_UNDEFINED 1720 ! write(200+mpi_enreg%me,*) "MPI_IDENT =",MPI_IDENT 1721 ! write(200+mpi_enreg%me,*) "MPI_CONGRUENT =",MPI_CONGRUENT 1722 ! write(200+mpi_enreg%me,*) "MPI_SIMILAR =",MPI_SIMILAR 1723 ! write(200+mpi_enreg%me,*) "MPI_UNEQUAL =",MPI_UNEQUAL 1724 ! write(200+mpi_enreg%me,*) "null_comm =",MPI_COMM_NULL 1725 ! write(200+mpi_enreg%me,*) "self_comm =",xmpi_comm_self 1726 ! write(200+mpi_enreg%me,*) "world_comm =",mpi_enreg%comm_world 1727 ! write(200+mpi_enreg%me,*) "empty_group =",MPI_GROUP_EMPTY 1728 ! write(200+mpi_enreg%me,*) "nimage =",mpi_enreg%my_nimage 1729 ! write(200+mpi_enreg%me,*) "nproc_per_image=",nproc_per_image 1730 ! call MPI_COMM_SIZE(mpi_enreg%comm_world,irank,ierr) 1731 ! write(200+mpi_enreg%me,*) "Size of world_comm =",irank 1732 ! call MPI_COMM_RANK(mpi_enreg%comm_world,irank,ierr) 1733 ! write(200+mpi_enreg%me,*) "My rank in world_comm =",irank 1734 ! if (option==1.or.option==-1) then 1735 ! write(200+mpi_enreg%me,*) "index_img=",mpi_enreg%my_imgtab(:) 1736 ! end if 1737 ! if (option==2.or.option==3.or.option==-1) then 1738 ! write(200+mpi_enreg%me,*) "nproc_cell =",mpi_enreg%nproc_cell 1739 ! write(200+mpi_enreg%me,*) "me_cell =",mpi_enreg%me_cell 1740 ! call xmpi_comm_group(mpi_enreg%comm_cell,group_cell,ierr) 1741 ! write(200+mpi_enreg%me,*) "group_cell =",group_cell 1742 ! write(200+mpi_enreg%me,*) "comm_cell =",mpi_enreg%comm_cell 1743 ! if (group_cell/=MPI_GROUP_EMPTY) then 1744 ! call MPI_GROUP_SIZE(group_cell,irank,ierr) 1745 ! write(200+mpi_enreg%me,*) "Size of group_cell =",irank 1746 ! call MPI_GROUP_RANK(group_cell,irank,ierr) 1747 ! write(200+mpi_enreg%me,*) "My rank in group_cell=",irank 1748 ! else 1749 ! write(200+mpi_enreg%me,*) "Size of group_cell =",0 1750 ! write(200+mpi_enreg%me,*) "My rank in group_cell=",-1 1751 ! end if 1752 ! if (mpi_enreg%comm_cell/=MPI_COMM_NULL) then 1753 ! call MPI_COMM_SIZE(mpi_enreg%comm_cell,irank,ierr) 1754 ! write(200+mpi_enreg%me,*) "Size of comm_cell =",irank 1755 ! call MPI_COMM_RANK(mpi_enreg%comm_cell,irank,ierr) 1756 ! write(200+mpi_enreg%me,*) "My rank in comm_cell=",irank 1757 ! call MPI_COMM_COMPARE(mpi_enreg%comm_world,mpi_enreg%comm_cell,irank,ierr) 1758 ! write(200+mpi_enreg%me,*) "Comparison world_comm/comm_cell=",irank 1759 ! call MPI_COMM_COMPARE(xmpi_comm_self,mpi_enreg%comm_cell,irank,ierr) 1760 ! write(200+mpi_enreg%me,*) "Comparison self_comm/comm_cell =",irank 1761 ! else 1762 ! write(200+mpi_enreg%me,*) "Size of comm_cell =",0 1763 ! write(200+mpi_enreg%me,*) "My rank in comm_cell=",-1 1764 ! write(200+mpi_enreg%me,*) "Comparison world_comm/comm_cell=",MPI_UNEQUAL 1765 ! write(200+mpi_enreg%me,*) "Comparison self_comm/comm_cell =",MPI_UNEQUAL 1766 ! end if 1767 ! end if 1768 ! if (option==3.or.option==-1) then 1769 ! write(200+mpi_enreg%me,*) "nproc_img =",mpi_enreg%nproc_img 1770 ! write(200+mpi_enreg%me,*) "me_img =",mpi_enreg%me_img 1771 ! write(200+mpi_enreg%me,*) "img_comm =",mpi_enreg%comm_img 1772 ! if (mpi_enreg%comm_img/=MPI_COMM_NULL) then 1773 ! call MPI_COMM_SIZE(mpi_enreg%comm_img,irank,ierr) 1774 ! write(200+mpi_enreg%me,*) "Size of img_comm =",irank 1775 ! call MPI_COMM_RANK(mpi_enreg%comm_img,irank,ierr) 1776 ! write(200+mpi_enreg%me,*) "My rank in img_comm=",irank 1777 ! call MPI_COMM_COMPARE(mpi_enreg%comm_world,mpi_enreg%comm_img,irank,ierr) 1778 ! write(200+mpi_enreg%me,*) "Comparison world_comm/img_comm=",irank 1779 ! call MPI_COMM_COMPARE(xmpi_comm_self,mpi_enreg%comm_img,irank,ierr) 1780 ! write(200+mpi_enreg%me,*) "Comparison self_comm/img_comm =",irank 1781 ! else 1782 ! write(200+mpi_enreg%me,*) "Size of img_comm =",0 1783 ! write(200+mpi_enreg%me,*) "My rank in img_comm=",-1 1784 ! write(200+mpi_enreg%me,*) "Comparison world_comm/img_comm=",MPI_UNEQUAL 1785 ! write(200+mpi_enreg%me,*) "Comparison self_comm/img_comm =",MPI_UNEQUAL 1786 ! end if 1787 ! write(200+mpi_enreg%me,*) "distrb_img=",mpi_enreg%distrb_img(:) 1788 ! end if 1789 ! write(200+mpi_enreg%me,*) 1790 ! call flush_unit(200+mpi_enreg%me) 1791 ! if (option==-1) stop 1792 ! end if 1793 1794 else 1795 1796 ! Do not activate flag for parallelization over images 1797 mpi_enreg%paral_img=0 1798 ! # of images treated by current proc 1799 if (dtset%optdriver==RUNL_GSTATE) then 1800 mpi_enreg%my_nimage=dtset%nimage 1801 else 1802 mpi_enreg%my_nimage=1 1803 end if 1804 ! Indexes of images treated by current proc 1805 if (option==1.or.option==-1) then 1806 ABI_MALLOC(mpi_enreg%my_imgtab,(mpi_enreg%my_nimage)) 1807 mpi_enreg%my_imgtab=(/(irank,irank=1,mpi_enreg%my_nimage)/) 1808 end if 1809 ! Communicator over all images 1810 if (option==2.or.option==3.or.option==-1) then 1811 ! Communicator for one image 1812 mpi_enreg%nproc_cell=mpi_enreg%nproc 1813 mpi_enreg%me_cell=mpi_enreg%me 1814 end if 1815 if (option==3.or.option==-1) then 1816 ! Communicator over all images 1817 mpi_enreg%nproc_img=1 1818 mpi_enreg%comm_img=xmpi_comm_self 1819 mpi_enreg%me_img=0 1820 ABI_MALLOC(mpi_enreg%distrb_img,(dtset%nimage)) 1821 mpi_enreg%distrb_img(:)=0 1822 end if 1823 end if 1824 1825 DBG_EXIT("COLL") 1826 1827 end subroutine initmpi_img
ABINIT/initmpi_pert [ Functions ]
NAME
initmpi_pert
FUNCTION
Creates group for Parallelization over Perturbations.
INPUTS
dtset <type(dataset_type)>=all input variables in this dataset
OUTPUT
SIDE EFFECTS
mpi_enreg=information about MPI parallelization
SOURCE
1894 subroutine initmpi_pert(dtset,mpi_enreg) 1895 1896 !Arguments ------------------------------------ 1897 !scalars 1898 type(MPI_type),intent(inout) :: mpi_enreg 1899 type(dataset_type),intent(in) :: dtset 1900 1901 !Local variables------------------------------- 1902 !scalars 1903 integer:: iprocmin,irank,npert,nproc_per_cell,nrank,numproc 1904 integer,allocatable :: ranks(:) 1905 !character(len=500) :: msg 1906 !arrays 1907 integer,pointer :: nkpt_rbz(:) 1908 real(dp),pointer :: nband_rbz(:,:) 1909 1910 ! *********************************************************************** 1911 1912 if (mpi_enreg%me_pert<0) then 1913 ABI_ERROR('Error in MPI distribution! Change your proc(s) distribution or use autoparal>0.') 1914 end if 1915 1916 call dtset%get_npert_rbz(nband_rbz, nkpt_rbz, npert) 1917 1918 if (dtset%nppert>=1) then 1919 if (mpi_enreg%comm_cell/=mpi_enreg%comm_world) then 1920 call xmpi_comm_free(mpi_enreg%comm_cell) 1921 end if 1922 mpi_enreg%comm_cell=mpi_enreg%comm_world 1923 1924 ! These values will be properly set in set_pert_comm 1925 mpi_enreg%me_cell=mpi_enreg%me 1926 mpi_enreg%nproc_cell=mpi_enreg%nproc 1927 1928 if (mpi_enreg%me>=0) then 1929 nproc_per_cell=mpi_enreg%nproc/dtset%nppert 1930 ABI_MALLOC(ranks,(dtset%nppert)) 1931 iprocmin=mod(mpi_enreg%me,nproc_per_cell) 1932 ranks=(/((iprocmin+(irank-1)*nproc_per_cell),irank=1,dtset%nppert)/) 1933 mpi_enreg%comm_pert=xmpi_subcomm(mpi_enreg%comm_world,dtset%nppert,ranks) 1934 ABI_FREE(ranks) 1935 mpi_enreg%me_pert=xmpi_comm_rank(mpi_enreg%comm_pert) 1936 mpi_enreg%nproc_pert=dtset%nppert 1937 if (iprocmin==0.and.mpi_enreg%me_pert==0.and.mpi_enreg%me/=0) then 1938 ABI_BUG('Error on me_pert!') 1939 end if 1940 ! Define mpi_enreg%distrb_pert 1941 ABI_MALLOC(mpi_enreg%distrb_pert,(npert)) 1942 nrank=0 1943 do irank=1,npert 1944 nrank=nrank+1 1945 mpi_enreg%distrb_pert(irank)=mod(nrank,dtset%nppert)-1 1946 if (mpi_enreg%distrb_pert(irank)==-1) mpi_enreg%distrb_pert(irank)=dtset%nppert-1 1947 end do 1948 ! Make sure that subrank 0 is working on the last perturbation 1949 ! Swap the ranks if necessary 1950 numproc=mpi_enreg%distrb_pert(npert) 1951 if(numproc/=0) then 1952 do irank=1,npert 1953 if (mpi_enreg%distrb_pert(irank)==numproc) mpi_enreg%distrb_pert(irank)=-2 1954 if (mpi_enreg%distrb_pert(irank)==0) mpi_enreg%distrb_pert(irank)=-3 1955 end do 1956 do irank=1,npert 1957 if (mpi_enreg%distrb_pert(irank)==-2) mpi_enreg%distrb_pert(irank)=0 1958 if (mpi_enreg%distrb_pert(irank)==-3) mpi_enreg%distrb_pert(irank)=numproc 1959 end do 1960 end if 1961 ! Communicator over one cell 1962 ABI_MALLOC(ranks,(nproc_per_cell)) 1963 iprocmin=(mpi_enreg%me/nproc_per_cell)*nproc_per_cell 1964 ranks=(/((iprocmin+irank-1),irank=1,nproc_per_cell)/) 1965 mpi_enreg%comm_cell_pert=xmpi_subcomm(mpi_enreg%comm_world,nproc_per_cell,ranks) 1966 ABI_FREE(ranks) 1967 end if 1968 1969 else !nppert<=1 1970 mpi_enreg%nproc_pert=1 1971 mpi_enreg%comm_pert=xmpi_comm_self 1972 mpi_enreg%me_pert=0 1973 ABI_MALLOC(mpi_enreg%distrb_pert,(npert)) 1974 mpi_enreg%distrb_pert(:)=0 1975 end if 1976 1977 ABI_FREE(nband_rbz) 1978 ABI_FREE(nkpt_rbz) 1979 1980 end subroutine initmpi_pert
ABINIT/initmpi_seq [ Functions ]
NAME
initmpi_seq
FUNCTION
Initializes the MPI information for sequential use.
INPUTS
OUTPUT
mpi_enreg=information about MPI parallelization
SOURCE
934 subroutine initmpi_seq(mpi_enreg) 935 936 !Arguments ------------------------------------ 937 type(MPI_type),intent(out) :: mpi_enreg 938 939 ! *********************************************************************** 940 941 DBG_ENTER("COLL") 942 943 !Set default seq values for scalars 944 mpi_enreg%bandpp=1 945 mpi_enreg%me=0 946 mpi_enreg%me_band=0 947 mpi_enreg%me_cell=0 948 mpi_enreg%me_fft=0 949 mpi_enreg%me_g0=1 950 mpi_enreg%me_g0_fft=1 951 mpi_enreg%me_img=0 952 mpi_enreg%me_hf=0 953 mpi_enreg%me_kpt=0 954 mpi_enreg%me_pert=0 955 mpi_enreg%me_spinor=0 956 mpi_enreg%me_wvl=0 957 mpi_enreg%my_natom=0 ! Should be natom 958 mpi_enreg%my_isppoltab=0 ! Should be (1,0) if nsppol=1 or (1,1) if nsppol=2 959 mpi_enreg%ngfft3_ionic=1 960 mpi_enreg%my_nimage=1 961 mpi_enreg%nproc=1 962 mpi_enreg%nproc_atom=1 963 mpi_enreg%nproc_band=1 964 mpi_enreg%nproc_cell=1 965 mpi_enreg%nproc_fft=1 966 mpi_enreg%nproc_img=1 967 mpi_enreg%nproc_hf=1 968 mpi_enreg%nproc_spkpt=1 969 mpi_enreg%nproc_pert=1 970 mpi_enreg%nproc_spinor=1 971 mpi_enreg%nproc_wvl=1 972 mpi_enreg%paralbd=0 973 mpi_enreg%paral_img=0 974 mpi_enreg%paral_hf=0 975 mpi_enreg%paral_kgb=0 976 mpi_enreg%paral_pert=0 977 mpi_enreg%paral_spinor=0 978 mpi_enreg%pw_unbal_thresh=-1._dp 979 980 !Set default seq values for communicators 981 mpi_enreg%comm_world = xmpi_world 982 mpi_enreg%comm_atom = xmpi_comm_self 983 mpi_enreg%comm_band = xmpi_comm_self 984 mpi_enreg%comm_bandspinor = xmpi_comm_self 985 mpi_enreg%comm_bandfft = xmpi_comm_self 986 mpi_enreg%comm_bandspinorfft = xmpi_comm_self 987 mpi_enreg%comm_cell = xmpi_comm_self 988 mpi_enreg%comm_cell_pert = xmpi_comm_self 989 mpi_enreg%comm_fft = xmpi_comm_self 990 mpi_enreg%comm_hf = xmpi_comm_self 991 mpi_enreg%comm_img = xmpi_comm_self 992 mpi_enreg%comm_kpt = xmpi_comm_self 993 mpi_enreg%comm_kptband = xmpi_comm_self 994 mpi_enreg%comm_pert = xmpi_comm_self 995 mpi_enreg%comm_spinor = xmpi_comm_self 996 mpi_enreg%comm_spinorfft = xmpi_comm_self 997 mpi_enreg%comm_wvl = xmpi_comm_self 998 999 !Nullify all pointers 1000 call nullify_mpi_enreg(mpi_enreg) 1001 1002 !Allocate and nullify distribfft datastructure 1003 ! This is not good since distribfft is not initialized here (even with 0s). 1004 ! It can be dangerous if use with no care (Valgrind might complain) 1005 ABI_MALLOC(mpi_enreg%distribfft,) 1006 1007 DBG_EXIT("COLL") 1008 1009 end subroutine initmpi_seq
ABINIT/initmpi_world [ Functions ]
NAME
initmpi_world
FUNCTION
%comm_world is redifined for the number of processors on which ABINIT is launched
INPUTS
OUTPUT
SOURCE
883 subroutine initmpi_world(mpi_enreg,nproc) 884 885 !Arguments ------------------------------------ 886 integer, intent(in)::nproc 887 type(MPI_type),intent(inout) :: mpi_enreg 888 889 !Local variables------------------------------- 890 !scalars 891 integer :: ii 892 !arrays 893 integer,allocatable :: ranks(:) 894 895 ! *********************************************************************** 896 897 DBG_ENTER("COLL") 898 899 if(nproc==mpi_enreg%nproc) return 900 901 ABI_MALLOC(ranks,(0:nproc-1)) 902 ranks(0:nproc-1)=(/((ii),ii=0,nproc-1)/) 903 mpi_enreg%comm_world=xmpi_subcomm(xmpi_world,nproc,ranks) 904 ABI_FREE(ranks) 905 906 if(mpi_enreg%me<nproc) then 907 mpi_enreg%me=xmpi_comm_rank(mpi_enreg%comm_world) 908 mpi_enreg%nproc=xmpi_comm_size(mpi_enreg%comm_world) 909 call abi_io_redirect(new_io_comm=mpi_enreg%comm_world) 910 call libpaw_write_comm_set(mpi_enreg%comm_world) 911 else 912 mpi_enreg%me=-1 913 end if 914 915 DBG_EXIT("COLL") 916 917 end subroutine initmpi_world
ABINIT/m_mpinfo [ Modules ]
NAME
m_mpinfo
FUNCTION
COPYRIGHT
Copyright (C) 2008-2024 ABINIT group (MT, GG, XG, FJ, AR, MB, CMartins) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt .
TODO
Change the name of the datatype: (MPI_|mpi_) is a reserved keyword and should not be used in client code!
SOURCE
19 #if defined HAVE_CONFIG_H 20 #include "config.h" 21 #endif 22 23 #include "abi_common.h" 24 25 MODULE m_mpinfo 26 27 use defs_basis 28 use m_errors 29 use m_abicore 30 #if defined HAVE_MPI2 31 use mpi 32 #endif 33 use m_xmpi 34 use m_sort 35 use m_distribfft 36 use m_dtset 37 38 use defs_abitypes, only : MPI_type 39 use m_fstrings, only : sjoin, ltoa 40 use m_io_tools, only : file_exists, open_file 41 use m_libpaw_tools, only : libpaw_write_comm_set 42 use m_paral_atom, only : get_my_natom, get_my_atmtab 43 44 implicit none 45 46 private 47 48 #if defined HAVE_MPI1 49 include 'mpif.h' 50 #endif 51 52 public :: init_mpi_enreg ! Initialise a mpi_enreg structure with dataset independent values. 53 public :: nullify_mpi_enreg ! nullify a mpi_enreg datastructure 54 public :: destroy_mpi_enreg ! Free memory 55 public :: copy_mpi_enreg ! Copy a mpi_enreg datastructure into another. 56 public :: set_mpi_enreg_fft ! Set the content of a MPI datastructure in order to call fourwf/fourdp 57 public :: unset_mpi_enreg_fft ! Unset the content of a MPI datastructure used to call fourwf/fourdp 58 public :: ptabs_fourdp ! Return *pointers* to the internal MPI-FFT tables used in fourdp 59 public :: ptabs_fourwf ! Return *pointers* to the internal MPI-FFT tables used in fourwf 60 public :: mpi_distrib_is_ok ! Check if a MPI datastructure contains number of processors 61 ! compatible (in terms of efficiency) with the number of spins/kpts/bands 62 public :: proc_distrb_cycle ! Test a condition to cycle 63 public :: proc_distrb_nband ! Return number of bands present on this cpu 64 public :: proc_distrb_cycle_bands ! Return array of logicals for bands to cycle at this k and spin 65 public :: proc_distrb_band ! Return array of me indices for bands at this k and spin 66 67 public :: initmpi_seq ! Initializes the MPI information for sequential use. 68 public :: initmpi_world ! %comm_world is redifined for the number of processors on which ABINIT is launched 69 70 public :: initmpi_atom ! Initializes the mpi information for parallelism over atoms (PAW). 71 public :: clnmpi_atom ! Cleans-up the mpi information for the parallelism over atoms (PAW). 72 73 public :: initmpi_grid ! Initializes the MPI information for the grid 74 public :: clnmpi_grid ! Cleans-up the mpi information for parallelism over grid (kpt/band/fft). 75 76 public :: initmpi_img ! Initializes the mpi information for parallelism over images of the cell (npimage>1). 77 public :: clnmpi_img ! Cleans-up the mpi information for parallelism over images of the cell (npimage>1). 78 79 public :: initmpi_pert ! Creates group for Parallelization over Perturbations. 80 public :: clnmpi_pert ! Cleans-up the mpi information for parallelization over perturbations. 81 82 public :: initmpi_band ! Initializes the mpi information for band parallelism (paralbd=1). 83 84 ! Helper functions. 85 public :: pre_gather 86 public :: pre_scatter 87 public :: iwrite_fftdatar ! Select the subset of processors that will write density/potential files. 88 89 public :: distrb2 ! Creates the tabs of repartition of processors for sharing the jobs on k-points, spins and bands. 90 public :: distrb2_hf ! Creates the tabs of repartition for Hartree-Fock calculations.
ABINIT/proc_distrb_band [ Functions ]
NAME
proc_distrb_band
FUNCTION
return `rank_band` array with the rank of the processor in comm_band treating `band`
INPUTS
SOURCE
840 subroutine proc_distrb_band(rank_band,distrib,ikpt,isppol,nband,me_band,me_kpt,comm_band) 841 842 !Arguments ------------------------------------ 843 !scalars 844 integer,intent(in) :: nband, ikpt, isppol 845 integer,intent(in) :: me_band,me_kpt,comm_band 846 integer,allocatable,intent(in) :: distrib(:,:,:) 847 integer,intent(out) :: rank_band(nband) 848 849 integer :: ierr, iband 850 851 ! ************************************************************************* 852 853 rank_band = 0 854 855 if (allocated(distrib)) then 856 do iband=1, nband 857 ! is this (k, band, spin) on current proc? 858 if (distrib(ikpt,iband,isppol)/=me_kpt) cycle 859 ! if so save rank in band subcommunicator 860 rank_band(iband) = me_band+1 861 end do 862 call xmpi_sum(rank_band,comm_band,ierr) 863 end if 864 865 rank_band = rank_band-1 866 867 end subroutine proc_distrb_band
ABINIT/proc_distrb_cycle [ Functions ]
NAME
proc_distrb_cycle
FUNCTION
test a condition to cycle over bands and k-points which do not belong to the present processor if return value is true, you can cycle over the given range of bands for this k and sppol
INPUTS
SOURCE
693 function proc_distrb_cycle(distrb,ikpt,iband1,iband2,isppol,me) 694 695 !Arguments ------------------------------------ 696 !scalars 697 integer,intent(in) :: ikpt,iband1,iband2,isppol,me 698 integer,allocatable,intent(in) :: distrb(:,:,:) 699 logical :: proc_distrb_cycle 700 701 ! ************************************************************************* 702 703 proc_distrb_cycle=.false. 704 if (allocated(distrb)) then 705 if (isppol==-1) then 706 ! in this condition, if one of the distrb is for me, then the minval will be == 0, so it returns false 707 proc_distrb_cycle=(minval(abs(distrb(ikpt,iband1:iband2,:)-me))/=0) 708 else 709 proc_distrb_cycle=(minval(abs(distrb(ikpt,iband1:iband2,isppol)-me))/=0) 710 end if 711 end if 712 713 end function proc_distrb_cycle
ABINIT/proc_distrb_cycle_bands [ Functions ]
NAME
proc_distrb_cycle_bands
FUNCTION
return vector of logicals for each band being on present proc
INPUTS
SOURCE
763 subroutine proc_distrb_cycle_bands(cycle_bands,distrb,ikpt,isppol,me) 764 765 !Arguments ------------------------------------ 766 !scalars 767 integer,intent(in) :: ikpt,isppol,me 768 integer,allocatable,intent(in) :: distrb(:,:,:) 769 logical,allocatable,intent(out) :: cycle_bands(:) 770 character(len=500) :: msg 771 772 ! ************************************************************************* 773 774 ABI_REMALLOC (cycle_bands, (size(distrb, 2))) 775 cycle_bands=.false. 776 if (allocated(distrb)) then 777 if (isppol==-1) then 778 ! TODO : should raise error here - the output rank will be all wrong 779 ! could return an OR of the two spin channels, if appropriate 780 cycle_bands=(distrb(ikpt,:,1)/=me) 781 write (msg, "(a)") " for the moment proc_distrb_cycle_bands does not handle the 'any spin' option nsppol -1" 782 ABI_ERROR(msg) 783 else 784 cycle_bands=(distrb(ikpt,:,isppol)/=me) 785 end if 786 end if 787 788 end subroutine proc_distrb_cycle_bands
ABINIT/proc_distrb_kptband [ Functions ]
NAME
proc_distrb_kptband
FUNCTION
return vector of processor indices for each band, within the full kpt communicator
INPUTS
SOURCE
802 subroutine proc_distrb_kptband(kpt_band_procs,distrb,ikpt,isppol) 803 804 !Arguments ------------------------------------ 805 !scalars 806 integer,intent(in) :: ikpt,isppol 807 integer,allocatable,intent(in) :: distrb(:,:,:) 808 integer,allocatable,intent(out) :: kpt_band_procs(:) 809 character(len=500) :: msg 810 811 ! ************************************************************************* 812 813 ABI_REMALLOC(kpt_band_procs, (size(distrb, 2))) 814 kpt_band_procs=-1 815 if (allocated(distrb)) then 816 if (isppol==-1) then 817 ! TODO : should raise error here - the output rank will be all wrong for isppol 2! 818 kpt_band_procs=distrb(ikpt,:,1) 819 write (msg, "(a)") " for the moment proc_distrb_kptband does not handle the 'any spin' option nsppol -1" 820 ABI_ERROR(msg) 821 else 822 kpt_band_procs=distrb(ikpt,:,isppol) 823 end if 824 end if 825 826 end subroutine proc_distrb_kptband
ABINIT/proc_distrb_nband [ Functions ]
NAME
proc_distrb_nband
FUNCTION
return number of bands at this k and spin which are on current proc "me" NB: could replace proc_distrb_cycle with "proc_distrb_nband > 0" isppol -1 means for all spin channels, I think...
INPUTS
SOURCE
729 function proc_distrb_nband(distrb,ikpt,nband_k,isppol,me) 730 731 !Arguments ------------------------------------ 732 !scalars 733 integer,intent(in) :: ikpt,isppol,me,nband_k 734 integer,allocatable,intent(in) :: distrb(:,:,:) 735 integer :: proc_distrb_nband 736 737 ! ************************************************************************* 738 739 proc_distrb_nband=0 740 if (allocated(distrb)) then 741 if (isppol==-1) then 742 !TODO: check this is used correctly : in nsppol 2 case you could end up with 2*nband 743 proc_distrb_nband=(count(distrb(ikpt,1:nband_k,:)==me)) 744 else 745 proc_distrb_nband=(count(distrb(ikpt,1:nband_k,isppol)==me)) 746 end if 747 end if 748 749 end function proc_distrb_nband
m_mpinfo/clnmpi_atom [ Functions ]
[ Top ] [ m_mpinfo ] [ Functions ]
NAME
clnmpi_atom
FUNCTION
Cleans-up the mpi information for the parallelism over atoms (PAW).
SIDE EFFECTS
mpi_enreg=information about MPI parallelization
SOURCE
1134 subroutine clnmpi_atom(mpi_enreg) 1135 1136 !Arguments ------------------------------------ 1137 type(MPI_type), intent(inout) :: mpi_enreg 1138 1139 ! *********************************************************************** 1140 1141 DBG_ENTER("COLL") 1142 1143 if (xmpi_paral==0) return 1144 1145 if (mpi_enreg%comm_atom/=mpi_enreg%comm_world) then 1146 call xmpi_comm_free(mpi_enreg%comm_atom) 1147 mpi_enreg%comm_atom=xmpi_comm_null 1148 end if 1149 1150 if(associated(mpi_enreg%my_atmtab)) then 1151 ABI_FREE(mpi_enreg%my_atmtab) 1152 end if 1153 1154 mpi_enreg%nproc_atom=1 1155 mpi_enreg%my_natom=0 ! should be natom 1156 1157 DBG_EXIT("COLL") 1158 1159 end subroutine clnmpi_atom
m_mpinfo/clnmpi_grid [ Functions ]
[ Top ] [ m_mpinfo ] [ Functions ]
NAME
clnmpi_grid
FUNCTION
Cleans-up the mpi information for parallelism over grid (kpt/band/fft).
SIDE EFFECTS
mpi_enreg=information about MPI parallelization
SOURCE
1455 subroutine clnmpi_grid(mpi_enreg) 1456 1457 !Arguments ------------------------------------ 1458 type(MPI_type), intent(inout) :: mpi_enreg 1459 1460 ! *********************************************************************** 1461 1462 DBG_ENTER("COLL") 1463 1464 if (xmpi_paral==0) return 1465 1466 if (mpi_enreg%comm_bandspinorfft/=mpi_enreg%comm_world) then 1467 call xmpi_comm_free(mpi_enreg%comm_bandspinorfft) 1468 mpi_enreg%comm_bandspinorfft=xmpi_comm_null 1469 end if 1470 1471 if (mpi_enreg%comm_bandfft/=mpi_enreg%comm_world) then 1472 call xmpi_comm_free(mpi_enreg%comm_bandfft) 1473 mpi_enreg%comm_bandfft=xmpi_comm_null 1474 end if 1475 1476 if (mpi_enreg%comm_spinorfft/=mpi_enreg%comm_world) then 1477 call xmpi_comm_free(mpi_enreg%comm_spinorfft) 1478 mpi_enreg%comm_spinorfft=xmpi_comm_null 1479 end if 1480 1481 if (mpi_enreg%comm_bandspinor/=mpi_enreg%comm_world) then 1482 call xmpi_comm_free(mpi_enreg%comm_bandspinor) 1483 mpi_enreg%comm_bandspinor=xmpi_comm_null 1484 end if 1485 1486 if (mpi_enreg%comm_kptband/=mpi_enreg%comm_world) then 1487 call xmpi_comm_free(mpi_enreg%comm_kptband) 1488 mpi_enreg%comm_kptband=xmpi_comm_null 1489 end if 1490 1491 if (mpi_enreg%comm_fft/=mpi_enreg%comm_world) then 1492 call xmpi_comm_free(mpi_enreg%comm_fft) 1493 mpi_enreg%comm_fft=xmpi_comm_null 1494 end if 1495 1496 if (mpi_enreg%comm_band/=mpi_enreg%comm_world) then 1497 call xmpi_comm_free(mpi_enreg%comm_band) 1498 mpi_enreg%comm_band=xmpi_comm_null 1499 end if 1500 1501 if (mpi_enreg%comm_spinor/=mpi_enreg%comm_world) then 1502 call xmpi_comm_free(mpi_enreg%comm_spinor) 1503 mpi_enreg%comm_spinor=xmpi_comm_null 1504 end if 1505 1506 if (mpi_enreg%comm_kpt/=mpi_enreg%comm_world) then 1507 call xmpi_comm_free(mpi_enreg%comm_kpt) 1508 mpi_enreg%comm_kpt=xmpi_comm_null 1509 end if 1510 1511 DBG_EXIT("COLL") 1512 1513 end subroutine clnmpi_grid
m_mpinfo/clnmpi_img [ Functions ]
[ Top ] [ m_mpinfo ] [ Functions ]
NAME
clnmpi_img
FUNCTION
Cleans-up the mpi information for parallelism over images of the cell (npimage>1).
SOURCE
1841 subroutine clnmpi_img(mpi_enreg) 1842 1843 !Arguments ------------------------------------ 1844 type(MPI_type), intent(inout) :: mpi_enreg 1845 1846 ! *********************************************************************** 1847 1848 DBG_ENTER("COLL") 1849 1850 if (xmpi_paral==0) return 1851 1852 if (mpi_enreg%comm_cell/=mpi_enreg%comm_world) then 1853 call xmpi_comm_free(mpi_enreg%comm_cell) 1854 mpi_enreg%comm_cell=xmpi_comm_null 1855 end if 1856 1857 if (mpi_enreg%comm_img/=mpi_enreg%comm_world) then 1858 call xmpi_comm_free(mpi_enreg%comm_img) 1859 mpi_enreg%comm_img=xmpi_comm_null 1860 end if 1861 1862 ABI_SFREE(mpi_enreg%my_imgtab) 1863 ABI_SFREE(mpi_enreg%distrb_img) 1864 1865 mpi_enreg%paral_img=0 1866 mpi_enreg%my_nimage=1 1867 mpi_enreg%me_img=0 1868 mpi_enreg%me_cell=0 1869 mpi_enreg%nproc_img=1 1870 mpi_enreg%nproc_cell=1 1871 1872 DBG_EXIT("COLL") 1873 1874 end subroutine clnmpi_img
m_mpinfo/clnmpi_pert [ Functions ]
[ Top ] [ m_mpinfo ] [ Functions ]
NAME
clnmpi_pert
FUNCTION
Cleans-up the mpi information for parallelization over perturbations.
INPUTS
SOURCE
1996 subroutine clnmpi_pert(mpi_enreg) 1997 1998 !Arguments ------------------------------------ 1999 type(MPI_type),intent(inout) :: mpi_enreg 2000 2001 ! *********************************************************************** 2002 2003 DBG_ENTER("COLL") 2004 2005 if (xmpi_paral==0) return 2006 2007 if(mpi_enreg%paral_pert == 1) then 2008 2009 ! Reset communicators 2010 if (mpi_enreg%comm_pert/=mpi_enreg%comm_world) then 2011 call xmpi_comm_free(mpi_enreg%comm_pert) 2012 mpi_enreg%comm_pert=xmpi_comm_null 2013 end if 2014 2015 ABI_SFREE(mpi_enreg%distrb_pert) 2016 2017 mpi_enreg%me_pert=0 2018 mpi_enreg%me_cell=0 2019 mpi_enreg%nproc_pert=1 2020 mpi_enreg%nproc_cell=1 2021 end if 2022 2023 DBG_EXIT("COLL") 2024 2025 end subroutine clnmpi_pert
m_mpinfo/copy_mpi_enreg [ Functions ]
[ Top ] [ m_mpinfo ] [ Functions ]
NAME
copy_mpi_enreg
FUNCTION
Copy a mpi_enreg datastructure into another
INPUTS
MPI_enreg1<MPI_type>=input mpi_enreg datastructure
OUTPUT
MPI_enreg2<MPI_type>=output mpi_enreg datastructure
SOURCE
238 subroutine copy_mpi_enreg(MPI_enreg1,MPI_enreg2) 239 240 !Arguments ------------------------------------ 241 !scalars 242 type(MPI_type),intent(in) :: mpi_enreg1 243 type(MPI_type),intent(out) :: MPI_enreg2 244 245 !Local variables------------------------------- 246 !scalars 247 integer :: sz1,sz2,sz3 248 249 ! ********************************************************************* 250 251 !scalars 252 mpi_enreg2%comm_world=mpi_enreg1%comm_world 253 mpi_enreg2%me=mpi_enreg1%me 254 mpi_enreg2%nproc=mpi_enreg1%nproc 255 mpi_enreg2%paral_spinor=mpi_enreg1%paral_spinor 256 mpi_enreg2%paralbd=mpi_enreg1%paralbd 257 mpi_enreg2%me_fft=mpi_enreg1%me_fft 258 mpi_enreg2%me_band=mpi_enreg1%me_band 259 mpi_enreg2%nproc_fft=mpi_enreg1%nproc_fft 260 mpi_enreg2%paral_kgb=mpi_enreg1%paral_kgb 261 mpi_enreg2%me_g0=mpi_enreg1%me_g0 262 mpi_enreg2%me_g0_fft=mpi_enreg1%me_g0_fft 263 mpi_enreg2%paral_pert=mpi_enreg1%paral_pert 264 mpi_enreg2%me_pert=mpi_enreg1%me_pert 265 mpi_enreg2%nproc_pert=mpi_enreg1%nproc_pert 266 mpi_enreg2%comm_pert=mpi_enreg1%comm_pert 267 mpi_enreg2%comm_bandfft=mpi_enreg1%comm_bandfft 268 mpi_enreg2%comm_band=mpi_enreg1%comm_band 269 mpi_enreg2%comm_fft=mpi_enreg1%comm_fft 270 mpi_enreg2%nproc_band=mpi_enreg1%nproc_band 271 mpi_enreg2%comm_bandspinorfft=mpi_enreg1%comm_bandspinorfft 272 mpi_enreg2%comm_kpt=mpi_enreg1%comm_kpt 273 mpi_enreg2%me_kpt=mpi_enreg1%me_kpt 274 mpi_enreg2%nproc_spkpt=mpi_enreg1%nproc_spkpt 275 mpi_enreg2%my_isppoltab=mpi_enreg1%my_isppoltab 276 mpi_enreg2%my_natom=mpi_enreg1%my_natom 277 mpi_enreg2%comm_atom=mpi_enreg1%comm_atom 278 mpi_enreg2%nproc_atom=mpi_enreg1%nproc_atom 279 mpi_enreg2%comm_kptband=mpi_enreg1%comm_kptband 280 mpi_enreg2%bandpp=mpi_enreg1%bandpp 281 mpi_enreg2%paral_img=mpi_enreg1%paral_img 282 mpi_enreg2%comm_img=mpi_enreg1%comm_img 283 mpi_enreg2%me_img=mpi_enreg1%me_img 284 mpi_enreg2%nproc_img=mpi_enreg1%nproc_img 285 mpi_enreg2%comm_cell=mpi_enreg1%comm_cell 286 mpi_enreg2%comm_cell_pert=mpi_enreg1%comm_cell_pert 287 mpi_enreg2%me_cell=mpi_enreg1%me_cell 288 mpi_enreg2%nproc_cell=mpi_enreg1%nproc_cell 289 mpi_enreg2%nproc_spinor=mpi_enreg1%nproc_spinor 290 mpi_enreg2%me_spinor=mpi_enreg1%me_spinor 291 mpi_enreg2%comm_spinorfft=mpi_enreg1%comm_spinorfft 292 mpi_enreg2%me_wvl =mpi_enreg1%me_wvl 293 mpi_enreg2%nproc_wvl =mpi_enreg1%nproc_wvl 294 mpi_enreg2%comm_wvl =mpi_enreg1%comm_wvl 295 mpi_enreg2%me_hf =mpi_enreg1%me_hf 296 mpi_enreg2%nproc_hf =mpi_enreg1%nproc_hf 297 mpi_enreg2%comm_hf =mpi_enreg1%comm_hf 298 mpi_enreg2%paral_hf=mpi_enreg1%paral_hf 299 mpi_enreg2%pw_unbal_thresh=mpi_enreg1%pw_unbal_thresh 300 301 !pointers 302 if (associated(mpi_enreg1%distribfft)) then 303 if (.not.associated(mpi_enreg2%distribfft)) then 304 ABI_MALLOC(mpi_enreg2%distribfft,) 305 end if 306 call copy_distribfft(mpi_enreg1%distribfft,mpi_enreg2%distribfft) 307 end if 308 309 if (allocated(mpi_enreg1%proc_distrb)) then 310 sz1=size(mpi_enreg1%proc_distrb,1) 311 sz2=size(mpi_enreg1%proc_distrb,2) 312 sz3=size(mpi_enreg1%proc_distrb,3) 313 ABI_MALLOC(mpi_enreg2%proc_distrb,(sz1,sz2,sz3)) 314 mpi_enreg2%proc_distrb=mpi_enreg1%proc_distrb 315 end if 316 if (allocated(mpi_enreg1%kptdstrb)) then 317 sz1=size(mpi_enreg1%kptdstrb,1) 318 sz2=size(mpi_enreg1%kptdstrb,2) 319 sz3=size(mpi_enreg1%kptdstrb,3) 320 ABI_MALLOC(mpi_enreg2%kptdstrb,(sz1,sz2,sz3)) 321 mpi_enreg2%kptdstrb=mpi_enreg1%kptdstrb 322 end if 323 if (allocated(mpi_enreg1%kpt_loc2fbz_sp)) then 324 sz1=size(mpi_enreg1%kpt_loc2fbz_sp,1)-1 325 sz2=size(mpi_enreg1%kpt_loc2fbz_sp,2) 326 sz3=size(mpi_enreg1%kpt_loc2fbz_sp,3) 327 ABI_MALLOC(mpi_enreg2%kpt_loc2fbz_sp,(0:sz1,1:sz2,1:sz3)) 328 mpi_enreg2%kpt_loc2fbz_sp=mpi_enreg1%kpt_loc2fbz_sp 329 end if 330 if (allocated(mpi_enreg1%kpt_loc2ibz_sp)) then 331 sz1=size(mpi_enreg1%kpt_loc2ibz_sp,1)-1 332 sz2=size(mpi_enreg1%kpt_loc2ibz_sp,2) 333 sz3=size(mpi_enreg1%kpt_loc2ibz_sp,3) 334 ABI_MALLOC(mpi_enreg2%kpt_loc2ibz_sp,(0:sz1,1:sz2,1:sz3)) 335 mpi_enreg2%kpt_loc2ibz_sp=mpi_enreg1%kpt_loc2ibz_sp 336 end if 337 if (allocated(mpi_enreg1%mkmem)) then 338 ABI_MALLOC(mpi_enreg2%mkmem,(0:size(mpi_enreg1%mkmem,1)-1)) 339 mpi_enreg2%mkmem=mpi_enreg1%mkmem 340 end if 341 if (associated(mpi_enreg1%my_atmtab)) then 342 ABI_MALLOC(mpi_enreg2%my_atmtab,(size(mpi_enreg1%my_atmtab))) 343 mpi_enreg2%my_atmtab=mpi_enreg1%my_atmtab 344 else 345 nullify(mpi_enreg2%my_atmtab) 346 end if 347 if (allocated(mpi_enreg1%my_kgtab)) then 348 sz1=size(mpi_enreg1%my_kgtab,1) 349 sz2=size(mpi_enreg1%my_kgtab,2) 350 ABI_MALLOC(mpi_enreg2%my_kgtab,(sz1,sz2)) 351 mpi_enreg2%my_kgtab=mpi_enreg1%my_kgtab 352 end if 353 if (allocated(mpi_enreg1%distrb_pert)) then 354 ABI_MALLOC(mpi_enreg2%distrb_pert,(size(mpi_enreg1%distrb_pert))) 355 mpi_enreg2%distrb_pert=mpi_enreg1%distrb_pert 356 end if 357 if (allocated(mpi_enreg1%distrb_img)) then 358 ABI_MALLOC(mpi_enreg2%distrb_img,(size(mpi_enreg1%distrb_img))) 359 mpi_enreg2%distrb_img=mpi_enreg1%distrb_img 360 end if 361 if (allocated(mpi_enreg1%my_imgtab)) then 362 ABI_MALLOC(mpi_enreg2%my_imgtab,(size(mpi_enreg1%my_imgtab))) 363 mpi_enreg2%my_imgtab=mpi_enreg1%my_imgtab 364 end if 365 if (allocated(mpi_enreg1%distrb_hf)) then 366 sz1=size(mpi_enreg1%distrb_hf,1) 367 sz2=size(mpi_enreg1%distrb_hf,2) 368 sz3=size(mpi_enreg1%distrb_hf,3) 369 ABI_MALLOC(mpi_enreg2%distrb_hf,(sz1,sz2,sz3)) 370 mpi_enreg2%distrb_hf=mpi_enreg1%distrb_hf 371 end if 372 373 !Optional pointers 374 if (allocated(mpi_enreg1%my_kpttab)) then 375 ABI_MALLOC(mpi_enreg2%my_kpttab,(size(mpi_enreg1%my_kpttab))) 376 mpi_enreg2%my_kpttab=mpi_enreg1%my_kpttab 377 end if 378 379 !Do not copy wavelet pointers, just associate. 380 mpi_enreg2%nscatterarr => mpi_enreg1%nscatterarr 381 mpi_enreg2%ngatherarr => mpi_enreg1%ngatherarr 382 383 end subroutine copy_mpi_enreg
m_mpinfo/destroy_mpi_enreg [ Functions ]
[ Top ] [ m_mpinfo ] [ Functions ]
NAME
destroy_mpi_enreg
FUNCTION
Destroy a mpi_enreg datastructure
SIDE EFFECTS
MPI_enreg<MPI_type>=Datatype gathering information on the parallelism.
SOURCE
186 subroutine destroy_mpi_enreg(MPI_enreg) 187 188 !Arguments ------------------------------------ 189 !scalars 190 type(MPI_type),intent(inout) :: MPI_enreg 191 192 ! ********************************************************************* 193 194 if (associated(mpi_enreg%distribfft)) then 195 call destroy_distribfft(mpi_enreg%distribfft) 196 ABI_FREE(mpi_enreg%distribfft) 197 nullify(mpi_enreg%distribfft) 198 end if 199 200 ABI_SFREE(mpi_enreg%proc_distrb) 201 ABI_SFREE(mpi_enreg%kptdstrb) 202 ABI_SFREE(mpi_enreg%kpt_loc2fbz_sp) 203 ABI_SFREE(mpi_enreg%kpt_loc2ibz_sp) 204 ABI_SFREE(mpi_enreg%mkmem) 205 ABI_SFREE(mpi_enreg%my_kpttab) 206 if (associated(mpi_enreg%my_atmtab)) then 207 ABI_FREE(mpi_enreg%my_atmtab) 208 nullify(mpi_enreg%my_atmtab) 209 end if 210 ABI_SFREE(mpi_enreg%distrb_pert) 211 ABI_SFREE(mpi_enreg%distrb_img) 212 ABI_SFREE(mpi_enreg%my_imgtab) 213 ABI_SFREE(mpi_enreg%my_kgtab) 214 ABI_SFREE(mpi_enreg%distrb_hf) 215 216 !Do not deallocate wavelet denspot distribution arrays, they are handled by BigDFT. 217 218 end subroutine destroy_mpi_enreg
m_mpinfo/init_mpi_enreg [ Functions ]
[ Top ] [ m_mpinfo ] [ Functions ]
NAME
init_mpi_enreg
FUNCTION
Initialise a mpi_enreg structure with dataset independent values. Other values of mpi_enreg are dataset dependent, and should NOT be initialized inside abinit.F90 . XG 071118 : At present several other values are initialized temporarily inside invars1.F90, FROM THE DTSET VALUES. In order to releave the present constraint of having mpi_enreg equal for all datasets, they should be reinitialized from the dtset values inside invars2m.F90 (where there is a loop over datasets, and finally, reinitialized from the dataset values inside each big routine called by driver, according to the kind of parallelisation that is needed there. One should have one init_mpi_dtset routine (or another name) per big routine (well, there is also the problem of TDDFT ...). Also, one should have a clean_mpi_dtset called at the end of each big routine, as well as invars1.F90 or invars2m.F90 .
INPUTS
SIDE EFFECTS
MPI_enreg<MPI_type>=All pointer set to null().
SOURCE
122 subroutine init_mpi_enreg(mpi_enreg) 123 124 !Arguments ------------------------------------ 125 !scalars 126 type(MPI_type),intent(inout) :: MPI_enreg 127 128 ! ********************************************************************* 129 130 !Default for sequential use 131 call initmpi_seq(mpi_enreg) 132 !Initialize MPI 133 #if defined HAVE_MPI 134 mpi_enreg%comm_world=xmpi_world 135 mpi_enreg%me = xmpi_comm_rank(xmpi_world) 136 mpi_enreg%nproc = xmpi_comm_size(xmpi_world) 137 #endif 138 139 end subroutine init_mpi_enreg
m_mpinfo/iwrite_fftdatar [ Functions ]
[ Top ] [ m_mpinfo ] [ Functions ]
NAME
iwrite_fftdatar
FUNCTION
This function selects the subset of processors that should write density/potential Return True if the processors should do IO.
INPUTS
mpi_enreg<MPI_type>=Datatype gathering information on the parallelism
SOURCE
2247 logical function iwrite_fftdatar(mpi_enreg) result(ans) 2248 2249 !Arguments ------------------------------------ 2250 !scalars 2251 type(MPI_type),intent(in) :: mpi_enreg 2252 2253 ! ********************************************************************* 2254 2255 ans = (xmpi_paral==0 .or. & ! No MPI 2256 (mpi_enreg%paral_kgb==0 .and. mpi_enreg%me==0) .or. & ! paral_kgb=0 does not use MPI-FFT and cartesian communicators. 2257 (mpi_enreg%paral_kgb==1 .and. mpi_enreg%me_band==0 .and. & ! select procs in one FFT communicator. 2258 mpi_enreg%me_kpt==0 .and. mpi_enreg%me_spinor==0) .or. & 2259 (mpi_enreg%paral_pert==1 .and. mpi_enreg%me_cell==0) & ! Group master in perturbation communicator. 2260 ) 2261 2262 end function iwrite_fftdatar
m_mpinfo/mpi_distrib_is_ok [ Functions ]
[ Top ] [ m_mpinfo ] [ Functions ]
NAME
mpi_distrib_is_ok
FUNCTION
Check if a MPI datastructure contains number of processors compatible (in terms of efficiency) with the number of spins/k-points/bands
INPUTS
MPI_enreg<MPI_type>=Datatype gathering information on the parallelism nband=number of bands nkpt=number of k-points nptk_current_proc=number of k-points handled by current MPI process nsppol= number of spins (1 or 2)
OUTPUT
mpi_distrib_is_ok (current function)=TRUE if the current MPI distribution is optimal FALSE otherwise [msg]= -optional- warning message to be printed out
SOURCE
641 logical function mpi_distrib_is_ok(MPI_enreg,nband,nkpt,nkpt_current_proc,nsppol,msg) 642 643 !Arguments ------------------------------------ 644 !scalars 645 integer,intent(in) :: nband,nkpt,nkpt_current_proc,nsppol 646 type(MPI_type),intent(in) :: MPI_enreg 647 character(len=*),optional,intent(out) :: msg 648 649 ! ********************************************************************* 650 651 mpi_distrib_is_ok=.true. 652 653 if (MPI_enreg%paralbd==0) then 654 if (MPI_enreg%nproc_spkpt-floor(nsppol*nkpt*one/nkpt_current_proc)>=nkpt_current_proc) then 655 mpi_distrib_is_ok=.false. 656 if (present(msg)) then 657 write(msg,'(a,i0,4a,i0,3a)') & 658 'Your number of spins*k-points (=',nsppol*nkpt,') ',& 659 'will not distribute correctly',ch10, & 660 'with the current number of processors (=',MPI_enreg%nproc_spkpt,').',ch10,& 661 'You will leave some empty.' 662 end if 663 end if 664 else 665 if (mod(nband,max(1,MPI_enreg%nproc_spkpt/(nsppol*nkpt)))/=0) then 666 mpi_distrib_is_ok=.false. 667 if (present(msg)) then 668 write(msg,'(a,i0,2a,i0,4a,i0,7a)')& 669 'Your number of spins*k-points (=',nsppol*nkpt,') ',& 670 'and bands (=',nband,') ',& 671 'will not distribute correctly',ch10,& 672 'with the current number of processors (=',MPI_enreg%nproc_spkpt,').',ch10,& 673 'You will leave some empty.' 674 end if 675 end if 676 end if 677 678 end function mpi_distrib_is_ok
m_mpinfo/nullify_mpi_enreg [ Functions ]
[ Top ] [ m_mpinfo ] [ Functions ]
NAME
nullify_mpi_enreg
FUNCTION
nullify a mpi_enreg datastructure
SIDE EFFECTS
MPI_enreg<MPI_type>=All pointer set to null().
SOURCE
156 subroutine nullify_mpi_enreg(MPI_enreg) 157 158 !Arguments ------------------------------------ 159 !scalars 160 type(MPI_type),intent(inout) :: MPI_enreg 161 162 ! ********************************************************************* 163 164 nullify(mpi_enreg%nscatterarr) 165 nullify(mpi_enreg%ngatherarr) 166 nullify(mpi_enreg%my_atmtab) 167 nullify(mpi_enreg%distribfft) 168 169 end subroutine nullify_mpi_enreg
m_mpinfo/pre_gather [ Functions ]
[ Top ] [ m_mpinfo ] [ Functions ]
NAME
pre_gather
FUNCTION
Gathers data from FFT processors.
INPUTS
n1,n2,n3= FFT grid dimensions n4= n3/mpi_enreg%nproc_fft array= data to gather among procs
OUTPUT
None
SIDE EFFECTS
array_allgather= gathered data
SOURCE
2178 subroutine pre_gather(array,array_allgather,n1,n2,n3,n4,mpi_enreg) 2179 2180 !Arguments ------------------------------------ 2181 integer,intent(in) :: n1,n2,n3,n4 2182 real(dp),intent(in) :: array(n1,n2,n4,1) 2183 real(dp),intent(inout) :: array_allgather(n1,n2,n3,1) 2184 type(mpi_type),intent(in) :: mpi_enreg 2185 2186 !Local variables------------------------------- 2187 integer :: ier 2188 2189 ! ********************************************************************* 2190 2191 !Gather the array on all procs 2192 call xmpi_allgather(array,n1*n2*n3/mpi_enreg%nproc_fft,array_allgather,mpi_enreg%comm_fft,ier) 2193 2194 end subroutine pre_gather
m_mpinfo/pre_scatter [ Functions ]
[ Top ] [ m_mpinfo ] [ Functions ]
NAME
pre_scatter
FUNCTION
Scatters data among FFT processors.
INPUTS
n1,n2,n3= FFT grid dimensions n4= n3/mpi_enreg%nproc_fft array_allgather= data to scatter among FFT procs
OUTPUT
array= scattered data
SOURCE
2217 subroutine pre_scatter(array,array_allgather,n1,n2,n3,n4,mpi_enreg) 2218 2219 !Arguments ------------------------------------ 2220 integer,intent(in) :: n1,n2,n3,n4 2221 real(dp),intent(out) :: array(n1,n2,n4,1) 2222 real(dp),intent(in) :: array_allgather(n1,n2,n3,1) 2223 type(mpi_type),intent(in) :: mpi_enreg 2224 2225 ! ********************************************************************* 2226 2227 !Perform the reverse operation 2228 array(:,:,:,:) = & 2229 & array_allgather(:,:,n3/mpi_enreg%nproc_fft*mpi_enreg%me_fft+1:n3/mpi_enreg%nproc_fft*(mpi_enreg%me_fft+1),:) 2230 2231 end subroutine pre_scatter
m_mpinfo/ptabs_fourdp [ Functions ]
[ Top ] [ m_mpinfo ] [ Functions ]
NAME
ptabs_fourdp
FUNCTION
Returns pointers to the tables used for the MPI FFT of densities and potentials (fourdp routine).
NOTES
1) These pointers are references to the internal tables stored in MPI_enreg hence *** DO NOT DEALLOCATE THE POINTERS YOU HAVE RECEIVED! *** 2) Client code should declare the pointers with the attribute ABI_CONTIGUOUS (this macro expands to F2008 CONTIGUOUS if the compiler supports it)
INPUTS
MPI_enreg<MPI_type>=Datatype gathering information on the parallelism. n2,n3=Number of FFT divisions along y and z
OUTPUT
fftn2_distrib(:)= rank of the processor which own fft planes in 2nd dimension for fourdp ffti2_local(:) = local i2 indices in fourdp fftn3_distrib(:) = rank of the processor which own fft planes in 3rd dimension for fourdp ffti3_local(:) = local i3 indices in fourdp
SOURCE
498 subroutine ptabs_fourdp(MPI_enreg,n2,n3,fftn2_distrib,ffti2_local,fftn3_distrib,ffti3_local) 499 500 !Arguments ------------------------------------ 501 !scalars 502 integer,intent(in) :: n2,n3 503 type(MPI_type),intent(in) :: MPI_enreg 504 integer, ABI_CONTIGUOUS pointer :: fftn2_distrib(:),ffti2_local(:) 505 integer, ABI_CONTIGUOUS pointer :: fftn3_distrib(:),ffti3_local(:) 506 507 !Local variables------------------------------- 508 !scalars 509 logical :: grid_found 510 511 ! ********************************************************************* 512 513 grid_found=.false. 514 515 ! Get the distrib associated with this fft_grid => for i2 and i3 planes 516 grid_found=.false. 517 if (n2== mpi_enreg%distribfft%n2_coarse) then 518 if( n3 == size(mpi_enreg%distribfft%tab_fftdp3_distrib) )then 519 fftn2_distrib => mpi_enreg%distribfft%tab_fftdp2_distrib 520 ffti2_local => mpi_enreg%distribfft%tab_fftdp2_local 521 fftn3_distrib => mpi_enreg%distribfft%tab_fftdp3_distrib 522 ffti3_local => mpi_enreg%distribfft%tab_fftdp3_local 523 grid_found=.true. 524 end if 525 end if 526 527 if((n2 == mpi_enreg%distribfft%n2_fine).and.(.not.(grid_found))) then 528 if( n3 == size(mpi_enreg%distribfft%tab_fftdp3dg_distrib) )then 529 fftn2_distrib => mpi_enreg%distribfft%tab_fftdp2dg_distrib 530 ffti2_local => mpi_enreg%distribfft%tab_fftdp2dg_local 531 fftn3_distrib => mpi_enreg%distribfft%tab_fftdp3dg_distrib 532 ffti3_local => mpi_enreg%distribfft%tab_fftdp3dg_local 533 grid_found=.true. 534 end if 535 end if 536 537 if (.not.grid_found) then 538 ABI_BUG(sjoin("Unable to find an allocated distrib for this fft grid with n2, n3 = ", ltoa([n2, n3]))) 539 end if 540 541 end subroutine ptabs_fourdp
m_mpinfo/ptabs_fourwf [ Functions ]
[ Top ] [ m_mpinfo ] [ Functions ]
NAME
ptabs_fourwf
FUNCTION
Returns pointers to the tables used for the MPI FFT of the wavefunctions (fourwf routine).
NOTES
1) These pointers are references to the internal tables stored in MPI_enreg hence *** DO NOT DEALLOCATE THE POINTERS YOU HAVE RECEIVED! *** 2) Client code should declare the pointers with the attribute ABI_CONTIGUOUS (this macro expands to F2008 CONTIGUOUS if the compiler supports it)
INPUTS
MPI_enreg<MPI_type>=Datatype gathering information on the parallelism. n2,n3=Number of FFT divisions along y and z
OUTPUT
fftn2_distrib(:)= rank of the processors which own fft planes in 2nd dimension for fourwf ffti2_local(:) = local i2 indices in fourwf fftn3_distrib(:) = rank of the processors which own fft planes in 3rd dimension for fourwf ffti3_local(:) = local i3 indices in fourwf
SOURCE
572 subroutine ptabs_fourwf(MPI_enreg,n2,n3,fftn2_distrib,ffti2_local,fftn3_distrib,ffti3_local) 573 574 !Arguments ------------------------------------ 575 !scalars 576 integer,intent(in) :: n2,n3 577 type(MPI_type),intent(in) :: MPI_enreg 578 integer, ABI_CONTIGUOUS pointer :: fftn2_distrib(:),ffti2_local(:) 579 integer, ABI_CONTIGUOUS pointer :: fftn3_distrib(:),ffti3_local(:) 580 581 !Local variables------------------------------- 582 !scalars 583 logical :: grid_found 584 585 ! ********************************************************************* 586 587 grid_found=.false. 588 589 ! Get the distrib associated with this fft_grid => for i2 and i3 planes 590 if (n2 == mpi_enreg%distribfft%n2_coarse) then 591 if (n3 == size(mpi_enreg%distribfft%tab_fftdp3_distrib))then 592 fftn2_distrib => mpi_enreg%distribfft%tab_fftwf2_distrib 593 ffti2_local => mpi_enreg%distribfft%tab_fftwf2_local 594 fftn3_distrib => mpi_enreg%distribfft%tab_fftdp3_distrib 595 ffti3_local => mpi_enreg%distribfft%tab_fftdp3_local 596 grid_found=.true. 597 end if 598 end if 599 600 if((n2 == mpi_enreg%distribfft%n2_fine).and.(.not.(grid_found))) then 601 if (n3 == size(mpi_enreg%distribfft%tab_fftdp3dg_distrib) )then 602 fftn2_distrib => mpi_enreg%distribfft%tab_fftwf2dg_distrib 603 ffti2_local => mpi_enreg%distribfft%tab_fftwf2dg_local 604 fftn3_distrib => mpi_enreg%distribfft%tab_fftdp3dg_distrib 605 ffti3_local => mpi_enreg%distribfft%tab_fftdp3dg_local 606 grid_found=.true. 607 end if 608 end if 609 610 if(.not. grid_found) then 611 ABI_BUG(sjoin("Unable to find an allocated distrib for this fft grid", ltoa([n2, n3]))) 612 end if 613 614 end subroutine ptabs_fourwf
m_mpinfo/set_mpi_enreg_fft [ Functions ]
[ Top ] [ m_mpinfo ] [ Functions ]
NAME
set_mpi_enreg_fft
FUNCTION
Set the content of a MPI datastructure in order to call fourwf/fourdp (in view of a wrapper for these routines)
INPUTS
me_g0=1 if the current process treat the g=0 plane-wave comm_fft= MPI communicator over FFT components paral_kgb= flag used to activate "band-FFT" parallelism
SIDE EFFECTS
MPI_enreg<MPI_type>=FFT pointer/flags intialized
SOURCE
406 subroutine set_mpi_enreg_fft(MPI_enreg,comm_fft,distribfft,me_g0,paral_kgb) 407 408 !Arguments ------------------------------------ 409 !scalars 410 integer,intent(in) :: me_g0,comm_fft,paral_kgb 411 type(distribfft_type),intent(in),target :: distribfft 412 type(MPI_type),intent(inout) :: MPI_enreg 413 414 ! ********************************************************************* 415 416 mpi_enreg%comm_fft=comm_fft 417 mpi_enreg%paral_kgb=paral_kgb 418 mpi_enreg%me_g0=me_g0 419 mpi_enreg%nproc_fft=xmpi_comm_size(comm_fft) 420 mpi_enreg%me_fft=xmpi_comm_rank(comm_fft) 421 mpi_enreg%me_g0_fft=0 422 if (mpi_enreg%me_fft==0) then 423 mpi_enreg%me_g0_fft=1 424 end if 425 if (associated(mpi_enreg%distribfft)) then 426 call destroy_distribfft(mpi_enreg%distribfft) 427 ABI_FREE(mpi_enreg%distribfft) 428 end if 429 mpi_enreg%distribfft => distribfft 430 431 end subroutine set_mpi_enreg_fft
m_mpinfo/unset_mpi_enreg_fft [ Functions ]
[ Top ] [ m_mpinfo ] [ Functions ]
NAME
unset_mpi_enreg_fft
FUNCTION
Unset the content of a MPI datastructure used to call fourwf/fourdp (in view of a wrapper for these routines)
INPUTS
SIDE EFFECTS
MPI_enreg<MPI_type>=FFT pointer/flags intialized
SOURCE
451 subroutine unset_mpi_enreg_fft(MPI_enreg) 452 453 !Arguments ------------------------------------ 454 !scalars 455 type(MPI_type),intent(inout) :: MPI_enreg 456 457 ! ********************************************************************* 458 459 mpi_enreg%me_g0=1 460 mpi_enreg%comm_fft=xmpi_comm_self 461 mpi_enreg%nproc_fft=1 462 mpi_enreg%me_fft=0 463 mpi_enreg%me_g0_fft=1 464 mpi_enreg%paral_kgb=0 465 nullify(mpi_enreg%distribfft) 466 467 end subroutine unset_mpi_enreg_fft