The Elk Code
 
Loading...
Searching...
No Matches
sfacinit.f90
Go to the documentation of this file.
1
2! Copyright (C) 2010 A. I. Baranov and F. Wagner.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine sfacinit
7use modmain
8use modpw
9implicit none
10! local variables
11logical trhonorm0
12integer ik,ist,is,ias
13! allocatable arrays
14real(8), allocatable :: occcr0(:,:)
15! initialise universal variables
16call init0
17call init1
18! read density and potentials from file
19call readstate
20! use existing density if wsfac is default
21if ((wsfac(1) <= -1.d6).and.(wsfac(2) >= 1.d6)) goto 10
22! make a copy of the core state occupation numbers
23allocate(occcr0(nstspmax,natmtot))
24occcr0(:,:)=occcr(:,:)
25! zero the core state occupation numbers for eigenvalues not in energy window
26do ias=1,natmtot
27 is=idxis(ias)
28 do ist=1,nstsp(is)
29 if (spcore(ist,is)) then
30 if ((evalcr(ist,ias) < wsfac(1)).or.(evalcr(ist,ias) > wsfac(2))) then
31 occcr(ist,ias)=0.d0
32 end if
33 end if
34 end do
35end do
36! generate the core wavefunctions and densities
37call gencore
38! restore the core state occupation numbers
39occcr(:,:)=occcr0(:,:)
40deallocate(occcr0)
41! find the new linearisation energies
42call linengy
43! generate the APW radial functions
44call genapwfr
45! generate the local-orbital radial functions
46call genlofr
47do ik=1,nkpt
48! get the eigenvalues and occupation numbers from file
49 call getevalsv(filext,ik,vkl(:,ik),evalsv(:,ik))
50 call getoccsv(filext,ik,vkl(:,ik),occsv(:,ik))
51! zero occupation numbers for eigenvalues not in energy window
52 do ist=1,nstsv
53 if ((evalsv(ist,ik) < wsfac(1)).or.(evalsv(ist,ik) > wsfac(2))) then
54 occsv(ist,ik)=0.d0
55 end if
56 end do
57end do
58! computed density should not be normalised
59trhonorm0=trhonorm
60trhonorm=.false.
61! generate the density and magnetisation
62call rhomag
63trhonorm=trhonorm0
6410 continue
65! generate the H-vectors
66call genhvec
67end subroutine
68
subroutine genapwfr
Definition genapwfr.f90:10
subroutine gencore
Definition gencore.f90:10
subroutine genhvec
Definition genhvec.f90:7
subroutine genlofr
Definition genlofr.f90:10
subroutine getevalsv(fext, ikp, vpl, evalsv_)
Definition getevalsv.f90:7
subroutine getoccsv(fext, ikp, vpl, occsvp)
Definition getoccsv.f90:7
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
subroutine linengy
Definition linengy.f90:10
real(8), dimension(:,:), allocatable evalcr
Definition modmain.f90:934
logical trhonorm
Definition modmain.f90:618
character(256) filext
Definition modmain.f90:1300
logical, dimension(maxstsp, maxspecies) spcore
Definition modmain.f90:127
integer nstspmax
Definition modmain.f90:115
real(8), dimension(2) wsfac
Definition modmain.f90:1105
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer, dimension(maxspecies) nstsp
Definition modmain.f90:113
integer nkpt
Definition modmain.f90:461
integer natmtot
Definition modmain.f90:40
integer nstsv
Definition modmain.f90:886
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
real(8), dimension(:,:), allocatable occsv
Definition modmain.f90:902
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
real(8), dimension(:,:), allocatable occcr
Definition modmain.f90:932
Definition modpw.f90:6
subroutine readstate
Definition readstate.f90:10
subroutine rhomag
Definition rhomag.f90:7
subroutine sfacinit
Definition sfacinit.f90:7