The Elk Code
 
Loading...
Searching...
No Matches
gwdmat.f90
Go to the documentation of this file.
1
2! Copyright (C) 2018 J. K. Dewhurst, S. Sharma 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 gwdmat
7use modmain
8use modgw
9use modmpi
10use modomp
11implicit none
12! local variables
13integer ik,lp,nthd
14! initialise universal variables
15call init0
16call init1
17call init3
18! read Fermi energy from file
19call readfermi
20! get the eigenvalues from file
21do ik=1,nkpt
22 call getevalsv(filext,ik,vkl(:,ik),evalsv(:,ik))
23end do
24! determine the GW Fermi energy
25call gwefermi
26! compute the GW density matrices and write the natural orbitals and occupation
27! numbers to EVECSV.OUT and OCCSV.OUT, respectively
28call holdthd(nkpt/np_mpi,nthd)
29!$OMP PARALLEL DO DEFAULT(SHARED) &
30!$OMP SCHEDULE(DYNAMIC) &
31!$OMP NUM_THREADS(nthd)
32do ik=1,nkpt
33! distribute among MPI processes
34 if (mod(ik-1,np_mpi) /= lp_mpi) cycle
35 call gwdmatk(ik)
36end do
37!$OMP END PARALLEL DO
38call freethd(nthd)
39! broadcast occupation number array to every MPI process
40if (np_mpi > 1) then
41 do ik=1,nkpt
42 lp=mod(ik-1,np_mpi)
43 call mpi_bcast(occsv(:,ik),nstsv,mpi_double_precision,lp,mpicom,ierror)
44 end do
45end if
46if (mp_mpi) then
47! write the occupation numbers to file
48 do ik=1,nkpt
49 call putoccsv(filext,ik,occsv(:,ik))
50 end do
51 write(*,*)
52 write(*,'("Info(gwdmat):")')
53 write(*,'(" GW density matrices determined for each k-point")')
54 write(*,*)
55 write(*,'(" Natural orbitals and occupation numbers written to")')
56 write(*,'(" EVECSV.OUT and OCCSV.OUT, respectively")')
57end if
58end subroutine
59
subroutine getevalsv(fext, ikp, vpl, evalsv_)
Definition getevalsv.f90:7
subroutine gwdmat
Definition gwdmat.f90:7
subroutine gwdmatk(ik)
Definition gwdmatk.f90:7
subroutine gwefermi
Definition gwefermi.f90:7
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
subroutine init3
Definition init3.f90:7
Definition modgw.f90:6
character(256) filext
Definition modmain.f90:1300
integer nkpt
Definition modmain.f90:461
integer nstsv
Definition modmain.f90:886
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
real(8), dimension(:,:), allocatable occsv
Definition modmain.f90:902
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
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 holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine putoccsv(fext, ik, occsvp)
Definition putoccsv.f90:7
subroutine readfermi
Definition readfermi.f90:10