The Elk Code
symrvfir.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: symrvfir
8 ! !INTERFACE:
9 subroutine symrvfir(tspin,tnc,ngridg_,ngtot_,ngvec_,nfgrz_,igfft_,igrzf_,ld, &
10  rvfir)
11 ! !USES:
12 use modmain
13 ! !INPUT/OUTPUT PARAMETERS:
14 ! tspin : .true. if spin rotations should be used (in,logical)
15 ! tnc : .true. if the vector field is non-collinear, otherwise it is
16 ! collinear along the z-axis (in,logical)
17 ! ngridg_ : G-vector grid sizes (in,integer(3))
18 ! ngtot_ : total number of G-vectors (in,integer)
19 ! ngvec_ : number of G-vectors within cut-off (in,integer)
20 ! nfgrz_ : number of FFT elements for real-complex transforms (in,integer)
21 ! igfft_ : map from G-vector index to FFT array (in,integer(ngvec_))
22 ! igrzf_ : map from real-complex FFT index to G-point index
23 ! (in,integer(nfgrz_))
24 ! ld : leading dimension (in,integer)
25 ! rvfir : real interstitial vector function (inout,real(ld,*))
26 ! !DESCRIPTION:
27 ! Symmetrises a real interstitial vector function. See routines {\tt symrvf}
28 ! and {\tt symrfir} for details.
29 !
30 ! !REVISION HISTORY:
31 ! Created July 2007 (JKD)
32 !EOP
33 !BOC
34 implicit none
35 ! arguments
36 logical, intent(in) :: tspin,tnc
37 integer, intent(in) :: ngridg_(3),ngtot_,ngvec_,nfgrz_
38 integer, intent(in) :: igfft_(ngvec_),igrzf_(nfgrz_),ld
39 real(8), intent(inout) :: rvfir(ld,*)
40 ! local variables
41 logical tv0
42 integer isym,lspl,lspn,sym(3,3)
43 integer nd,ig,jg,ifg,jfg
44 integer i1,i2,i3,j1,j2,j3,i
45 real(8) sc(3,3),v1,v2,v3,t1
46 complex(8) z0,z1,z2,z3
47 ! allocatable arrays
48 complex(8), allocatable :: zfft1(:,:),zfft2(:,:)
49 ! dimension of the vector field
50 if (tnc) then
51  nd=3
52 else
53  nd=1
54 end if
55 allocate(zfft1(ngtot_,nd),zfft2(nfgrz_,nd))
56 ! Fourier transform vector function to G-space
57 do i=1,nd
58  zfft1(1:ngtot_,i)=rvfir(1:ngtot_,i)
59  call zfftifc(3,ngridg_,-1,zfft1(:,i))
60 end do
61 zfft2(1:nfgrz_,1:nd)=0.d0
62 ! loop over crystal symmetries
63 do isym=1,nsymcrys
64 ! zero translation vector flag
65  tv0=tv0symc(isym)
66 ! translation vector in Cartesian coordinates
67  if (.not.tv0) then
68  v1=vtcsymc(1,isym)
69  v2=vtcsymc(2,isym)
70  v3=vtcsymc(3,isym)
71  end if
72 ! index to spatial rotation lattice symmetry
73  lspl=lsplsymc(isym)
74  sym(1:3,1:3)=symlat(1:3,1:3,lspl)
75  if (tspin) then
76 ! global spin proper rotation in Cartesian coordinates
77  lspn=lspnsymc(isym)
78  sc(1:3,1:3)=symlatd(lspn)*symlatc(1:3,1:3,lspn)
79  else
80 ! set spin rotation equal to spatial rotation
81  lspn=lspl
82  sc(1:3,1:3)=symlatc(1:3,1:3,lspl)
83  end if
84  do ifg=1,nfgrz_
85  ig=igrzf_(ifg)
86  if (ig > ngvec_) cycle
87 ! multiply the transpose of the inverse symmetry matrix with the G-vector
88  if (lspl == 1) then
89  jg=ig
90  else
91  i1=ivg(1,ig); i2=ivg(2,ig); i3=ivg(3,ig)
92  j1=sym(1,1)*i1+sym(2,1)*i2+sym(3,1)*i3
93  j2=sym(1,2)*i1+sym(2,2)*i2+sym(3,2)*i3
94  j3=sym(1,3)*i1+sym(2,3)*i2+sym(3,3)*i3
95  jg=ivgig(j1,j2,j3)
96  end if
97  jfg=igfft_(jg)
98 ! translation, spatial rotation and global spin rotation
99  if (tv0) then
100 ! zero translation vector
101  if (lspn == 1) then
102 ! global spin symmetry is the identity
103  zfft2(ifg,1:nd)=zfft2(ifg,1:nd)+zfft1(jfg,1:nd)
104  else
105  if (tnc) then
106 ! non-collinear case
107  z1=zfft1(jfg,1); z2=zfft1(jfg,2); z3=zfft1(jfg,3)
108  zfft2(ifg,1)=zfft2(ifg,1)+sc(1,1)*z1+sc(1,2)*z2+sc(1,3)*z3
109  zfft2(ifg,2)=zfft2(ifg,2)+sc(2,1)*z1+sc(2,2)*z2+sc(2,3)*z3
110  zfft2(ifg,3)=zfft2(ifg,3)+sc(3,1)*z1+sc(3,2)*z2+sc(3,3)*z3
111  else
112 ! collinear case
113  zfft2(ifg,1)=zfft2(ifg,1)+sc(3,3)*zfft1(jfg,1)
114  end if
115  end if
116  else
117 ! complex phase factor for translation
118  t1=vgc(1,jg)*v1+vgc(2,jg)*v2+vgc(3,jg)*v3
119  z0=cmplx(cos(t1),-sin(t1),8)
120  if (lspn == 1) then
121  zfft2(ifg,1:nd)=zfft2(ifg,1:nd)+z0*zfft1(jfg,1:nd)
122  else
123  if (tnc) then
124  z1=zfft1(jfg,1); z2=zfft1(jfg,2); z3=zfft1(jfg,3)
125  zfft2(ifg,1)=zfft2(ifg,1)+z0*(sc(1,1)*z1+sc(1,2)*z2+sc(1,3)*z3)
126  zfft2(ifg,2)=zfft2(ifg,2)+z0*(sc(2,1)*z1+sc(2,2)*z2+sc(2,3)*z3)
127  zfft2(ifg,3)=zfft2(ifg,3)+z0*(sc(3,1)*z1+sc(3,2)*z2+sc(3,3)*z3)
128  else
129  zfft2(ifg,1)=zfft2(ifg,1)+sc(3,3)*z0*zfft1(jfg,1)
130  end if
131  end if
132  end if
133  end do
134 end do
135 ! Fourier transform to real-space and normalise
136 t1=1.d0/dble(nsymcrys)
137 do i=1,nd
138  call rzfftifc(3,ngridg_,1,rvfir(:,i),zfft2(:,i))
139  rvfir(1:ngtot_,i)=t1*rvfir(1:ngtot_,i)
140 end do
141 deallocate(zfft1,zfft2)
142 end subroutine
143 !EOC
144 
integer, dimension(maxsymcrys) lspnsymc
Definition: modmain.f90:366
subroutine symrvfir(tspin, tnc, ngridg_, ngtot_, ngvec_, nfgrz_, igfft_, igrzf_, ld, rvfir)
Definition: symrvfir.f90:11
integer, dimension(:,:,:), allocatable ivgig
Definition: modmain.f90:402
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
real(8), dimension(3, 3, 48) symlatc
Definition: modmain.f90:350
integer, dimension(48) symlatd
Definition: modmain.f90:346
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