The Elk Code
 
Loading...
Searching...
No Matches
genjr.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 genjr
7use modmain
8use modtddft
9implicit none
10! local variables
11integer is,ias,np,l
12real(8) ca,t1
13! generate the paramagnetic current
14call genjpr
15! add the diamagnetic term if a time-dependent A-field is present
16if (tafieldt) then
17! coupling constant of the external A-field (-1/c)
18 ca=-1.d0/solsc
19! muffin-tin part
20 do l=1,3
21 t1=ca*afieldt(l,itimes)
22 do ias=1,natmtot
23 is=idxis(ias)
24 np=npmt(is)
25 jrmt(1:np,ias,l)=jrmt(1:np,ias,l)+t1*(rhomt(1:np,ias)-rhosmt(1:np,ias,l))
26 end do
27 end do
28! interstitial part
29 do l=1,3
30 t1=ca*afieldt(l,itimes)
31 jrir(1:ngtot,l)=jrir(1:ngtot,l)+t1*(rhoir(1:ngtot)-rhosir(1:ngtot,l))
32 end do
33end if
34end subroutine
35
subroutine genjpr
Definition genjpr.f90:7
subroutine genjr
Definition genjr.f90:7
integer ngtot
Definition modmain.f90:390
real(8), dimension(:), pointer, contiguous rhoir
Definition modmain.f90:614
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer natmtot
Definition modmain.f90:40
real(8), dimension(:,:,:), allocatable jrmt
Definition modmain.f90:622
integer, dimension(maxspecies) npmt
Definition modmain.f90:213
real(8) solsc
Definition modmain.f90:1252
real(8), dimension(:,:), pointer, contiguous rhomt
Definition modmain.f90:614
real(8), dimension(:,:), allocatable jrir
Definition modmain.f90:622
logical tafieldt
Definition modtddft.f90:56
real(8), dimension(:,:,:), allocatable rhosmt
Definition modtddft.f90:82
real(8), dimension(:,:), allocatable afieldt
Definition modtddft.f90:58
real(8), dimension(:,:), allocatable rhosir
Definition modtddft.f90:82
integer itimes
Definition modtddft.f90:46