The Elk Code
 
Loading...
Searching...
No Matches
efieldmt.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
6!BOP
7! !ROUTINE: efieldmt
8! !INTERFACE:
9subroutine efieldmt
10! !USES:
11use modmain
12! !DESCRIPTION:
13! Calculates the average electric field in each muffin-tin from the gradient
14! of the Coulomb potential:
15! \begin{align*}
16! {\bf E}_{\alpha}&\equiv-\frac{3}{4\pi R_{\alpha}^3}
17! \int_{{\rm MT}_\alpha} \nabla V_{\rm C}({\bf r})\,d^3r \\
18! &=-\frac{3}{4\pi R_{\alpha}^3}\int_{{\rm MT}_\alpha}
19! V_{\rm C}({\bf r})\,\hat{\bf n}\,dS,
20! \end{align*}
21! where $R_{\alpha}$ is the radius of muffin-tin $\alpha$.
22!
23! !REVISION HISTORY:
24! Created April 2024 (JKD)
25!EOP
26!BOC
27implicit none
28! local variables
29integer is,ias,nr,nri,i
30real(8) t1,t2
31! automatic arrays
32real(8) grfmt(npmtmax,3)
33! external functions
34real(8), external :: rfmtint
35do ias=1,natmtot
36 is=idxis(ias)
37 nr=nrmt(is)
38 nri=nrmti(is)
39 call gradrfmt(nr,nri,rlmt(:,-1,is),wcrmt(:,:,is),vclmt(:,ias),npmtmax,grfmt)
40 t1=-(fourpi/3.d0)*rmt(is)**3
41 do i=1,3
42 t2=rfmtint(nr,nri,wr2mt(:,is),grfmt(:,i))
43 efcmt(i,ias)=t2/t1
44 end do
45end do
46end subroutine
47!EOC
48
subroutine efieldmt
Definition efieldmt.f90:10
subroutine gradrfmt(nr, nri, ri, wcr, rfmt, ld, grfmt)
Definition gradrfmt.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
real(8), dimension(maxspecies) rmt
Definition modmain.f90:162
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
real(8), parameter fourpi
Definition modmain.f90:1231
integer natmtot
Definition modmain.f90:40
real(8), dimension(:,:), allocatable vclmt
Definition modmain.f90:624
real(8), dimension(:,:), allocatable wr2mt
Definition modmain.f90:183
real(8), dimension(:,:), allocatable efcmt
Definition modmain.f90:316
real(8), dimension(:,:,:), allocatable rlmt
Definition modmain.f90:179
pure real(8) function rfmtint(nr, nri, wr, rfmt)
Definition rfmtint.f90:7