The Elk Code
 
Loading...
Searching...
No Matches
gradrhomt.f90
Go to the documentation of this file.
1
2! Copyright (C) 2012 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine gradrhomt
7use modmain
8use modphonon
9implicit none
10! local variables
11integer nr,nri,np
12! allocatable arrays
13complex(8), allocatable :: zfmt(:),grhomt(:,:)
14allocate(zfmt(npmtmax),grhomt(npmtmax,3))
15! add gradient contribution from rigid shift of muffin-tin
16nr=nrmt(isph)
17nri=nrmti(isph)
18np=npmt(isph)
19! convert the density to complex spherical harmonic expansion
20call rtozfmt(nr,nri,rhomt(:,iasph),zfmt)
21! compute the gradient
22call gradzfmt(nr,nri,rlmt(:,-1,isph),wcrmt(:,:,isph),zfmt,npmtmax,grhomt)
23! store density derivative for displaced atom
24if (tlast) drhomt(1:np,natmtot+1)=drhomt(1:np,iasph)
25! subtract from the density derivative
26drhomt(1:np,iasph)=drhomt(1:np,iasph)-grhomt(1:np,ipph)
27deallocate(zfmt,grhomt)
28end subroutine
29
subroutine gradrhomt
Definition gradrhomt.f90:7
subroutine gradzfmt(nr, nri, ri, wcr, zfmt, ld, gzfmt)
Definition gradzfmt.f90:10
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
real(8), dimension(:,:,:), allocatable wcrmt
Definition modmain.f90:187
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
integer natmtot
Definition modmain.f90:40
logical tlast
Definition modmain.f90:1050
integer npmtmax
Definition modmain.f90:216
integer, dimension(maxspecies) npmt
Definition modmain.f90:213
real(8), dimension(:,:), pointer, contiguous rhomt
Definition modmain.f90:614
real(8), dimension(:,:,:), allocatable rlmt
Definition modmain.f90:179
integer ipph
Definition modphonon.f90:15
integer iasph
Definition modphonon.f90:15
integer isph
Definition modphonon.f90:15
complex(8), dimension(:,:), allocatable drhomt
Definition modphonon.f90:98
pure subroutine rtozfmt(nr, nri, rfmt, zfmt)
Definition rtozfmt.f90:7