1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146
|
!
! Copyright (C) 2006-2010 Quantum ESPRESSO group
! This file is distributed under the terms of the
! GNU General Public License. See the file `License'
! in the root directory of the present distribution,
! or http://www.gnu.org/copyleft/gpl.txt .
!
!
!
!----------------------------------------------------------------------
! FFT base Module.
! Written by Carlo Cavazzoni
!----------------------------------------------------------------------
!
!=----------------------------------------------------------------------=!
MODULE fft_base
!=----------------------------------------------------------------------=!
USE kinds, ONLY: DP
USE parallel_include
USE fft_types, ONLY: fft_dlay_descriptor
IMPLICIT NONE
! ... data structure containing all information
! ... about fft data distribution for a given
! ... potential grid, and its wave functions sub-grid.
TYPE ( fft_dlay_descriptor ) :: dfftp ! descriptor for dense grid
! Dimensions of the 3D real and reciprocal space FFT grid
! relative to the charge density and potential ("dense" grid)
TYPE ( fft_dlay_descriptor ) :: dffts ! descriptor for smooth grid
! Dimensions of the 3D real and reciprocal space
! FFT grid relative to the smooth part of the charge density
! (may differ from the full charge density grid for USPP )
TYPE ( fft_dlay_descriptor ) :: dfftb ! descriptor for box grids
! Dimensions of the 3D real and reciprocal space
! FFT grid relative to the "small box" computation
! of the atomic augmentation part of the
! charge density used in USPP (to speed up CPV iterations)
SAVE
PRIVATE
PUBLIC :: fft_scatter, grid_gather, grid_scatter
PUBLIC :: dfftp, dffts, dfftb, fft_dlay_descriptor
PUBLIC :: cgather_sym, cgather_smooth, cgather_custom
PUBLIC :: cscatter_sym, cscatter_smooth, cscatter_custom
PUBLIC :: gather_smooth, scatter_smooth
PUBLIC :: tg_gather, tg_cgather
!=----------------------------------------------------------------------=!
CONTAINS
!=----------------------------------------------------------------------=!
!
!
!
! ALLTOALL based SCATTER, should be better on network
! with a defined topology, like on bluegene and cray machine
!
!-----------------------------------------------------------------------
SUBROUTINE fft_scatter ( dfft, f_in, nr3x, nxx_, f_aux, ncp_, npp_, isgn, use_tg )
!-----------------------------------------------------------------------
!
! transpose the fft grid across nodes
! a) From columns to planes (isgn > 0)
!
! "columns" (or "pencil") representation:
! processor "me" has ncp_(me) contiguous columns along z
! Each column has nr3x elements for a fft of order nr3
! nr3x can be =nr3+1 in order to reduce memory conflicts.
!
! The transpose take places in two steps:
! 1) on each processor the columns are divided into slices along z
! that are stored contiguously. On processor "me", slices for
! processor "proc" are npp_(proc)*ncp_(me) big
! 2) all processors communicate to exchange slices
! (all columns with z in the slice belonging to "me"
! must be received, all the others must be sent to "proc")
! Finally one gets the "planes" representation:
! processor "me" has npp_(me) complete xy planes
!
! b) From planes to columns (isgn < 0)
!
! Quite the same in the opposite direction
!
! The output is overwritten on f_in ; f_aux is used as work space
!
! If optional argument "use_tg" is true the subroutines performs
! the trasposition using the Task Groups distribution
!
#ifdef __MPI
USE parallel_include
#endif
USE kinds, ONLY : DP
IMPLICIT NONE
TYPE (fft_dlay_descriptor), INTENT(in) :: dfft
INTEGER, INTENT(in) :: nr3x, nxx_, isgn, ncp_ (:), npp_ (:)
COMPLEX (DP), INTENT(inout) :: f_in (nxx_), f_aux (nxx_)
LOGICAL, OPTIONAL, INTENT(in) :: use_tg
#ifdef __MPI
INTEGER :: dest, from, k, offset, proc, ierr, me, nprocp, gproc, gcomm, i, kdest, kfrom
INTEGER :: me_p, nppx, mc, j, npp, nnp, ii, it, ip, ioff, sendsiz, ncpx
!
LOGICAL :: use_tg_
#if defined __HPM
! CALL f_hpmstart( 10, 'scatter' )
#endif
!
! Task Groups
use_tg_ = .false.
IF( present( use_tg ) ) use_tg_ = use_tg
me = dfft%mype + 1
!
IF( use_tg_ ) THEN
! This is the number of procs. in the plane-wave group
nprocp = dfft%npgrp
ELSE
nprocp = dfft%nproc
ENDIF
!
CALL start_clock ('fft_scatter')
!
ncpx = 0
nppx = 0
IF( use_tg_ ) THEN
DO proc = 1, nprocp
gproc = dfft%nplist( proc ) + 1
ncpx = max( ncpx, ncp_ ( gproc ) )
nppx = max( nppx, npp_ ( gproc ) )
ENDDO
ELSE
DO proc = 1, nprocp
ncpx = max( ncpx, ncp_ ( proc ) )
nppx = max( nppx, npp_ ( proc ) )
ENDDO
IF ( dfft%nproc == 1 ) THEN
nppx = dfft%nr3x
END IF
ENDIF
sendsiz = ncpx * nppx
!
ierr = 0
IF (isgn.gt.0) THEN
IF (nprocp==1) GO TO 10
!
! "forward" scatter from columns to planes
!
! step one: store contiguously the slices
!
offset = 1
DO proc = 1, nprocp
from = offset
IF( use_tg_ ) THEN
gproc = dfft%nplist(proc)+1
ELSE
gproc = proc
ENDIF
dest = 1 + ( proc - 1 ) * sendsiz
!
DO k = 1, ncp_ (me)
kdest = dest + (k - 1) * nppx - 1
kfrom = from + (k - 1) * nr3x - 1
DO i = 1, npp_ ( gproc )
f_aux ( kdest + i ) = f_in ( kfrom + i )
ENDDO
ENDDO
offset = offset + npp_ ( gproc )
ENDDO
!
! maybe useless; ensures that no garbage is present in the output
!
f_in = 0.0_DP
!
! step two: communication
!
IF( use_tg_ ) THEN
gcomm = dfft%pgrp_comm
ELSE
gcomm = dfft%comm
ENDIF
! CALL mpi_barrier (gcomm, ierr) ! why barrier? for buggy openmpi over ib
CALL mpi_alltoall (f_aux(1), sendsiz, MPI_DOUBLE_COMPLEX, f_in(1), sendsiz, MPI_DOUBLE_COMPLEX, gcomm, ierr)
IF( abs(ierr) /= 0 ) CALL errore ('fft_scatter', 'info<>0', abs(ierr) )
!
10 CONTINUE
!
f_aux = (0.d0, 0.d0)
!
IF( isgn == 1 ) THEN
!!$omp parallel default(none) private(ip,ioff,i,mc,it,j) shared(dfft,nppx,sendsiz,me,f_in,f_aux)
!!$omp do
DO ip = 1, dfft%nproc
ioff = dfft%iss( ip )
DO i = 1, dfft%nsp( ip )
mc = dfft%ismap( i + ioff )
it = ( i - 1 ) * nppx + ( ip - 1 ) * sendsiz
DO j = 1, dfft%npp( me )
f_aux( mc + ( j - 1 ) * dfft%nnp ) = f_in( j + it )
ENDDO
ENDDO
ENDDO
!!$omp end do
!!$omp end parallel
ELSE
IF( use_tg_ ) THEN
npp = dfft%tg_npp( me )
nnp = dfft%nr1x * dfft%nr2x
ELSE
npp = dfft%npp( me )
nnp = dfft%nnp
ENDIF
!
!!$omp parallel default(none) private(ip,ioff,i,mc,it,j,gproc,ii) shared(dfft,nppx,npp,nnp,sendsiz,use_tg_,f_in,f_aux)
!!$omp do
DO ip = 1, dfft%nproc
IF( use_tg_ ) THEN
gproc = ( ip - 1 ) / dfft%nogrp + 1
IF( MOD( ip - 1, dfft%nogrp ) == 0 ) ii = 0
ELSE
gproc = ip
ii = 0
ENDIF
!
ioff = dfft%iss( ip )
!
DO i = 1, dfft%nsw( ip )
!
mc = dfft%ismap( i + ioff )
!
it = ii * nppx + ( gproc - 1 ) * sendsiz
!
DO j = 1, npp
f_aux( mc + ( j - 1 ) * nnp ) = f_in( j + it )
ENDDO
!
ii = ii + 1
!
ENDDO
!
ENDDO
!!$omp end do
!!$omp end parallel
END IF
ELSE
!
! "backward" scatter from planes to columns
!
IF( isgn == -1 ) THEN
DO ip = 1, dfft%nproc
ioff = dfft%iss( ip )
DO i = 1, dfft%nsp( ip )
mc = dfft%ismap( i + ioff )
it = ( i - 1 ) * nppx + ( ip - 1 ) * sendsiz
DO j = 1, dfft%npp( me )
f_in( j + it ) = f_aux( mc + ( j - 1 ) * dfft%nnp )
ENDDO
ENDDO
ENDDO
ELSE
IF( use_tg_ ) THEN
npp = dfft%tg_npp( me )
nnp = dfft%nr1x * dfft%nr2x
ELSE
npp = dfft%npp( me )
nnp = dfft%nnp
ENDIF
DO ip = 1, dfft%nproc
IF( use_tg_ ) THEN
gproc = ( ip - 1 ) / dfft%nogrp + 1
IF( MOD( ip - 1, dfft%nogrp ) == 0 ) ii = 0
ELSE
gproc = ip
ii = 0
ENDIF
!
ioff = dfft%iss( ip )
!
DO i = 1, dfft%nsw( ip )
!
mc = dfft%ismap( i + ioff )
!
it = ii * nppx + ( gproc - 1 ) * sendsiz
!
DO j = 1, npp
f_in( j + it ) = f_aux( mc + ( j - 1 ) * nnp )
ENDDO
!
ii = ii + 1
!
ENDDO
ENDDO
END IF
IF( nprocp == 1 ) GO TO 20
!
! step two: communication
!
IF( use_tg_ ) THEN
gcomm = dfft%pgrp_comm
ELSE
gcomm = dfft%comm
ENDIF
! CALL mpi_barrier (gcomm, ierr) ! why barrier? for buggy openmpi over ib
CALL mpi_alltoall (f_in(1), sendsiz, MPI_DOUBLE_COMPLEX, f_aux(1), sendsiz, MPI_DOUBLE_COMPLEX, gcomm, ierr)
IF( abs(ierr) /= 0 ) CALL errore ('fft_scatter', 'info<>0', abs(ierr) )
!
! step one: store contiguously the columns
!
f_in = 0.0_DP
!
offset = 1
DO proc = 1, nprocp
from = offset
IF( use_tg_ ) THEN
gproc = dfft%nplist(proc)+1
ELSE
gproc = proc
ENDIF
dest = 1 + ( proc - 1 ) * sendsiz
!
DO k = 1, ncp_ (me)
kdest = dest + (k - 1) * nppx - 1
kfrom = from + (k - 1) * nr3x - 1
DO i = 1, npp_ ( gproc )
f_in ( kfrom + i ) = f_aux ( kdest + i )
ENDDO
ENDDO
offset = offset + npp_ ( gproc )
ENDDO
20 CONTINUE
ENDIF
CALL stop_clock ('fft_scatter')
#endif
#if defined __HPM
! CALL f_hpmstop( 10 )
#endif
RETURN
END SUBROUTINE fft_scatter
!----------------------------------------------------------------------------
SUBROUTINE grid_gather( f_in, f_out )
!----------------------------------------------------------------------------
!
! ... gathers nproc distributed data on the first processor of every pool
!
! ... REAL*8 f_in = distributed variable (nxx)
! ... REAL*8 f_out = gathered variable (nr1x*nr2x*nr3x)
!
USE kinds, ONLY : DP
USE parallel_include
!
IMPLICIT NONE
!
REAL(DP) :: f_in( : ), f_out( : )
!
#if defined (__MPI)
!
INTEGER :: proc, info
INTEGER :: displs(0:dfftp%nproc-1), recvcount(0:dfftp%nproc-1)
!
IF( size( f_in ) < dfftp%nnr ) &
CALL errore( ' grid_gather ', ' f_in too small ', dfftp%nnr - size( f_in ) )
!
CALL start_clock( 'gather' )
!
DO proc = 0, ( dfftp%nproc - 1 )
!
recvcount(proc) = dfftp%nnp * dfftp%npp(proc+1)
!
IF ( proc == 0 ) THEN
!
displs(proc) = 0
!
ELSE
!
displs(proc) = displs(proc-1) + recvcount(proc-1)
!
ENDIF
!
ENDDO
!
info = size( f_out ) - displs( dfftp%nproc - 1 ) - recvcount( dfftp%nproc - 1 )
!
IF( info < 0 ) &
CALL errore( ' grid_gather ', ' f_out too small ', -info )
!
info = 0
!
CALL MPI_GATHERV( f_in, recvcount(dfftp%mype), MPI_DOUBLE_PRECISION, f_out, &
recvcount, displs, MPI_DOUBLE_PRECISION, dfftp%root, &
dfftp%comm, info )
!
CALL errore( 'grid_gather', 'info<>0', info )
!
CALL stop_clock( 'gather' )
!
#else
CALL errore('grid_gather', 'do not use in serial execution', 1)
#endif
!
RETURN
!
END SUBROUTINE grid_gather
!----------------------------------------------------------------------------
SUBROUTINE grid_scatter( f_in, f_out )
!----------------------------------------------------------------------------
!
! ... scatters data from the first processor of every pool
!
! ... REAL*8 f_in = gathered variable (nr1x*nr2x*nr3x)
! ... REAL*8 f_out = distributed variable (nxx)
!
USE kinds, ONLY : DP
USE parallel_include
!
IMPLICIT NONE
!
REAL(DP) :: f_in( : ), f_out( : )
!
#if defined (__MPI)
!
INTEGER :: proc, info
INTEGER :: displs(0:dfftp%nproc-1), sendcount(0:dfftp%nproc-1)
!
IF( size( f_out ) < dfftp%nnr ) &
CALL errore( ' grid_scatter ', ' f_out too small ', dfftp%nnr - size( f_in ) )
!
CALL start_clock( 'scatter' )
!
DO proc = 0, ( dfftp%nproc - 1 )
!
sendcount(proc) = dfftp%nnp * dfftp%npp(proc+1)
!
IF ( proc == 0 ) THEN
!
displs(proc) = 0
!
ELSE
!
displs(proc) = displs(proc-1) + sendcount(proc-1)
!
ENDIF
!
ENDDO
!
info = size( f_in ) - displs( dfftp%nproc - 1 ) - sendcount( dfftp%nproc - 1 )
!
IF( info < 0 ) &
CALL errore( ' grid_scatter ', ' f_in too small ', -info )
!
info = 0
!
CALL MPI_SCATTERV( f_in, sendcount, displs, MPI_DOUBLE_PRECISION, &
f_out, sendcount(dfftp%mype), MPI_DOUBLE_PRECISION, &
dfftp%root, dfftp%comm, info )
!
CALL errore( 'grid_scatter', 'info<>0', info )
!
IF ( sendcount(dfftp%mype) /= dfftp%nnr ) f_out(sendcount(dfftp%mype)+1:dfftp%nnr) = 0.D0
!
CALL stop_clock( 'scatter' )
!
#else
CALL errore('grid_scatter', 'do not use in serial execution', 1)
#endif
!
RETURN
!
END SUBROUTINE grid_scatter
!
! ... "gather"-like subroutines
!
!-----------------------------------------------------------------------
SUBROUTINE cgather_sym( f_in, f_out )
!-----------------------------------------------------------------------
!
! ... gather complex data for symmetrization (in phonon code)
! ... COMPLEX*16 f_in = distributed variable (nrxx)
! ... COMPLEX*16 f_out = gathered variable (nr1x*nr2x*nr3x)
!
USE mp, ONLY : mp_barrier
USE parallel_include
!
IMPLICIT NONE
!
COMPLEX(DP) :: f_in( : ), f_out(:)
!
#if defined (__MPI)
!
INTEGER :: proc, info
INTEGER :: displs(0:dfftp%nproc-1), recvcount(0:dfftp%nproc-1)
!
!
CALL start_clock( 'cgather' )
!
DO proc = 0, ( dfftp%nproc - 1 )
!
recvcount(proc) = 2 * dfftp%nnp * dfftp%npp(proc+1)
!
IF ( proc == 0 ) THEN
!
displs(proc) = 0
!
ELSE
!
displs(proc) = displs(proc-1) + recvcount(proc-1)
!
ENDIF
!
ENDDO
!
CALL mp_barrier( dfftp%comm )
!
CALL MPI_ALLGATHERV( f_in, recvcount(dfftp%mype), MPI_DOUBLE_PRECISION, &
f_out, recvcount, displs, MPI_DOUBLE_PRECISION, &
dfftp%comm, info )
!
CALL errore( 'cgather_sym', 'info<>0', info )
!
! CALL mp_barrier( dfftp%comm )
!
CALL stop_clock( 'cgather' )
!
#else
CALL errore('cgather_sym', 'do not use in serial execution', 1)
#endif
!
RETURN
!
END SUBROUTINE cgather_sym
!
!----------------------------------------------------------------------------
SUBROUTINE cgather_smooth ( f_in, f_out )
!----------------------------------------------------------------------------
!
! ... gathers data on the smooth AND complex fft grid
!
! ... gathers nproc distributed data on the first processor of every pool
!
! ... COMPLEX*16 f_in = distributed variable ( dffts%nnr )
! ... COMPLEX*16 f_out = gathered variable (nr1sx*nr2sx*nr3sx)
!
USE mp, ONLY : mp_barrier
USE kinds, ONLY : DP
USE parallel_include
!
IMPLICIT NONE
!
COMPLEX(DP) :: f_in(:), f_out(:)
!
#if defined (__MPI)
!
INTEGER :: proc, info
INTEGER :: displs(0:dfftp%nproc-1), recvcount(0:dfftp%nproc-1)
!
!
CALL start_clock( 'gather' )
!
DO proc = 0, ( dfftp%nproc - 1 )
!
recvcount(proc) = 2 * dffts%nnp * dffts%npp(proc+1)
!
IF ( proc == 0 ) THEN
!
displs(proc) = 0
!
ELSE
!
displs(proc) = displs(proc-1) + recvcount(proc-1)
!
ENDIF
!
ENDDO
!
CALL mp_barrier( dfftp%comm )
!
CALL MPI_GATHERV( f_in, recvcount(dfftp%mype), MPI_DOUBLE_PRECISION, f_out, &
recvcount, displs, MPI_DOUBLE_PRECISION, dfftp%root, &
dfftp%comm, info )
!
CALL errore( 'cgather_smooth', 'info<>0', info )
!
CALL stop_clock( 'gather' )
!
#else
CALL errore('cgather_smooth', 'do not use in serial execution', 1)
#endif
!
RETURN
!
END SUBROUTINE cgather_smooth
!
!----------------------------------------------------------------------------
SUBROUTINE cgather_custom ( f_in, f_out, dfftt )
!----------------------------------------------------------------------------
!
! ... gathers data on the custom AND complex fft grid
!
! ... gathers nproc distributed data on the first processor of every pool
!
! ... COMPLEX*16 f_in = distributed variable ( dfftt%nnr )
! ... COMPLEX*16 f_out = gathered variable (nr1sx*nr2sx*nr3sx)
!
USE mp, ONLY : mp_barrier
USE kinds, ONLY : DP
USE parallel_include
!
IMPLICIT NONE
!
COMPLEX(DP) :: f_in(:), f_out(:)
TYPE ( fft_dlay_descriptor ), INTENT(IN) :: dfftt
!
#if defined (__MPI)
!
INTEGER :: proc, info
INTEGER :: displs(0:dfftp%nproc-1), recvcount(0:dfftp%nproc-1)
!
!
CALL start_clock( 'gather' )
!
DO proc = 0, ( dfftp%nproc - 1 )
!
recvcount(proc) = 2 * dfftt%nnp * dfftt%npp(proc+1)
!
IF ( proc == 0 ) THEN
!
displs(proc) = 0
!
ELSE
!
displs(proc) = displs(proc-1) + recvcount(proc-1)
!
ENDIF
!
ENDDO
!
CALL mp_barrier( dfftp%comm )
!
CALL MPI_GATHERV( f_in, recvcount(dfftp%mype), MPI_DOUBLE_PRECISION, f_out, &
recvcount, displs, MPI_DOUBLE_PRECISION, dfftp%root, &
dfftp%comm, info )
!
CALL errore( 'cgather_custom', 'info<>0', info )
!
CALL stop_clock( 'gather' )
!
#else
CALL errore('cgather_custom', 'do not use in serial execution', 1)
#endif
!
RETURN
!
END SUBROUTINE cgather_custom
!
! ... "scatter"-like subroutines
!
!----------------------------------------------------------------------------
SUBROUTINE cscatter_sym( f_in, f_out )
!----------------------------------------------------------------------------
!
! ... scatters data from the first processor of every pool
!
! ... COMPLEX*16 f_in = gathered variable (nr1x*nr2x*nr3x)
! ... COMPLEX*16 f_out = distributed variable (nxx)
!
USE mp, ONLY : mp_barrier
USE kinds, ONLY : DP
USE parallel_include
!
IMPLICIT NONE
!
COMPLEX(DP) :: f_in(:), f_out(:)
!
#if defined (__MPI)
!
INTEGER :: proc, info
INTEGER :: displs(0:dfftp%nproc-1), sendcount(0:dfftp%nproc-1)
!
!
CALL start_clock( 'cscatter_sym' )
!
DO proc = 0, ( dfftp%nproc - 1 )
!
sendcount(proc) = 2 * dfftp%nnp * dfftp%npp(proc+1)
!
IF ( proc == 0 ) THEN
!
displs(proc) = 0
!
ELSE
!
displs(proc) = displs(proc-1) + sendcount(proc-1)
!
ENDIF
!
ENDDO
!
CALL mp_barrier( dfftp%comm )
!
CALL MPI_SCATTERV( f_in, sendcount, displs, MPI_DOUBLE_PRECISION, &
f_out, sendcount(dfftp%mype), MPI_DOUBLE_PRECISION, &
dfftp%root, dfftp%comm, info )
!
CALL errore( 'cscatter_sym', 'info<>0', info )
!
IF ( sendcount(dfftp%mype) /= dfftp%nnr ) f_out(sendcount(dfftp%mype)+1: dfftp%nnr ) = 0.D0
!
CALL stop_clock( 'cscatter_sym' )
!
#else
CALL errore('cscatter_sym', 'do not use in serial execution', 1)
#endif
!
RETURN
!
END SUBROUTINE cscatter_sym
!
!----------------------------------------------------------------------------
SUBROUTINE cscatter_smooth( f_in, f_out )
!----------------------------------------------------------------------------
!
! ... scatters data on the smooth AND complex fft grid
! ... scatters data from the first processor of every pool
!
! ... COMPLEX*16 f_in = gathered variable (nr1sx*nr2sx*nr3sx)
! ... COMPLEX*16 f_out = distributed variable ( dffts%nnr)
!
USE mp, ONLY : mp_barrier
USE kinds, ONLY : DP
USE parallel_include
!
IMPLICIT NONE
!
COMPLEX(DP) :: f_in(:), f_out(:)
!
#if defined (__MPI)
!
INTEGER :: proc, info
INTEGER :: displs(0:dfftp%nproc-1), sendcount(0:dfftp%nproc-1)
!
!
CALL start_clock( 'scatter' )
!
DO proc = 0, ( dfftp%nproc - 1 )
!
sendcount(proc) = 2 * dffts%nnp * dffts%npp(proc+1)
!
IF ( proc == 0 ) THEN
!
displs(proc) = 0
!
ELSE
!
displs(proc) = displs(proc-1) + sendcount(proc-1)
!
ENDIF
!
ENDDO
!
CALL mp_barrier( dfftp%comm )
!
CALL MPI_SCATTERV( f_in, sendcount, displs, MPI_DOUBLE_PRECISION, &
f_out, sendcount(dfftp%mype), MPI_DOUBLE_PRECISION, &
dfftp%root, dfftp%comm, info )
!
CALL errore( 'scatter', 'info<>0', info )
!
IF ( sendcount(dfftp%mype) /= dffts%nnr ) f_out(sendcount(dfftp%mype)+1: dffts%nnr ) = 0.D0
!
CALL stop_clock( 'scatter' )
!
#else
CALL errore('cscatter_smooth', 'do not use in serial execution', 1)
#endif
!
RETURN
!
END SUBROUTINE cscatter_smooth
!
!----------------------------------------------------------------------------
SUBROUTINE cscatter_custom( f_in, f_out, dfftt )
!----------------------------------------------------------------------------
!
! ... scatters data on the custom AND complex fft grid
! ... scatters data from the first processor of every pool
!
! ... COMPLEX*16 f_in = gathered variable (nr1sx*nr2sx*nr3sx)
! ... COMPLEX*16 f_out = distributed variable ( dfftt%nnr)
!
USE mp, ONLY : mp_barrier
USE kinds, ONLY : DP
USE parallel_include
!
IMPLICIT NONE
!
COMPLEX(DP) :: f_in(:), f_out(:)
TYPE ( fft_dlay_descriptor ), INTENT(IN) :: dfftt
!
#if defined (__MPI)
!
INTEGER :: proc, info
INTEGER :: displs(0:dfftp%nproc-1), sendcount(0:dfftp%nproc-1)
!
!
CALL start_clock( 'scatter' )
!
DO proc = 0, ( dfftp%nproc - 1 )
!
sendcount(proc) = 2 * dfftt%nnp * dfftt%npp(proc+1)
!
IF ( proc == 0 ) THEN
!
displs(proc) = 0
!
ELSE
!
displs(proc) = displs(proc-1) + sendcount(proc-1)
!
ENDIF
!
ENDDO
!
CALL mp_barrier( dfftp%comm )
!
CALL MPI_SCATTERV( f_in, sendcount, displs, MPI_DOUBLE_PRECISION, &
f_out, sendcount(dfftp%mype), MPI_DOUBLE_PRECISION, &
dfftp%root, dfftp%comm, info )
!
CALL errore( 'scatter', 'info<>0', info )
!
IF ( sendcount(dfftp%mype) /= dfftt%nnr ) f_out(sendcount(dfftp%mype)+1: dfftt%nnr ) = 0.D0
!
CALL stop_clock( 'scatter' )
!
#else
CALL errore('cscatter_custom', 'do not use in serial execution', 1)
#endif
!
RETURN
!
END SUBROUTINE cscatter_custom
!
!----------------------------------------------------------------------------
SUBROUTINE gather_smooth ( f_in, f_out )
!----------------------------------------------------------------------------
!
! ... gathers data on the smooth AND real fft grid
!
! ... gathers nproc distributed data on the first processor of every pool
!
! ... REAL*8 f_in = distributed variable ( dffts%nnr )
! ... REAL*8 f_out = gathered variable (nr1sx*nr2sx*nr3sx)
!
USE mp, ONLY : mp_barrier
USE kinds, ONLY : DP
USE parallel_include
!
IMPLICIT NONE
!
REAL(DP) :: f_in(:), f_out(:)
!
#if defined (__MPI)
!
INTEGER :: proc, info
INTEGER :: displs(0:dffts%nproc-1), recvcount(0:dffts%nproc-1)
!
!
CALL start_clock( 'gather' )
!
DO proc = 0, ( dffts%nproc - 1 )
!
recvcount(proc) = dffts%nnp * dffts%npp(proc+1)
!
IF ( proc == 0 ) THEN
!
displs(proc) = 0
!
ELSE
!
displs(proc) = displs(proc-1) + recvcount(proc-1)
!
ENDIF
!
ENDDO
!
CALL mp_barrier( dffts%comm )
!
CALL MPI_GATHERV( f_in, recvcount(dffts%mype), MPI_DOUBLE_PRECISION, f_out, &
recvcount, displs, MPI_DOUBLE_PRECISION, dffts%root, &
dffts%comm, info )
!
CALL errore( 'gather', 'info<>0', info )
!
CALL stop_clock( 'gather' )
!
#else
CALL errore('gather_smooth', 'do not use in serial execution', 1)
#endif
!
RETURN
!
END SUBROUTINE gather_smooth
!
!----------------------------------------------------------------------------
SUBROUTINE scatter_smooth( f_in, f_out )
!----------------------------------------------------------------------------
!
! ... scatters data on the smooth AND real fft grid
! ... scatters data from the first processor of every pool
!
! ... REAL*8 f_in = gathered variable (nr1sx*nr2sx*nr3sx)
! ... REAL*8 f_out = distributed variable ( dffts%nnr)
!
USE mp, ONLY : mp_barrier
USE kinds, ONLY : DP
USE parallel_include
!
IMPLICIT NONE
!
REAL(DP) :: f_in(:), f_out(:)
!
#if defined (__MPI)
!
INTEGER :: proc, info
INTEGER :: displs(0:dffts%nproc-1), sendcount(0:dffts%nproc-1)
!
!
CALL start_clock( 'scatter' )
!
DO proc = 0, ( dffts%nproc - 1 )
!
sendcount(proc) = dffts%nnp * dffts%npp(proc+1)
!
IF ( proc == 0 ) THEN
!
displs(proc) = 0
!
ELSE
!
displs(proc) = displs(proc-1) + sendcount(proc-1)
!
ENDIF
!
ENDDO
!
CALL mp_barrier( dffts%comm )
!
CALL MPI_SCATTERV( f_in, sendcount, displs, MPI_DOUBLE_PRECISION, &
f_out, sendcount(dffts%mype), MPI_DOUBLE_PRECISION, &
dffts%root, dffts%comm, info )
!
CALL errore( 'scatter', 'info<>0', info )
!
IF ( sendcount(dffts%mype) /= dffts%nnr ) f_out(sendcount(dffts%mype)+1: dffts%nnr ) = 0.D0
!
CALL stop_clock( 'scatter' )
!
#else
CALL errore('scatter_smooth', 'do not use in serial execution', 1)
#endif
!
RETURN
!
END SUBROUTINE scatter_smooth
!
SUBROUTINE tg_gather( dffts, v, tg_v )
!
USE parallel_include
!
USE fft_types, ONLY : fft_dlay_descriptor
! T.G.
! NOGRP: Number of processors per orbital task group
IMPLICIT NONE
TYPE(fft_dlay_descriptor), INTENT(in) :: dffts
REAL(DP) :: v(:)
REAL(DP) :: tg_v(:)
INTEGER :: nsiz, i, ierr, nsiz_tg
INTEGER :: recv_cnt( dffts%nogrp ), recv_displ( dffts%nogrp )
nsiz_tg = dffts%tg_nnr * dffts%nogrp
IF( size( tg_v ) < nsiz_tg ) &
CALL errore( ' tg_gather ', ' tg_v too small ', ( nsiz_tg - size( tg_v ) ) )
nsiz = dffts%npp( dffts%mype+1 ) * dffts%nr1x * dffts%nr2x
IF( size( v ) < nsiz ) &
CALL errore( ' tg_gather ', ' v too small ', ( nsiz - size( v ) ) )
!
! The potential in v is distributed across all processors
! We need to redistribute it so that it is completely contained in the
! processors of an orbital TASK-GROUP
!
recv_cnt(1) = dffts%npp( dffts%nolist(1) + 1 ) * dffts%nr1x * dffts%nr2x
recv_displ(1) = 0
DO i = 2, dffts%nogrp
recv_cnt(i) = dffts%npp( dffts%nolist(i) + 1 ) * dffts%nr1x * dffts%nr2x
recv_displ(i) = recv_displ(i-1) + recv_cnt(i-1)
ENDDO
! clean only elements that will not be overwritten
!
DO i = recv_displ(dffts%nogrp) + recv_cnt( dffts%nogrp ) + 1, size( tg_v )
tg_v( i ) = 0.0d0
ENDDO
#if defined (__MPI)
CALL MPI_Allgatherv( v(1), nsiz, MPI_DOUBLE_PRECISION, &
tg_v(1), recv_cnt, recv_displ, MPI_DOUBLE_PRECISION, dffts%ogrp_comm, IERR)
IF( ierr /= 0 ) &
CALL errore( ' tg_gather ', ' MPI_Allgatherv ', abs( ierr ) )
#endif
END SUBROUTINE tg_gather
!
! Complex version of previous routine
!
SUBROUTINE tg_cgather( dffts, v, tg_v )
!
USE parallel_include
!
USE fft_types, ONLY : fft_dlay_descriptor
! T.G.
! NOGRP: Number of processors per orbital task group
IMPLICIT NONE
TYPE(fft_dlay_descriptor), INTENT(in) :: dffts
COMPLEX(DP) :: v(:)
COMPLEX(DP) :: tg_v(:)
INTEGER :: nsiz, i, ierr, nsiz_tg
INTEGER :: recv_cnt( dffts%nogrp ), recv_displ( dffts%nogrp )
nsiz_tg = dffts%tg_nnr * dffts%nogrp
IF( size( tg_v ) < nsiz_tg ) &
CALL errore( ' tg_gather ', ' tg_v too small ', ( nsiz_tg - size( tg_v ) ) )
nsiz = dffts%npp( dffts%mype+1 ) * dffts%nr1x * dffts%nr2x
IF( size( v ) < nsiz ) &
CALL errore( ' tg_gather ', ' v too small ', ( nsiz - size( v ) ) )
!
! The potential in v is distributed across all processors
! We need to redistribute it so that it is completely contained in the
! processors of an orbital TASK-GROUP
!
recv_cnt(1) = dffts%npp( dffts%nolist(1) + 1 ) * dffts%nr1x * dffts%nr2x
recv_displ(1) = 0
DO i = 2, dffts%nogrp
recv_cnt(i) = dffts%npp( dffts%nolist(i) + 1 ) * dffts%nr1x * dffts%nr2x
recv_displ(i) = recv_displ(i-1) + recv_cnt(i-1)
ENDDO
! clean only elements that will not be overwritten
!
DO i = recv_displ(dffts%nogrp) + recv_cnt( dffts%nogrp ) + 1, size( tg_v )
tg_v( i ) = (0.0d0,0.0d0)
ENDDO
!
! The quantities are complex, multiply the cunters by 2 and gather
! real numbers
!
nsiz = 2 * nsiz
recv_cnt = 2 * recv_cnt
recv_displ = 2 * recv_displ
#if defined (__MPI)
CALL MPI_Allgatherv( v(1), nsiz, MPI_DOUBLE_PRECISION, &
tg_v(1), recv_cnt, recv_displ, MPI_DOUBLE_PRECISION, dffts%ogrp_comm, IERR)
IF( ierr /= 0 ) &
CALL errore( ' tg_cgather ', ' MPI_Allgatherv ', abs( ierr ) )
#endif
END SUBROUTINE tg_cgather
!=----------------------------------------------------------------------=!
END MODULE fft_base
!=----------------------------------------------------------------------=!
|