The Elk Code
 
Loading...
Searching...
No Matches
genjpr.f90
Go to the documentation of this file.
1
2! Copyright (C) 2018 S. Sharma, J. K. Dewhurst 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 genjpr
7use modmain
8use modmpi
9use modomp
10implicit none
11! local variables
12integer ik,is,ias,n,i,nthd
13! allocatable arrays
14real(8), allocatable :: jrmt_(:,:,:),jrir_(:,:)
15! current density cannot be computed if wavefunctions do not exist
16if (iscl < 1) then
17 jrmt(:,:,:)=0.d0
18 jrir(:,:)=0.d0
19 return
20end if
21! set the current density to zero
22allocate(jrmt_(npcmtmax,natmtot,3),jrir_(ngtc,3))
23jrmt_(:,:,:)=0.d0
24jrir_(:,:)=0.d0
25call holdthd(nkpt/np_mpi,nthd)
26!$OMP PARALLEL DO DEFAULT(SHARED) &
27!$OMP REDUCTION(+:jrmt_,jrir_) &
28!$OMP SCHEDULE(STATIC) &
29!$OMP NUM_THREADS(nthd)
30do ik=1,nkpt
31! distribute among MPI processes
32 if (mod(ik-1,np_mpi) /= lp_mpi) cycle
33 call genjprk(ik,jrmt_,jrir_)
34end do
35!$OMP END PARALLEL DO
36call freethd(nthd)
37! add current densities from each process and redistribute
38if (np_mpi > 1) then
40 call mpi_allreduce(mpi_in_place,jrmt_,n,mpi_double_precision,mpi_sum,mpicom, &
41 ierror)
42 n=ngtc*3
43 call mpi_allreduce(mpi_in_place,jrir_,n,mpi_double_precision,mpi_sum,mpicom, &
44 ierror)
45end if
46! copy current density to global arrays
47do i=1,3
48 do ias=1,natmtot
49 is=idxis(ias)
50 jrmt(1:npcmt(is),ias,i)=jrmt_(1:npcmt(is),ias,i)
51 end do
52end do
53jrir(1:ngtc,1:3)=jrir_(1:ngtc,1:3)
54deallocate(jrmt_,jrir_)
55! convert muffin-tin current density to spherical harmonics
56call holdthd(natmtot,nthd)
57!$OMP PARALLEL DEFAULT(SHARED) &
58!$OMP PRIVATE(ias,is,i) &
59!$OMP NUM_THREADS(nthd)
60do i=1,3
61!$OMP DO SCHEDULE(DYNAMIC)
62 do ias=1,natmtot
63 is=idxis(ias)
64 call rfshtip(nrcmt(is),nrcmti(is),jrmt(:,ias,i))
65 end do
66!$OMP END DO NOWAIT
67end do
68!$OMP END PARALLEL
69call freethd(nthd)
70! symmetrise the current density
71call symrvf(.false.,.true.,nrcmt,nrcmti,npcmt,ngdgc,ngtc,ngvc,nfgrzc,igfc, &
73! convert muffin-tin and interstitial current density from coarse to fine grids
74call holdthd(6,nthd)
75!$OMP PARALLEL DEFAULT(SHARED) &
76!$OMP NUM_THREADS(nthd)
77!$OMP DO SCHEDULE(DYNAMIC)
78do i=1,3
79 call rfmtctof(jrmt(:,:,i))
80end do
81!$OMP END DO NOWAIT
82!$OMP DO SCHEDULE(DYNAMIC)
83do i=1,3
84 call rfirctof(jrir(:,i),jrir(:,i))
85end do
86!$OMP END DO
87!$OMP END PARALLEL
88call freethd(nthd)
89end subroutine
90
subroutine genjpr
Definition genjpr.f90:7
subroutine genjprk(ik, jrmt_, jrir_)
Definition genjprk.f90:7
integer nfgrzc
Definition modmain.f90:414
integer ngtot
Definition modmain.f90:390
integer, dimension(3) ngdgc
Definition modmain.f90:388
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
integer, dimension(:), allocatable igfc
Definition modmain.f90:410
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer ngtc
Definition modmain.f90:392
integer nkpt
Definition modmain.f90:461
integer natmtot
Definition modmain.f90:40
integer npcmtmax
Definition modmain.f90:216
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
integer npmtmax
Definition modmain.f90:216
real(8), dimension(:,:,:), allocatable jrmt
Definition modmain.f90:622
integer, dimension(:), allocatable igrzfc
Definition modmain.f90:418
integer ngvc
Definition modmain.f90:398
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
integer iscl
Definition modmain.f90:1048
real(8), dimension(:,:), allocatable jrir
Definition modmain.f90:622
integer lp_mpi
Definition modmpi.f90:15
integer ierror
Definition modmpi.f90:19
integer mpicom
Definition modmpi.f90:11
integer np_mpi
Definition modmpi.f90:13
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine rfirctof(rfirc, rfir)
Definition rfirctof.f90:10
subroutine rfmtctof(rfmt)
Definition rfmtctof.f90:10
subroutine rfshtip(nr, nri, rfmt)
Definition rfshtip.f90:7
subroutine symrvf(tspin, tnc, nrmt_, nrmti_, npmt_, ngridg_, ngtot_, ngvec_, nfgrz_, igfft_, igrzf_, ld1, rvfmt, ld2, rvfir)
Definition symrvf.f90:11