The Elk Code
 
Loading...
Searching...
No Matches
writeefield.f90
Go to the documentation of this file.
1
2! Copyright (C) 2024 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 writeefield(fnum)
7use modmain
8implicit none
9! arguments
10integer, intent(in) :: fnum
11! local variables
12integer is,ia,ias,i
13real(8) efc(3),t1
14if (natmtot == 0) return
15! determine the average electric field in each muffin-tin
16call efieldmt
17! write the electric fields to file
18write(fnum,*)
19write(fnum,'("Average electric field in each muffin-tin")')
20do is=1,nspecies
21 write(fnum,'(" species : ",I4," (",A,")")') is,trim(spsymb(is))
22 do ia=1,natoms(is)
23 ias=idxas(ia,is)
24 write(fnum,'(" atom ",I4,T30,": ",3G18.10)') ia,efcmt(:,ias)
25 end do
26end do
27! compute the average electric field
28do i=1,3
29 efc(i)=sum(efcmt(i,1:natmtot))/dble(natmtot)
30end do
31write(fnum,*)
32write(fnum,'("Average of muffin-tin electric fields :")')
33write(fnum,'(3G18.10)') efc
34t1=norm2(efc(1:3))
35write(fnum,'(" magnitude : ",G18.10)') t1
36write(fnum,'(" volts/nanometer : ",G18.10)') t1*ef_si/1.d9
37end subroutine
38
subroutine efieldmt
Definition efieldmt.f90:10
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
integer natmtot
Definition modmain.f90:40
character(64), dimension(maxspecies) spsymb
Definition modmain.f90:78
integer nspecies
Definition modmain.f90:34
real(8), dimension(:,:), allocatable efcmt
Definition modmain.f90:316
real(8), parameter ef_si
Definition modmain.f90:1272
subroutine writeefield(fnum)