The Elk Code
symrfir.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2007 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: symrfir
8 ! !INTERFACE:
9 subroutine symrfir(ngridg_,ngtot_,ngvec_,nfgrz_,igfft_,igrzf_,rfir)
10 ! !USES:
11 use modmain
12 ! !INPUT/OUTPUT PARAMETERS:
13 ! ngridg_ : G-vector grid sizes (in,integer(3))
14 ! ngtot_ : total number of G-vectors (in,integer)
15 ! ngvec_ : number of G-vectors within cut-off (in,integer)
16 ! nfgrz_ : number of FFT elements for real-complex transforms (in,integer)
17 ! igfft_ : map from G-vector index to FFT array (in,integer(ngvec_))
18 ! igrzf_ : map from real-complex FFT index to G-point index
19 ! (in,integer(nfgrz_))
20 ! rfir : real intersitial function (inout,real(ngtot_))
21 ! !DESCRIPTION:
22 ! Symmetrises a real scalar interstitial function. The function is first
23 ! Fourier transformed to $G$-space, and then averaged over each symmetry by
24 ! rotating the Fourier coefficients and multiplying them by a phase factor
25 ! corresponding to the symmetry translation.
26 !
27 ! !REVISION HISTORY:
28 ! Created July 2007 (JKD)
29 !EOP
30 !BOC
31 implicit none
32 ! arguments
33 integer, intent(in) :: ngridg_(3),ngtot_,ngvec_,nfgrz_
34 integer, intent(in) :: igfft_(ngvec_),igrzf_(nfgrz_)
35 real(8), intent(inout) :: rfir(ngtot_)
36 ! local variables
37 logical tv0
38 integer isym,lspl,sym(3,3)
39 integer ig,jg,ifg,jfg
40 integer i1,i2,i3,j1,j2,j3
41 real(8) v1,v2,v3,t1
42 ! allocatable arrays
43 complex(8), allocatable :: zfft1(:),zfft2(:)
44 allocate(zfft1(ngtot_),zfft2(nfgrz_))
45 ! Fourier transform function to G-space
46 zfft1(1:ngtot_)=rfir(1:ngtot_)
47 call zfftifc(3,ngridg_,-1,zfft1)
48 zfft2(1:nfgrz_)=0.d0
49 ! loop over crystal symmetries
50 do isym=1,nsymcrys
51 ! zero translation vector flag
52  tv0=tv0symc(isym)
53 ! translation vector in Cartesian coordinates
54  if (.not.tv0) then
55  v1=vtcsymc(1,isym)
56  v2=vtcsymc(2,isym)
57  v3=vtcsymc(3,isym)
58  end if
59 ! index to spatial rotation lattice symmetry
60  lspl=lsplsymc(isym)
61  sym(1:3,1:3)=symlat(1:3,1:3,lspl)
62  do ifg=1,nfgrz_
63  ig=igrzf_(ifg)
64  if (ig > ngvec_) cycle
65 ! multiply the transpose of the inverse symmetry matrix with the G-vector
66  if (lspl == 1) then
67  jg=ig
68  else
69  i1=ivg(1,ig); i2=ivg(2,ig); i3=ivg(3,ig)
70  j1=sym(1,1)*i1+sym(2,1)*i2+sym(3,1)*i3
71  j2=sym(1,2)*i1+sym(2,2)*i2+sym(3,2)*i3
72  j3=sym(1,3)*i1+sym(2,3)*i2+sym(3,3)*i3
73  jg=ivgig(j1,j2,j3)
74  end if
75  jfg=igfft_(jg)
76 ! translation and spatial rotation
77  if (tv0) then
78 ! zero translation vector
79  zfft2(ifg)=zfft2(ifg)+zfft1(jfg)
80  else
81 ! complex phase factor for translation
82  t1=vgc(1,jg)*v1+vgc(2,jg)*v2+vgc(3,jg)*v3
83  zfft2(ifg)=zfft2(ifg)+zfft1(jfg)*cmplx(cos(t1),-sin(t1),8)
84  end if
85  end do
86 end do
87 ! Fourier transform to real-space and normalise
88 call rzfftifc(3,ngridg_,1,rfir,zfft2)
89 t1=1.d0/dble(nsymcrys)
90 rfir(1:ngtot_)=t1*rfir(1:ngtot_)
91 deallocate(zfft1,zfft2)
92 end subroutine
93 !EOC
94 
integer, dimension(:,:,:), allocatable ivgig
Definition: modmain.f90:402
subroutine symrfir(ngridg_, ngtot_, ngvec_, nfgrz_, igfft_, igrzf_, rfir)
Definition: symrfir.f90:10
integer nsymcrys
Definition: modmain.f90:358
integer, dimension(3, 3, 48) symlat
Definition: modmain.f90:344
logical, dimension(maxsymcrys) tv0symc
Definition: modmain.f90:362
real(8), dimension(:,:), allocatable vgc
Definition: modmain.f90:420
subroutine zfftifc(nd, n, sgn, z)
Definition: zfftifc_fftw.f90:7
integer, dimension(maxsymcrys) lsplsymc
Definition: modmain.f90:364
integer, dimension(:,:), allocatable ivg
Definition: modmain.f90:400
subroutine rzfftifc(nd, n, sgn, r, z)
real(8), dimension(3, maxsymcrys) vtcsymc
Definition: modmain.f90:360