The Elk Code
 
Loading...
Searching...
No Matches
gencfun.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: gencfun
8! !INTERFACE:
9subroutine gencfun
10! !USES:
11use modmain
12! !DESCRIPTION:
13! Generates the smooth characteristic function. This is the function which is
14! 0 within the muffin-tins and 1 in the intersitial region and is constructed
15! from radial step function form factors with $G<G_{\rm max}$. The form
16! factors are given by
17! $$ \tilde{\Theta}_i(G)=\begin{cases}
18! \frac{4\pi R_i^3}{3 \Omega} & G=0 \\
19! \frac{4\pi R_i^3}{\Omega}\frac{j_1(GR_i)}{GR_i} & 0<G\le G_{\rm max} \\
20! 0 & G>G_{\rm max}\end{cases} $$
21! where $R_i$ is the muffin-tin radius of the $i$th species and $\Omega$ is
22! the unit cell volume. Therefore the characteristic function in $G$-space is
23! $$ \tilde{\Theta}({\bf G})=\delta_{G,0}-\sum_{ij}\exp(-i{\bf G}\cdot
24! {\bf r}_{ij})\tilde{\Theta}_i(G), $$
25! where ${\bf r}_{ij}$ is the position of the $j$th atom of the $i$th species.
26!
27! !REVISION HISTORY:
28! Created January 2003 (JKD)
29!EOP
30!BOC
31implicit none
32! local variables
33integer is,ia,ig,ifg
34real(8) v1,v2,v3,t1
35! allocatable arrays
36complex(8), allocatable :: zfft(:)
37allocate(zfft(ngtot))
38zfft(1)=1.d0
39zfft(2:ngtot)=0.d0
40! begin loop over species
41do is=1,nspecies
42! loop over atoms
43 do ia=1,natoms(is)
44 v1=atposc(1,ia,is); v2=atposc(2,ia,is); v3=atposc(3,ia,is)
45 do ig=1,ngtot
46 ifg=igfft(ig)
47! structure factor
48 t1=vgc(1,ig)*v1+vgc(2,ig)*v2+vgc(3,ig)*v3
49! add to characteristic function in G-space
50 zfft(ifg)=zfft(ifg)-ffacg(ig,is)*cmplx(cos(t1),-sin(t1),8)
51 end do
52 end do
53end do
54! allocate global characteristic function arrays
55if (allocated(cfunig)) deallocate(cfunig)
56allocate(cfunig(ngvec))
57if (allocated(cfunir)) deallocate(cfunir)
58allocate(cfunir(ngtot))
59cfunig(1:ngvec)=zfft(igfft(1:ngvec))
60! Fourier transform to real-space
61call zfftifc(3,ngridg,1,zfft)
62cfunir(1:ngtot)=dble(zfft(1:ngtot))
63deallocate(zfft)
64end subroutine
65!EOC
66
subroutine gencfun
Definition gencfun.f90:10
integer ngtot
Definition modmain.f90:390
integer, dimension(3) ngridg
Definition modmain.f90:386
real(8), dimension(3, maxatoms, maxspecies) atposc
Definition modmain.f90:54
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
integer ngvec
Definition modmain.f90:396
real(8), dimension(:,:), allocatable ffacg
Definition modmain.f90:432
complex(8), dimension(:), allocatable cfunig
Definition modmain.f90:434
real(8), dimension(:), allocatable cfunir
Definition modmain.f90:436
integer, dimension(:), allocatable igfft
Definition modmain.f90:406
real(8), dimension(:,:), allocatable vgc
Definition modmain.f90:420
integer nspecies
Definition modmain.f90:34
subroutine zfftifc(nd, n, sgn, z)