The Elk Code
 
Loading...
Searching...
No Matches
ggair_sp_1.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: ggair_sp_1
8! !INTERFACE:
9subroutine ggair_sp_1(rhoup,rhodn,grho,gup,gdn,g2up,g2dn,g3rho,g3up,g3dn)
10! !INPUT/OUTPUT PARAMETERS:
11! rhoup : spin-up density (in,real(ngtot))
12! rhodn : spin-down density (in,real(ngtot))
13! grho : |grad rho| (out,real(ngtot))
14! gup : |grad rhoup| (out,real(ngtot))
15! gdn : |grad rhodn| (out,real(ngtot))
16! g2up : grad^2 rhoup (out,real(ngtot))
17! g2dn : grad^2 rhodn (out,real(ngtot))
18! g3rho : (grad rho).(grad |grad rho|) (out,real(ngtot))
19! g3up : (grad rhoup).(grad |grad rhoup|) (out,real(ngtot))
20! g3dn : (grad rhodn).(grad |grad rhodn|) (out,real(ngtot))
21! !DESCRIPTION:
22! Computes $|\nabla\rho|$, $|\nabla\rho^{\uparrow}|$,
23! $|\nabla\rho^{\downarrow}|$, $\nabla^2\rho^{\uparrow}$,
24! $\nabla^2\rho^{\downarrow}$, $\nabla\rho\cdot(\nabla|\nabla\rho|)$,
25! $\nabla\rho^{\uparrow}\cdot(\nabla|\nabla\rho^{\uparrow}|)$ and
26! $\nabla\rho^{\downarrow}\cdot(\nabla|\nabla\rho^{\downarrow}|)$ for the
27! interstitial charge density, as required by the generalised gradient
28! approximation functionals of type 1 for spin-polarised densities. See
29! routines {\tt potxc} and {\tt modxcifc}.
30!
31! !REVISION HISTORY:
32! Created October 2004 (JKD)
33! Simplified and improved, October 2009 (JKD)
34!EOP
35!BOC
36use modmain
37implicit none
38! arguments
39real(8), intent(in) :: rhoup(ngtot),rhodn(ngtot)
40real(8), intent(out) :: grho(ngtot),gup(ngtot),gdn(ngtot)
41real(8), intent(out) :: g2up(ngtot),g2dn(ngtot)
42real(8), intent(out) :: g3rho(ngtot),g3up(ngtot),g3dn(ngtot)
43! local variables
44integer ig,ifg,i
45! allocatable arrays
46real(8), allocatable :: gvup(:,:),gvdn(:,:),rfir(:)
47complex(8), allocatable :: zfft1(:),zfft2(:)
48allocate(gvup(ngtot,3),gvdn(ngtot,3),rfir(ngtot))
49allocate(zfft1(nfgrz),zfft2(nfgrz))
50!----------------!
51! rho up !
52!----------------!
53call rzfftifc(3,ngridg,-1,rhoup,zfft1)
54! |grad rhoup|
55do i=1,3
56 do ifg=1,nfgrz
57 ig=igrzf(ifg)
58 if (ig <= ngvc) then
59 zfft2(ifg)=vgc(i,ig)*zi*zfft1(ifg)
60 else
61 zfft2(ifg)=0.d0
62 end if
63 end do
64 call rzfftifc(3,ngridg,1,gvup(:,i),zfft2)
65end do
66gup(:)=sqrt(gvup(:,1)**2+gvup(:,2)**2+gvup(:,3)**2)
67! grad^2 rhoup
68do ifg=1,nfgrz
69 ig=igrzf(ifg)
70 if (ig <= ngvc) then
71 zfft2(ifg)=-(gc(ig)**2)*zfft1(ifg)
72 else
73 zfft2(ifg)=0.d0
74 end if
75end do
76call rzfftifc(3,ngridg,1,g2up,zfft2)
77! (grad rhoup).(grad |grad rhoup|)
78call rzfftifc(3,ngridg,-1,gup,zfft1)
79g3up(:)=0.d0
80do i=1,3
81 do ifg=1,nfgrz
82 ig=igrzf(ifg)
83 zfft2(ifg)=vgc(i,ig)*zi*zfft1(ifg)
84 end do
85 call rzfftifc(3,ngridg,1,rfir,zfft2)
86 g3up(:)=g3up(:)+gvup(:,i)*rfir(:)
87end do
88!------------------!
89! rho down !
90!------------------!
91call rzfftifc(3,ngridg,-1,rhodn,zfft1)
92! |grad rhodn|
93do i=1,3
94 do ifg=1,nfgrz
95 ig=igrzf(ifg)
96 if (ig <= ngvc) then
97 zfft2(ifg)=vgc(i,ig)*zi*zfft1(ifg)
98 else
99 zfft2(ifg)=0.d0
100 end if
101 end do
102 call rzfftifc(3,ngridg,1,gvdn(:,i),zfft2)
103end do
104gdn(:)=sqrt(gvdn(:,1)**2+gvdn(:,2)**2+gvdn(:,3)**2)
105! grad^2 rhodn
106do ifg=1,nfgrz
107 ig=igrzf(ifg)
108 if (ig <= ngvc) then
109 zfft2(ifg)=-(gc(ig)**2)*zfft1(ifg)
110 else
111 zfft2(ifg)=0.d0
112 end if
113end do
114call rzfftifc(3,ngridg,1,g2dn,zfft2)
115! (grad rhodn).(grad |grad rhodn|)
116call rzfftifc(3,ngridg,-1,gdn,zfft1)
117g3dn(:)=0.d0
118do i=1,3
119 do ifg=1,nfgrz
120 ig=igrzf(ifg)
121 zfft2(ifg)=vgc(i,ig)*zi*zfft1(ifg)
122 end do
123 call rzfftifc(3,ngridg,1,rfir,zfft2)
124 g3dn(:)=g3dn(:)+gvdn(:,i)*rfir(:)
125end do
126!-------------!
127! rho !
128!-------------!
129! |grad rho|
130grho(:)=sqrt((gvup(:,1)+gvdn(:,1))**2 &
131 +(gvup(:,2)+gvdn(:,2))**2 &
132 +(gvup(:,3)+gvdn(:,3))**2)
133! (grad rho).(grad |grad rho|)
134call rzfftifc(3,ngridg,-1,grho,zfft1)
135g3rho(:)=0.d0
136do i=1,3
137 do ifg=1,nfgrz
138 ig=igrzf(ifg)
139 zfft2(ifg)=vgc(i,ig)*zi*zfft1(ifg)
140 end do
141 call rzfftifc(3,ngridg,1,rfir,zfft2)
142 g3rho(:)=g3rho(:)+(gvup(:,i)+gvdn(:,i))*rfir(:)
143end do
144deallocate(gvup,gvdn,rfir,zfft1,zfft2)
145end subroutine
146!EOC
147
subroutine ggair_sp_1(rhoup, rhodn, grho, gup, gdn, g2up, g2dn, g3rho, g3up, g3dn)
integer, dimension(3) ngridg
Definition modmain.f90:386
integer nfgrz
Definition modmain.f90:412
integer, dimension(:), allocatable igrzf
Definition modmain.f90:416
complex(8), parameter zi
Definition modmain.f90:1239
real(8), dimension(:,:), allocatable vgc
Definition modmain.f90:420
integer ngvc
Definition modmain.f90:398
real(8), dimension(:), allocatable gc
Definition modmain.f90:422
subroutine rzfftifc(nd, n, sgn, r, z)