The Elk Code
gensfacgp.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 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: gensfacgp
8 ! !INTERFACE:
9 pure subroutine gensfacgp(ngp,vgpc,ld,sfacgp)
10 ! !USES:
11 use modmain
12 ! !INPUT/OUTPUT PARAMETERS:
13 ! ngp : number of G+p-vectors (in,integer)
14 ! vgpc : G+p-vectors in Cartesian coordinates (in,real(3,*))
15 ! ld : leading dimension (in,integer)
16 ! sfacgp : structure factors of G+p-vectors (out,complex(ld,natmtot))
17 ! !DESCRIPTION:
18 ! Generates the atomic structure factors for a set of ${\bf G+p}$-vectors:
19 ! $$ S_{\alpha}({\bf G+p})=\exp(i({\bf G+p})\cdot{\bf r}_{\alpha}), $$
20 ! where ${\bf r}_{\alpha}$ is the position of atom $\alpha$.
21 !
22 ! !REVISION HISTORY:
23 ! Created January 2003 (JKD)
24 !EOP
25 !BOC
26 implicit none
27 ! arguments
28 integer, intent(in) :: ngp
29 real(8), intent(in) :: vgpc(3,ngp)
30 integer, intent(in) :: ld
31 complex(8), intent(out) :: sfacgp(ld,natmtot)
32 ! local variables
33 integer is,ia,ias,igp
34 real(8) v1,v2,v3,t1
35 do ias=1,natmtot
36  is=idxis(ias)
37  ia=idxia(ias)
38  v1=atposc(1,ia,is); v2=atposc(2,ia,is); v3=atposc(3,ia,is)
39 !$OMP SIMD PRIVATE(t1) SIMDLEN(8)
40  do igp=1,ngp
41  t1=vgpc(1,igp)*v1+vgpc(2,igp)*v2+vgpc(3,igp)*v3
42  sfacgp(igp,ias)=cmplx(cos(t1),sin(t1),8)
43  end do
44 end do
45 end subroutine
46 !EOC
47 
pure subroutine gensfacgp(ngp, vgpc, ld, sfacgp)
Definition: gensfacgp.f90:10
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer, dimension(maxatoms *maxspecies) idxia
Definition: modmain.f90:45
integer natmtot
Definition: modmain.f90:40
real(8), dimension(3, maxatoms, maxspecies) atposc
Definition: modmain.f90:54