The Elk Code
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 
6 subroutine curlrvf(rvfmt,rvfir,curlmt,curlir)
7 use modmain
8 use modomp
9 implicit none
10 ! arguments
11 real(8), intent(in) :: rvfmt(npmtmax,natmtot,3),rvfir(ngtot,3)
12 real(8), intent(out) :: curlmt(npmtmax,natmtot,3),curlir(ngtot,3)
13 ! local variables
14 integer is,ias,np,i,nthd
15 ! allocatable arrays
16 real(8), allocatable :: grfmt(:,:,:,:),grfir(:,:,:)
17 allocate(grfmt(npmtmax,natmtot,3,3),grfir(ngtot,3,3))
18 ! compute the gradients
19 call 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)
24 do i=1,3
25  call gradrf(rvfmt(:,:,i),rvfir(:,i),grfmt(:,:,:,i),grfir(:,:,i))
26 end do
27 !$OMP END DO
28 ! determine the muffin-tin and interstitial curl
29 !$OMP SECTIONS
30 !$OMP SECTION
31 do 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)
37 end do
38 !$OMP SECTION
39 curlir(1:ngtot,1)=grfir(1:ngtot,2,3)-grfir(1:ngtot,3,2)
40 curlir(1:ngtot,2)=grfir(1:ngtot,3,1)-grfir(1:ngtot,1,3)
41 curlir(1:ngtot,3)=grfir(1:ngtot,1,2)-grfir(1:ngtot,2,1)
42 !$OMP END SECTIONS
43 !$OMP END PARALLEL
44 call freethd(nthd)
45 deallocate(grfmt,grfir)
46 end subroutine
47 
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
Definition: modomp.f90:6
subroutine gradrf(rfmt, rfir, grfmt, grfir)
Definition: gradrf.f90:7
subroutine curlrvf(rvfmt, rvfir, curlmt, curlir)
Definition: curlrvf.f90:7
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78