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,n,i,j
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
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,4,natmtot))
56  i=0
57  do ias=1,natmtot
58  is=idxis(ias)
59  n=0
60  do j=1,4
61  l=lprojw90(j,is)
62  if ((l < 0).or.(l > 3)) exit
63  n=n+1
64  do m=-l,l
65  do ispn=1,nspinor
66  i=i+1
67  if (m == -l) prjwn(ispn,j,ias)=i
68  end do
69  end do
70  end do
71  nprojw90(is)=n
72  end do
73  num_wann=i
74  if (num_wann == 0) then
75  write(*,*)
76  write(*,'("Error(init0): no projectors specified")')
77  write(*,*)
78  stop
79  end if
80 else if (num_wann < 1) then
82  num_wann=max(num_wann,1)
83 end if
84 if (num_wann > num_bands) then
85  write(*,*)
86  write(*,'("Error(initw90): num_wann > num_bands :",2(X,I0))') num_wann, &
87  num_bands
88  write(*,*)
89  stop
90 end if
91 ! read density and potentials from file
92 call readstate
93 ! find the new linearisation energies
94 call linengy
95 ! generate the APW and local-orbital radial functions and integrals
96 call genapwlofr
97 ! read in the second-variational eigenvalues
98 do ik=1,nkpt
99  call getevalsv(filext,ik,vkl(:,ik),evalsv(:,ik))
100 end do
101 end subroutine
102 !EOC
103 
integer num_bands
Definition: modw90.f90:22
character(256) filext
Definition: modmain.f90:1288
integer, dimension(maxspecies) nprojw90
Definition: modw90.f90:30
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:910
integer nkpt
Definition: modmain.f90:461
logical projw90
Definition: modw90.f90:26
subroutine getevalsv(fext, ikp, vpl, evalsv_)
Definition: getevalsv.f90:7
integer num_wann
Definition: modw90.f90:20
subroutine genapwlofr
Definition: genapwlofr.f90:7
subroutine linengy
Definition: linengy.f90:10
subroutine initw90
Definition: initw90.f90:10
integer nstsv
Definition: modmain.f90:880
Definition: modw90.f90:6
integer, dimension(:), allocatable idxw90
Definition: modw90.f90:24
integer nspinor
Definition: modmain.f90:267
integer, dimension(:,:,:), allocatable prjwn
Definition: modw90.f90:32
subroutine init1
Definition: init1.f90:10
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
subroutine readstate
Definition: readstate.f90:10
subroutine init0
Definition: init0.f90:10
integer natmtot
Definition: modmain.f90:40
integer, dimension(4, maxspecies) lprojw90
Definition: modw90.f90:28