File: symrfir.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (82 lines) | stat: -rw-r--r-- 2,331 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82

! Copyright (C) 2007 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

!BOP
! !ROUTINE: symrfir
! !INTERFACE:
subroutine symrfir(rfir)
! !USES:
use modmain
! !INPUT/OUTPUT PARAMETERS:
!   rfir : real intersitial function (inout,real(ngtot))
! !DESCRIPTION:
!   Symmetrises a real scalar interstitial function. The function is first
!   Fourier transformed to $G$-space, and then averaged over each symmetry by
!   rotating the Fourier coefficients and multiplying them by a phase factor
!   corresponding to the symmetry translation.
!
! !REVISION HISTORY:
!   Created July 2007 (JKD)
!EOP
!BOC
implicit none
! arguments
real(8), intent(inout) :: rfir(ngtot)
! local variables
logical tvz
integer isym,lspl,ilspl,sym(3,3)
integer iv(3),jv(3),ig,ifg,jfg
real(8) vtc(3),t1
complex(8) z1
! allocatable arrays
complex(8), allocatable :: zfft1(:),zfft2(:)
allocate(zfft1(ngtot),zfft2(ngtot))
! Fourier transform function to G-space
zfft1(:)=rfir(:)
call zfftifc(3,ngridg,-1,zfft1)
zfft2(:)=0.d0
! loop over crystal symmetries
do isym=1,nsymcrys
! translation in Cartesian coordinates
  vtc(:)=vtcsymc(:,isym)
! index to lattice symmetry of spatial rotation
  lspl=lsplsymc(isym)
! inverse rotation required for rotation of G-vectors
  ilspl=isymlat(lspl)
  sym(:,:)=symlat(:,:,ilspl)
! zero translation vector flag
  tvz=tv0symc(isym)
  do ig=1,ngvec
    ifg=igfft(ig)
! multiply the transpose of the inverse symmetry matrix with the G-vector
    if (lspl.eq.1) then
      jfg=ifg
    else
      iv(:)=ivg(:,ig)
      jv(1)=sym(1,1)*iv(1)+sym(2,1)*iv(2)+sym(3,1)*iv(3)
      jv(2)=sym(1,2)*iv(1)+sym(2,2)*iv(2)+sym(3,2)*iv(3)
      jv(3)=sym(1,3)*iv(1)+sym(2,3)*iv(2)+sym(3,3)*iv(3)
      jfg=igfft(ivgig(jv(1),jv(2),jv(3)))
    end if
    if (tvz) then
! zero translation vector
      zfft2(jfg)=zfft2(jfg)+zfft1(ifg)
    else
! complex phase factor for translation
      t1=-(vgc(1,ig)*vtc(1)+vgc(2,ig)*vtc(2)+vgc(3,ig)*vtc(3))
      z1=cmplx(cos(t1),sin(t1),8)
      zfft2(jfg)=zfft2(jfg)+z1*zfft1(ifg)
    end if
  end do
end do
! Fourier transform to real-space and normalise
call zfftifc(3,ngridg,1,zfft2)
t1=1.d0/dble(nsymcrys)
rfir(:)=t1*dble(zfft2(:))
deallocate(zfft1,zfft2)
return
end subroutine
!EOC