TABLE OF CONTENTS
ABINIT/zerosym [ Functions ]
NAME
zerosym
FUNCTION
Symmetrize an array on the FFT grid by vanishing some term on the boundaries.
COPYRIGHT
Copyright (C) 1998-2017 ABINIT group (GZ, 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 . For the initials of contributors, see ~abinit/doc/developers/contributors.txt.
INPUTS
cplex= if 1, input array is REAL, if 2, input array is COMPLEX mpi_enreg=informations about MPI parallelization n1,n2,n3=FFT dimensions nfft=n1*n2*n3 ig1,ig2,ig3=optional arguments= indexes of unbalanced g-vectors to cancel if not present, ig1=1+n1/2, ig2=1+n2/2, ig3=1+n3/2 for even n1,n2,n3 if igj=-1, nothing is done in direction j
OUTPUT
(see side effects)
SIDE EFFECTS
array(cplex,n1*n2*n3)=complex array to be symetrized
PARENTS
atm2fft,dfpt_atm2fft,dfpt_dyfro,forces,hartre,initro,m_fock,pawmknhat pawmknhat_psipsi,prcref,prcref_PMA,stress,transgrid
CHILDREN
SOURCE
38 #if defined HAVE_CONFIG_H 39 #include "config.h" 40 #endif 41 42 #include "abi_common.h" 43 44 subroutine zerosym(array,cplex,n1,n2,n3,& 45 & ig1,ig2,ig3,comm_fft,distribfft) ! Optional arguments 46 47 use defs_basis 48 use m_profiling_abi 49 use m_distribfft 50 use m_errors 51 use m_xmpi 52 53 !This section has been created automatically by the script Abilint (TD). 54 !Do not modify the following lines by hand. 55 #undef ABI_FUNC 56 #define ABI_FUNC 'zerosym' 57 !End of the abilint section 58 59 implicit none 60 61 !Arguments ------------------------------------ 62 !scalars 63 integer,intent(in) :: cplex,n1,n2,n3 64 integer,optional,intent(in) :: ig1,ig2,ig3,comm_fft 65 type(distribfft_type),intent(in),target,optional :: distribfft 66 !arrays 67 real(dp),intent(inout) :: array(cplex,n1*n2*n3) 68 69 !Local variables------------------------------- 70 !scalars 71 integer :: i1,i2,i3,ifft,ifft_proc,index,j,j1,j2,j3,me_fft,nd2 72 integer :: nproc_fft,n1sel,nn12,n2sel,n3sel,r2 73 !arrays 74 integer, ABI_CONTIGUOUS pointer :: fftn2_distrib(:),ffti2_local(:) 75 76 ! ********************************************************************** 77 78 DBG_ENTER("COLL") 79 80 me_fft=0;nproc_fft=1 81 if (present(comm_fft)) then 82 me_fft=xmpi_comm_rank(comm_fft) 83 nproc_fft=xmpi_comm_size(comm_fft) 84 end if 85 nd2=(n2-1)/nproc_fft+1 86 nn12=n1*n2 87 88 !Get the distrib associated with this fft_grid 89 if (present(distribfft)) then 90 if (n2== distribfft%n2_coarse) then 91 fftn2_distrib => distribfft%tab_fftdp2_distrib 92 ffti2_local => distribfft%tab_fftdp2_local 93 else if(n2 == distribfft%n2_fine) then 94 fftn2_distrib => distribfft%tab_fftdp2dg_distrib 95 ffti2_local => distribfft%tab_fftdp2dg_local 96 else 97 MSG_BUG("Unable to find an allocated distrib for this fft grid") 98 end if 99 else 100 ABI_ALLOCATE(fftn2_distrib,(n2)) 101 ABI_ALLOCATE(ffti2_local,(n2)) 102 fftn2_distrib=0;ffti2_local=(/(i2,i2=1,n2)/) 103 end if 104 105 if (present(ig1)) then 106 n1sel=ig1 107 else if (mod(n1,2)==0) then 108 n1sel=1+n1/2 109 else 110 n1sel=-1 111 end if 112 if (present(ig2)) then 113 n2sel=ig2 114 else if (mod(n2,2)==0) then 115 n2sel=1+n2/2 116 else 117 n2sel=-1 118 end if 119 if (present(ig3)) then 120 n3sel=ig3 121 else if (mod(n3,2)==0) then 122 n3sel=1+n3/2 123 else 124 n3sel=-1 125 end if 126 127 if (n1sel>0) then 128 index=n1sel-nn12-n1 129 do i3=1,n3 130 index=index+nn12;ifft=index 131 do i2=1,n2 132 ifft=ifft+n1 133 if (nproc_fft>1) then 134 ! MPIWF: consider ifft only if it is treated by the current proc and compute its adress 135 j=ifft-1;j1=modulo(j,n1);j2=modulo(j/n1,n2);j3=j/(n1*n2) !;r2=modulo(j2,nd2) 136 if(fftn2_distrib(j2+1)==me_fft) then ! MPIWF this ifft is to be treated by me_fft 137 r2= ffti2_local(j2+1) - 1 138 ifft_proc=n1*(nd2*j3+r2)+j1+1 !this is ifft in the current proc 139 array(:,ifft_proc)=zero 140 end if 141 else 142 array(:,ifft)=zero 143 end if 144 end do 145 end do 146 end if 147 148 if (n2sel>0) then 149 index=n1*n2sel-nn12-n1 150 do i3=1,n3 151 index=index+nn12;ifft=index 152 do i1=1,n1 153 ifft=ifft+1 154 if (nproc_fft>1) then 155 ! MPIWF: consider ifft only if it is treated by the current proc and compute its adress 156 j=ifft-1;j1=modulo(j,n1);j2=modulo(j/n1,n2);j3=j/(n1*n2); 157 if(fftn2_distrib(j2+1)==me_fft) then ! MPIWF this ifft is to be treated by me_fft 158 r2= ffti2_local(j2+1) - 1 159 ifft_proc=n1*(nd2*j3+r2)+j1+1 !this is ifft in the current proc 160 array(:,ifft_proc)=zero 161 end if 162 else 163 array(:,ifft)=zero 164 end if 165 end do 166 end do 167 end if 168 169 if (n3sel>0) then 170 index=nn12*n3sel-nn12-n1 171 do i2=1,n2 172 index=index+n1;ifft=index 173 do i1=1,n1 174 ifft=ifft+1 175 if (nproc_fft>1) then 176 ! MPIWF: consider ifft only if it is treated by the current proc and compute its adress 177 j=ifft-1;j1=modulo(j,n1);j2=modulo(j/n1,n2);j3=j/(n1*n2) 178 if(fftn2_distrib(j2+1)==me_fft) then ! MPIWF this ifft is to be treated by me_fft 179 r2= ffti2_local(j2+1) - 1 180 ifft_proc=n1*(nd2*j3+r2)+j1+1 !this is ifft in the current proc 181 array(:,ifft_proc)=zero 182 end if 183 else 184 array(:,ifft)=zero 185 end if 186 end do 187 end do 188 end if 189 190 if (.not.present(distribfft)) then 191 ABI_DEALLOCATE(fftn2_distrib) 192 ABI_DEALLOCATE(ffti2_local) 193 end if 194 195 DBG_EXIT("COLL") 196 197 end subroutine zerosym