TABLE OF CONTENTS
- ABINIT/bandfft_kpt_mpi_recv
- ABINIT/bandfft_kpt_mpi_send
- ABINIT/m_bandfft_kpt
- ABINIT/prep_bandfft_tabs
- m_bandfft_kpt/bandfft_kpt_copy
- m_bandfft_kpt/bandfft_kpt_destroy
- m_bandfft_kpt/bandfft_kpt_destroy_array
- m_bandfft_kpt/bandfft_kpt_get_ikpt
- m_bandfft_kpt/bandfft_kpt_init1
- m_bandfft_kpt/bandfft_kpt_init2
- m_bandfft_kpt/bandfft_kpt_reset
- m_bandfft_kpt/bandfft_kpt_restoretabs
- m_bandfft_kpt/bandfft_kpt_savetabs
- m_bandfft_kpt/bandfft_kpt_set_ikpt
- m_bandfft_kpt/bandfft_kpt_type
ABINIT/bandfft_kpt_mpi_recv [ Functions ]
NAME
bandfft_kpt_mpi_recv
FUNCTION
Receive a bandfft_kpt_type inside a MPI communicator.
INPUTS
sender=ID of the sender in spaceComm tag=message tag spaceComm=MPI Communicator
OUTPUT
ierr=Error status output=# of on proc. sender
PARENTS
posdoppler
CHILDREN
SOURCE
1524 subroutine bandfft_kpt_mpi_recv(output,sender,tag,spaceComm,ierr) 1525 1526 1527 !This section has been created automatically by the script Abilint (TD). 1528 !Do not modify the following lines by hand. 1529 #undef ABI_FUNC 1530 #define ABI_FUNC 'bandfft_kpt_mpi_recv' 1531 !End of the abilint section 1532 1533 implicit none 1534 1535 !Arguments ------------------------------------ 1536 !scalars 1537 integer,intent(in) :: sender,spaceComm,tag 1538 integer,intent(out) :: ierr 1539 !arrays 1540 type(bandfft_kpt_type),intent(out) :: output 1541 1542 !Local variables------------------------------- 1543 !scalars 1544 integer :: ipck,nsize,size_dp,size_int 1545 integer :: size1_kg_k_gather=0,size2_kg_k_gather=0,size_recvcounts=0 1546 integer :: size_sendcounts=0,size_rdispls=0,size_sdispls=0,size1_gbound=0,size2_gbound 1547 integer :: size1_kg_k_gather_sym=0,size2_kg_k_gather_sym=0,size_rdispls_sym=0 1548 integer :: size_sdispls_sym=0,size_recvcounts_sym=0,size_recvcounts_sym_tot=0 1549 integer :: size_sendcounts_sym=0,size_sendcounts_sym_all=0,size_tab_proc=0 1550 integer :: size_indices_pw_fft=0,size_sendcount_fft=0,size_senddisp_fft=0 1551 integer :: size_recvcount_fft=0,size_recvdisp_fft=0,size1_kg_k_fft=0,size2_kg_k_fft=0 1552 integer :: size1_ffnl_gather=0,size2_ffnl_gather=0,size3_ffnl_gather=0,size4_ffnl_gather=0 1553 integer :: size_kinpw_gather=0,size1_ph3d_gather=0,size2_ph3d_gather=0,size3_ph3d_gather=0 1554 integer :: size1_kpg_k_gather=0,size2_kpg_k_gather=0,sz1,sz2,sz3,sz4 1555 !arrays 1556 integer,allocatable :: buffer_int(:) 1557 real(dp),allocatable :: buffer_dp(:) 1558 1559 ! ************************************************************************* 1560 1561 if (xmpi_comm_size(spaceComm)<=1) return 1562 1563 !=== Receive flags and array sizes ==== 1564 nsize=45 1565 ABI_ALLOCATE(buffer_int,(nsize)) 1566 call xmpi_recv(buffer_int,sender,3*tag-2,spaceComm,ierr) 1567 output%flag1_is_allocated=buffer_int(1) 1568 output%flag2_is_allocated=buffer_int(2) 1569 output%flag3_is_allocated=buffer_int(3) 1570 output%have_to_reequilibrate=(buffer_int(4)/=0) 1571 output%istwf_k=buffer_int(5) 1572 output%npw_tot=buffer_int(6) 1573 output%npw_fft=buffer_int(7) 1574 output%ndatarecv=buffer_int(8) 1575 output%idatarecv0=buffer_int(9) 1576 output%ndatarecv_tot=buffer_int(10) 1577 output%ndatasend_sym=buffer_int(11) 1578 size_recvcounts=buffer_int(12) 1579 size_sendcounts=buffer_int(13) 1580 size_rdispls=buffer_int(14) 1581 size_sdispls=buffer_int(15) 1582 size1_gbound=buffer_int(16) 1583 size2_gbound=buffer_int(17) 1584 size1_ffnl_gather=buffer_int(18) 1585 size2_ffnl_gather=buffer_int(19) 1586 size3_ffnl_gather=buffer_int(20) 1587 size4_ffnl_gather=buffer_int(21) 1588 size_kinpw_gather=buffer_int(22) 1589 size1_ph3d_gather=buffer_int(23) 1590 size2_ph3d_gather=buffer_int(24) 1591 size3_ph3d_gather=buffer_int(25) 1592 size1_kpg_k_gather=buffer_int(26) 1593 size2_kpg_k_gather=buffer_int(27) 1594 size_indices_pw_fft=buffer_int(28) 1595 size_sendcount_fft=buffer_int(29) 1596 size_senddisp_fft=buffer_int(30) 1597 size_recvcount_fft=buffer_int(31) 1598 size_recvdisp_fft=buffer_int(32) 1599 size1_kg_k_fft=buffer_int(33) 1600 size2_kg_k_fft=buffer_int(34) 1601 size1_kg_k_gather=buffer_int(35) 1602 size2_kg_k_gather=buffer_int(36) 1603 size1_kg_k_gather_sym=buffer_int(37) 1604 size2_kg_k_gather_sym=buffer_int(38) 1605 size_rdispls_sym=buffer_int(39) 1606 size_sdispls_sym=buffer_int(40) 1607 size_recvcounts_sym=buffer_int(41) 1608 size_recvcounts_sym_tot=buffer_int(42) 1609 size_sendcounts_sym=buffer_int(43) 1610 size_sendcounts_sym_all=buffer_int(44) 1611 size_tab_proc=buffer_int(45) 1612 ABI_DEALLOCATE(buffer_int) 1613 1614 !=== Compute amount of transmitted data ==== 1615 size_int=size1_kg_k_gather*size2_kg_k_gather & 1616 & +size_recvcounts+size_sendcounts+size_rdispls+size_sdispls & 1617 & +size1_gbound*size2_gbound & 1618 & +size1_kg_k_gather_sym*size2_kg_k_gather_sym & 1619 & +size_rdispls_sym+size_sdispls_sym & 1620 & +size_recvcounts_sym+size_recvcounts_sym_tot & 1621 & +size_sendcounts_sym+size_sendcounts_sym_all+size_tab_proc & 1622 & +size_indices_pw_fft+size_sendcount_fft+size_senddisp_fft & 1623 & +size_recvcount_fft+size_recvdisp_fft & 1624 & +size1_kg_k_fft*size2_kg_k_fft 1625 size_dp=size1_ffnl_gather*size2_ffnl_gather*size3_ffnl_gather*size4_ffnl_gather & 1626 & +size_kinpw_gather & 1627 & +size1_ph3d_gather*size2_ph3d_gather*size3_ph3d_gather & 1628 & +size1_kpg_k_gather*size2_kpg_k_gather 1629 1630 !=== Receive and unpack integer data ==== 1631 if (size_int>0) then 1632 ipck=0 1633 ABI_ALLOCATE(buffer_int,(size_int)) 1634 call xmpi_recv(buffer_int,sender,3*tag-1,spaceComm,ierr) 1635 if (allocated(output%kg_k_gather)) then 1636 ABI_DEALLOCATE(output%kg_k_gather) 1637 end if 1638 if (size1_kg_k_gather*size2_kg_k_gather>0) then 1639 nsize=size1_kg_k_gather*size2_kg_k_gather 1640 sz1=size1_kg_k_gather;sz2=size2_kg_k_gather 1641 ABI_ALLOCATE(output%kg_k_gather,(sz1,sz2)) 1642 output%kg_k_gather(:,:)=reshape(buffer_int(ipck+1:ipck+nsize),(/sz1,sz2/)) 1643 ipck=ipck+nsize 1644 end if 1645 if (allocated(output%recvcounts)) then 1646 ABI_DEALLOCATE(output%recvcounts) 1647 end if 1648 if (size_recvcounts>0) then 1649 ABI_ALLOCATE(output%recvcounts,(size_recvcounts)) 1650 output%recvcounts(:)=buffer_int(ipck+1:ipck+size_recvcounts) 1651 ipck=ipck+size_recvcounts 1652 end if 1653 if (allocated(output%sendcounts)) then 1654 ABI_DEALLOCATE(output%sendcounts) 1655 end if 1656 if (size_sendcounts>0) then 1657 ABI_ALLOCATE(output%sendcounts,(size_sendcounts)) 1658 output%sendcounts(:)=buffer_int(ipck+1:ipck+size_sendcounts) 1659 ipck=ipck+size_sendcounts 1660 end if 1661 if (allocated(output%rdispls)) then 1662 ABI_DEALLOCATE(output%rdispls) 1663 end if 1664 if (size_rdispls>0) then 1665 ABI_ALLOCATE(output%rdispls,(size_rdispls)) 1666 output%rdispls(:)=buffer_int(ipck+1:ipck+size_rdispls) 1667 ipck=ipck+size_rdispls 1668 end if 1669 if (allocated(output%sdispls)) then 1670 ABI_DEALLOCATE(output%sdispls) 1671 end if 1672 if (size_sdispls>0) then 1673 ABI_ALLOCATE(output%sdispls,(size_sdispls)) 1674 output%sdispls(:)=buffer_int(ipck+1:ipck+size_sdispls) 1675 ipck=ipck+size_sdispls 1676 end if 1677 if (allocated(output%gbound)) then 1678 ABI_DEALLOCATE(output%gbound) 1679 end if 1680 if (size1_gbound*size2_gbound>0) then 1681 nsize=size1_gbound*size2_gbound 1682 sz1=size1_gbound;sz2=size2_gbound 1683 ABI_ALLOCATE(output%gbound,(sz1,sz2)) 1684 output%gbound(:,:)=reshape(buffer_int(ipck+1:ipck+nsize),(/sz1,sz2/)) 1685 ipck=ipck+nsize 1686 end if 1687 if (allocated(output%kg_k_gather_sym)) then 1688 ABI_DEALLOCATE(output%kg_k_gather_sym) 1689 end if 1690 if (size1_kg_k_gather_sym*size2_kg_k_gather_sym>0) then 1691 nsize=size1_kg_k_gather_sym*size2_kg_k_gather_sym 1692 sz1=size1_kg_k_gather_sym;sz2=size2_kg_k_gather_sym 1693 ABI_ALLOCATE(output%kg_k_gather_sym,(sz1,sz2)) 1694 output%kg_k_gather_sym(:,:)=reshape(buffer_int(ipck+1:ipck+nsize),(/sz1,sz2/)) 1695 ipck=ipck+nsize 1696 end if 1697 if (allocated(output%rdispls_sym)) then 1698 ABI_DEALLOCATE(output%rdispls_sym) 1699 end if 1700 if (size_rdispls_sym>0) then 1701 ABI_ALLOCATE(output%rdispls_sym,(size_rdispls_sym)) 1702 output%rdispls_sym(:)=buffer_int(ipck+1:ipck+size_rdispls_sym) 1703 ipck=ipck+size_rdispls_sym 1704 end if 1705 if (allocated(output%sdispls_sym)) then 1706 ABI_DEALLOCATE(output%sdispls_sym) 1707 end if 1708 if (size_sdispls_sym>0) then 1709 ABI_ALLOCATE(output%sdispls_sym,(size_sdispls_sym)) 1710 output%sdispls_sym(:)=buffer_int(ipck+1:ipck+size_sdispls_sym) 1711 ipck=ipck+size_sdispls_sym 1712 end if 1713 if (allocated(output%recvcounts_sym)) then 1714 ABI_DEALLOCATE(output%recvcounts_sym) 1715 end if 1716 if (size_recvcounts_sym>0) then 1717 ABI_ALLOCATE(output%recvcounts_sym,(size_recvcounts_sym)) 1718 output%recvcounts_sym(:)=buffer_int(ipck+1:ipck+size_recvcounts_sym) 1719 ipck=ipck+size_recvcounts_sym 1720 end if 1721 if (allocated(output%recvcounts_sym_tot)) then 1722 ABI_DEALLOCATE(output%recvcounts_sym_tot) 1723 end if 1724 if (size_recvcounts_sym_tot>0) then 1725 ABI_ALLOCATE(output%recvcounts_sym_tot,(size_recvcounts_sym_tot)) 1726 output%recvcounts_sym_tot(:)=buffer_int(ipck+1:ipck+size_recvcounts_sym_tot) 1727 ipck=ipck+size_recvcounts_sym_tot 1728 end if 1729 if (allocated(output%sendcounts_sym)) then 1730 ABI_DEALLOCATE(output%sendcounts_sym) 1731 end if 1732 if (size_sendcounts_sym>0) then 1733 ABI_ALLOCATE(output%sendcounts_sym,(size_sendcounts_sym)) 1734 output%sendcounts_sym(:)=buffer_int(ipck+1:ipck+size_sendcounts_sym) 1735 ipck=ipck+size_sendcounts_sym 1736 end if 1737 if (allocated(output%sendcounts_sym_all)) then 1738 ABI_DEALLOCATE(output%sendcounts_sym_all) 1739 end if 1740 if (size_sendcounts_sym_all>0) then 1741 ABI_ALLOCATE(output%sendcounts_sym_all,(size_sendcounts_sym_all)) 1742 output%sendcounts_sym_all(:)=buffer_int(ipck+1:ipck+size_sendcounts_sym_all) 1743 ipck=ipck+size_sendcounts_sym_all 1744 end if 1745 if (allocated(output%tab_proc)) then 1746 ABI_DEALLOCATE(output%tab_proc) 1747 end if 1748 if (size_tab_proc>0) then 1749 ABI_ALLOCATE(output%tab_proc,(size_tab_proc)) 1750 output%tab_proc(:)=buffer_int(ipck+1:ipck+size_tab_proc) 1751 ipck=ipck+size_tab_proc 1752 end if 1753 if (allocated(output%indices_pw_fft)) then 1754 ABI_DEALLOCATE(output%indices_pw_fft) 1755 end if 1756 if (size_indices_pw_fft>0) then 1757 ABI_ALLOCATE(output%indices_pw_fft,(size_indices_pw_fft)) 1758 output%indices_pw_fft(:)=buffer_int(ipck+1:ipck+size_indices_pw_fft) 1759 ipck=ipck+size_indices_pw_fft 1760 end if 1761 if (allocated(output%sendcount_fft)) then 1762 ABI_DEALLOCATE(output%sendcount_fft) 1763 end if 1764 if (size_sendcount_fft>0) then 1765 ABI_ALLOCATE(output%sendcount_fft,(size_sendcount_fft)) 1766 output%sendcount_fft(:)=buffer_int(ipck+1:ipck+size_sendcount_fft) 1767 ipck=ipck+size_sendcount_fft 1768 end if 1769 if (allocated(output%senddisp_fft)) then 1770 ABI_DEALLOCATE(output%senddisp_fft) 1771 end if 1772 if (size_senddisp_fft>0) then 1773 ABI_ALLOCATE(output%senddisp_fft,(size_senddisp_fft)) 1774 output%senddisp_fft(:)=buffer_int(ipck+1:ipck+size_senddisp_fft) 1775 ipck=ipck+size_senddisp_fft 1776 end if 1777 if (allocated(output%recvcount_fft)) then 1778 ABI_DEALLOCATE(output%recvcount_fft) 1779 end if 1780 if (size_recvcount_fft>0) then 1781 ABI_ALLOCATE(output%recvcount_fft,(size_recvcount_fft)) 1782 output%recvcount_fft(:)=buffer_int(ipck+1:ipck+size_recvcount_fft) 1783 ipck=ipck+size_recvcount_fft 1784 end if 1785 if (allocated(output%recvdisp_fft)) then 1786 ABI_DEALLOCATE(output%recvdisp_fft) 1787 end if 1788 if (size_recvdisp_fft>0) then 1789 ABI_ALLOCATE(output%recvdisp_fft,(size_recvdisp_fft)) 1790 output%recvdisp_fft(:)=buffer_int(ipck+1:ipck+size_recvdisp_fft) 1791 ipck=ipck+size_recvdisp_fft 1792 end if 1793 if (allocated(output%kg_k_fft)) then 1794 ABI_DEALLOCATE(output%kg_k_fft) 1795 end if 1796 if (size1_kg_k_fft*size2_kg_k_fft>0) then 1797 nsize=size1_kg_k_fft*size2_kg_k_fft 1798 sz1=size1_kg_k_fft;sz2=size2_kg_k_fft 1799 ABI_ALLOCATE(output%kg_k_fft,(sz1,sz2)) 1800 output%kg_k_fft(:,:)=reshape(buffer_int(ipck+1:ipck+nsize),(/sz1,sz2/)) 1801 ipck=ipck+nsize 1802 end if 1803 ABI_DEALLOCATE(buffer_int) 1804 end if 1805 1806 !=== Receive and unpack real data ==== 1807 if (size_dp>0) then 1808 ipck=0 1809 ABI_ALLOCATE(buffer_dp,(size_dp)) 1810 call xmpi_recv(buffer_dp,sender,3*tag,spaceComm,ierr) 1811 if (allocated(output%ffnl_gather)) then 1812 ABI_DEALLOCATE(output%ffnl_gather) 1813 end if 1814 if (size1_ffnl_gather*size2_ffnl_gather*size3_ffnl_gather*size4_ffnl_gather>0) then 1815 nsize=size1_ffnl_gather*size2_ffnl_gather*size3_ffnl_gather*size4_ffnl_gather 1816 sz1=size1_ffnl_gather;sz2=size2_ffnl_gather;sz3=size3_ffnl_gather;sz4=size4_ffnl_gather 1817 ABI_ALLOCATE(output%ffnl_gather,(sz1,sz2,sz3,sz4)) 1818 output%ffnl_gather(:,:,:,:)=reshape(buffer_dp(ipck+1:ipck+nsize),(/sz1,sz2,sz3,sz4/)) 1819 ipck=ipck+nsize 1820 end if 1821 if (allocated(output%kinpw_gather)) then 1822 ABI_DEALLOCATE(output%kinpw_gather) 1823 end if 1824 if (size_kinpw_gather>0) then 1825 ABI_ALLOCATE(output%kinpw_gather,(size_kinpw_gather)) 1826 output%kinpw_gather(:)=buffer_dp(ipck+1:ipck+size_kinpw_gather) 1827 ipck=ipck+size_kinpw_gather 1828 end if 1829 if (allocated(output%ph3d_gather)) then 1830 ABI_DEALLOCATE(output%ph3d_gather) 1831 end if 1832 if (size1_ph3d_gather*size2_ph3d_gather*size3_ph3d_gather>0) then 1833 nsize=size1_ph3d_gather*size2_ph3d_gather*size3_ph3d_gather 1834 sz1=size1_ph3d_gather;sz2=size2_ph3d_gather;sz3=size3_ph3d_gather 1835 ABI_ALLOCATE(output%ph3d_gather,(sz1,sz2,sz3)) 1836 output%ph3d_gather(:,:,:)=reshape(buffer_dp(ipck+1:ipck+nsize),(/sz1,sz2,sz3/)) 1837 ipck=ipck+nsize 1838 end if 1839 if (allocated(output%kpg_k_gather)) then 1840 ABI_DEALLOCATE(output%kpg_k_gather) 1841 end if 1842 if (size1_kpg_k_gather*size2_kpg_k_gather>0) then 1843 nsize=size1_kpg_k_gather*size2_kpg_k_gather 1844 sz1=size1_kpg_k_gather;sz2=size2_kpg_k_gather 1845 ABI_ALLOCATE(output%kpg_k_gather,(sz1,sz2)) 1846 output%kpg_k_gather(:,:)=reshape(buffer_dp(ipck+1:ipck+nsize),(/sz1,sz2/)) 1847 ipck=ipck+nsize 1848 end if 1849 ABI_DEALLOCATE(buffer_dp) 1850 end if 1851 1852 end subroutine bandfft_kpt_mpi_recv
ABINIT/bandfft_kpt_mpi_send [ Functions ]
NAME
bandfft_kpt_mpi_send
FUNCTION
Send a bandfft_kpt_type inside a MPI communicator
INPUTS
input=The datatype to be transmitted receiver=ID of the receiver in spaceComm tag=message tag spaceComm=MPI Communicator [profile]=(character, optional, default=full) Define the part of the datastructure to be sent; possible values: 'full'= the entire datastructure 'fourwf'= only arrays needed by the fourwf routine
OUTPUT
ierr=Error status
PARENTS
posdoppler
CHILDREN
SOURCE
1220 subroutine bandfft_kpt_mpi_send(input,receiver,tag,spaceComm,ierr,profile) 1221 1222 1223 !This section has been created automatically by the script Abilint (TD). 1224 !Do not modify the following lines by hand. 1225 #undef ABI_FUNC 1226 #define ABI_FUNC 'bandfft_kpt_mpi_send' 1227 !End of the abilint section 1228 1229 implicit none 1230 1231 !Arguments ------------------------------------ 1232 !scalars 1233 integer,intent(in) :: receiver,spaceComm,tag 1234 integer,intent(out) :: ierr 1235 character(len=*),optional,intent(in) :: profile 1236 !arrays 1237 type(bandfft_kpt_type),intent(in) :: input 1238 1239 !Local variables------------------------------- 1240 !scalars 1241 integer :: ipck,nsize,size_dp,size_int 1242 integer :: size1_kg_k_gather=0,size2_kg_k_gather=0,size_recvcounts=0 1243 integer :: size_sendcounts=0,size_rdispls=0,size_sdispls=0,size1_gbound=0,size2_gbound 1244 integer :: size1_kg_k_gather_sym=0,size2_kg_k_gather_sym=0,size_rdispls_sym=0 1245 integer :: size_sdispls_sym=0,size_recvcounts_sym=0,size_recvcounts_sym_tot=0 1246 integer :: size_sendcounts_sym=0,size_sendcounts_sym_all=0,size_tab_proc=0 1247 integer :: size_indices_pw_fft=0,size_sendcount_fft=0,size_senddisp_fft=0 1248 integer :: size_recvcount_fft=0,size_recvdisp_fft=0,size1_kg_k_fft=0,size2_kg_k_fft=0 1249 integer :: size1_ffnl_gather=0,size2_ffnl_gather=0,size3_ffnl_gather=0,size4_ffnl_gather=0 1250 integer :: size_kinpw_gather=0,size1_ph3d_gather=0,size2_ph3d_gather=0,size3_ph3d_gather=0 1251 integer :: size1_kpg_k_gather=0,size2_kpg_k_gather=0 1252 logical :: fourwf,full 1253 !arrays 1254 integer,allocatable :: buffer_int(:) 1255 real(dp),allocatable :: buffer_dp(:) 1256 1257 ! ************************************************************************* 1258 1259 if (xmpi_comm_size(spaceComm)<=1) return 1260 1261 fourwf=.false.;if (present(profile)) fourwf=(trim(profile)=='fourwf') 1262 full=(.not.fourwf) 1263 1264 !=== Store sizes ==== 1265 if (fourwf.or.full) then 1266 if (allocated(input%kg_k_gather)) size1_kg_k_gather=size(input%kg_k_gather,1) 1267 if (allocated(input%kg_k_gather)) size2_kg_k_gather=size(input%kg_k_gather,2) 1268 if (input%flag1_is_allocated==1) then 1269 if (allocated(input%recvcounts)) size_recvcounts=size(input%recvcounts) 1270 if (allocated(input%sendcounts)) size_sendcounts=size(input%sendcounts) 1271 if (allocated(input%rdispls)) size_rdispls=size(input%rdispls) 1272 if (allocated(input%sdispls)) size_sdispls=size(input%sdispls) 1273 if (allocated(input%gbound)) size1_gbound=size(input%gbound,1) 1274 if (allocated(input%gbound)) size2_gbound=size(input%gbound,2) 1275 end if 1276 if (input%flag3_is_allocated==1.and.input%istwf_k>1) then 1277 if (allocated(input%kg_k_gather_sym)) size1_kg_k_gather_sym=size(input%kg_k_gather_sym,1) 1278 if (allocated(input%kg_k_gather_sym)) size2_kg_k_gather_sym=size(input%kg_k_gather_sym,2) 1279 if (allocated(input%rdispls_sym)) size_rdispls_sym=size(input%rdispls_sym) 1280 if (allocated(input%sdispls_sym)) size_sdispls_sym=size(input%sdispls_sym) 1281 if (allocated(input%recvcounts_sym)) size_recvcounts_sym=size(input%recvcounts_sym) 1282 if (allocated(input%recvcounts_sym_tot)) size_recvcounts_sym_tot=size(input%recvcounts_sym_tot) 1283 if (allocated(input%sendcounts_sym)) size_sendcounts_sym=size(input%sendcounts_sym) 1284 if (allocated(input%sendcounts_sym_all)) size_sendcounts_sym_all=size(input%sendcounts_sym_all) 1285 if (allocated(input%tab_proc)) size_tab_proc=size(input%tab_proc) 1286 end if 1287 if (input%have_to_reequilibrate) then 1288 if (allocated(input%indices_pw_fft)) size_indices_pw_fft=size(input%indices_pw_fft) 1289 if (allocated(input%sendcount_fft)) size_sendcount_fft=size(input%sendcount_fft) 1290 if (allocated(input%senddisp_fft)) size_senddisp_fft=size(input%senddisp_fft) 1291 if (allocated(input%recvcount_fft)) size_recvcount_fft=size(input%recvcount_fft) 1292 if (allocated(input%recvdisp_fft)) size_recvdisp_fft=size(input%recvdisp_fft) 1293 if (allocated(input%kg_k_fft)) size1_kg_k_fft=size(input%kg_k_fft,1) 1294 if (allocated(input%kg_k_fft)) size2_kg_k_fft=size(input%kg_k_fft,2) 1295 end if 1296 end if 1297 if (input%flag2_is_allocated==1.and.full) then 1298 if (allocated(input%ffnl_gather)) size1_ffnl_gather=size(input%ffnl_gather,1) 1299 if (allocated(input%ffnl_gather)) size2_ffnl_gather=size(input%ffnl_gather,2) 1300 if (allocated(input%ffnl_gather)) size3_ffnl_gather=size(input%ffnl_gather,3) 1301 if (allocated(input%ffnl_gather)) size4_ffnl_gather=size(input%ffnl_gather,4) 1302 if (allocated(input%kinpw_gather)) size_kinpw_gather=size(input%kinpw_gather) 1303 if (allocated(input%ph3d_gather)) size1_ph3d_gather=size(input%ph3d_gather,1) 1304 if (allocated(input%ph3d_gather)) size2_ph3d_gather=size(input%ph3d_gather,2) 1305 if (allocated(input%ph3d_gather)) size3_ph3d_gather=size(input%ph3d_gather,3) 1306 if (allocated(input%kpg_k_gather)) size1_kpg_k_gather=size(input%kpg_k_gather,1) 1307 if (allocated(input%kpg_k_gather)) size2_kpg_k_gather=size(input%kpg_k_gather,2) 1308 end if 1309 1310 !=== Compute amount of data to transmit ==== 1311 size_int=size1_kg_k_gather*size2_kg_k_gather & 1312 & +size_recvcounts+size_sendcounts+size_rdispls+size_sdispls & 1313 & +size1_gbound*size2_gbound & 1314 & +size1_kg_k_gather_sym*size2_kg_k_gather_sym & 1315 & +size_rdispls_sym+size_sdispls_sym & 1316 & +size_recvcounts_sym+size_recvcounts_sym_tot & 1317 & +size_sendcounts_sym+size_sendcounts_sym_all+size_tab_proc & 1318 & +size_indices_pw_fft+size_sendcount_fft+size_senddisp_fft & 1319 & +size_recvcount_fft+size_recvdisp_fft & 1320 & +size1_kg_k_fft*size2_kg_k_fft 1321 size_dp=size1_ffnl_gather*size2_ffnl_gather*size3_ffnl_gather*size4_ffnl_gather & 1322 & +size_kinpw_gather & 1323 & +size1_ph3d_gather*size2_ph3d_gather*size3_ph3d_gather & 1324 & +size1_kpg_k_gather*size2_kpg_k_gather 1325 1326 !=== Send flags and array sizes ==== 1327 nsize=45 1328 ABI_ALLOCATE(buffer_int,(nsize)) 1329 buffer_int(1)=input%flag1_is_allocated 1330 buffer_int(2)=input%flag2_is_allocated 1331 buffer_int(3)=input%flag3_is_allocated 1332 buffer_int(4)=0;if (input%have_to_reequilibrate) buffer_int(4)=1 1333 buffer_int(5)=input%istwf_k 1334 buffer_int(6)=input%npw_tot 1335 buffer_int(7)=input%npw_fft 1336 buffer_int(8)=input%ndatarecv 1337 buffer_int(9)=input%idatarecv0 1338 buffer_int(10)=input%ndatarecv_tot 1339 buffer_int(11)=input%ndatasend_sym 1340 buffer_int(12)=size_recvcounts 1341 buffer_int(13)=size_sendcounts 1342 buffer_int(14)=size_rdispls 1343 buffer_int(15)=size_sdispls 1344 buffer_int(16)=size1_gbound 1345 buffer_int(17)=size2_gbound 1346 buffer_int(18)=size1_ffnl_gather 1347 buffer_int(19)=size2_ffnl_gather 1348 buffer_int(20)=size3_ffnl_gather 1349 buffer_int(21)=size4_ffnl_gather 1350 buffer_int(22)=size_kinpw_gather 1351 buffer_int(23)=size1_ph3d_gather 1352 buffer_int(24)=size2_ph3d_gather 1353 buffer_int(25)=size3_ph3d_gather 1354 buffer_int(26)=size1_kpg_k_gather 1355 buffer_int(27)=size2_kpg_k_gather 1356 buffer_int(28)=size_indices_pw_fft 1357 buffer_int(29)=size_sendcount_fft 1358 buffer_int(30)=size_senddisp_fft 1359 buffer_int(31)=size_recvcount_fft 1360 buffer_int(32)=size_recvdisp_fft 1361 buffer_int(33)=size1_kg_k_fft 1362 buffer_int(34)=size2_kg_k_fft 1363 buffer_int(35)=size1_kg_k_gather 1364 buffer_int(36)=size2_kg_k_gather 1365 buffer_int(37)=size1_kg_k_gather_sym 1366 buffer_int(38)=size2_kg_k_gather_sym 1367 buffer_int(39)=size_rdispls_sym 1368 buffer_int(40)=size_sdispls_sym 1369 buffer_int(41)=size_recvcounts_sym 1370 buffer_int(42)=size_recvcounts_sym_tot 1371 buffer_int(43)=size_sendcounts_sym 1372 buffer_int(44)=size_sendcounts_sym_all 1373 buffer_int(45)=size_tab_proc 1374 call xmpi_send(buffer_int,receiver,3*tag-2,spaceComm,ierr) 1375 ABI_DEALLOCATE(buffer_int) 1376 1377 !=== Pack and send integer data ==== 1378 if (size_int>0) then 1379 ipck=0 1380 ABI_ALLOCATE(buffer_int,(size_int)) 1381 if (size1_kg_k_gather*size2_kg_k_gather>0) then 1382 nsize=size1_kg_k_gather*size2_kg_k_gather 1383 buffer_int(ipck+1:ipck+nsize)=reshape(input%kg_k_gather(:,:),(/nsize/)) 1384 ipck=ipck+nsize 1385 end if 1386 if (size_recvcounts>0) then 1387 buffer_int(ipck+1:ipck+size_recvcounts)=input%recvcounts(:) 1388 ipck=ipck+size_recvcounts 1389 end if 1390 if (size_sendcounts>0) then 1391 buffer_int(ipck+1:ipck+size_sendcounts)=input%sendcounts(:) 1392 ipck=ipck+size_sendcounts 1393 end if 1394 if (size_rdispls>0) then 1395 buffer_int(ipck+1:ipck+size_rdispls)=input%rdispls(:) 1396 ipck=ipck+size_rdispls 1397 end if 1398 if (size_sdispls>0) then 1399 buffer_int(ipck+1:ipck+size_sdispls)=input%sdispls(:) 1400 ipck=ipck+size_sdispls 1401 end if 1402 if (size1_gbound*size2_gbound>0) then 1403 nsize=size1_gbound*size2_gbound 1404 buffer_int(ipck+1:ipck+nsize)=reshape(input%gbound(:,:),(/nsize/)) 1405 ipck=ipck+nsize 1406 end if 1407 if (size1_kg_k_gather_sym*size2_kg_k_gather_sym>0) then 1408 nsize=size1_kg_k_gather_sym*size2_kg_k_gather_sym 1409 buffer_int(ipck+1:ipck+nsize)=reshape(input%kg_k_gather_sym(:,:),(/nsize/)) 1410 ipck=ipck+nsize 1411 end if 1412 if (size_rdispls_sym>0) then 1413 buffer_int(ipck+1:ipck+size_rdispls_sym)=input%rdispls_sym(:) 1414 ipck=ipck+size_rdispls_sym 1415 end if 1416 if (size_sdispls_sym>0) then 1417 buffer_int(ipck+1:ipck+size_sdispls_sym)=input%sdispls_sym(:) 1418 ipck=ipck+size_sdispls_sym 1419 end if 1420 if (size_recvcounts_sym>0) then 1421 buffer_int(ipck+1:ipck+size_recvcounts_sym)=input%recvcounts_sym(:) 1422 ipck=ipck+size_recvcounts_sym 1423 end if 1424 if (size_recvcounts_sym_tot>0) then 1425 buffer_int(ipck+1:ipck+size_recvcounts_sym_tot)=input%recvcounts_sym_tot(:) 1426 ipck=ipck+size_recvcounts_sym_tot 1427 end if 1428 if (size_sendcounts_sym>0) then 1429 buffer_int(ipck+1:ipck+size_sendcounts_sym)=input%sendcounts_sym(:) 1430 ipck=ipck+size_sendcounts_sym 1431 end if 1432 if (size_sendcounts_sym_all>0) then 1433 buffer_int(ipck+1:ipck+size_sendcounts_sym_all)=input%sendcounts_sym_all(:) 1434 ipck=ipck+size_sendcounts_sym_all 1435 end if 1436 if (size_tab_proc>0) then 1437 buffer_int(ipck+1:ipck+size_tab_proc)=input%tab_proc(:) 1438 ipck=ipck+size_tab_proc 1439 end if 1440 if (size_indices_pw_fft>0) then 1441 buffer_int(ipck+1:ipck+size_indices_pw_fft)=input%indices_pw_fft(:) 1442 ipck=ipck+size_indices_pw_fft 1443 end if 1444 if (size_sendcount_fft>0) then 1445 buffer_int(ipck+1:ipck+size_sendcount_fft)=input%sendcount_fft(:) 1446 ipck=ipck+size_sendcount_fft 1447 end if 1448 if (size_senddisp_fft>0) then 1449 buffer_int(ipck+1:ipck+size_senddisp_fft)=input%senddisp_fft(:) 1450 ipck=ipck+size_senddisp_fft 1451 end if 1452 if (size_recvcount_fft>0) then 1453 buffer_int(ipck+1:ipck+size_recvcount_fft)=input%recvcount_fft(:) 1454 ipck=ipck+size_recvcount_fft 1455 end if 1456 if (size_recvdisp_fft>0) then 1457 buffer_int(ipck+1:ipck+size_recvdisp_fft)=input%recvdisp_fft(:) 1458 ipck=ipck+size_recvdisp_fft 1459 end if 1460 if (size1_kg_k_fft*size2_kg_k_fft>0) then 1461 nsize=size1_kg_k_fft*size2_kg_k_fft 1462 buffer_int(ipck+1:ipck+nsize)=reshape(input%kg_k_fft(:,:),(/nsize/)) 1463 ipck=ipck+nsize 1464 end if 1465 call xmpi_send(buffer_int,receiver,3*tag-1,spaceComm,ierr) 1466 ABI_DEALLOCATE(buffer_int) 1467 end if 1468 1469 !=== Pack and send real data ==== 1470 if (size_dp>0) then 1471 ipck=0 1472 ABI_ALLOCATE(buffer_dp,(size_dp)) 1473 if (size1_ffnl_gather*size2_ffnl_gather*size3_ffnl_gather*size4_ffnl_gather>0) then 1474 nsize=size1_ffnl_gather*size2_ffnl_gather*size3_ffnl_gather*size4_ffnl_gather 1475 buffer_dp(ipck+1:ipck+nsize)=reshape(input%ffnl_gather(:,:,:,:),(/nsize/)) 1476 ipck=ipck+nsize 1477 end if 1478 if (size_kinpw_gather>0) then 1479 buffer_dp(ipck+1:ipck+size_kinpw_gather)=input%kinpw_gather(:) 1480 ipck=ipck+size_kinpw_gather 1481 end if 1482 if (size1_ph3d_gather*size2_ph3d_gather*size3_ph3d_gather>0) then 1483 nsize=size1_ph3d_gather*size2_ph3d_gather*size3_ph3d_gather 1484 buffer_dp(ipck+1:ipck+nsize)=reshape(input%ph3d_gather(:,:,:),(/nsize/)) 1485 ipck=ipck+nsize 1486 end if 1487 if (size1_kpg_k_gather*size2_kpg_k_gather>0) then 1488 nsize=size1_kpg_k_gather*size2_kpg_k_gather 1489 buffer_dp(ipck+1:ipck+nsize)=reshape(input%kpg_k_gather(:,:),(/nsize/)) 1490 ipck=ipck+nsize 1491 end if 1492 call xmpi_send(buffer_dp,receiver,3*tag,spaceComm,ierr) 1493 ABI_DEALLOCATE(buffer_dp) 1494 end if 1495 1496 end subroutine bandfft_kpt_mpi_send
ABINIT/m_bandfft_kpt [ Modules ]
NAME
m_bandfft_kpt
FUNCTION
This module provides the definition of the bandfft_kpt_type used for kgb parallelization.
COPYRIGHT
Copyright (C) 2011-2018 ABINIT group (FJ, FB, MT) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt .
PARENTS
CHILDREN
SOURCE
22 #if defined HAVE_CONFIG_H 23 #include "config.h" 24 #endif 25 26 #include "abi_common.h" 27 28 MODULE m_bandfft_kpt 29 30 use defs_basis 31 use defs_abitypes 32 use m_abicore 33 use m_errors 34 use m_xmpi 35 36 use m_time, only : timab 37 use m_kg, only : mkkpg 38 use m_fftcore, only : sphereboundary 39 use m_mpinfo, only : proc_distrb_cycle 40 use m_hamiltonian, only : gs_hamiltonian_type 41 42 implicit none 43 44 private 45 public :: bandfft_kpt_init1 46 public :: bandfft_kpt_init2 47 public :: bandfft_kpt_reset 48 public :: bandfft_kpt_destroy 49 public :: bandfft_kpt_destroy_array 50 public :: bandfft_kpt_copy 51 public :: bandfft_kpt_mpi_send 52 public :: bandfft_kpt_mpi_recv 53 public :: bandfft_kpt_savetabs 54 public :: bandfft_kpt_restoretabs 55 public :: bandfft_kpt_set_ikpt 56 public :: bandfft_kpt_get_ikpt 57 public :: prep_bandfft_tabs
ABINIT/prep_bandfft_tabs [ Functions ]
NAME
prep_bandfft_tabs
FUNCTION
This routine transpose various tabs needed in bandfft parallelization
INPUTS
ikpt=index of the k-point gs_hamk <type(gs_hamiltonian_type)>=all data for the Hamiltonian at k mkmem=number of k points which can fit in memory; set to 0 if use disk
SIDE EFFECTS
bandfft_kpt tabs (defined in m_bandfft_kpt module)
PARENTS
energy,fock2ACE,forstrnps,vtorho
CHILDREN
bandfft_kpt_init2,bandfft_kpt_set_ikpt,mkkpg,timab,xmpi_allgatherv
SOURCE
2173 subroutine prep_bandfft_tabs(gs_hamk,ikpt,mkmem,mpi_enreg) 2174 2175 2176 !This section has been created automatically by the script Abilint (TD). 2177 !Do not modify the following lines by hand. 2178 #undef ABI_FUNC 2179 #define ABI_FUNC 'prep_bandfft_tabs' 2180 !End of the abilint section 2181 2182 implicit none 2183 2184 !Arguments ------------------------------- 2185 integer,intent(in) :: ikpt,mkmem 2186 type(gs_hamiltonian_type),intent(inout) :: gs_hamk 2187 type(MPI_type),intent(inout) :: mpi_enreg 2188 2189 !Local variables------------------------------- 2190 integer :: dimffnl,ierr,ikpt_this_proc,ipw,lmnmax,matblk,ndatarecv,nkpg,npw_k,ntypat,spaceComm 2191 logical :: tabs_allocated 2192 real(dp) :: tsec(2) 2193 character(len=500) :: message 2194 integer, allocatable :: recvcounts(:),rdispls(:) 2195 integer, allocatable :: recvcountsloc(:),rdisplsloc(:) 2196 real(dp),allocatable :: ffnl_gather(:,:,:,:),ffnl_little(:,:,:,:),ffnl_little_gather(:,:,:,:) 2197 real(dp),allocatable :: kinpw_gather(:),kpg_k_gather(:,:) 2198 real(dp),allocatable :: ph3d_gather(:,:,:),ph3d_little(:,:,:),ph3d_little_gather(:,:,:) 2199 2200 ! ********************************************************************* 2201 2202 call timab(575,1,tsec) 2203 2204 ikpt_this_proc=mpi_enreg%my_kpttab(ikpt) 2205 tabs_allocated=((bandfft_kpt(ikpt_this_proc)%flag1_is_allocated==1).and.& 2206 & (ikpt_this_proc <= mkmem).and.(ikpt_this_proc/=0)) 2207 2208 if (.not.tabs_allocated) then 2209 message = ' the bandfft tabs are not allocated !' 2210 MSG_BUG(message) 2211 end if 2212 2213 ntypat=gs_hamk%ntypat 2214 lmnmax=gs_hamk%lmnmax 2215 matblk=gs_hamk%matblk 2216 npw_k=gs_hamk%npw_k 2217 dimffnl=size(gs_hamk%ffnl_k,2) 2218 nkpg=size(gs_hamk%kpg_k,2) 2219 2220 ABI_ALLOCATE(rdispls,(mpi_enreg%nproc_band)) 2221 ABI_ALLOCATE(recvcounts,(mpi_enreg%nproc_band)) 2222 ABI_ALLOCATE(rdisplsloc,(mpi_enreg%nproc_band)) 2223 ABI_ALLOCATE(recvcountsloc,(mpi_enreg%nproc_band)) 2224 2225 spaceComm =mpi_enreg%comm_band 2226 ndatarecv =bandfft_kpt(ikpt_this_proc)%ndatarecv 2227 rdispls(:) =bandfft_kpt(ikpt_this_proc)%rdispls(:) 2228 recvcounts(:)=bandfft_kpt(ikpt_this_proc)%recvcounts(:) 2229 2230 !---- Process FFNL 2231 if (associated(gs_hamk%ffnl_k)) then 2232 ABI_ALLOCATE(ffnl_gather,(ndatarecv,dimffnl,lmnmax,ntypat)) 2233 ABI_ALLOCATE(ffnl_little,(dimffnl,lmnmax,ntypat,npw_k)) 2234 ABI_ALLOCATE(ffnl_little_gather,(dimffnl,lmnmax,ntypat,ndatarecv)) 2235 do ipw=1,npw_k 2236 ffnl_little(:,:,:,ipw)=gs_hamk%ffnl_k(ipw,:,:,:) 2237 end do 2238 recvcountsloc(:)=recvcounts(:)*dimffnl*lmnmax*ntypat 2239 rdisplsloc(:)=rdispls(:)*dimffnl*lmnmax*ntypat 2240 call xmpi_allgatherv(ffnl_little,npw_k*dimffnl*lmnmax*ntypat,ffnl_little_gather,& 2241 & recvcountsloc(:),rdisplsloc,spaceComm,ierr) 2242 do ipw=1,ndatarecv 2243 ffnl_gather(ipw,:,:,:)=ffnl_little_gather(:,:,:,ipw) 2244 end do 2245 ABI_DEALLOCATE(ffnl_little) 2246 ABI_DEALLOCATE(ffnl_little_gather) 2247 else 2248 ABI_ALLOCATE(ffnl_gather,(0,0,0,0)) 2249 end if 2250 2251 !---- Process PH3D 2252 if (associated(gs_hamk%ph3d_k)) then 2253 ABI_ALLOCATE(ph3d_gather,(2,ndatarecv,matblk)) 2254 ABI_ALLOCATE(ph3d_little,(2,matblk,npw_k)) 2255 ABI_ALLOCATE(ph3d_little_gather,(2,matblk,ndatarecv)) 2256 recvcountsloc(:)=recvcounts(:)*2*matblk 2257 rdisplsloc(:)=rdispls(:)*2*matblk 2258 do ipw=1,npw_k 2259 ph3d_little(:,:,ipw)=gs_hamk%ph3d_k(:,ipw,:) 2260 end do 2261 call xmpi_allgatherv(ph3d_little,npw_k*2*matblk,ph3d_little_gather,& 2262 & recvcountsloc(:),rdisplsloc,spaceComm,ierr) 2263 do ipw=1,ndatarecv 2264 ph3d_gather(:,ipw,:)=ph3d_little_gather(:,:,ipw) 2265 end do 2266 ABI_DEALLOCATE(ph3d_little_gather) 2267 ABI_DEALLOCATE(ph3d_little) 2268 else 2269 ABI_ALLOCATE(ph3d_gather,(0,0,0)) 2270 end if 2271 2272 !---- Process KPG_K 2273 if (associated(gs_hamk%kpg_k)) then 2274 ABI_ALLOCATE(kpg_k_gather,(ndatarecv,nkpg)) 2275 if (nkpg>0) then 2276 call mkkpg(bandfft_kpt(ikpt_this_proc)%kg_k_gather,kpg_k_gather,gs_hamk%kpt_k,nkpg,ndatarecv) 2277 end if 2278 else 2279 ABI_ALLOCATE(kpg_k_gather,(0,0)) 2280 end if 2281 2282 !---- Process KINPW 2283 if (associated(gs_hamk%kinpw_k)) then 2284 ABI_ALLOCATE(kinpw_gather,(ndatarecv)) 2285 recvcountsloc(:)=recvcounts(:) 2286 rdisplsloc(:)=rdispls(:) 2287 call xmpi_allgatherv(gs_hamk%kinpw_k,npw_k,kinpw_gather,recvcountsloc(:),rdisplsloc,spaceComm,ierr) 2288 else 2289 ABI_ALLOCATE(kinpw_gather,(0)) 2290 end if 2291 2292 ABI_DEALLOCATE(recvcounts) 2293 ABI_DEALLOCATE(rdispls) 2294 ABI_DEALLOCATE(recvcountsloc) 2295 ABI_DEALLOCATE(rdisplsloc) 2296 2297 call bandfft_kpt_init2(bandfft_kpt,dimffnl,ffnl_gather,ikpt_this_proc,kinpw_gather,kpg_k_gather,& 2298 & lmnmax,matblk,mkmem,ndatarecv,nkpg,ntypat,ph3d_gather) 2299 2300 ABI_DEALLOCATE(ffnl_gather) 2301 ABI_DEALLOCATE(ph3d_gather) 2302 ABI_DEALLOCATE(kinpw_gather) 2303 ABI_DEALLOCATE(kpg_k_gather) 2304 2305 !---- Store current kpt index 2306 call bandfft_kpt_set_ikpt(ikpt,mpi_enreg) 2307 2308 call timab(575,2,tsec) 2309 2310 end subroutine prep_bandfft_tabs
m_bandfft_kpt/bandfft_kpt_copy [ Functions ]
[ Top ] [ m_bandfft_kpt ] [ Functions ]
NAME
bandfft_kpt_copy
FUNCTION
Copy a bandfft_kpt datastructure into another
INPUTS
bandfft_kpt_in=<type(bandfft_kpt_type)>=input bandfft_kpt datastructure
OUTPUT
bandfft_kpt_out=<type(bandfft_kpt_type)>=output bandfft_kpt datastructure
PARENTS
CHILDREN
SOURCE
1000 subroutine bandfft_kpt_copy(bandfft_kpt_in,bandfft_kpt_out,mpi_enreg1,opt_bandfft) 1001 1002 1003 !This section has been created automatically by the script Abilint (TD). 1004 !Do not modify the following lines by hand. 1005 #undef ABI_FUNC 1006 #define ABI_FUNC 'bandfft_kpt_copy' 1007 !End of the abilint section 1008 1009 implicit none 1010 1011 !Arguments ------------------------------------ 1012 !scalars 1013 integer,intent(in) :: opt_bandfft 1014 type(bandfft_kpt_type),pointer :: bandfft_kpt_in(:) 1015 type(bandfft_kpt_type),pointer :: bandfft_kpt_out(:) 1016 type(MPI_type),intent(inout) :: mpi_enreg1 1017 1018 !Local variables------------------------------- 1019 !scalars 1020 integer :: ikpt,isppol,jkpt,sz1,sz2,sz3,sz4 1021 1022 ! ********************************************************************* 1023 1024 !Optional pointers 1025 if (opt_bandfft==0) then 1026 nullify(bandfft_kpt_out) 1027 else if (opt_bandfft==1) then 1028 if (associated(bandfft_kpt_in)) then 1029 ABI_DATATYPE_ALLOCATE(bandfft_kpt_out,(size(bandfft_kpt_in))) 1030 do isppol=1,size(mpi_enreg1%proc_distrb,3) 1031 do ikpt=1,size(mpi_enreg1%proc_distrb,1) 1032 sz1=size(mpi_enreg1%proc_distrb,2) 1033 if(proc_distrb_cycle(mpi_enreg1%proc_distrb,ikpt,1,sz1,isppol,mpi_enreg1%me_kpt)) then 1034 cycle 1035 end if 1036 jkpt=mpi_enreg1%my_kpttab(ikpt) 1037 ! if (allocated(bandfft_kpt_in(jkpt)%ind_kg_mpi_to_seq)) then 1038 ! sz1=size(bandfft_kpt_in(jkpt)%ind_kg_mpi_to_seq) 1039 ! ABI_ALLOCATE(bandfft_kpt_out(jkpt)%ind_kg_mpi_to_seq,(sz1)) 1040 ! bandfft_kpt_out(jkpt)%ind_kg_mpi_to_seq= & 1041 ! & bandfft_kpt_in(jkpt)%ind_kg_mpi_to_seq 1042 ! end if 1043 if (allocated(bandfft_kpt_in(jkpt)%kg_k_gather)) then 1044 sz1=size(bandfft_kpt_in(jkpt)%kg_k_gather,1) 1045 sz2=size(bandfft_kpt_in(jkpt)%kg_k_gather,2) 1046 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%kg_k_gather,(sz1,sz2)) 1047 bandfft_kpt_out(jkpt)%kg_k_gather= & 1048 & bandfft_kpt_in(jkpt)%kg_k_gather 1049 end if 1050 bandfft_kpt_out(jkpt)%flag1_is_allocated=bandfft_kpt_in(jkpt)%flag1_is_allocated 1051 if (allocated(bandfft_kpt_in(jkpt)%gbound)) then 1052 sz1=size(bandfft_kpt_in(jkpt)%gbound,1) 1053 sz2=size(bandfft_kpt_in(jkpt)%gbound,2) 1054 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%gbound,(sz1,sz2)) 1055 bandfft_kpt_out(jkpt)%gbound= & 1056 & bandfft_kpt_in(jkpt)%gbound 1057 end if 1058 if (allocated(bandfft_kpt_in(jkpt)%recvcounts)) then 1059 sz1=size(bandfft_kpt_in(jkpt)%recvcounts) 1060 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%recvcounts,(sz1)) 1061 bandfft_kpt_out(jkpt)%recvcounts= & 1062 & bandfft_kpt_in(jkpt)%recvcounts 1063 end if 1064 if (allocated(bandfft_kpt_in(jkpt)%sendcounts)) then 1065 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%sendcounts,(size(bandfft_kpt_in(jkpt)%sendcounts))) 1066 bandfft_kpt_out(jkpt)%sendcounts= & 1067 & bandfft_kpt_in(jkpt)%sendcounts 1068 end if 1069 if (allocated(bandfft_kpt_in(jkpt)%rdispls)) then 1070 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%rdispls,(size(bandfft_kpt_in(jkpt)%rdispls))) 1071 bandfft_kpt_out(jkpt)%rdispls= & 1072 & bandfft_kpt_in(jkpt)%rdispls 1073 end if 1074 if (allocated(bandfft_kpt_in(jkpt)%sdispls)) then 1075 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%sdispls,(size(bandfft_kpt_in(jkpt)%sdispls))) 1076 bandfft_kpt_out(jkpt)%sdispls= & 1077 & bandfft_kpt_in(jkpt)%sdispls 1078 end if 1079 bandfft_kpt_out(jkpt)%flag2_is_allocated=bandfft_kpt_in(jkpt)%flag2_is_allocated 1080 if (allocated(bandfft_kpt_in(jkpt)%ffnl_gather)) then 1081 sz1=size(bandfft_kpt_in(jkpt)%ffnl_gather,1) 1082 sz2=size(bandfft_kpt_in(jkpt)%ffnl_gather,2) 1083 sz3=size(bandfft_kpt_in(jkpt)%ffnl_gather,3) 1084 sz4=size(bandfft_kpt_in(jkpt)%ffnl_gather,4) 1085 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%ffnl_gather,(sz1,sz2,sz3,sz4)) 1086 bandfft_kpt_out(jkpt)%ffnl_gather= & 1087 & bandfft_kpt_in(jkpt)%ffnl_gather 1088 end if 1089 if (allocated(bandfft_kpt_in(jkpt)%kinpw_gather)) then 1090 sz1=size(bandfft_kpt_in(jkpt)%kinpw_gather) 1091 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%kinpw_gather,(sz1)) 1092 bandfft_kpt_out(jkpt)%kinpw_gather= & 1093 & bandfft_kpt_in(jkpt)%kinpw_gather 1094 end if 1095 if (allocated(bandfft_kpt_in(jkpt)%ph3d_gather)) then 1096 sz1=size(bandfft_kpt_in(jkpt)%ph3d_gather,1) 1097 sz2=size(bandfft_kpt_in(jkpt)%ph3d_gather,2) 1098 sz3=size(bandfft_kpt_in(jkpt)%ph3d_gather,3) 1099 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%ph3d_gather,(sz1,sz2,sz3)) 1100 bandfft_kpt_out(jkpt)%ph3d_gather= & 1101 & bandfft_kpt_in(jkpt)%ph3d_gather 1102 end if 1103 if (allocated(bandfft_kpt_in(jkpt)%kpg_k_gather)) then 1104 sz1=size(bandfft_kpt_in(jkpt)%kpg_k_gather,1) 1105 sz2=size(bandfft_kpt_in(jkpt)%kpg_k_gather,2) 1106 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%kpg_k_gather,(sz1,sz2)) 1107 bandfft_kpt_out(jkpt)%kpg_k_gather= & 1108 & bandfft_kpt_in(jkpt)%kpg_k_gather 1109 end if 1110 bandfft_kpt_out(jkpt)%flag3_is_allocated=bandfft_kpt_in(jkpt)%flag3_is_allocated 1111 if (allocated(bandfft_kpt_in(jkpt)%kg_k_gather_sym)) then 1112 sz1=size(bandfft_kpt_in(jkpt)%kg_k_gather_sym,1) 1113 sz2=size(bandfft_kpt_in(jkpt)%kg_k_gather_sym,2) 1114 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%kg_k_gather_sym,(sz1,sz2)) 1115 bandfft_kpt_out(jkpt)%kg_k_gather_sym= & 1116 & bandfft_kpt_in(jkpt)%kg_k_gather_sym 1117 end if 1118 if (allocated(bandfft_kpt_in(jkpt)%rdispls_sym)) then 1119 sz1=size(bandfft_kpt_in(jkpt)%rdispls_sym) 1120 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%rdispls_sym,(sz1)) 1121 bandfft_kpt_out(jkpt)%rdispls_sym= & 1122 & bandfft_kpt_in(jkpt)%rdispls_sym 1123 end if 1124 if (allocated(bandfft_kpt_in(jkpt)%recvcounts_sym)) then 1125 sz1=size(bandfft_kpt_in(jkpt)%recvcounts_sym) 1126 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%recvcounts_sym,(sz1)) 1127 bandfft_kpt_out(jkpt)%recvcounts_sym= & 1128 & bandfft_kpt_in(jkpt)%recvcounts_sym 1129 end if 1130 if (allocated(bandfft_kpt_in(jkpt)%recvcounts_sym_tot)) then 1131 sz1=size(bandfft_kpt_in(jkpt)%recvcounts_sym_tot) 1132 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%recvcounts_sym_tot,(sz1)) 1133 bandfft_kpt_out(jkpt)%recvcounts_sym_tot= & 1134 & bandfft_kpt_in(jkpt)%recvcounts_sym_tot 1135 end if 1136 if (allocated(bandfft_kpt_in(jkpt)%sdispls_sym)) then 1137 sz1=size(bandfft_kpt_in(jkpt)%sdispls_sym) 1138 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%sdispls_sym,(sz1)) 1139 bandfft_kpt_out(jkpt)%sdispls_sym= & 1140 & bandfft_kpt_in(jkpt)%sdispls_sym 1141 end if 1142 if (allocated(bandfft_kpt_in(jkpt)%sendcounts_sym)) then 1143 sz1=size(bandfft_kpt_in(jkpt)%sendcounts_sym) 1144 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%sendcounts_sym,(sz1)) 1145 bandfft_kpt_out(jkpt)%sendcounts_sym= & 1146 & bandfft_kpt_in(jkpt)%sendcounts_sym 1147 end if 1148 if (allocated(bandfft_kpt_in(jkpt)%sendcounts_sym_all)) then 1149 sz1=size(bandfft_kpt_in(jkpt)%sendcounts_sym_all) 1150 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%sendcounts_sym_all,(sz1)) 1151 bandfft_kpt_out(jkpt)%sendcounts_sym_all= & 1152 & bandfft_kpt_in(jkpt)%sendcounts_sym_all 1153 end if 1154 if (allocated(bandfft_kpt_in(jkpt)%tab_proc)) then 1155 sz1=size(bandfft_kpt_in(jkpt)%tab_proc) 1156 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%tab_proc,(sz1)) 1157 bandfft_kpt_out(jkpt)%tab_proc= & 1158 bandfft_kpt_in(jkpt)%tab_proc 1159 end if 1160 if (bandfft_kpt_in(jkpt)%have_to_reequilibrate) then 1161 bandfft_kpt_out(jkpt)%have_to_reequilibrate = .true. 1162 bandfft_kpt_out(jkpt)%npw_fft = bandfft_kpt_in(jkpt)%npw_fft 1163 sz1=size(bandfft_kpt_in(jkpt)%indices_pw_fft) 1164 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%indices_pw_fft,(sz1)) 1165 bandfft_kpt_out(jkpt)%indices_pw_fft=bandfft_kpt_in(jkpt)%indices_pw_fft 1166 sz1=size(bandfft_kpt_in(jkpt)%sendcount_fft) 1167 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%sendcount_fft,(sz1)) 1168 bandfft_kpt_out(jkpt)%sendcount_fft=bandfft_kpt_in(jkpt)%sendcount_fft 1169 sz1=size(bandfft_kpt_in(jkpt)%senddisp_fft) 1170 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%senddisp_fft,(sz1)) 1171 bandfft_kpt_out(jkpt)%senddisp_fft=bandfft_kpt_in(jkpt)%senddisp_fft 1172 sz1=size(bandfft_kpt_in(jkpt)%recvcount_fft) 1173 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%recvcount_fft,(sz1)) 1174 bandfft_kpt_out(jkpt)%recvcount_fft=bandfft_kpt_in(jkpt)%recvcount_fft 1175 sz1=size(bandfft_kpt_in(jkpt)%recvdisp_fft) 1176 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%recvdisp_fft,(sz1)) 1177 bandfft_kpt_out(jkpt)%recvdisp_fft=bandfft_kpt_in(jkpt)%recvdisp_fft 1178 sz1=size(bandfft_kpt_in(jkpt)%kg_k_fft,1) 1179 sz2=size(bandfft_kpt_in(jkpt)%kg_k_fft,2) 1180 ABI_ALLOCATE(bandfft_kpt_out(jkpt)%kg_k_fft,(sz1,sz2)) 1181 bandfft_kpt_out(jkpt)%kg_k_fft=bandfft_kpt_in(jkpt)%kg_k_fft 1182 end if 1183 end do 1184 end do 1185 end if 1186 end if 1187 1188 end subroutine bandfft_kpt_copy
m_bandfft_kpt/bandfft_kpt_destroy [ Functions ]
[ Top ] [ m_bandfft_kpt ] [ Functions ]
NAME
bandfft_kpt_destroy
FUNCTION
Destroy a bandfft_kpt datastructure
INPUTS
OUTPUT
SIDE EFFECTS
bandfft_kpt_in=the datastructure to destroy
PARENTS
m_bandfft_kpt,posdoppler
CHILDREN
SOURCE
802 subroutine bandfft_kpt_destroy(bandfft_kpt_in) 803 804 805 !This section has been created automatically by the script Abilint (TD). 806 !Do not modify the following lines by hand. 807 #undef ABI_FUNC 808 #define ABI_FUNC 'bandfft_kpt_destroy' 809 !End of the abilint section 810 811 implicit none 812 813 !Arguments ------------------------------------ 814 type(bandfft_kpt_type) :: bandfft_kpt_in 815 !Local variables------------------------------- 816 817 ! *********************************************************************** 818 819 bandfft_kpt_in%flag1_is_allocated=0 820 bandfft_kpt_in%flag2_is_allocated=0 821 bandfft_kpt_in%flag3_is_allocated=0 822 bandfft_kpt_in%have_to_reequilibrate=.false. 823 824 if (allocated(bandfft_kpt_in%kg_k_gather)) then 825 ABI_DEALLOCATE(bandfft_kpt_in%kg_k_gather) 826 end if 827 if (allocated(bandfft_kpt_in%gbound)) then 828 ABI_DEALLOCATE(bandfft_kpt_in%gbound) 829 end if 830 if (allocated(bandfft_kpt_in%recvcounts)) then 831 ABI_DEALLOCATE(bandfft_kpt_in%recvcounts) 832 end if 833 if (allocated(bandfft_kpt_in%sendcounts)) then 834 ABI_DEALLOCATE(bandfft_kpt_in%sendcounts) 835 end if 836 if (allocated(bandfft_kpt_in%rdispls)) then 837 ABI_DEALLOCATE(bandfft_kpt_in%rdispls) 838 end if 839 if (allocated(bandfft_kpt_in%sdispls)) then 840 ABI_DEALLOCATE(bandfft_kpt_in%sdispls) 841 end if 842 if (allocated(bandfft_kpt_in%ffnl_gather)) then 843 ABI_DEALLOCATE(bandfft_kpt_in%ffnl_gather) 844 end if 845 if (allocated(bandfft_kpt_in%kinpw_gather)) then 846 ABI_DEALLOCATE(bandfft_kpt_in%kinpw_gather) 847 end if 848 if (allocated(bandfft_kpt_in%kpg_k_gather)) then 849 ABI_DEALLOCATE(bandfft_kpt_in%kpg_k_gather) 850 end if 851 if (allocated(bandfft_kpt_in%ph3d_gather)) then 852 ABI_DEALLOCATE(bandfft_kpt_in%ph3d_gather) 853 end if 854 if (allocated(bandfft_kpt_in%kg_k_gather_sym)) then 855 ABI_DEALLOCATE(bandfft_kpt_in%kg_k_gather_sym) 856 end if 857 if (allocated(bandfft_kpt_in%rdispls_sym)) then 858 ABI_DEALLOCATE(bandfft_kpt_in%rdispls_sym) 859 end if 860 if (allocated(bandfft_kpt_in%recvcounts_sym)) then 861 ABI_DEALLOCATE(bandfft_kpt_in%recvcounts_sym) 862 end if 863 if (allocated(bandfft_kpt_in%recvcounts_sym_tot)) then 864 ABI_DEALLOCATE(bandfft_kpt_in%recvcounts_sym_tot) 865 end if 866 if (allocated(bandfft_kpt_in%sdispls_sym)) then 867 ABI_DEALLOCATE(bandfft_kpt_in%sdispls_sym) 868 end if 869 if (allocated(bandfft_kpt_in%sendcounts_sym)) then 870 ABI_DEALLOCATE(bandfft_kpt_in%sendcounts_sym) 871 end if 872 if (allocated(bandfft_kpt_in%sendcounts_sym_all)) then 873 ABI_DEALLOCATE(bandfft_kpt_in%sendcounts_sym_all) 874 end if 875 if (allocated(bandfft_kpt_in%tab_proc)) then 876 ABI_DEALLOCATE(bandfft_kpt_in%tab_proc) 877 end if 878 if (allocated(bandfft_kpt_in%indices_pw_fft)) then 879 ABI_DEALLOCATE(bandfft_kpt_in%indices_pw_fft) 880 end if 881 if (allocated(bandfft_kpt_in%sendcount_fft)) then 882 ABI_DEALLOCATE(bandfft_kpt_in%sendcount_fft) 883 end if 884 if (allocated(bandfft_kpt_in%senddisp_fft)) then 885 ABI_DEALLOCATE(bandfft_kpt_in%senddisp_fft) 886 end if 887 if (allocated(bandfft_kpt_in%recvcount_fft)) then 888 ABI_DEALLOCATE(bandfft_kpt_in%recvcount_fft) 889 end if 890 if (allocated(bandfft_kpt_in%recvdisp_fft)) then 891 ABI_DEALLOCATE(bandfft_kpt_in%recvdisp_fft) 892 end if 893 if (allocated(bandfft_kpt_in%kg_k_fft)) then 894 ABI_DEALLOCATE(bandfft_kpt_in%kg_k_fft) 895 end if 896 897 end subroutine bandfft_kpt_destroy
m_bandfft_kpt/bandfft_kpt_destroy_array [ Functions ]
[ Top ] [ m_bandfft_kpt ] [ Functions ]
NAME
bandfft_kpt_destroy_array
FUNCTION
Clean and destroy an array of bandfft_kpt datastructures
INPUTS
OUTPUT
SIDE EFFECTS
bandfft_kpt_in(:)=the array of datastructure to destroy
PARENTS
gstate,gwls_hamiltonian
CHILDREN
SOURCE
924 subroutine bandfft_kpt_destroy_array(bandfft_kpt_in,mpi_enreg) 925 926 927 !This section has been created automatically by the script Abilint (TD). 928 !Do not modify the following lines by hand. 929 #undef ABI_FUNC 930 #define ABI_FUNC 'bandfft_kpt_destroy_array' 931 !End of the abilint section 932 933 implicit none 934 935 !Arguments ------------------------------------ 936 type(bandfft_kpt_type),pointer :: bandfft_kpt_in(:) 937 type(MPI_type), intent(inout) :: mpi_enreg 938 939 !Local variables------------------------------- 940 integer :: ikpt_this_proc,isppol,ikpt,mkmem,nband,nkpt,nsppol 941 character(len=500) :: msg 942 943 ! *********************************************************************** 944 945 if (xmpi_paral==0) return 946 947 if (associated(bandfft_kpt_in)) then 948 mkmem =size(bandfft_kpt_in) 949 nkpt=size(mpi_enreg%proc_distrb,1) 950 nband=size(mpi_enreg%proc_distrb,2) 951 nsppol=size(mpi_enreg%proc_distrb,3) 952 if (nsppol==0.or.nkpt==0) then 953 msg=' mpi_enreg%proc_distrb should be allocated !' 954 MSG_BUG(msg) 955 end if 956 nkpt=size(mpi_enreg%my_kpttab) 957 if (nkpt==0) then 958 msg=' mpi_enreg%my_kpttab should be allocated !' 959 MSG_BUG(msg) 960 end if 961 do isppol=1,nsppol 962 do ikpt=1,nkpt 963 if(proc_distrb_cycle(mpi_enreg%proc_distrb,ikpt,1,nband,isppol,mpi_enreg%me_kpt)) cycle 964 ikpt_this_proc=mpi_enreg%my_kpttab(ikpt) 965 if ((ikpt_this_proc>mkmem) .or.(ikpt_this_proc<=0)) then 966 msg=' The bandfft tab cannot be deallocated !' 967 MSG_BUG(msg) 968 end if 969 call bandfft_kpt_destroy(bandfft_kpt_in(ikpt_this_proc)) 970 end do 971 end do 972 ABI_DATATYPE_DEALLOCATE(bandfft_kpt_in) 973 nullify(bandfft_kpt_in) 974 end if 975 976 end subroutine bandfft_kpt_destroy_array
m_bandfft_kpt/bandfft_kpt_get_ikpt [ Functions ]
[ Top ] [ m_bandfft_kpt ] [ Functions ]
NAME
bandfft_kpt_get_ikpt
FUNCTION
Get the value of bandfft_kpt_current_ikpt variable, i.e. current index of bandfft_kpt loaded in memory INPUT
OUTPUT
bandfft_kpt_get_ikpt= current index of bandfft_kpt
PARENTS
apply_invovl,chebfi,getghc,make_invovl,prep_fourwf prep_getghc,prep_nonlop
CHILDREN
SOURCE
2127 function bandfft_kpt_get_ikpt() 2128 2129 2130 !This section has been created automatically by the script Abilint (TD). 2131 !Do not modify the following lines by hand. 2132 #undef ABI_FUNC 2133 #define ABI_FUNC 'bandfft_kpt_get_ikpt' 2134 !End of the abilint section 2135 2136 implicit none 2137 2138 !Arguments ------------------------------- 2139 integer :: bandfft_kpt_get_ikpt 2140 !Local variables------------------------------- 2141 2142 ! ********************************************************************* 2143 2144 bandfft_kpt_get_ikpt=bandfft_kpt_current_ikpt 2145 2146 end function bandfft_kpt_get_ikpt
m_bandfft_kpt/bandfft_kpt_init1 [ Functions ]
[ Top ] [ m_bandfft_kpt ] [ Functions ]
NAME
bandfft_kpt_init1
FUNCTION
Init all (or part of) scalars and pointers in a bandfft_kpt datastructure
INPUTS
istwfk(nkpt) = input option parameter that describes the storage of wfs kg(3,mpw*mkmem) = dimensionless coords of G vecs in basis sphere at k point mgfft = maximum single fft dimension (IN) mkmem = number of k points which can fit in memory; set to 0 if use disk mpi_enreg = information about MPI parallelization mpw = maximum number of planewaves as dimensioned in calling routine nband(nkpt*nsppol) = number of bands at each k point nkpt = number of k points npwarr(nkpt) = array holding npw for each k point, taking into account the effect of istwfk, and the spreading over processors nsppol = 1 for unpolarized, 2 for polarized
OUTPUT
within the bandfft_kpt data_type : Initialize and compute gbound = sphere boundary info idatarecv0 = position of the planewave coordinates (0,0,0) istwf_k = input option parameter that describes the storage of wfs kg_k_gather = planewave coordinates (of the processor + sended by other processors band) kg_k_gather_sym = planewave coordinates (kg_k_gather + opposited planewave coordinates sended by the processors fft) ndatarecv = total number of values received by the processor and sended by the other processors band ndatasend_sym = number of sended values to the processors fft to create opposited planewave coordinates ndatarecv_tot = total number of received values by the processor (ndatarecv + number of received opposited planewave coordinates) recvcounts = number of values received by the processor from each processor band recvcounts_sym = number of values received by the processor from each processor fft recvcounts_sym_tot = number of values received by each processor from the other processors fft rdispls = positions of values received by the processor from each processor band rdispls_sym = positions of values received by the processor from each processor fft sendcounts = number of values sended by the processor to each processor band sendcounts_sym = number of values sended by the processor to each processor fft sendcounts_sym_all = number of values sended by each processor to the other processors fft sdispls = postions of values sended by the processor to each processor band sdispls_sym = postions of values sended by the processor to each processor fft tab_proc = positions of opposited planewave coordinates in the list of the processors fft
SIDE EFFECTS
bandfft_kpt_in=<type(bandfft_kpt)>=bandfft_kpt datastructure
PARENTS
gstate
CHILDREN
SOURCE
204 subroutine bandfft_kpt_init1(bandfft_kpt_in,istwfk,kg,mgfft,mkmem,mpi_enreg,mpw,nband,nkpt,npwarr,nsppol) 205 206 207 !This section has been created automatically by the script Abilint (TD). 208 !Do not modify the following lines by hand. 209 #undef ABI_FUNC 210 #define ABI_FUNC 'bandfft_kpt_init1' 211 !End of the abilint section 212 213 implicit none 214 215 !Arguments ------------------------------------ 216 !scalars 217 integer,intent(in) :: mgfft,mkmem,mpw,nkpt,nsppol 218 type(bandfft_kpt_type),pointer :: bandfft_kpt_in(:) 219 type(MPI_type),intent(inout) :: mpi_enreg 220 !arrays 221 integer,intent(in) :: istwfk(nkpt),nband(nkpt*nsppol) 222 integer,intent(in) :: kg(3,mpw*mkmem),npwarr(nkpt) 223 224 !Local variables------------------------------- 225 !scalars 226 integer :: comm_band,comm_fft,idatarecv,idatarecv0,ierr,ikg,ikpt,ikpt_this_proc,iproc,isppol,istwf_k,itest,jsendloc 227 integer :: me_fft,me_kpt,n2,nband_k,ndatarecv,ndatarecv_tot,ndatasend_sym,nproc_band,nproc_fft,npw_fft,npw_k,npw_tot 228 logical :: reequilibrate_allocated 229 character(len=500) :: message 230 !arrays 231 integer,allocatable :: buff_kg(:,:),gbound(:,:),indices_pw_fft(:),kg_k(:,:),kg_k_fft(:,:),kg_k_gather(:,:) 232 integer,allocatable :: kg_k_gather_all(:,:),kg_k_gather_send(:,:),kg_k_gather_sym(:,:) 233 integer,allocatable :: npw_per_proc(:) 234 integer,allocatable :: rdispls(:),rdispls_all(:),rdispls_sym(:),rdispls_sym_loc(:) 235 integer,allocatable :: recvcounts(:),recvcounts_sym(:),recvcounts_fft(:),recvcounts_sym_loc(:),recvcounts_sym_tot(:) 236 integer,allocatable :: recvdisp_fft(:),sdispls(:),sdispls_sym_loc(:),sdispls_sym(:) 237 integer,allocatable :: sendcounts(:),sendcounts_fft(:),sendcounts_sym(:),sendcounts_sym_all(:),sendcounts_sym_loc(:) 238 integer,allocatable :: senddisp_fft(:),sum_kg(:),tab_proc(:) 239 240 ! ********************************************************************* 241 242 If(mpi_enreg%paral_kgb/=1) then 243 nullify(bandfft_kpt_in) 244 return 245 end if 246 247 !--------------------------------------------- 248 !Initialisation 249 !--------------------------------------------- 250 nproc_fft = mpi_enreg%nproc_fft 251 nproc_band = mpi_enreg%nproc_band 252 253 me_fft = mpi_enreg%me_fft 254 me_kpt = mpi_enreg%me_kpt 255 256 comm_band = mpi_enreg%comm_band 257 comm_fft = mpi_enreg%comm_fft 258 259 !============================================================================= 260 !Compute and store various tabs in bandfft_kpt(ikpt) data_struc 261 !These ones will be used in following subroutines: 262 !vtorho, mkrho, prep_nonlop, prep_fourwf, prep_getghc... 263 !============================================================================= 264 265 ABI_DATATYPE_ALLOCATE(bandfft_kpt_in,(mkmem)) 266 267 ABI_ALLOCATE(sdispls ,(nproc_band)) 268 ABI_ALLOCATE(sendcounts ,(nproc_band)) 269 ABI_ALLOCATE(rdispls ,(nproc_band)) 270 ABI_ALLOCATE(recvcounts ,(nproc_band)) 271 ABI_ALLOCATE(sendcounts_fft,(nproc_fft)) 272 ABI_ALLOCATE(senddisp_fft ,(nproc_fft)) 273 ABI_ALLOCATE(recvcounts_fft,(nproc_fft)) 274 ABI_ALLOCATE(recvdisp_fft ,(nproc_fft)) 275 276 bandfft_kpt_in(:)%flag1_is_allocated=0 277 bandfft_kpt_in(:)%flag2_is_allocated=0 278 bandfft_kpt_in(:)%flag3_is_allocated=0 279 280 do isppol=1,nsppol 281 ikg=0 282 do ikpt=1,nkpt 283 npw_k=npwarr(ikpt) 284 istwf_k=istwfk(ikpt) 285 nband_k=nband(ikpt+(isppol-1)*nkpt) 286 if(proc_distrb_cycle(mpi_enreg%proc_distrb,ikpt,1,nband_k,isppol,me_kpt))cycle 287 ikpt_this_proc=mpi_enreg%my_kpttab(ikpt) 288 if ((ikpt_this_proc > mkmem).or.(ikpt_this_proc==0)) then 289 message = ' this bandfft tab is not allocated !' 290 MSG_ERROR(message) 291 end if 292 if (bandfft_kpt_in(ikpt_this_proc)%flag1_is_allocated==0) then 293 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%gbound ,(2*mgfft+8,2)) 294 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%recvcounts,(nproc_band)) 295 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%sendcounts,(nproc_band)) 296 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%rdispls ,(nproc_band)) 297 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%sdispls ,(nproc_band)) 298 bandfft_kpt_in(ikpt_this_proc)%flag1_is_allocated=1 299 end if 300 301 ! Initialize various quantities for the reequilibration step 302 reequilibrate_allocated=(isppol==2.and.bandfft_kpt_in(ikpt_this_proc)%have_to_reequilibrate) 303 bandfft_kpt_in(ikpt_this_proc)%have_to_reequilibrate = .false. 304 bandfft_kpt_in(ikpt_this_proc)%npw_fft = 0 305 306 call xmpi_allgather(npw_k,recvcounts,comm_band,ierr) 307 rdispls(1)=0 308 do iproc=2,nproc_band 309 rdispls(iproc)=rdispls(iproc-1)+recvcounts(iproc-1) 310 end do 311 ndatarecv=rdispls(nproc_band)+recvcounts(nproc_band) 312 313 ABI_ALLOCATE(kg_k_gather,(3,ndatarecv)) 314 ABI_ALLOCATE(kg_k,(3,npw_k)) 315 kg_k(:,1:npw_k)=kg(:,1+ikg:npw_k+ikg) 316 call xmpi_allgatherv(kg_k,3*npw_k,kg_k_gather,3*recvcounts(:),3*rdispls(:),comm_band,ierr) 317 318 sendcounts(:)=npw_k*mpi_enreg%bandpp 319 do iproc=1,nproc_band 320 sdispls(iproc)=(iproc-1)*npw_k*mpi_enreg%bandpp 321 end do 322 323 ! ============================================================================ 324 ! Here we compute gbound, as well for istwf_k=1 as for istwf_k=2 and store it 325 ! ============================================================================ 326 !MG: Why don't we compute gbound locally by just computing the full G-sphere from ecut 327 ABI_ALLOCATE(npw_per_proc,(nproc_fft)) 328 ABI_ALLOCATE(rdispls_all,(nproc_fft)) 329 ABI_ALLOCATE(gbound,(2*mgfft+8,2)) 330 if (mgfft>0) gbound(:,:)=0 331 if (istwf_k==1) then 332 call xmpi_allgather(ndatarecv,npw_per_proc,mpi_enreg%comm_fft,ierr) 333 rdispls_all(1)=0 334 do iproc=2,nproc_fft 335 rdispls_all(iproc)=rdispls_all(iproc-1)+npw_per_proc(iproc-1) 336 end do 337 npw_tot=rdispls_all(nproc_fft)+npw_per_proc(nproc_fft) 338 ABI_ALLOCATE(kg_k_gather_all,(3,npw_tot)) 339 call xmpi_allgatherv(kg_k_gather,& 340 & 3*ndatarecv,kg_k_gather_all,3*npw_per_proc(:),3*rdispls_all,mpi_enreg%comm_fft,ierr) 341 if (mgfft>0) then 342 call sphereboundary(gbound,istwf_k,kg_k_gather_all,mgfft,npw_tot) 343 end if 344 345 else if (istwf_k==2) then 346 347 ! ============================================================================ 348 ! In this case, we have to add the opposite values in the kg_k_gather tab 349 ! before computing gbound 350 ! ============================================================================ 351 352 ! Allocation 353 ABI_ALLOCATE(tab_proc ,(ndatarecv)) 354 ABI_ALLOCATE(sendcounts_sym ,(nproc_fft)) 355 ABI_ALLOCATE(sendcounts_sym_all,(nproc_fft*nproc_fft)) 356 ABI_ALLOCATE(sdispls_sym ,(nproc_fft)) 357 ABI_ALLOCATE(recvcounts_sym ,(nproc_fft)) 358 ABI_ALLOCATE(recvcounts_sym_tot,(nproc_fft)) 359 ABI_ALLOCATE(rdispls_sym ,(nproc_fft)) 360 ABI_ALLOCATE(sendcounts_sym_loc,(nproc_fft)) 361 ABI_ALLOCATE(sdispls_sym_loc ,(nproc_fft)) 362 ABI_ALLOCATE(recvcounts_sym_loc,(nproc_fft)) 363 ABI_ALLOCATE(rdispls_sym_loc ,(nproc_fft)) 364 365 ! Initialisation 366 tab_proc(:) = 0 367 sendcounts_sym(:) = 0 368 sendcounts_sym_all(:) = 0 369 sdispls_sym(:) = 0 370 recvcounts_sym(:) = 0 371 recvcounts_sym_tot(:) = 0 372 373 ! Localisation of kg_k==[0 0 0] 374 ABI_ALLOCATE(sum_kg,(ndatarecv)) 375 idatarecv0 = -1 376 ndatasend_sym = ndatarecv 377 sum_kg=sum(abs(kg_k_gather),1) 378 if (count(sum_kg==0)/=0) then 379 do idatarecv=1,ndatarecv 380 if (sum_kg(idatarecv)==0) idatarecv0=idatarecv 381 end do 382 ndatasend_sym = ndatarecv-1 383 end if 384 385 ! Localisation of the processor where the vector -k2 is 386 n2 = mpi_enreg%distribfft%n2_coarse 387 do idatarecv=1,ndatarecv 388 if (idatarecv/=idatarecv0) then 389 tab_proc(idatarecv) = mpi_enreg%distribfft%tab_fftwf2_distrib(modulo(-kg_k_gather(2,idatarecv),n2) + 1) 390 else 391 tab_proc(idatarecv) = -1 392 end if 393 end do 394 395 ! Number of values send by the processor to the others 396 do iproc=1,nproc_fft 397 sendcounts_sym(iproc) = count(tab_proc(:)==(iproc-1)) 398 end do 399 400 ! Save sendcounts_sym for each processor in sendcounts_sym_all 401 ! knowed by all processors of comm_fft 402 rdispls_sym(1)=0 403 do iproc=2,nproc_fft 404 rdispls_sym(iproc)= nproc_fft*(iproc-1) 405 end do 406 recvcounts_sym(:)=nproc_fft 407 call xmpi_allgatherv(sendcounts_sym(:),nproc_fft,& 408 & sendcounts_sym_all(:),recvcounts_sym,rdispls_sym,comm_fft,ierr) 409 410 ! Calculation of the dimension of kg_k_gather_sym for each processor 411 ! recvcounts_sym_tot is knowed by all processors of comm_fft 412 call xmpi_sum(sendcounts_sym,recvcounts_sym_tot,nproc_fft,comm_fft,ierr) 413 414 ! Dimension of kg_k_gather_sym 415 ndatarecv_tot = ndatarecv+recvcounts_sym_tot(me_fft+1) 416 417 ! Intialize kg_k_gather_sym 418 ABI_ALLOCATE(kg_k_gather_sym,(3,ndatarecv_tot)) 419 kg_k_gather_sym(:,:)=0 420 kg_k_gather_sym(:,1:ndatarecv) = kg_k_gather(:,:) 421 422 ! Allocation and initialisation 423 ABI_ALLOCATE(kg_k_gather_send,(3,ndatasend_sym)) 424 kg_k_gather_send(:,:)=0 425 426 ! The values are sorted in blocks 427 jsendloc=0 428 do iproc=1,nproc_fft 429 430 ! Position of the beginning of the block 431 sdispls_sym(iproc)=jsendloc 432 433 ! Creation of the blocks 434 do idatarecv=1,ndatarecv 435 if (tab_proc(idatarecv)==(iproc-1)) then 436 jsendloc=jsendloc+1 437 kg_k_gather_send(:,jsendloc) = -kg_k_gather(:,idatarecv) 438 end if 439 end do 440 end do 441 442 ! Position of received data 443 rdispls_sym(1)= ndatarecv 444 recvcounts_sym(1)= sendcounts_sym_all((me_fft+1)) 445 do iproc=2,nproc_fft 446 rdispls_sym(iproc) = rdispls_sym(iproc-1) + & 447 sendcounts_sym_all((me_fft+1)+(iproc-2)*nproc_fft) 448 recvcounts_sym(iproc) = sendcounts_sym_all((me_fft+1)+(iproc-1)*nproc_fft) 449 end do 450 451 ! Exchange of kg_k 452 sendcounts_sym_loc = sendcounts_sym*3 453 sdispls_sym_loc = sdispls_sym *3 454 recvcounts_sym_loc = recvcounts_sym*3 455 rdispls_sym_loc = rdispls_sym *3 456 call xmpi_alltoallv(kg_k_gather_send(:,:),sendcounts_sym_loc,sdispls_sym_loc,& 457 & kg_k_gather_sym(:,:) ,recvcounts_sym_loc,rdispls_sym_loc,comm_fft,ierr) 458 459 ! Store the following data in the bandfft_kpt_in data_struc 460 ikpt_this_proc=mpi_enreg%my_kpttab(ikpt) 461 if (bandfft_kpt_in(ikpt_this_proc)%flag3_is_allocated==0) then 462 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%kg_k_gather_sym,(3,ndatarecv_tot)) 463 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%rdispls_sym,(nproc_fft)) 464 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%recvcounts_sym,(nproc_fft)) 465 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%recvcounts_sym_tot,(nproc_fft)) 466 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%sdispls_sym,(nproc_fft)) 467 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%sendcounts_sym,(nproc_fft)) 468 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%sendcounts_sym_all,(nproc_fft*nproc_fft)) 469 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%tab_proc,(ndatarecv)) 470 bandfft_kpt_in(ikpt_this_proc)%flag3_is_allocated=1 471 end if 472 473 bandfft_kpt_in(ikpt_this_proc)%idatarecv0 =idatarecv0 474 bandfft_kpt_in(ikpt_this_proc)%ndatarecv_tot =ndatarecv_tot 475 bandfft_kpt_in(ikpt_this_proc)%ndatasend_sym =ndatasend_sym 476 bandfft_kpt_in(ikpt_this_proc)%kg_k_gather_sym(:,:) =kg_k_gather_sym(:,:) 477 bandfft_kpt_in(ikpt_this_proc)%rdispls_sym(:) =rdispls_sym(:) 478 bandfft_kpt_in(ikpt_this_proc)%recvcounts_sym(:) =recvcounts_sym(:) 479 bandfft_kpt_in(ikpt_this_proc)%recvcounts_sym_tot(:)=recvcounts_sym_tot(:) 480 bandfft_kpt_in(ikpt_this_proc)%sdispls_sym(:) =sdispls_sym(:) 481 bandfft_kpt_in(ikpt_this_proc)%sendcounts_sym(:) =sendcounts_sym(:) 482 bandfft_kpt_in(ikpt_this_proc)%sendcounts_sym_all(:)=sendcounts_sym_all(:) 483 bandfft_kpt_in(ikpt_this_proc)%tab_proc(:) =tab_proc(:) 484 485 ABI_DEALLOCATE(tab_proc) 486 ABI_DEALLOCATE(sendcounts_sym) 487 ABI_DEALLOCATE(sendcounts_sym_all) 488 ABI_DEALLOCATE(sdispls_sym) 489 ABI_DEALLOCATE(recvcounts_sym) 490 ABI_DEALLOCATE(recvcounts_sym_tot) 491 ABI_DEALLOCATE(rdispls_sym) 492 ABI_DEALLOCATE(kg_k_gather_sym) 493 ABI_DEALLOCATE(sendcounts_sym_loc) 494 ABI_DEALLOCATE(recvcounts_sym_loc) 495 ABI_DEALLOCATE(sdispls_sym_loc) 496 ABI_DEALLOCATE(rdispls_sym_loc) 497 ABI_DEALLOCATE(kg_k_gather_send) 498 ABI_DEALLOCATE(sum_kg) 499 500 ! Then compute gbound 501 call xmpi_allgather(ndatarecv_tot,npw_per_proc,mpi_enreg%comm_fft,ierr) 502 rdispls_all(1)=0 503 do iproc=2,nproc_fft 504 rdispls_all(iproc)=rdispls_all(iproc-1)+npw_per_proc(iproc-1) 505 end do 506 npw_tot=rdispls_all(nproc_fft)+npw_per_proc(nproc_fft) 507 ABI_ALLOCATE(kg_k_gather_all,(3,npw_tot)) 508 call xmpi_allgatherv(bandfft_kpt_in(ikpt_this_proc)%kg_k_gather_sym,& 509 & 3*ndatarecv_tot,kg_k_gather_all,3*npw_per_proc(:),3*rdispls_all,mpi_enreg%comm_fft,ierr) 510 if (mgfft>0) then 511 call sphereboundary(gbound,istwf_k,kg_k_gather_all,mgfft,npw_tot) 512 end if 513 514 ! Only calculations with istwfk=1 or 2 515 else 516 write(message, '(a,i0,a)' )' the value istwfk=',istwf_k,' is not allowed in case of bandfft parallelization!' 517 MSG_BUG(message) 518 end if 519 ABI_DEALLOCATE(kg_k_gather_all) 520 ABI_DEALLOCATE(npw_per_proc) 521 ABI_DEALLOCATE(rdispls_all) 522 ! ============================================================================ 523 ! End of gbound 524 ! ============================================================================ 525 526 ! Check if there is k_G vectors have been redistributed (after unbalancing dectecion) 527 ! Note that FFT load balancing is directly related to unbalancing detection 528 ! made in kpgsph routine. 529 itest=0 ; n2 = mpi_enreg%distribfft%n2_coarse 530 do idatarecv=1,ndatarecv 531 iproc=mpi_enreg%distribfft%tab_fftwf2_distrib(modulo(kg_k_gather(2,idatarecv),n2)+1) 532 if (iproc/=me_fft) itest=itest+1 533 end do 534 call xmpi_sum(itest,mpi_enreg%comm_fft,ierr) 535 if (itest>0) then 536 write(message, '(a,i4,3a)' ) & 537 & 'There is a load unbalancing for the FFT parallelization (kpt',ikpt,').',ch10,& 538 & 'Plane-wave components will be redistributed before each FFT!' 539 MSG_COMMENT(message) 540 end if 541 bandfft_kpt_in(ikpt_this_proc)%have_to_reequilibrate=(itest>0) 542 543 ! If yes, store relevant data 544 if(bandfft_kpt_in(ikpt_this_proc)%have_to_reequilibrate) then 545 n2 = mpi_enreg%distribfft%n2_coarse 546 sendcounts_fft(:) = 0 547 recvcounts_fft(:) = 0 548 senddisp_fft(:) = 0 549 recvdisp_fft(:) = 0 550 do idatarecv = 1 ,ndatarecv 551 iproc = mpi_enreg%distribfft%tab_fftwf2_distrib(modulo(kg_k_gather(2,idatarecv),n2) + 1) 552 sendcounts_fft(iproc + 1) = sendcounts_fft(iproc + 1) + 1 553 end do 554 call xmpi_alltoall(sendcounts_fft,1,recvcounts_fft,1,mpi_enreg%comm_fft,ierr) 555 do iproc =2, nproc_fft 556 senddisp_fft(iproc) = senddisp_fft(iproc - 1) + sendcounts_fft(iproc-1) 557 recvdisp_fft(iproc) = recvdisp_fft(iproc - 1) + recvcounts_fft(iproc-1) 558 end do 559 npw_fft = recvdisp_fft(nproc_fft) + recvcounts_fft(nproc_fft) ! nb plane wave for fourwf call 560 ABI_ALLOCATE(buff_kg,(3,ndatarecv)) ! for sorting kg_k 561 ABI_ALLOCATE(kg_k_fft,(3,npw_fft)) 562 ABI_ALLOCATE(indices_pw_fft,(ndatarecv)) 563 !filling of sorted send buffers 564 sendcounts_fft(:) = 0 565 do idatarecv = 1 ,ndatarecv 566 iproc = mpi_enreg%distribfft%tab_fftwf2_distrib(modulo(kg_k_gather(2,idatarecv),n2) + 1) 567 sendcounts_fft(iproc + 1) = sendcounts_fft(iproc + 1) + 1 568 indices_pw_fft(idatarecv) = senddisp_fft(iproc+1) + sendcounts_fft(iproc+1) 569 buff_kg(1:3,indices_pw_fft(idatarecv)) = kg_k_gather(1:3,idatarecv) 570 end do 571 572 call xmpi_alltoallv(buff_kg, 3*sendcounts_fft, 3*senddisp_fft, & 573 & kg_k_fft,3*recvcounts_fft, 3*recvdisp_fft, mpi_enreg%comm_fft,ierr) 574 575 if (.not.reequilibrate_allocated) then 576 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%kg_k_fft, (3,npw_fft) ) 577 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%indices_pw_fft, (ndatarecv)) 578 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%sendcount_fft, (nproc_fft)) 579 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%senddisp_fft, (nproc_fft)) 580 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%recvcount_fft, (nproc_fft)) 581 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%recvdisp_fft, (nproc_fft)) 582 end if 583 bandfft_kpt_in(ikpt_this_proc)%npw_fft = npw_fft 584 bandfft_kpt_in(ikpt_this_proc)%indices_pw_fft(:) = indices_pw_fft(:) 585 bandfft_kpt_in(ikpt_this_proc)%sendcount_fft (:) = sendcounts_fft (:) 586 bandfft_kpt_in(ikpt_this_proc)%senddisp_fft(:) = senddisp_fft(:) 587 bandfft_kpt_in(ikpt_this_proc)%recvcount_fft(:) = recvcounts_fft(:) 588 bandfft_kpt_in(ikpt_this_proc)%recvdisp_fft(:) = recvdisp_fft(:) 589 bandfft_kpt_in(ikpt_this_proc)%kg_k_fft(:,:) = kg_k_fft(:,:) 590 ABI_DEALLOCATE(buff_kg ) 591 ABI_DEALLOCATE(kg_k_fft) 592 ABI_DEALLOCATE(indices_pw_fft) 593 end if 594 595 ! Tabs which are common to istwf_k=1 and 2 596 if (.not. allocated(bandfft_kpt_in(ikpt_this_proc)%kg_k_gather)) then 597 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%kg_k_gather,(3,ndatarecv)) 598 end if 599 bandfft_kpt_in(ikpt_this_proc)%recvcounts(:) =recvcounts(:) 600 bandfft_kpt_in(ikpt_this_proc)%sendcounts(:) =sendcounts(:) 601 bandfft_kpt_in(ikpt_this_proc)%rdispls(:) =rdispls(:) 602 bandfft_kpt_in(ikpt_this_proc)%sdispls(:) =sdispls(:) 603 bandfft_kpt_in(ikpt_this_proc)%gbound(:,:) =gbound(:,:) 604 bandfft_kpt_in(ikpt_this_proc)%kg_k_gather(:,:)=kg_k_gather(:,:) 605 bandfft_kpt_in(ikpt_this_proc)%ndatarecv =ndatarecv 606 bandfft_kpt_in(ikpt_this_proc)%istwf_k =istwf_k 607 bandfft_kpt_in(ikpt_this_proc)%npw_tot =npw_tot 608 ABI_DEALLOCATE(kg_k_gather) 609 ABI_DEALLOCATE(kg_k) 610 ABI_DEALLOCATE(gbound) 611 612 ikg=ikg+npw_k 613 end do 614 end do 615 ABI_DEALLOCATE(recvcounts) 616 ABI_DEALLOCATE(sendcounts) 617 ABI_DEALLOCATE(rdispls) 618 ABI_DEALLOCATE(sdispls) 619 ABI_DEALLOCATE(sendcounts_fft) 620 ABI_DEALLOCATE(senddisp_fft) 621 ABI_DEALLOCATE(recvcounts_fft) 622 ABI_DEALLOCATE(recvdisp_fft) 623 624 !============================================================================= 625 !End of computation and storage of the bandfft_kpt(ikpt) data_struc 626 !============================================================================= 627 628 end subroutine bandfft_kpt_init1
m_bandfft_kpt/bandfft_kpt_init2 [ Functions ]
[ Top ] [ m_bandfft_kpt ] [ Functions ]
NAME
bandfft_kpt_init2
FUNCTION
Init part of scalars and pointers in a bandfft_kpt datastructure
INPUTS
dimffnl=second dimension of ffnl (1+number of derivatives) ffnl_gather(ndatarecv,dimffnl,lmnmax,ntypat)=nonlocal form factors on basis sphere. kinpw_gather(:)=(modified) kinetic energy for each plane wave (Hartree) kpg_k_gather(ndatarecv,nkpg)=k+G vector for a given k point lmnmax=if useylm=1, max number of (l,m,n) comp. over all type of psps =if useylm=0, max number of (l,n) comp. over all type of psps matblk=dimension of the array ph3d mkmem =number of k points which can fit in memory; set to 0 if use disk ndatarecv= dimension of the arrays nkpg=second dimension of kpg_k (0 if useylm=0) ntypat=number of types of atoms in unit cell. ph3d_gather(2,ndatarecv,matblk)=3-dim structure factors, for each atom and plane wave.
OUTPUT
SIDE EFFECTS
PARENTS
prep_bandfft_tabs
CHILDREN
SOURCE
665 subroutine bandfft_kpt_init2(bandfft_kpt_in,dimffnl,ffnl_gather,ikpt_this_proc,kinpw_gather,& 666 & kpg_k_gather,lmnmax,matblk,mkmem,ndatarecv,nkpg,ntypat,ph3d_gather) 667 668 669 !This section has been created automatically by the script Abilint (TD). 670 !Do not modify the following lines by hand. 671 #undef ABI_FUNC 672 #define ABI_FUNC 'bandfft_kpt_init2' 673 !End of the abilint section 674 675 implicit none 676 677 !Arguments ------------------------------- 678 integer, intent(in) :: dimffnl,ikpt_this_proc,lmnmax,matblk,mkmem,ndatarecv,nkpg,ntypat 679 !Local variables------------------------------- 680 real(dp),intent(in) :: ffnl_gather(:,:,:,:),kinpw_gather(:) 681 real(dp),intent(in) :: kpg_k_gather(:,:),ph3d_gather(:,:,:) 682 type(bandfft_kpt_type), intent(inout) :: bandfft_kpt_in (mkmem) 683 684 ! ********************************************************************* 685 686 if (allocated(bandfft_kpt_in(ikpt_this_proc)%ffnl_gather)) then 687 ABI_DEALLOCATE(bandfft_kpt_in(ikpt_this_proc)%ffnl_gather) 688 end if 689 if (size(ffnl_gather)>0) then 690 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%ffnl_gather,(ndatarecv,dimffnl,lmnmax,ntypat)) 691 bandfft_kpt_in(ikpt_this_proc)%ffnl_gather(:,:,:,:)=ffnl_gather(:,:,:,:) 692 else 693 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%ffnl_gather,(0,0,0,0)) 694 end if 695 696 if (allocated(bandfft_kpt_in(ikpt_this_proc)%ph3d_gather)) then 697 ABI_DEALLOCATE(bandfft_kpt_in(ikpt_this_proc)%ph3d_gather) 698 end if 699 if (size(ph3d_gather)>0) then 700 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%ph3d_gather,(2,ndatarecv,matblk)) 701 bandfft_kpt_in(ikpt_this_proc)%ph3d_gather(:,:,:) =ph3d_gather(:,:,:) 702 else 703 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%ph3d_gather,(0,0,0)) 704 end if 705 706 if (allocated(bandfft_kpt_in(ikpt_this_proc)%kpg_k_gather)) then 707 ABI_DEALLOCATE(bandfft_kpt_in(ikpt_this_proc)%kpg_k_gather) 708 end if 709 if (size(kpg_k_gather)>0) then 710 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%kpg_k_gather,(ndatarecv,nkpg)) 711 bandfft_kpt_in(ikpt_this_proc)%kpg_k_gather(:,:) =kpg_k_gather(:,:) 712 else 713 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%kpg_k_gather,(0,0)) 714 end if 715 716 if (allocated(bandfft_kpt_in(ikpt_this_proc)%kinpw_gather)) then 717 ABI_DEALLOCATE(bandfft_kpt_in(ikpt_this_proc)%kinpw_gather) 718 end if 719 if (size(kinpw_gather)>0) then 720 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%kinpw_gather,(ndatarecv)) 721 bandfft_kpt_in(ikpt_this_proc)%kinpw_gather(:) =kinpw_gather(:) 722 else 723 ABI_ALLOCATE(bandfft_kpt_in(ikpt_this_proc)%kinpw_gather,(0)) 724 end if 725 726 bandfft_kpt_in(ikpt_this_proc)%flag2_is_allocated=1 727 728 end subroutine bandfft_kpt_init2
m_bandfft_kpt/bandfft_kpt_reset [ Functions ]
[ Top ] [ m_bandfft_kpt ] [ Functions ]
NAME
bandfft_kpt_reset
FUNCTION
Reset flags a bandfft_kpt datastructure
INPUTS
OUTPUT
SIDE EFFECTS
bandfft_kpt_in=the datastructure to nullify
PARENTS
posdoppler
CHILDREN
SOURCE
754 subroutine bandfft_kpt_reset(bandfft_kpt_in) 755 756 757 !This section has been created automatically by the script Abilint (TD). 758 !Do not modify the following lines by hand. 759 #undef ABI_FUNC 760 #define ABI_FUNC 'bandfft_kpt_reset' 761 !End of the abilint section 762 763 implicit none 764 765 !Arguments ------------------------------------ 766 type(bandfft_kpt_type) :: bandfft_kpt_in 767 !Local variables------------------------------- 768 769 ! *********************************************************************** 770 771 bandfft_kpt_in%flag1_is_allocated=0 772 bandfft_kpt_in%flag2_is_allocated=0 773 bandfft_kpt_in%flag3_is_allocated=0 774 bandfft_kpt_in%have_to_reequilibrate=.false. 775 776 end subroutine bandfft_kpt_reset
m_bandfft_kpt/bandfft_kpt_restoretabs [ Functions ]
[ Top ] [ m_bandfft_kpt ] [ Functions ]
NAME
bandfft_kpt_restoretabs
FUNCTION
Restore some (kpt-dependent) tabs into a bandfft_kpt datastructure INPUT
OUTPUT
=== Arrays to be eventually restored ===== ffnl(:,:,:,:)=nonlocal form factors on basis sphere ph3d(:,:,:)=3-dim structure factors, for each atom and plane wave kpg(:,:)=k+G vector for a given k point kinpw(:)=kinetic energy for each plane wave (Hartree)
SIDE EFFECTS
bandfft_kpt_out=<type(bandfft_kpt)>=bandfft_kpt datastructure
PARENTS
energy,fock2ACE,forstrnps
CHILDREN
SOURCE
1977 subroutine bandfft_kpt_restoretabs(bandfft_kpt_out,ffnl,ph3d,kpg,kinpw) 1978 1979 1980 !This section has been created automatically by the script Abilint (TD). 1981 !Do not modify the following lines by hand. 1982 #undef ABI_FUNC 1983 #define ABI_FUNC 'bandfft_kpt_restoretabs' 1984 !End of the abilint section 1985 1986 implicit none 1987 1988 !Arguments ------------------------------- 1989 type(bandfft_kpt_type),intent(inout) :: bandfft_kpt_out 1990 real(dp),intent(inout),allocatable,optional :: ffnl(:,:,:,:),ph3d(:,:,:),kpg(:,:),kinpw(:) 1991 !Local variables------------------------------- 1992 integer :: is1,is2,is3,is4 1993 1994 ! ********************************************************************* 1995 1996 if (present(ffnl)) then 1997 if (allocated(bandfft_kpt_out%ffnl_gather)) then 1998 ABI_DEALLOCATE(bandfft_kpt_out%ffnl_gather) 1999 end if 2000 if (allocated(ffnl)) then 2001 is1=size(ffnl,1) 2002 is2=size(ffnl,2) 2003 is3=size(ffnl,3) 2004 is4=size(ffnl,4) 2005 ABI_ALLOCATE(bandfft_kpt_out%ffnl_gather,(is1,is2,is3,is4)) 2006 bandfft_kpt_out%ffnl_gather(:,:,:,:)=ffnl(:,:,:,:) 2007 ABI_DEALLOCATE(ffnl) 2008 bandfft_kpt_out%flag2_is_allocated=1 2009 end if 2010 end if 2011 if (present(ph3d)) then 2012 if (allocated(bandfft_kpt_out%ph3d_gather)) then 2013 ABI_DEALLOCATE(bandfft_kpt_out%ph3d_gather) 2014 end if 2015 if (allocated(ph3d)) then 2016 is1=size(ph3d,1) 2017 is2=size(ph3d,2) 2018 is3=size(ph3d,3) 2019 ABI_ALLOCATE(bandfft_kpt_out%ph3d_gather,(is1,is2,is3)) 2020 bandfft_kpt_out%ph3d_gather(:,:,:)=ph3d(:,:,:) 2021 ABI_DEALLOCATE(ph3d) 2022 bandfft_kpt_out%flag2_is_allocated=1 2023 end if 2024 end if 2025 if (present(kpg)) then 2026 if (allocated(bandfft_kpt_out%kpg_k_gather)) then 2027 ABI_DEALLOCATE(bandfft_kpt_out%kpg_k_gather) 2028 end if 2029 if (allocated(kpg)) then 2030 is1=size(kpg,1) 2031 is2=size(kpg,2) 2032 ABI_ALLOCATE(bandfft_kpt_out%kpg_k_gather,(is1,is2)) 2033 bandfft_kpt_out%kpg_k_gather(:,:)=kpg(:,:) 2034 ABI_DEALLOCATE(kpg) 2035 bandfft_kpt_out%flag2_is_allocated=1 2036 end if 2037 end if 2038 if (present(kinpw)) then 2039 if (allocated(bandfft_kpt_out%kinpw_gather)) then 2040 ABI_DEALLOCATE(bandfft_kpt_out%kinpw_gather) 2041 end if 2042 if (allocated(kinpw)) then 2043 is1=size(kinpw,1) 2044 ABI_ALLOCATE(bandfft_kpt_out%kinpw_gather,(is1)) 2045 bandfft_kpt_out%kinpw_gather(:)=kinpw(:) 2046 ABI_DEALLOCATE(kinpw) 2047 end if 2048 end if 2049 2050 end subroutine bandfft_kpt_restoretabs
m_bandfft_kpt/bandfft_kpt_savetabs [ Functions ]
[ Top ] [ m_bandfft_kpt ] [ Functions ]
NAME
bandfft_kpt_savetabs
FUNCTION
Save some (kpt-dependent) tabs from a bandfft_kpt datastructure
INPUTS
bandfft_kpt_in=<type(bandfft_kpt)>=bandfft_kpt datastructure
OUTPUT
ffnl(:,:,:,:)=nonlocal form factors on basis sphere ph3d(:,:,:)=3-dim structure factors, for each atom and plane wave kpg(:,:)=k+G vector for a given k point kinpw(:)=kinetic energy for each plane wave (Hartree)
PARENTS
energy,fock2ACE,forstrnps
CHILDREN
SOURCE
1880 subroutine bandfft_kpt_savetabs(bandfft_kpt_in,ffnl,ph3d,kpg,kinpw) 1881 1882 1883 !This section has been created automatically by the script Abilint (TD). 1884 !Do not modify the following lines by hand. 1885 #undef ABI_FUNC 1886 #define ABI_FUNC 'bandfft_kpt_savetabs' 1887 !End of the abilint section 1888 1889 implicit none 1890 1891 !Arguments ------------------------------- 1892 type(bandfft_kpt_type),intent(inout) :: bandfft_kpt_in 1893 real(dp),intent(inout),allocatable,optional :: ffnl(:,:,:,:),ph3d(:,:,:),kpg(:,:),kinpw(:) 1894 !Local variables------------------------------- 1895 integer :: is1,is2,is3,is4 1896 1897 ! ********************************************************************* 1898 1899 if (present(ffnl)) then 1900 if (allocated(ffnl)) then 1901 ABI_DEALLOCATE(ffnl) 1902 end if 1903 if (allocated(bandfft_kpt_in%ffnl_gather)) then 1904 is1=size(bandfft_kpt_in%ffnl_gather,1) 1905 is2=size(bandfft_kpt_in%ffnl_gather,2) 1906 is3=size(bandfft_kpt_in%ffnl_gather,3) 1907 is4=size(bandfft_kpt_in%ffnl_gather,4) 1908 ABI_ALLOCATE(ffnl,(is1,is2,is3,is4)) 1909 ffnl(:,:,:,:)=bandfft_kpt_in%ffnl_gather(:,:,:,:) 1910 end if 1911 end if 1912 if (present(ph3d)) then 1913 if (allocated(ph3d)) then 1914 ABI_DEALLOCATE(ph3d) 1915 end if 1916 if (allocated(bandfft_kpt_in%ph3d_gather)) then 1917 is1=size(bandfft_kpt_in%ph3d_gather,1) 1918 is2=size(bandfft_kpt_in%ph3d_gather,2) 1919 is3=size(bandfft_kpt_in%ph3d_gather,3) 1920 ABI_ALLOCATE(ph3d,(is1,is2,is3)) 1921 ph3d(:,:,:)=bandfft_kpt_in%ph3d_gather(:,:,:) 1922 end if 1923 end if 1924 if (present(kpg)) then 1925 if (allocated(kpg)) then 1926 ABI_DEALLOCATE(kpg) 1927 end if 1928 if (allocated(bandfft_kpt_in%kpg_k_gather)) then 1929 is1=size(bandfft_kpt_in%kpg_k_gather,1) 1930 is2=size(bandfft_kpt_in%kpg_k_gather,2) 1931 ABI_ALLOCATE(kpg,(is1,is2)) 1932 kpg(:,:)=bandfft_kpt_in%kpg_k_gather(:,:) 1933 end if 1934 end if 1935 if (present(kinpw)) then 1936 if (allocated(kinpw)) then 1937 ABI_DEALLOCATE(kinpw) 1938 end if 1939 if (allocated(bandfft_kpt_in%kinpw_gather)) then 1940 is1=size(bandfft_kpt_in%kinpw_gather,1) 1941 ABI_ALLOCATE(kinpw,(is1)) 1942 kinpw(:)=bandfft_kpt_in%kinpw_gather(:) 1943 end if 1944 end if 1945 1946 end subroutine bandfft_kpt_savetabs
m_bandfft_kpt/bandfft_kpt_set_ikpt [ Functions ]
[ Top ] [ m_bandfft_kpt ] [ Functions ]
NAME
bandfft_kpt_set_ikpt
FUNCTION
Set the value of bandfft_kpt_current_ikpt variable, i.e. current index of bandfft_kpt loaded in memory INPUT ikpt=index of k-point (in the global array dtset%kpt) mpi_enreg= information about MPI parallelization
OUTPUT
bandfft_kpt_current_ikpt value changed
PARENTS
mkrho,prep_bandfft_tabs,vtorho
CHILDREN
SOURCE
2077 subroutine bandfft_kpt_set_ikpt(ikpt,mpi_enreg) 2078 2079 2080 !This section has been created automatically by the script Abilint (TD). 2081 !Do not modify the following lines by hand. 2082 #undef ABI_FUNC 2083 #define ABI_FUNC 'bandfft_kpt_set_ikpt' 2084 !End of the abilint section 2085 2086 implicit none 2087 2088 !Arguments ------------------------------- 2089 integer,intent(in) :: ikpt 2090 type(MPI_type),intent(inout) :: mpi_enreg 2091 !Local variables------------------------------- 2092 2093 ! ********************************************************************* 2094 2095 if (ikpt>0) then 2096 bandfft_kpt_current_ikpt=mpi_enreg%my_kpttab(ikpt) 2097 else 2098 bandfft_kpt_current_ikpt=-1 2099 end if 2100 2101 end subroutine bandfft_kpt_set_ikpt
m_bandfft_kpt/bandfft_kpt_type [ Types ]
[ Top ] [ m_bandfft_kpt ] [ Types ]
NAME
bandfft_kpt_type
FUNCTION
The bandfft_kpt_type structured datatype gather different information about the triple band-fft-kpt parallelisation : tabs which are distributed over all the three dimensions and stored during the calculation, dimensions of messages exchange during the calculations... i.e.: all the information which were spread over the entire code before and recomputed at each iline, istep or itime STEP with a large probability to make a mistake.
SOURCE
75 type, public :: bandfft_kpt_type 76 77 ! WARNING : if you modify this datatype, please check whether there might be creation/destruction/copy routines, 78 ! declared in another part of ABINIT, that might need to take into account your modification. 79 80 integer :: flag1_is_allocated ! determine if the following data are allocated or not 81 integer :: npw_tot ! array holding the total number of plane waves for each k point 82 integer :: ndatarecv ! total number of values received by the processor and sent 83 ! by the other processors band 84 integer, allocatable :: kg_k_gather(:,:) ! planewave coordinates 85 ! (of the processor + sent by other processors band) 86 integer, allocatable :: recvcounts(:) ! number of values received by the processor from each processor band 87 integer, allocatable :: sendcounts(:) ! number of values sent by the processor to each processor band 88 integer, allocatable :: rdispls (:) ! positions of values received by the processor from each processor band 89 integer, allocatable :: sdispls (:) ! postions of values sent by the processor to each processor band 90 integer, allocatable :: gbound(:,:) ! sphere boundary info: gbound(2*mgfft+8,2) 91 92 integer :: flag2_is_allocated ! determine if the following data are allocated or not 93 real(dp), allocatable :: ffnl_gather(:,:,:,:) ! ffnl tab (of the processor + sent by other processors band) 94 real(dp), allocatable :: kinpw_gather(:) ! kinpw tab (of the processor + sent by other processors band) 95 real(dp), allocatable :: ph3d_gather(:,:,:) ! ph3d tab (of the processor + sent by other processors band) 96 real(dp), allocatable :: kpg_k_gather(:,:) ! kpg_k tab (of the processor + sent by other processors band) 97 98 99 integer :: flag3_is_allocated ! determine if the following data are allocated or not 100 integer :: istwf_k ! input option parameter that describes the storage of wfs 101 integer :: idatarecv0 ! position of the planewave coordinates (0,0,0) 102 integer :: ndatarecv_tot ! total number of received values by the processor 103 ! (ndatarecv + number of received opposited planewave coordinates) 104 integer :: ndatasend_sym ! number of sent values to the processors fft to create opposited 105 ! planewave coordinates 106 integer, allocatable :: kg_k_gather_sym(:,:) ! planewave coordinates 107 ! (kg_k_gather + opposited planewave coordinates sent by the processors fft) 108 integer, allocatable :: rdispls_sym(:) ! positions of values received by the processor from each processor fft 109 integer, allocatable :: recvcounts_sym(:) ! number of values received by the processor from each processor fft 110 integer, allocatable :: recvcounts_sym_tot(:) ! number of values received by each processor from the other processors fft 111 integer, allocatable :: sdispls_sym(:) ! postions of values sent by the processor to each processor fft 112 integer, allocatable :: sendcounts_sym(:) ! number of values sent by the processor to each processor fft 113 integer, allocatable :: sendcounts_sym_all(:) ! number of values sent by each processor to the other processors fft 114 integer, allocatable :: tab_proc(:) ! positions of opposited planewave coordinates in the list of the processors fft 115 116 logical :: have_to_reequilibrate ! indicates weather we will have to reequilibrate and allocate all these stuff 117 integer :: npw_fft ! Number of plane waves during fft step 118 integer, allocatable :: indices_pw_fft(:) ! Indices for sorting pw like pw 119 integer, allocatable :: sendcount_fft (:) ! Number of pw to send to others proc fft 120 integer, allocatable :: senddisp_fft(:) ! Positions for sending 121 integer, allocatable :: recvcount_fft(:) ! Number of pw to receive from others proc fft 122 integer, allocatable :: recvdisp_fft(:) ! Positions for receiving 123 integer, allocatable :: kg_k_fft(:,:) ! planewaves coordinates 124 125 end type bandfft_kpt_type