The Elk Code
genidxbse.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2010 S. Sharma, J. K. Dewhurst and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine genidxbse
7 use modmain
8 implicit none
9 integer ik,jk,a,ntop
10 integer ist,jst,i,j,k
11 ! automatic arrays
12 integer idx(nstsv)
13 ! check if the BSE extra valence or conduction states are in range
14 do i=1,nvxbse
15  ist=istxbse(i)
16  if ((ist < 1).or.(ist > nstsv)) then
17  write(*,*)
18  write(*,'("Error(genidxbse): extra valence state out of range : ",I8)') ist
19  write(*,*)
20  stop
21  end if
22 end do
23 do j=1,ncxbse
24  jst=jstxbse(j)
25  if ((jst < 1).or.(jst > nstsv)) then
26  write(*,*)
27  write(*,'("Error(genidxbse): extra conduction state out of range : ",I8)') &
28  jst
29  write(*,*)
30  stop
31  end if
32 end do
33 ! number of valence states for transitions
35 ! number of conduction states for transitions
37 if ((nvbse < 1).or.(ncbse < 1)) then
38  write(*,*)
39  write(*,'("Error(genidxbse): invalid number of valence or conduction &
40  &transition states : ",2I8)') nvbse,ncbse
41  write(*,*)
42  stop
43 end if
44 ! total number of transitions
46 ! block size in BSE matrix
48 ! BSE matrix size
49 if (bsefull) then
50  nmbse=2*nbbse
51 else
52  nmbse=nbbse
53 end if
54 ! allocate global BSE index arrays
55 if (allocated(istbse)) deallocate(istbse)
56 allocate(istbse(nvbse,nkptnr))
57 if (allocated(jstbse)) deallocate(jstbse)
58 allocate(jstbse(ncbse,nkptnr))
59 if (allocated(ijkbse)) deallocate(ijkbse)
60 allocate(ijkbse(nvbse,ncbse,nkptnr))
61 a=0
62 ! loop over non-reduced k-points
63 do ik=1,nkptnr
64 ! equivalent reduced k-point
65  jk=ivkik(ivk(1,ik),ivk(2,ik),ivk(3,ik))
66 ! index for sorting the eigenvalues into ascending order
67  call sortidx(nstsv,evalsv(:,jk),idx)
68 ! find the topmost occupied band
69  ntop=nstsv
70  do ist=nstsv,1,-1
71  if (evalsv(idx(ist),jk) < efermi) then
72  ntop=ist
73  exit
74  end if
75  end do
76  if ((ntop-nvbse0+1) < 1) then
77  write(*,*)
78  write(*,'("Error(genidxbse): not enough valence states, reduce nvbse")')
79  write(*,*)
80  stop
81  end if
82  if ((ntop+ncbse0) > nstsv) then
83  write(*,*)
84  write(*,'("Error(genidxbse): not enough conduction states")')
85  write(*,'(" reduce ncbse or increase nempty")')
86  write(*,*)
87  stop
88  end if
89 ! index from BSE valence states to second-variational state numbers
90  do i=1,nvbse0
91  istbse(i,ik)=idx(ntop-nvbse0+i)
92  end do
93 ! index from BSE conduction states to second-variational state numbers
94  do j=1,ncbse0
95  jstbse(j,ik)=idx(ntop+j)
96  end do
97 ! add extra states to the list
98  do i=1,nvxbse
99  ist=istxbse(i)
100  if (evalsv(ist,jk) > efermi) then
101  write(*,*)
102  write(*,'("Error(genidxbse): extra valence state above Fermi energy : ",&
103  &I6)') ist
104  write(*,'(" for k-point ",I8)') jk
105  write(*,*)
106  stop
107  end if
108  do k=1,nvbse0+i-1
109  if (ist == istbse(k,ik)) then
110  write(*,*)
111  write(*,'("Error(genidxbse): redundant extra valence state : ",I6)') ist
112  write(*,'(" for k-point ",I8)') jk
113  write(*,*)
114  stop
115  end if
116  end do
117  istbse(nvbse0+i,ik)=ist
118  end do
119  do j=1,ncxbse
120  jst=jstxbse(j)
121  if (evalsv(jst,jk) < efermi) then
122  write(*,*)
123  write(*,'("Error(genidxbse): extra conduction state below Fermi &
124  &energy : ",I6)') jst
125  write(*,'(" for k-point ",I8)') jk
126  write(*,*)
127  stop
128  end if
129  do k=1,ncbse0+j-1
130  if (jst == jstbse(k,ik)) then
131  write(*,*)
132  write(*,'("Error(genidxbse): redundant extra conduction state : ",&
133  &I6)') jst
134  write(*,'(" for k-point ",I8)') jk
135  write(*,*)
136  stop
137  end if
138  end do
139  jstbse(ncbse0+j,ik)=jst
140  end do
141 ! index from BSE valence-conduction pair and k-point to location in BSE matrix
142  do i=1,nvbse
143  do j=1,ncbse
144  a=a+1
145  ijkbse(i,j,ik)=a
146  end do
147  end do
148 ! end loop over non-reduced k-points
149 end do
150 end subroutine
151 
real(8) efermi
Definition: modmain.f90:907
integer ncxbse
Definition: modmain.f90:1185
integer, dimension(maxxbse) jstxbse
Definition: modmain.f90:1187
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:921
integer nvcbse
Definition: modmain.f90:1189
integer, dimension(:,:,:), allocatable ijkbse
Definition: modmain.f90:1199
integer nbbse
Definition: modmain.f90:1191
integer, dimension(maxxbse) istxbse
Definition: modmain.f90:1187
integer, dimension(:,:,:), allocatable ivkik
Definition: modmain.f90:467
integer nvxbse
Definition: modmain.f90:1185
integer nkptnr
Definition: modmain.f90:463
integer, dimension(:,:), allocatable jstbse
Definition: modmain.f90:1197
integer ncbse0
Definition: modmain.f90:1181
integer nvbse0
Definition: modmain.f90:1181
pure subroutine sortidx(n, x, idx)
Definition: sortidx.f90:10
integer nmbse
Definition: modmain.f90:1193
subroutine genidxbse
Definition: genidxbse.f90:7
logical bsefull
Definition: modmain.f90:1206
integer nvbse
Definition: modmain.f90:1179
integer ncbse
Definition: modmain.f90:1179
integer, dimension(:,:), allocatable ivk
Definition: modmain.f90:465
integer, dimension(:,:), allocatable istbse
Definition: modmain.f90:1195