The Elk Code
 
Loading...
Searching...
No Matches
curlrvf.f90
Go to the documentation of this file.
1
2! Copyright (C) 2018 T. Mueller, J. K. Dewhurst, S. Sharma 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 curlrvf(rvfmt,rvfir,curlmt,curlir)
7use modmain
8use modomp
9implicit none
10! arguments
11real(8), intent(in) :: rvfmt(npmtmax,natmtot,3),rvfir(ngtot,3)
12real(8), intent(out) :: curlmt(npmtmax,natmtot,3),curlir(ngtot,3)
13! local variables
14integer is,ias,np,i,nthd
15! allocatable arrays
16real(8), allocatable :: grfmt(:,:,:,:),grfir(:,:,:)
17allocate(grfmt(npmtmax,natmtot,3,3),grfir(ngtot,3,3))
18! compute the gradients
19call holdthd(3,nthd)
20!$OMP PARALLEL DEFAULT(SHARED) &
21!$OMP PRIVATE(i,ias,is,np) &
22!$OMP NUM_THREADS(nthd)
23!$OMP DO SCHEDULE(DYNAMIC)
24do i=1,3
25 call gradrf(rvfmt(:,:,i),rvfir(:,i),grfmt(:,:,:,i),grfir(:,:,i))
26end do
27!$OMP END DO
28! determine the muffin-tin and interstitial curl
29!$OMP SECTIONS
30!$OMP SECTION
31do ias=1,natmtot
32 is=idxis(ias)
33 np=npmt(is)
34 curlmt(1:np,ias,1)=grfmt(1:np,ias,2,3)-grfmt(1:np,ias,3,2)
35 curlmt(1:np,ias,2)=grfmt(1:np,ias,3,1)-grfmt(1:np,ias,1,3)
36 curlmt(1:np,ias,3)=grfmt(1:np,ias,1,2)-grfmt(1:np,ias,2,1)
37end do
38!$OMP SECTION
39curlir(1:ngtot,1)=grfir(1:ngtot,2,3)-grfir(1:ngtot,3,2)
40curlir(1:ngtot,2)=grfir(1:ngtot,3,1)-grfir(1:ngtot,1,3)
41curlir(1:ngtot,3)=grfir(1:ngtot,1,2)-grfir(1:ngtot,2,1)
42!$OMP END SECTIONS
43!$OMP END PARALLEL
44call freethd(nthd)
45deallocate(grfmt,grfir)
46end subroutine
47
subroutine curlrvf(rvfmt, rvfir, curlmt, curlir)
Definition curlrvf.f90:7
subroutine gradrf(rfmt, rfir, grfmt, grfir)
Definition gradrf.f90:7
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer, dimension(maxspecies) npmt
Definition modmain.f90:213
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106