The Elk Code
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:
9
subroutine
initw90
10
! !USES:
11
use
modmain
12
use
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
20
implicit none
21
! local variables
22
integer
ik,ist,ispn
23
integer
is,ias,l,m,i
24
! initialise universal variables
25
call
init0
26
call
init1
27
if
(
num_bands
>
nstsv
)
then
28
write
(*,*)
29
write
(*,
'("Error(initw90): num_bands > nstsv :",2(X,I0))'
)
num_bands
,
nstsv
30
write
(*,*)
31
stop
32
end if
33
! if num_bands is not positive then assume all states are used
34
if
(
num_bands
< 1)
then
35
if
(
allocated
(
idxw90
))
deallocate
(
idxw90
)
36
allocate
(
idxw90
(
nstsv
))
37
do
ist=1,
nstsv
38
idxw90
(ist)=ist
39
end do
40
num_bands
=
nstsv
41
end if
42
! check that each state index is in range
43
do
i=1,
num_bands
44
ist=
idxw90
(i)
45
if
((ist < 1).or.(ist >
nstsv
))
then
46
write
(*,*)
47
write
(*,
'("Error(initw90): state index out of range : ",I0)'
) ist
48
write
(*,*)
49
stop
50
end if
51
end do
52
if
(
projw90
)
then
53
! determine the map from the spin, l and atom projector to the Wannier function
54
if
(
allocated
(
prjwn
))
deallocate
(
prjwn
)
55
allocate
(
prjwn
(
nspinor
,0:3,
natmtot
))
56
i=0
57
do
ias=1,
natmtot
58
is=
idxis
(ias)
59
do
l=0,
lprojw90
(is)
60
do
m=-l,l
61
do
ispn=1,
nspinor
62
i=i+1
63
if
(m == -l)
prjwn
(ispn,l,ias)=i
64
end do
65
end do
66
end do
67
end do
68
num_wann
=i
69
else
if
(
num_wann
< 1)
then
70
num_wann
=
num_bands
+
num_wann
71
num_wann
=max(
num_wann
,1)
72
end if
73
if
(
num_wann
>
num_bands
)
then
74
write
(*,*)
75
write
(*,
'("Error(initw90): num_wann > num_bands :",2(X,I0))'
)
num_wann
, &
76
num_bands
77
write
(*,*)
78
stop
79
end if
80
! read density and potentials from file
81
call
readstate
82
! find the new linearisation energies
83
call
linengy
84
! generate the APW and local-orbital radial functions and integrals
85
call
genapwlofr
86
! read in the second-variational eigenvalues
87
do
ik=1,
nkpt
88
call
getevalsv
(
filext
,ik,
vkl
(:,ik),
evalsv
(:,ik))
89
end do
90
end subroutine
91
!EOC
92
modw90::num_bands
integer num_bands
Definition:
modw90.f90:22
modmain::filext
character(256) filext
Definition:
modmain.f90:1293
modmain::evalsv
real(8), dimension(:,:), allocatable evalsv
Definition:
modmain.f90:913
modmain::nkpt
integer nkpt
Definition:
modmain.f90:461
modw90::projw90
logical projw90
Definition:
modw90.f90:26
getevalsv
subroutine getevalsv(fext, ikp, vpl, evalsv_)
Definition:
getevalsv.f90:7
modw90::num_wann
integer num_wann
Definition:
modw90.f90:20
genapwlofr
subroutine genapwlofr
Definition:
genapwlofr.f90:7
linengy
subroutine linengy
Definition:
linengy.f90:10
initw90
subroutine initw90
Definition:
initw90.f90:10
modmain
Definition:
modmain.f90:6
modmain::nstsv
integer nstsv
Definition:
modmain.f90:883
modw90
Definition:
modw90.f90:6
modw90::idxw90
integer, dimension(:), allocatable idxw90
Definition:
modw90.f90:24
modmain::nspinor
integer nspinor
Definition:
modmain.f90:267
modw90::prjwn
integer, dimension(:,:,:), allocatable prjwn
Definition:
modw90.f90:30
init1
subroutine init1
Definition:
init1.f90:10
modmain::vkl
real(8), dimension(:,:), allocatable vkl
Definition:
modmain.f90:471
modmain::idxis
integer, dimension(maxatoms *maxspecies) idxis
Definition:
modmain.f90:44
readstate
subroutine readstate
Definition:
readstate.f90:10
init0
subroutine init0
Definition:
init0.f90:10
modmain::natmtot
integer natmtot
Definition:
modmain.f90:40
modw90::lprojw90
integer, dimension(maxspecies) lprojw90
Definition:
modw90.f90:28
initw90.f90
Generated by
1.8.14