The Elk Code
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 
6 subroutine writeemd
7 use modmain
8 use modpw
9 use modmpi
10 use modomp
11 use moddelf
12 implicit none
13 ! local variables
14 integer ik,ihk,recl
15 integer ist,ispn,nthd
16 real(8) sm
17 complex(8) z1
18 ! allocatable arrays
19 real(8), allocatable :: emd(:)
20 complex(8), allocatable :: wfpw(:,:,:)
21 if (spinsprl) then
22  write(*,*)
23  write(*,'("Error(writeemd): electron momentum density not available for &
24  &spin-spirals")')
25  write(*,*)
26  stop
27 end if
28 ! initialise universal variables
29 call init0
30 call init1
31 call init4
32 ! read density and potentials from file
33 call readstate
34 ! find the new linearisation energies
35 call linengy
36 ! generate the APW radial functions
37 call genapwfr
38 ! generate the local-orbital radial functions
39 call genlofr
40 ! get the occupation numbers from file
41 do ik=1,nkpt
42  call getoccsv(filext,ik,vkl(:,ik),occsv(:,ik))
43 end do
44 ! delete existing EMD.OUT
45 if (mp_mpi) call delfile('EMD.OUT')
46 ! synchronise MPI processes
47 call mpi_barrier(mpicom,ierror)
48 allocate(emd(nhkmax))
49 inquire(iolength=recl) vkl(:,1),nhkmax,emd
50 deallocate(emd)
51 open(250,file='EMD.OUT',form='UNFORMATTED',access='DIRECT',recl=recl)
52 ! loop over k-points
53 call 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)
58 allocate(emd(nhkmax),wfpw(nhkmax,nspinor,nstsv))
59 !$OMP DO SCHEDULE(DYNAMIC)
60 do 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)
85 end do
86 !$OMP END DO
87 deallocate(emd,wfpw)
88 !$OMP END PARALLEL
89 call freethd(nthd)
90 close(250)
91 ! synchronise MPI processes
92 call mpi_barrier(mpicom,ierror)
93 if (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")')
97 end if
98 end subroutine
99 
real(8), dimension(:,:,:,:), allocatable vhkc
Definition: modpw.f90:48
character(256) filext
Definition: modmain.f90:1301
logical mp_mpi
Definition: modmpi.f90:17
integer nkpt
Definition: modmain.f90:461
subroutine genlofr
Definition: genlofr.f90:10
logical spinsprl
Definition: modmain.f90:283
Definition: modomp.f90:6
subroutine writeemd
Definition: writeemd.f90:7
integer np_mpi
Definition: modmpi.f90:13
complex(8), dimension(:,:,:,:), allocatable sfacgk
Definition: modmain.f90:509
subroutine linengy
Definition: linengy.f90:10
integer, dimension(:,:), allocatable nhk
Definition: modpw.f90:40
integer nstsv
Definition: modmain.f90:889
integer, dimension(:,:), allocatable ngk
Definition: modmain.f90:497
subroutine genwfpw(vpl, ngp, igpig, vgpl, vgpc, gpc, sfacgp, nhp, vhpc, hpc, sfachp, wfpw)
Definition: genwfpw.f90:7
real(8), dimension(:,:,:,:), allocatable vgkl
Definition: modmain.f90:503
real(8), dimension(:,:), allocatable occsv
Definition: modmain.f90:905
integer nspinor
Definition: modmain.f90:267
real(8), dimension(:,:,:,:), allocatable vgkc
Definition: modmain.f90:505
subroutine init1
Definition: init1.f90:10
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
real(8), dimension(:,:,:), allocatable hkc
Definition: modpw.f90:50
subroutine delfile(fname)
Definition: moddelf.f90:15
subroutine genapwfr
Definition: genapwfr.f90:10
Definition: modmpi.f90:6
subroutine getoccsv(fext, ikp, vpl, occsvp)
Definition: getoccsv.f90:7
real(8), dimension(:,:,:), allocatable gkc
Definition: modmain.f90:507
subroutine init4
Definition: init4.f90:7
Definition: modpw.f90:6
subroutine readstate
Definition: readstate.f90:10
integer lp_mpi
Definition: modmpi.f90:15
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
subroutine init0
Definition: init0.f90:10
integer, dimension(:,:,:), allocatable igkig
Definition: modmain.f90:501
complex(8), dimension(:,:,:,:), allocatable sfachk
Definition: modpw.f90:52
integer nhkmax
Definition: modpw.f90:42
integer mpicom
Definition: modmpi.f90:11
integer ierror
Definition: modmpi.f90:19