The Elk Code
 
Loading...
Searching...
No Matches
writeemd.f90
Go to the documentation of this file.
1
2! Copyright (C) 2012 S. Dugdale, D. Ernsting and J. K. Dewhurst.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine writeemd
7use modmain
8use modpw
9use modmpi
10use modomp
11use moddelf
12implicit none
13! local variables
14integer ik,ihk,recl
15integer ist,ispn,nthd
16real(8) sm
17complex(8) z1
18! allocatable arrays
19real(8), allocatable :: emd(:)
20complex(8), allocatable :: wfpw(:,:,:)
21if (spinsprl) then
22 write(*,*)
23 write(*,'("Error(writeemd): electron momentum density not available for &
24 &spin-spirals")')
25 write(*,*)
26 stop
27end if
28! initialise universal variables
29call init0
30call init1
31call init4
32! read density and potentials from file
33call readstate
34! find the new linearisation energies
35call linengy
36! generate the APW radial functions
37call genapwfr
38! generate the local-orbital radial functions
39call genlofr
40! get the occupation numbers from file
41do ik=1,nkpt
42 call getoccsv(filext,ik,vkl(:,ik),occsv(:,ik))
43end do
44! delete existing EMD.OUT
45if (mp_mpi) call delfile('EMD.OUT')
46! synchronise MPI processes
47call mpi_barrier(mpicom,ierror)
48allocate(emd(nhkmax))
49inquire(iolength=recl) vkl(:,1),nhkmax,emd
50deallocate(emd)
51open(250,file='EMD.OUT',form='UNFORMATTED',access='DIRECT',recl=recl)
52! loop over k-points
53call holdthd(nkpt/np_mpi,nthd)
54!$OMP PARALLEL DEFAULT(SHARED) &
55!$OMP PRIVATE(emd,wfpw,ihk,sm) &
56!$OMP PRIVATE(ist,ispn,z1) &
57!$OMP NUM_THREADS(nthd)
58allocate(emd(nhkmax),wfpw(nhkmax,nspinor,nstsv))
59!$OMP DO SCHEDULE(DYNAMIC)
60do ik=1,nkpt
61! distribute among MPI processes
62 if (mod(ik-1,np_mpi) /= lp_mpi) cycle
63!$OMP CRITICAL(writeemd_)
64 write(*,'("Info(writeemd): ",I6," of ",I6," k-points")') ik,nkpt
65!$OMP END CRITICAL(writeemd_)
66! Fourier transform the wavefunctions
67 call genwfpw(vkl(:,ik),ngk(1,ik),igkig(:,1,ik),vgkl(:,:,1,ik), &
68 vgkc(:,:,1,ik),gkc(:,1,ik),sfacgk(:,:,1,ik),nhk(1,ik),vhkc(:,:,1,ik), &
69 hkc(:,1,ik),sfachk(:,:,1,ik),wfpw)
70! loop over all H+k-vectors
71 do ihk=1,nhk(1,ik)
72! sum over occupied states and spins
73 sm=0.d0
74 do ist=1,nstsv
75 do ispn=1,nspinor
76 z1=wfpw(ihk,ispn,ist)
77 sm=sm+occsv(ist,ik)*(z1%re**2+z1%im**2)
78 end do
79 end do
80 emd(ihk)=sm
81 end do
82!$OMP CRITICAL(u250)
83 write(250,rec=ik) vkl(:,ik),nhk(1,ik),emd
84!$OMP END CRITICAL(u250)
85end do
86!$OMP END DO
87deallocate(emd,wfpw)
88!$OMP END PARALLEL
89call freethd(nthd)
90close(250)
91! synchronise MPI processes
92call mpi_barrier(mpicom,ierror)
93if (mp_mpi) then
94 write(*,*)
95 write(*,'("Info(writeemd): electron momentum density written to EMD.OUT")')
96 write(*,'(" for all H+k-vectors up to |H+k| < hkmax")')
97end if
98end subroutine
99
subroutine genapwfr
Definition genapwfr.f90:10
subroutine genlofr
Definition genlofr.f90:10
subroutine genwfpw(vpl, ngp, igpig, vgpl, vgpc, gpc, sfacgp, nhp, vhpc, hpc, sfachp, wfpw)
Definition genwfpw.f90:7
subroutine getoccsv(fext, ikp, vpl, occsvp)
Definition getoccsv.f90:7
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
subroutine init4
Definition init4.f90:7
subroutine linengy
Definition linengy.f90:10
subroutine delfile(fname)
Definition moddelf.f90:15
real(8), dimension(:,:,:,:), allocatable vgkc
Definition modmain.f90:505
real(8), dimension(:,:,:), allocatable gkc
Definition modmain.f90:507
integer nspinor
Definition modmain.f90:267
character(256) filext
Definition modmain.f90:1300
integer, dimension(:,:), allocatable ngk
Definition modmain.f90:497
integer, dimension(:,:,:), allocatable igkig
Definition modmain.f90:501
integer nkpt
Definition modmain.f90:461
real(8), dimension(:,:,:,:), allocatable vgkl
Definition modmain.f90:503
integer nstsv
Definition modmain.f90:886
logical spinsprl
Definition modmain.f90:283
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
complex(8), dimension(:,:,:,:), allocatable sfacgk
Definition modmain.f90:509
real(8), dimension(:,:), allocatable occsv
Definition modmain.f90:902
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
Definition modpw.f90:6
real(8), dimension(:,:,:), allocatable hkc
Definition modpw.f90:50
real(8), dimension(:,:,:,:), allocatable vhkc
Definition modpw.f90:48
integer nhkmax
Definition modpw.f90:42
integer, dimension(:,:), allocatable nhk
Definition modpw.f90:40
complex(8), dimension(:,:,:,:), allocatable sfachk
Definition modpw.f90:52
subroutine readstate
Definition readstate.f90:10
subroutine writeemd
Definition writeemd.f90:7