The Elk Code
 
Loading...
Searching...
No Matches
writeefg.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2006 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: writeefg
8! !INTERFACE:
9subroutine writeefg
10! !USES:
11use modmain
12use modtest
13! !DESCRIPTION:
14! Computes the electric field gradient (EFG) tensor for each atom, $\alpha$,
15! and writes it to the file {\tt EFG.OUT} along with its eigenvalues. The EFG
16! is defined by
17! $$ V^{\alpha}_{ij}\equiv\left.\frac{\partial^2 V'_{\rm C}({\bf r})}
18! {\partial{\bf r}_i\partial{\bf r}_j}\right|_{{\bf r}={\bf r}_{\alpha}}, $$
19! where $V'_{\rm C}$ is the Coulomb potential with the $l=m=0$ component
20! removed in each muffin-tin. The derivatives are computed explicitly using
21! the routine {\tt gradrfmt}.
22!
23! !REVISION HISTORY:
24! Created May 2004 (JKD)
25! Fixed serious problem, November 2006 (JKD)
26!EOP
27!BOC
28implicit none
29! local variables
30integer, parameter :: lwork=10
31integer is,ia,ias
32integer nr,nri,ir
33integer np,i,j,info
34real(8) efg(3,3),a(3,3)
35real(8) w(3),work(lwork)
36! allocatable arrays
37real(8), allocatable :: rfmt(:),grfmt1(:,:),grfmt2(:,:)
38if (lmaxi < 2) then
39 write(*,*)
40 write(*,'("Error(writeefg): lmaxi too small for calculating the EFG : ",&
41 &I4)') lmaxi
42 write(*,'(" Run the ground-state calculation again with lmaxi >= 2")')
43 write(*,*)
44 stop
45end if
46! initialise universal variables
47call init0
48! read density and potentials from file
49call readstate
50! allocate local arrays
51allocate(rfmt(npmtmax),grfmt1(npmtmax,3),grfmt2(npmtmax,3))
52open(50,file='EFG.OUT',form='FORMATTED')
53write(50,*)
54write(50,'("(electric field gradient tensor is in Cartesian coordinates)")')
55do is=1,nspecies
56 nr=nrmt(is)
57 nri=nrmti(is)
58 np=npmt(is)
59 do ia=1,natoms(is)
60 ias=idxas(ia,is)
61 write(50,*)
62 write(50,*)
63 write(50,'("Species : ",I4," (",A,"), atom : ",I4)') is,trim(spsymb(is)),ia
64! remove the l = m = 0 part of the potential
65 rfmt(1:np)=vclmt(1:np,ias)
66 i=1
67 do ir=1,nri
68 rfmt(i)=0.d0
69 i=i+lmmaxi
70 end do
71 do ir=nri+1,nr
72 rfmt(i)=0.d0
73 i=i+lmmaxo
74 end do
75! compute the gradient of the Coulomb potential
76 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),rfmt,npmtmax,grfmt1)
77 do i=1,3
78! compute the gradient of the gradient
79 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),grfmt1(:,i),npmtmax, &
80 grfmt2)
81 do j=1,3
82 efg(i,j)=grfmt2(1,j)*y00
83 end do
84 end do
85! symmetrise the EFG
86 do i=1,3
87 do j=i+1,3
88 efg(i,j)=0.5d0*(efg(i,j)+efg(j,i))
89 efg(j,i)=efg(i,j)
90 end do
91 end do
92 write(50,*)
93 write(50,'(" EFG tensor :")')
94 do i=1,3
95 write(50,'(3G18.10)') (efg(i,j),j=1,3)
96 end do
97 write(50,'(" trace : ",G18.10)') efg(1,1)+efg(2,2)+efg(3,3)
98! diagonalise the EFG
99 a(:,:)=efg(:,:)
100 call dsyev('N','U',3,a,3,w,work,lwork,info)
101 write(50,'(" eigenvalues :")')
102 write(50,'(3G18.10)') w
103 end do
104end do
105close(50)
106write(*,*)
107write(*,'("Info(writeefg): electric field gradient written to EFG.OUT")')
108deallocate(rfmt,grfmt1,grfmt2)
109! write EFG of last atom to test file
110call writetest(115,'electric field gradient',nv=9,tol=1.d-3,rva=efg)
111end subroutine
112!EOC
113
subroutine gradrfmt(nr, nri, ri, wcr, rfmt, ld, grfmt)
Definition gradrfmt.f90:10
subroutine init0
Definition init0.f90:10
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
real(8), parameter y00
Definition modmain.f90:1233
real(8), dimension(:,:,:), allocatable wcrmt
Definition modmain.f90:187
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
integer lmmaxi
Definition modmain.f90:207
character(64), dimension(maxspecies) spsymb
Definition modmain.f90:78
integer npmtmax
Definition modmain.f90:216
integer, dimension(maxspecies) npmt
Definition modmain.f90:213
real(8), dimension(:,:), allocatable vclmt
Definition modmain.f90:624
integer lmaxi
Definition modmain.f90:205
integer lmmaxo
Definition modmain.f90:203
integer nspecies
Definition modmain.f90:34
real(8), dimension(:,:,:), allocatable rlmt
Definition modmain.f90:179
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition modtest.f90:16
subroutine readstate
Definition readstate.f90:10
subroutine writeefg
Definition writeefg.f90:10