The Elk Code
gengqvec.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 gengqvec(iq)
7 use modmain
8 use modphonon
9 implicit none
10 ! arguments
11 integer, intent(in) :: iq
12 ! local variables
13 integer ig
14 do ig=1,ngvec
15 ! G+q-vector in Cartesian coordinates
16  vgqc(1:3,ig)=vgc(1:3,ig)+vqc(1:3,iq)
17 ! G+q-vector length
18  gqc(ig)=sqrt(vgqc(1,ig)**2+vgqc(2,ig)**2+vgqc(3,ig)**2)
19 ! spherical harmonics for G+q-vectors
20  call genylmv(.true.,lmaxo,vgqc(:,ig),ylmgq(:,ig))
21 end do
22 ! compute the spherical Bessel functions j_l(|G+q|R_mt)
24 ! structure factors for G+q
26 ! generate the smooth step function form factors for G+q
28 end subroutine
29 
pure subroutine gensfacgp(ngp, vgpc, ld, sfacgp)
Definition: gensfacgp.f90:10
real(8), dimension(:,:,:), allocatable jlgqrmt
Definition: modphonon.f90:68
integer lmaxo
Definition: modmain.f90:201
pure subroutine genylmv(t4pil, lmax, v, ylm)
Definition: genylmv.f90:10
real(8), dimension(:,:), allocatable vgc
Definition: modmain.f90:420
real(8), dimension(:,:), allocatable vqc
Definition: modmain.f90:547
integer ngvec
Definition: modmain.f90:396
real(8), dimension(:,:), allocatable vgqc
Definition: modphonon.f90:62
real(8), dimension(:,:), allocatable ffacgq
Definition: modphonon.f90:74
pure subroutine genffacgp(ngp, gpc, ld, ffacgp)
Definition: genffacgp.f90:10
complex(8), dimension(:,:), allocatable sfacgq
Definition: modphonon.f90:72
integer lnpsd
Definition: modmain.f90:628
subroutine genjlgprmt(lmax, ngp, gpc, ld, jlgprmt)
Definition: genjlgprmt.f90:10
real(8), dimension(:), allocatable gqc
Definition: modphonon.f90:64
subroutine gengqvec(iq)
Definition: gengqvec.f90:7
complex(8), dimension(:,:), allocatable ylmgq
Definition: modphonon.f90:70