File: gengqvec.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 (34 lines) | stat: -rw-r--r-- 926 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

! Copyright (C) 2014 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

subroutine gengqvec(iq)
use modmain
use modphonon
implicit none
! arguments
integer, intent(in) :: iq
! local variables
integer is,ig
real(8) tp(2)
! loop over G-vectors
do ig=1,ngtot
! G+q-vector in Cartesian coordinates
  vgqc(:,ig)=vgc(:,ig)+vqc(:,iq)
! G+q-vector length and (theta, phi) coordinates
  call sphcrd(vgqc(:,ig),gqc(ig),tp)
! spherical harmonics for G+q-vectors
  call genylm(lmaxo,tp,ylmgq(:,ig))
end do
! compute the spherical Bessel functions j_l(|G+q|R_mt)
call genjlgprmt(lnpsd,ngvec,gqc,ngvec,jlgqrmt)
! structure factors for G+q
call gensfacgp(ngvec,vgqc,ngvec,sfacgq)
! generate the smooth step function form factors for G+q
do is=1,nspecies
  call genffacgp(is,gqc,ffacgq(:,is))
end do
return
end subroutine