File: writeemd.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (105 lines) | stat: -rw-r--r-- 2,772 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105

! Copyright (C) 2012 S. Dugdale, D. Ernsting and J. K. Dewhurst.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

subroutine writeemd
use modmain
use modpw
use modmpi
use modomp
implicit none
! local variables
integer ik,ihk,recl
integer ist,ispn,nthd
real(8) sum,t1
complex(8) z1
! allocatable arrays
real(8), allocatable :: emd(:)
complex(8), allocatable :: wfpw(:,:,:)
if (spinsprl) then
  write(*,*)
  write(*,'("Error(writeemd): electron momentum density not available for &
   &spin-spirals")')
  write(*,*)
  stop
end if
! initialise universal variables
call init0
call init1
call init4
! read density and potentials from file
call readstate
! read Fermi energy from file
call readfermi
! find the new linearisation energies
call linengy
! generate the APW radial functions
call genapwfr
! generate the local-orbital radial functions
call genlofr
! get the occupancies from file
do ik=1,nkpt
  call getoccsv(filext,ik,vkl(:,ik),occsv(:,ik))
end do
! delete existing EMD.OUT
if (mp_mpi) then
  open(160,file='EMD.OUT')
  close(160,status='DELETE')
end if
! synchronise MPI processes
call mpi_barrier(mpicom,ierror)
allocate(emd(nhkmax))
inquire(iolength=recl) vkl(:,1),nhkmax,emd
deallocate(emd)
open(160,file='EMD.OUT',form='UNFORMATTED',access='DIRECT',recl=recl)
! loop over k-points
call omp_hold(nkpt/np_mpi,nthd)
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP PRIVATE(emd,wfpw,ihk,sum) &
!$OMP PRIVATE(ist,ispn,z1,t1) &
!$OMP NUM_THREADS(nthd)
!$OMP DO
do ik=1,nkpt
! distribute among MPI processes
  if (mod(ik-1,np_mpi).ne.lp_mpi) cycle
  allocate(emd(nhkmax),wfpw(nhkmax,nspinor,nstsv))
!$OMP CRITICAL(writeemd_)
  write(*,'("Info(writeemd): ",I6," of ",I6," k-points")') ik,nkpt
!$OMP END CRITICAL(writeemd_)
! Fourier transform the wavefunctions
  call genwfpw(vkl(:,ik),ngk(1,ik),igkig(:,1,ik),vgkl(:,:,1,ik), &
   vgkc(:,:,1,ik),gkc(:,1,ik),tpgkc(:,:,1,ik),sfacgk(:,:,1,ik),nhk(1,ik), &
   vhkc(:,:,1,ik),hkc(:,1,ik),tphkc(:,:,1,ik),sfachk(:,:,1,ik),wfpw)
! loop over all H+k-vectors
  do ihk=1,nhk(1,ik)
! sum over occupied states and spins
    sum=0.d0
    do ist=1,nstsv
      do ispn=1,nspinor
        z1=wfpw(ihk,ispn,ist)
        t1=dble(z1)**2+aimag(z1)**2
        sum=sum+occsv(ist,ik)*t1
      end do
    end do
    emd(ihk)=sum
  end do
!$OMP CRITICAL(u160)
  write(160,rec=ik) vkl(:,ik),nhk(1,ik),emd
!$OMP END CRITICAL(u160)
  deallocate(emd,wfpw)
end do
!$OMP END DO
!$OMP END PARALLEL
call omp_free(nthd)
close(160)
! synchronise MPI processes
call mpi_barrier(mpicom,ierror)
if (mp_mpi) then
  write(*,*)
  write(*,'("Info(writeemd): electron momentum density written to EMD.OUT")')
  write(*,'(" for all H+k-vectors up to |H+k| < hkmax")')
end if
return
end subroutine