The Elk Code
 
Loading...
Searching...
No Matches
writehmlbse.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 writehmlbse
7use modmain
8use modmpi
9! sets up the BSE matrix and writes it to file
10implicit none
11! local variables
12integer ik,jk,ist,jst
13integer a,b,i,j,m,n
14real(8) t1
15! initialise global variables
16call init0
17call init1
18call init2
19call init3
20! read density and potentials from file
21call readstate
22! find the new linearisation energies
23call linengy
24! generate the APW radial functions
25call genapwfr
26! generate the local-orbital radial functions
27call genlofr
28! get the eigenvalues and occupation numbers from file
29do ik=1,nkpt
30 call getevalsv(filext,ik,vkl(:,ik),evalsv(:,ik))
31 call getoccsv(filext,ik,vkl(:,ik),occsv(:,ik))
32end do
33! check if system is metallic
34t1=minval(abs(0.5d0-occsv(:,:)/occmax))
35if (abs(t1-0.5d0) > 0.01d0) then
36 write(*,*)
37 write(*,'("Warning(writehmlbse): system is metallic, the BSE may fail")')
38 write(*,'("Try using a different vkloff or reducing swidth")')
39end if
40! generate the BSE state index arrays
41call genidxbse
42if (allocated(hmlbse)) deallocate(hmlbse)
43allocate(hmlbse(nmbse,nmbse))
44! synchronise MPI processes
45call mpi_barrier(mpicom,ierror)
46if (mp_mpi) then
47 write(*,*)
48 write(*,'("Info(writehmlbse): setting up BSE Hamiltonian matrix")')
49end if
50! zero the BSE Hamiltonian
51hmlbse(:,:)=0.d0
52! compute diagonal matrix elements
53do ik=1,nkptnr
54! distribute among MPI processes
55 if (mod(ik-1,np_mpi) /= lp_mpi) cycle
56 jk=ivkik(ivk(1,ik),ivk(2,ik),ivk(3,ik))
57 do i=1,nvbse
58 ist=istbse(i,ik)
59 do j=1,ncbse
60 jst=jstbse(j,ik)
61 a=ijkbse(i,j,ik)
62 hmlbse(a,a)=evalsv(jst,jk)-evalsv(ist,jk)
63 if (bsefull) then
64 b=a+nbbse
65 hmlbse(b,b)=-hmlbse(a,a)
66 end if
67 end do
68 end do
69end do
70! add the exchange matrix elements
71if (hxbse) call hmlxbse
72! add the direct matrix elements
73if (hdbse) call hmldbse
74! add matrices from all processes and redistribute
75if (np_mpi > 1) then
76! ensure that the number of elements transmitted by MPI is not larger than the
77! maximum packet size (assumed to be 1024³ bytes)
78 m=67108864/nmbse
79 do b=1,nmbse,m
80 n=min(nmbse-b+1,m)
81 n=nmbse*n
82 call mpi_allreduce(mpi_in_place,hmlbse(1,b),n,mpi_double_complex,mpi_sum, &
84 end do
85end if
86! write the BSE matrix to HMLBSE.OUT
87if (mp_mpi) then
88 open(140,file='HMLBSE.OUT',form='UNFORMATTED',action='WRITE')
89 write(140) nmbse
90 write(140) hmlbse
91 close(140)
92 write(*,*)
93 write(*,'("Info(writehmlbse): BSE Hamiltonian matrix written to HMLBSE.OUT")')
94end if
95! deallocate global BSE arrays
96deallocate(istbse,jstbse,ijkbse,hmlbse)
97end subroutine
98
subroutine genapwfr
Definition genapwfr.f90:10
subroutine genidxbse
Definition genidxbse.f90:7
subroutine genlofr
Definition genlofr.f90:10
subroutine getevalsv(fext, ikp, vpl, evalsv_)
Definition getevalsv.f90:7
subroutine getoccsv(fext, ikp, vpl, occsvp)
Definition getoccsv.f90:7
subroutine hmldbse
Definition hmldbse.f90:7
subroutine hmlxbse
Definition hmlxbse.f90:7
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
subroutine init2
Definition init2.f90:7
subroutine init3
Definition init3.f90:7
subroutine linengy
Definition linengy.f90:10
logical bsefull
Definition modmain.f90:1203
logical hxbse
Definition modmain.f90:1206
integer nkptnr
Definition modmain.f90:463
integer, dimension(:,:), allocatable istbse
Definition modmain.f90:1192
character(256) filext
Definition modmain.f90:1300
integer nbbse
Definition modmain.f90:1188
integer, dimension(:,:), allocatable ivk
Definition modmain.f90:465
complex(8), dimension(:,:), allocatable hmlbse
Definition modmain.f90:1198
integer nkpt
Definition modmain.f90:461
integer, dimension(:,:,:), allocatable ijkbse
Definition modmain.f90:1196
integer, dimension(:,:), allocatable jstbse
Definition modmain.f90:1194
logical hdbse
Definition modmain.f90:1206
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
real(8) occmax
Definition modmain.f90:898
integer, dimension(:,:,:), allocatable ivkik
Definition modmain.f90:467
integer nmbse
Definition modmain.f90:1190
real(8), dimension(:,:), allocatable occsv
Definition modmain.f90:902
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
integer ncbse
Definition modmain.f90:1176
integer nvbse
Definition modmain.f90:1176
integer lp_mpi
Definition modmpi.f90:15
integer ierror
Definition modmpi.f90:19
integer mpicom
Definition modmpi.f90:11
integer np_mpi
Definition modmpi.f90:13
logical mp_mpi
Definition modmpi.f90:17
subroutine readstate
Definition readstate.f90:10
subroutine writehmlbse