The Elk Code
allatoms.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 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 !BOP
7 ! !ROUTINE: allatoms
8 ! !INTERFACE:
9 subroutine allatoms
10 ! !USES:
11 use modmain
12 use modxcifc
13 use modomp
14 ! !DESCRIPTION:
15 ! Solves the Kohn-Sham-Dirac equations for each atom type in the solid and
16 ! finds the self-consistent radial wavefunctions, eigenvalues, charge
17 ! densities and potentials. The atomic densities can then be used to
18 ! initialise the crystal densities, and the atomic self-consistent potentials
19 ! can be appended to the muffin-tin potentials to solve for the core states.
20 ! Note that, irrespective of the value of {\tt xctype}, exchange-correlation
21 ! functional type 3 is used. See also {\tt atoms}, {\tt rhoinit},
22 ! {\tt gencore} and {\tt modxcifc}.
23 !
24 ! !REVISION HISTORY:
25 ! Created September 2002 (JKD)
26 ! Modified for GGA, June 2007 (JKD)
27 !EOP
28 !BOC
29 implicit none
30 logical hybrid_
31 integer xcspin_,xcgrad_
32 integer is,nthd
33 real(8) hybridc_
34 character(264) xcdescr_
35 ! allocatable arrays
36 real(8), allocatable :: rwf(:,:,:)
37 ! allocate global species charge density and potential arrays
38 if (allocated(rhosp)) deallocate(rhosp)
39 allocate(rhosp(nrspmax,nspecies))
40 if (allocated(vrsp)) deallocate(vrsp)
41 allocate(vrsp(nrspmax,nspecies))
42 ! get the exchange-correlation functional data
43 call getxcdata(xctsp,xcdescr_,xcspin_,xcgrad_,hybrid_,hybridc_)
44 call holdthd(nspecies,nthd)
45 !$OMP PARALLEL DEFAULT(SHARED) &
46 !$OMP PRIVATE(rwf) &
47 !$OMP NUM_THREADS(nthd)
48 allocate(rwf(nrspmax,2,nstspmax))
49 !$OMP DO SCHEDULE(DYNAMIC)
50 do is=1,nspecies
51  call atom(solsc,ptnucl,spzn(is),nstsp(is),nsp(:,is),lsp(:,is),ksp(:,is), &
52  occsp(:,is),xctsp,xcgrad_,nrsp(is),rsp(:,is),evalsp(:,is),rhosp(:,is), &
53  vrsp(:,is),rwf)
54 end do
55 !$OMP END DO
56 deallocate(rwf)
57 !$OMP END PARALLEL
58 call freethd(nthd)
59 end subroutine
60 !EOC
61 
real(8), dimension(maxstsp, maxspecies) occsp
Definition: modmain.f90:133
real(8), dimension(:,:), allocatable rhosp
Definition: modmain.f90:137
subroutine atom(sol, ptnucl, zn, nst, n, l, k, occ, xctype, xcgrad, nr, r, eval, rho, vr, rwf)
Definition: atom.f90:10
integer, dimension(maxstsp, maxspecies) ksp
Definition: modmain.f90:125
integer, dimension(maxstsp, maxspecies) lsp
Definition: modmain.f90:123
integer, dimension(3) xctsp
Definition: modmain.f90:142
Definition: modomp.f90:6
integer nrspmax
Definition: modmain.f90:109
subroutine getxcdata(xctype, xcdescr, xcspin, xcgrad, hybrid, hybridc)
Definition: modxcifc.f90:411
real(8), dimension(maxstsp, maxspecies) evalsp
Definition: modmain.f90:131
subroutine allatoms
Definition: allatoms.f90:10
real(8) solsc
Definition: modmain.f90:1253
integer, dimension(maxspecies) nrsp
Definition: modmain.f90:107
integer, dimension(maxstsp, maxspecies) nsp
Definition: modmain.f90:121
logical ptnucl
Definition: modmain.f90:83
integer nspecies
Definition: modmain.f90:34
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer, dimension(maxspecies) nstsp
Definition: modmain.f90:113
real(8), dimension(maxspecies) spzn
Definition: modmain.f90:80
real(8), dimension(:,:), allocatable rsp
Definition: modmain.f90:135
real(8), dimension(:,:), allocatable vrsp
Definition: modmain.f90:139
integer nstspmax
Definition: modmain.f90:115