The Elk Code
 
Loading...
Searching...
No Matches
gwspecf.f90
Go to the documentation of this file.
1
2! Copyright (C) 2016 A. Davydov, A. Sanna, J. K. Dewhurst, S. Sharma and
3! E. K. U. Gross. This file is distributed under the terms of the GNU General
4! Public License. See the file COPYING for license details.
5
6subroutine gwspecf
7use modmain
8use modgw
9use modmpi
10use modomp
11use modtest
12implicit none
13! local variables
14integer ik,iw,nthd
15real(8) dw
16! allocatable arrays
17real(8), allocatable :: wr(:),sft(:),sf(:)
18complex(8), allocatable :: se(:,:,:)
19! initialise universal variables
20call init0
21call init1
22call init2
23call init3
24! read Fermi energy from file
25call readfermi
26! get the eigenvalues and occupation numbers from file
27do ik=1,nkpt
28 call getevalsv(filext,ik,vkl(:,ik),evalsv(:,ik))
29 call getoccsv(filext,ik,vkl(:,ik),occsv(:,ik))
30end do
31! real axis frequencies
32allocate(wr(nwplot))
33dw=(wplot(2)-wplot(1))/dble(nwplot)
34do iw=1,nwplot
35 wr(iw)=dw*dble(iw-1)+wplot(1)
36end do
37! allocate and zero the total spectral function
38allocate(sft(nwplot))
39sft(1:nwplot)=0.d0
40if (mp_mpi) write(*,*)
41! synchronise MPI processes
42call mpi_barrier(mpicom,ierror)
43! loop over reduced k-point set
44call holdthd(nkpt/np_mpi,nthd)
45!$OMP PARALLEL DEFAULT(SHARED) &
46!$OMP PRIVATE(sf,se) &
47!$OMP NUM_THREADS(nthd)
48allocate(sf(nwplot),se(nstsv,nstsv,0:nwfm))
49!$OMP DO SCHEDULE(DYNAMIC)
50do ik=1,nkpt
51! distribute among MPI processes
52 if (mod(ik-1,np_mpi) /= lp_mpi) cycle
53!$OMP CRITICAL(gwspecf_1)
54 write(*,'("Info(gwspecf): ",I6," of ",I6," k-points")') ik,nkpt
55!$OMP END CRITICAL(gwspecf_1)
56! get the self-energy at the fermionic frequencies from file
57 call getgwsefm(ik,se)
58! solve the Dyson equation on the real axis
59 call dysonr(ik,wr,se,sf)
60! write the spectral function to file
61!$OMP CRITICAL(gwspecf_2)
62 call writegwsf(ik,wr,sf)
63!$OMP END CRITICAL(gwspecf_2)
64! add to the total spectral function
65!$OMP CRITICAL(gwspecf_3)
66 sft(1:nwplot)=sft(1:nwplot)+wkpt(ik)*sf(1:nwplot)
67!$OMP END CRITICAL(gwspecf_3)
68end do
69!$OMP END DO
70deallocate(sf,se)
71!$OMP END PARALLEL
72call freethd(nthd)
73! add total spectral function from each process
74if (np_mpi > 1) then
75 call mpi_allreduce(mpi_in_place,sft,nwplot,mpi_double_precision,mpi_sum, &
77end if
78! write the total spectral function to file (MPI master process only)
79if (mp_mpi) then
80 open(50,file='GWTSF.OUT',form='FORMATTED')
81 do iw=1,nwplot
82 write(50,'(2G18.10)') wr(iw),sft(iw)
83 end do
84 close(50)
85 write(*,*)
86 write(*,'("Info(gw):")')
87 write(*,'(" GW spectral functions and Kohn-Sham eigenvalues written to &
88 &GWSF_Kkkkkkk.OUT")')
89 write(*,'(" for all k-points")')
90 write(*,*)
91 write(*,'(" Total GW spectral function written to GWTSF.OUT")')
92 write(*,*)
93 write(*,'(" Fermi energy for the Kohn-Sham eigenvalues is at zero in plots")')
94 write(*,'(" Fermi energy for the GW spectral function is undetermined")')
95 write(*,*)
96 write(*,'(" Spectral function units are states/Hartree/unit cell")')
97end if
98! write the total GW spectral function to test file
99call writetest(610,'total GW spectral function',nv=nwplot,tol=5.d-2,rva=sft)
100deallocate(wr,sft)
101end subroutine
102
subroutine dysonr(ik, wr, sem, sf)
Definition dysonr.f90:7
subroutine getevalsv(fext, ikp, vpl, evalsv_)
Definition getevalsv.f90:7
subroutine getgwsefm(ik, se)
Definition getgwsefm.f90:7
subroutine getoccsv(fext, ikp, vpl, occsvp)
Definition getoccsv.f90:7
subroutine gwspecf
Definition gwspecf.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
Definition modgw.f90:6
integer nwfm
Definition modgw.f90:19
real(8), dimension(:), allocatable wkpt
Definition modmain.f90:475
character(256) filext
Definition modmain.f90:1300
integer nwplot
Definition modmain.f90:1070
integer nkpt
Definition modmain.f90:461
real(8), dimension(2) wplot
Definition modmain.f90:1076
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 writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition modtest.f90:16
subroutine readfermi
Definition readfermi.f90:10
subroutine writegwsf(ik, wr, sf)
Definition writegwsf.f90:7