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