The Elk Code
 
Loading...
Searching...
No Matches
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:
9subroutine allatoms
10! !USES:
11use modmain
12use modxcifc
13use 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
29implicit none
30logical hybrid_
31integer xcspin_,xcgrad_
32integer is,nthd
33real(8) hybridc_
34character(264) xcdescr_
35! allocatable arrays
36real(8), allocatable :: rwf(:,:,:)
37! allocate global species charge density and potential arrays
38if (allocated(rhosp)) deallocate(rhosp)
39allocate(rhosp(nrspmax,nspecies))
40if (allocated(vrsp)) deallocate(vrsp)
41allocate(vrsp(nrspmax,nspecies))
42! get the exchange-correlation functional data
43call getxcdata(xctsp,xcdescr_,xcspin_,xcgrad_,hybrid_,hybridc_)
44call holdthd(nspecies,nthd)
45!$OMP PARALLEL DEFAULT(SHARED) &
46!$OMP PRIVATE(rwf) &
47!$OMP NUM_THREADS(nthd)
48allocate(rwf(nrspmax,2,nstspmax))
49!$OMP DO SCHEDULE(DYNAMIC)
50do 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)
54end do
55!$OMP END DO
56deallocate(rwf)
57!$OMP END PARALLEL
58call freethd(nthd)
59end subroutine
60!EOC
61
subroutine allatoms
Definition allatoms.f90:10
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) lsp
Definition modmain.f90:123
real(8), dimension(:,:), allocatable rhosp
Definition modmain.f90:137
integer, dimension(maxstsp, maxspecies) nsp
Definition modmain.f90:121
logical ptnucl
Definition modmain.f90:83
integer nstspmax
Definition modmain.f90:115
real(8), dimension(:,:), allocatable rsp
Definition modmain.f90:135
real(8), dimension(maxstsp, maxspecies) occsp
Definition modmain.f90:133
integer, dimension(maxspecies) nstsp
Definition modmain.f90:113
integer, dimension(3) xctsp
Definition modmain.f90:142
integer, dimension(maxstsp, maxspecies) ksp
Definition modmain.f90:125
real(8) solsc
Definition modmain.f90:1252
real(8), dimension(:,:), allocatable vrsp
Definition modmain.f90:139
real(8), dimension(maxstsp, maxspecies) evalsp
Definition modmain.f90:131
integer, dimension(maxspecies) nrsp
Definition modmain.f90:107
integer nrspmax
Definition modmain.f90:109
integer nspecies
Definition modmain.f90:34
real(8), dimension(maxspecies) spzn
Definition modmain.f90:80
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine getxcdata(xctype, xcdescr, xcspin, xcgrad, hybrid, hybridc)
Definition modxcifc.f90:411