The Elk Code
 
Loading...
Searching...
No Matches
writewfpw.f90
Go to the documentation of this file.
1
2! Copyright (C) 2010 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 writewfpw
7use modmain
8use modpw
9use modmpi
10use modomp
11use moddelf
12implicit none
13! local variables
14integer ik,recl,nthd
15! allocatable arrays
16complex(8), allocatable :: wfpw(:,:,:)
17! initialise global variables
18call init0
19call init1
20call init4
21! read density and potentials from file
22call readstate
23! find the new linearisation energies
24call linengy
25! generate the APW radial functions
26call genapwfr
27! generate the local-orbital radial functions
28call genlofr
29! delete existing WFPW.OUT
30if (mp_mpi) call delfile('WFPW.OUT')
31! synchronise MPI processes
32call mpi_barrier(mpicom,ierror)
33! determine the record length and open WFPW.OUT
34allocate(wfpw(nhkmax,nspinor,nstsv))
35inquire(iolength=recl) vkl(:,1),nhkmax,nspinor,nstsv,wfpw
36deallocate(wfpw)
37open(270,file='WFPW.OUT',form='UNFORMATTED',access='DIRECT',recl=recl)
38! begin parallel loop over k-points
39call holdthd(nkpt/np_mpi,nthd)
40!$OMP PARALLEL DEFAULT(SHARED) &
41!$OMP PRIVATE(wfpw) &
42!$OMP NUM_THREADS(nthd)
43allocate(wfpw(nhkmax,nspinor,nstsv))
44!$OMP DO SCHEDULE(DYNAMIC)
45do ik=1,nkpt
46! distribute among MPI processes
47 if (mod(ik-1,np_mpi) /= lp_mpi) cycle
48!$OMP CRITICAL(writewfpw_)
49 write(*,'("Info(writewfpw): ",I6," of ",I6," k-points")') ik,nkpt
50!$OMP END CRITICAL(writewfpw_)
51! generate the plane wave wavefunctions
52 call genwfpw(vkl(:,ik),ngk(:,ik),igkig(:,:,ik),vgkl(:,:,:,ik), &
53 vgkc(:,:,:,ik),gkc(:,:,ik),sfacgk(:,:,:,ik),nhk(:,ik),vhkc(:,:,:,ik), &
54 hkc(:,:,ik),sfachk(:,:,:,ik),wfpw)
55!$OMP CRITICAL(u270)
56 write(270,rec=ik) vkl(:,ik),nhkmax,nspinor,nstsv,wfpw
57!$OMP END CRITICAL(u270)
58end do
59!$OMP END DO
60deallocate(wfpw)
61!$OMP END PARALLEL
62call freethd(nthd)
63close(270)
64! synchronise MPI processes
65call mpi_barrier(mpicom,ierror)
66if (mp_mpi) then
67 write(*,*)
68 write(*,'("Info(writewfpw): plane wave wavefunctions written to WFPW.OUT")')
69 write(*,'(" for all H+k-vectors up to |H+k| < hkmax")')
70end if
71end subroutine
72
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 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
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
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
complex(8), dimension(:,:,:,:), allocatable sfacgk
Definition modmain.f90:509
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 writewfpw
Definition writewfpw.f90:7