The Elk Code
gencfrm.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2010 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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 gencfrm(wfmt11,wfmt12,wfir11,wfir12,wfmt21,wfmt22,wfir21,wfir22, &
7  crhomt,crhoir,cmagmt,cmagir)
8 use modmain
9 use modomp
10 implicit none
11 ! arguments
12 complex(4), intent(in) :: wfmt11(npcmtmax,natmtot),wfmt12(npcmtmax,natmtot)
13 complex(4), intent(in) :: wfir11(ngtot),wfir12(ngtot)
14 complex(4), intent(in) :: wfmt21(npcmtmax,natmtot),wfmt22(npcmtmax,natmtot)
15 complex(4), intent(in) :: wfir21(ngtot),wfir22(ngtot)
16 complex(4), intent(out) :: crhomt(npcmtmax,natmtot),crhoir(ngtot)
17 complex(4), intent(out) :: cmagmt(npcmtmax,natmtot,ndmag),cmagir(ngtot,ndmag)
18 ! local variables
19 integer ld,is,ias,nthd
20 ld=npcmtmax*natmtot
21 call holdthd(natmtot+1,nthd)
22 !$OMP PARALLEL DEFAULT(SHARED) &
23 !$OMP PRIVATE(ias,is) &
24 !$OMP NUM_THREADS(nthd)
25 !-------------------------!
26 ! muffin-tin part !
27 !-------------------------!
28 !$OMP DO SCHEDULE(DYNAMIC)
29 do ias=1,natmtot
30  is=idxis(ias)
31  call gencrm(npcmt(is),wfmt11(:,ias),wfmt12(:,ias),wfmt21(:,ias), &
32  wfmt22(:,ias),crhomt(:,ias),ld,cmagmt(:,ias,1))
33 end do
34 !$OMP END DO NOWAIT
35 !---------------------------!
36 ! interstitial part !
37 !---------------------------!
38 !$OMP SINGLE
39 call gencrm(ngtot,wfir11,wfir12,wfir21,wfir22,crhoir,ngtot,cmagir)
40 !$OMP END SINGLE
41 !$OMP END PARALLEL
42 call freethd(nthd)
43 end subroutine
44 
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
subroutine gencfrm(wfmt11, wfmt12, wfir11, wfir12, wfmt21, wfmt22, wfir21, wfir22, crhomt, crhoir, cmagmt, cmagir)
Definition: gencfrm.f90:8
Definition: modomp.f90:6
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
pure subroutine gencrm(n, wf11, wf12, wf21, wf22, crho, ld, cmag)
Definition: gencrm.f90:7
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78