The Elk Code
 
Loading...
Searching...
No Matches
initw90.f90
Go to the documentation of this file.
1
2! Copyright (C) 2017-18 Arsenii Gerasimov, Yaroslav Kvashnin and Lars Nordstrom.
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: initw90
8! !INTERFACE:
9subroutine initw90
10! !USES:
11use modmain
12use modw90
13! !DESCRIPTION:
14! Initialises global variables for the Wannier90 interface.
15!
16! !REVISION HISTORY:
17! Created November 2018 (Arsenii Gerasimov)
18!EOP
19!BOC
20implicit none
21! local variables
22integer ik,ist,i
23! initialise universal variables
24call init0
25call init1
26if (num_bands > nstsv) then
27 write(*,*)
28 write(*,'("Error(initw90): num_bands > nstsv : ",2I8)') num_bands,nstsv
29 write(*,*)
30 stop
31end if
32! if num_bands is not positive then assume all states are used
33if (num_bands < 1) then
34 if (allocated(idxw90)) deallocate(idxw90)
35 allocate(idxw90(nstsv))
36 do ist=1,nstsv
37 idxw90(ist)=ist
38 end do
40end if
41! check that each state index is in range
42do i=1,num_bands
43 ist=idxw90(i)
44 if ((ist < 1).or.(ist > nstsv)) then
45 write(*,*)
46 write(*,'("Error(initw90): state index out of range : ",I8)') ist
47 write(*,*)
48 stop
49 end if
50end do
51! set the number of Wannier functions
52if (num_wann < 1) then
54 num_wann=max(num_wann,1)
55end if
56! read density and potentials from file
57call readstate
58! find the new linearisation energies
59call linengy
60! generate the APW and local-orbital radial functions and integrals
61call genapwlofr
62! read in the second-variational eigenvalues
63do ik=1,nkpt
64 call getevalsv(filext,ik,vkl(:,ik),evalsv(:,ik))
65end do
66end subroutine
67!EOC
68
subroutine genapwlofr
Definition genapwlofr.f90:7
subroutine getevalsv(fext, ikp, vpl, evalsv_)
Definition getevalsv.f90:7
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
subroutine initw90
Definition initw90.f90:10
subroutine linengy
Definition linengy.f90:10
character(256) filext
Definition modmain.f90:1300
integer nkpt
Definition modmain.f90:461
integer nstsv
Definition modmain.f90:886
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
integer num_wann
Definition modw90.f90:18
integer, dimension(:), allocatable idxw90
Definition modw90.f90:22
integer num_bands
Definition modw90.f90:20
subroutine readstate
Definition readstate.f90:10