The Elk Code
 
Loading...
Searching...
No Matches
rhostatic.f90
Go to the documentation of this file.
1
2! Copyright (C) 2020 J. K. Dewhurst and S. Sharma.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine rhostatic
7use modmain
8use modtddft
9use modmpi
10implicit none
11! local variables
12integer is,ia,ias,np,i
13real(8) t0,t1
14! allocatable arrays
15real(8), allocatable :: jrmt0(:,:,:),jrir0(:,:)
16! external functions
17real(8), external :: rfmtint,rfint
18! store original parameters
20afieldc0(:)=afieldc(:)
23! initialise global variables
24call init0
25! only one self-consistent loop required
26maxscl=1
27! read potential from STATE.OUT
28trdstate=.true.
29! calculate current density with zero A-field
30tjr=.true.
31afieldc(:)=0.d0
32! switch off forces
33tforce=.false.
34! run the ground-state calculation
35call gndstate
36! store the current density
37allocate(jrmt0(npmtmax,natmtot,3),jrir0(ngtot,3))
38do i=1,3
39 do ias=1,natmtot
40 is=idxis(ias)
41 jrmt0(1:npmt(is),ias,i)=jrmt(1:npmt(is),ias,i)
42 end do
43 jrir0(1:ngtot,i)=jrir(1:ngtot,i)
44end do
45! allocate static density and charge global arrays
46if (allocated(rhosmt)) deallocate(rhosmt)
47allocate(rhosmt(npmtmax,natmtot,3))
48if (allocated(rhosir)) deallocate(rhosir)
49allocate(rhosir(ngtot,3))
50if (allocated(chgsmt)) deallocate(chgsmt)
51allocate(chgsmt(natmtot,3))
52! magnitude of applied external A-field
53t0=0.75d0*solsc
54t1=solsc/t0
55! loop over three directions
56do i=1,3
57 afieldc(:)=0.d0
58 afieldc(i)=t0
59! run the ground-state calculation
60 call gndstate
61! muffin-tin static density
62 do ias=1,natmtot
63 is=idxis(ias)
64 np=npmt(is)
65 rhosmt(1:np,ias,i)=rhomt(1:np,ias)+t1*(jrmt0(1:np,ias,i)-jrmt(1:np,ias,i))
66! compute the muffin-tin static charge
67 chgsmt(ias,i)=rfmtint(nrmt(is),nrmti(is),wr2mt(:,is),rhosmt(:,ias,i))
68 end do
69! interstitial static density
70 rhosir(1:ngtot,i)=rhoir(1:ngtot)+t1*(jrir0(1:ngtot,i)-jrir(1:ngtot,i))
71! compute the static charge
72 chgstot(i)=rfint(rhosmt(:,:,i),rhosir(:,i))
73end do
74if (mp_mpi) then
75! write static density and charges to binary file
76 open(100,file='RHOSTAT.OUT',form='UNFORMATTED',action='WRITE')
77 write(100) natmtot
78 write(100) npmtmax
79 write(100) ngtot
80 write(100) rhosmt,rhosir
81 write(100) chgsmt
82 write(100) chgstot
83 close(100)
84! write static charges to text file
85 open(50,file='CHGSTAT.OUT',form='FORMATTED',action='WRITE')
86 write(50,'("Muffin-tin static charges :")')
87 do is=1,nspecies
88 write(50,'(" species : ",I4," (",A,")")') is,trim(spsymb(is))
89 write(50,'(" core charge",T25,": ",G18.10)') chgcr(is)
90 do ia=1,natoms(is)
91 ias=idxas(ia,is)
92 write(50,'(" atom ",I4,T25,": ",3G18.10)') ia,chgsmt(ias,:)
93 end do
94 end do
95 write(50,*)
96 write(50,'("Total static charge",T25,": ",3G18.10)') chgstot(:)
97 close(50)
98end if
99deallocate(jrmt0,jrir0)
100! restore original input parameters
102afieldc(:)=afieldc0(:)
104tjr=tjr0
105! synchronise MPI processes
106call mpi_barrier(mpicom,ierror)
107end subroutine
108
subroutine gndstate
Definition gndstate.f90:10
subroutine init0
Definition init0.f90:10
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
logical tjr0
Definition modmain.f90:620
integer ngtot
Definition modmain.f90:390
real(8), dimension(3) afieldc0
Definition modmain.f90:325
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
real(8), dimension(:), pointer, contiguous rhoir
Definition modmain.f90:614
real(8), dimension(3) afieldc
Definition modmain.f90:325
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
real(8), dimension(maxspecies) chgcr
Definition modmain.f90:716
logical tjr
Definition modmain.f90:620
integer natmtot
Definition modmain.f90:40
character(64), dimension(maxspecies) spsymb
Definition modmain.f90:78
integer maxscl0
Definition modmain.f90:1046
logical trdstate
Definition modmain.f90:682
integer npmtmax
Definition modmain.f90:216
real(8), dimension(:,:,:), allocatable jrmt
Definition modmain.f90:622
integer, dimension(maxspecies) npmt
Definition modmain.f90:213
real(8) solsc
Definition modmain.f90:1252
integer maxscl
Definition modmain.f90:1046
real(8), dimension(:,:), allocatable wr2mt
Definition modmain.f90:183
logical tforce
Definition modmain.f90:986
real(8), dimension(:,:), pointer, contiguous rhomt
Definition modmain.f90:614
integer nspecies
Definition modmain.f90:34
logical tforce0
Definition modmain.f90:986
real(8), dimension(:,:), allocatable jrir
Definition modmain.f90:622
integer ierror
Definition modmpi.f90:19
integer mpicom
Definition modmpi.f90:11
logical mp_mpi
Definition modmpi.f90:17
real(8), dimension(:,:,:), allocatable rhosmt
Definition modtddft.f90:82
real(8), dimension(:,:), allocatable chgsmt
Definition modtddft.f90:86
real(8), dimension(:,:), allocatable rhosir
Definition modtddft.f90:82
real(8), dimension(3) chgstot
Definition modtddft.f90:84
real(8) function rfint(rfmt, rfir)
Definition rfint.f90:7
pure real(8) function rfmtint(nr, nri, wr, rfmt)
Definition rfmtint.f90:7
subroutine rhostatic
Definition rhostatic.f90:7