The Elk Code
writeevbse.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 writeevbse
7 use modmain
8 implicit none
9 ! local variables
10 integer ik,a
11 integer ios,nmbse_
12 ! allocatable arrays
13 complex(8), allocatable :: w(:)
14 ! initialise global variables
15 call init0
16 call init1
17 ! read Fermi energy from a file
18 call readefm
19 ! get the eigenvalues from file
20 do ik=1,nkpt
21  call getevalsv(filext,ik,vkl(:,ik),evalsv(:,ik))
22 end do
23 ! generate the BSE state index arrays
24 call genidxbse
25 ! allocate global BSE arrays
26 if (allocated(evalbse)) deallocate(evalbse)
27 allocate(evalbse(nmbse))
28 if (allocated(hmlbse)) deallocate(hmlbse)
29 allocate(hmlbse(nmbse,nmbse))
30 ! read in BSE Hamiltonian matrix
31 open(140,file='HMLBSE.OUT',form='UNFORMATTED',action='READ',status='OLD', &
32  iostat=ios)
33 if (ios /= 0) then
34  write(*,*)
35  write(*,'("Error(writeevbse): error opening HMLBSE.OUT")')
36  write(*,*)
37  stop
38 end if
39 read(140) nmbse_
40 if (nmbse /= nmbse_) then
41  write(*,*)
42  write(*,'("Error(writeevbse): differing nmbse")')
43  write(*,'(" current : ",I6)') nmbse
44  write(*,'(" HMLBSE.OUT : ",I6)') nmbse_
45  write(*,*)
46  stop
47 end if
48 read(140) hmlbse
49 close(140)
50 write(*,'("Info(writeevbse): diagonalising the BSE Hamiltonian matrix")')
51 if (bsefull) then
52 ! full non-Hermitian matrix
53  allocate(w(nmbse))
54  call eveqnzg(nmbse,nmbse,hmlbse,w)
55  evalbse(:)=dble(w(:))
56 else
57 ! Hermitian block only
59 end if
60 ! write the BSE eigenvectors and eigenvalues to file
61 open(140,file='EVBSE.OUT',form='UNFORMATTED',action='WRITE')
62 write(140) nmbse
63 write(140) evalbse
64 write(140) hmlbse
65 close(140)
66 ! write the BSE eigenvalues to file
67 open(50,file='EIGVAL_BSE.OUT',form='FORMATTED',action='WRITE')
68 write(50,'(I6," : nmbse")') nmbse
69 if (bsefull) then
70  do a=1,nmbse
71  write(50,'(I6,2G18.10)') a,dble(w(a)),aimag(w(a))
72  end do
73  deallocate(w)
74 else
75  do a=1,nmbse
76  write(50,'(I6,G18.10)') a,evalbse(a)
77  end do
78 end if
79 close(50)
80 write(*,*)
81 write(*,'("Info(writeevbse):")')
82 write(*,'(" BSE eigenvectors and eigenvalues written to EVBSE.OUT")')
83 write(*,'(" BSE eigenvalues written to EIGVAL_BSE.OUT")')
84 ! deallocate global BSE arrays
85 deallocate(evalbse,hmlbse)
86 end subroutine
87 
character(256) filext
Definition: modmain.f90:1301
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:921
integer nkpt
Definition: modmain.f90:461
subroutine getevalsv(fext, ikp, vpl, evalsv_)
Definition: getevalsv.f90:7
subroutine writeevbse
Definition: writeevbse.f90:7
subroutine eveqnzh(n, ld, a, w)
Definition: eveqnzh.f90:7
subroutine readefm
Definition: readefm.f90:10
integer nmbse
Definition: modmain.f90:1193
subroutine init1
Definition: init1.f90:10
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
subroutine genidxbse
Definition: genidxbse.f90:7
complex(8), dimension(:,:), allocatable hmlbse
Definition: modmain.f90:1201
logical bsefull
Definition: modmain.f90:1206
subroutine init0
Definition: init0.f90:10
real(8), dimension(:), allocatable evalbse
Definition: modmain.f90:1203
subroutine eveqnzg(n, ld, a, w)
Definition: eveqnzg.f90:7