The Elk Code
gengkqvec.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2014 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine gengkqvec(iq)
7 use modmain
8 use modphonon
9 implicit none
10 ! arguments
11 integer, intent(in) :: iq
12 ! local variables
13 integer ik,jspn
14 real(8) vl(3),vc(3)
15 ! loop over non-reduced k-point set
16 do ik=1,nkptnr
17 ! k+q-vectors in lattice and Cartesian coordinates
18  vkql(1:3,ik)=vkl(1:3,ik)+vql(1:3,iq)
19  vkqc(1:3,ik)=vkc(1:3,ik)+vqc(1:3,iq)
20  do jspn=1,nspnfv
21  vl(1:3)=vkql(1:3,ik)
22  vc(1:3)=vkqc(1:3,ik)
23 ! spin-spiral case
24  if (spinsprl) then
25  if (jspn == 1) then
26  vl(1:3)=vl(1:3)+0.5d0*vqlss(1:3)
27  vc(1:3)=vc(1:3)+0.5d0*vqcss(1:3)
28  else
29  vl(1:3)=vl(1:3)-0.5d0*vqlss(1:3)
30  vc(1:3)=vc(1:3)-0.5d0*vqcss(1:3)
31  end if
32  end if
33 ! generate the G+k+q-vectors
34  call gengkvec(ngvc,ivg,vgc,vl,vc,gkmax,ngkmax,ngkq(jspn,ik), &
35  igkqig(:,jspn,ik),vgkql(:,:,jspn,ik),vgkqc(:,:,jspn,ik),gkqc(:,jspn,ik))
36 ! generate structure factors for the G+k+q-vectors
37  call gensfacgp(ngkq(jspn,ik),vgkqc(:,:,jspn,ik),ngkmax,sfacgkq(:,:,jspn,ik))
38  end do
39 end do
40 end subroutine
41 
pure subroutine gensfacgp(ngp, vgpc, ld, sfacgp)
Definition: gensfacgp.f90:10
real(8), dimension(:,:,:), allocatable gkqc
Definition: modphonon.f90:90
integer ngkmax
Definition: modmain.f90:499
logical spinsprl
Definition: modmain.f90:283
real(8), dimension(:,:), allocatable vkqc
Definition: modphonon.f90:56
real(8), dimension(3) vqlss
Definition: modmain.f90:293
integer nkptnr
Definition: modmain.f90:463
integer ngvc
Definition: modmain.f90:398
subroutine gengkqvec(iq)
Definition: gengkqvec.f90:7
real(8), dimension(:,:), allocatable vkc
Definition: modmain.f90:473
real(8), dimension(:,:), allocatable vgc
Definition: modmain.f90:420
real(8), dimension(:,:), allocatable vqc
Definition: modmain.f90:547
complex(8), dimension(:,:,:,:), allocatable sfacgkq
Definition: modphonon.f90:92
real(8), dimension(:,:,:,:), allocatable vgkql
Definition: modphonon.f90:88
real(8), dimension(:,:), allocatable vql
Definition: modmain.f90:545
integer, dimension(:,:), allocatable ivg
Definition: modmain.f90:400
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
real(8), dimension(3) vqcss
Definition: modmain.f90:295
pure subroutine gengkvec(ngv, ivg, vgc, vkl, vkc, gkmax, ngkmax, ngk, igkig, vgkl, vgkc, gkc)
Definition: gengkvec.f90:11
integer, dimension(:,:), allocatable ngkq
Definition: modphonon.f90:84
real(8) gkmax
Definition: modmain.f90:495
real(8), dimension(:,:,:,:), allocatable vgkqc
Definition: modphonon.f90:88
real(8), dimension(:,:), allocatable vkql
Definition: modphonon.f90:54
integer nspnfv
Definition: modmain.f90:289
integer, dimension(:,:,:), allocatable igkqig
Definition: modphonon.f90:86