The Elk Code
Loading...
Searching...
No Matches
writepmat.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3
! This file is distributed under the terms of the GNU General Public License.
4
! See the file COPYING for license details.
5
6
!BOP
7
! !ROUTINE: writepmat
8
! !INTERFACE:
9
subroutine
writepmat
10
! !USES:
11
use
modmain
12
use
modmpi
13
use
modramdisk
14
! !DESCRIPTION:
15
! Calculates the momentum matrix elements using routine {\tt genpmat} and
16
! writes them to direct access file {\tt PMAT.OUT}.
17
!
18
! !REVISION HISTORY:
19
! Created November 2003 (Sharma)
20
!EOP
21
!BOC
22
implicit none
23
! ensure momentum matrix elements are written to disk
24
wrtdsk0
=
wrtdsk
25
wrtdsk
=.true.
26
! initialise universal variables
27
call
init0
28
call
init1
29
! read in the density and potentials from file
30
call
readstate
31
! find the new linearisation energies
32
call
linengy
33
! generate the APW radial functions
34
call
genapwfr
35
! generate the local-orbital radial functions
36
call
genlofr
37
! write the momentum matrix elements in the second-variational basis to file
38
call
genpmat
39
wrtdsk
=
wrtdsk0
40
if
(
mp_mpi
)
then
41
write
(*,*)
42
write
(*,
'("Info(writepmat):")'
)
43
write
(*,
'(" momentum matrix elements written to file PMAT.OUT")'
)
44
end if
45
end subroutine
46
!EOC
47
genapwfr
subroutine genapwfr
Definition
genapwfr.f90:10
genlofr
subroutine genlofr
Definition
genlofr.f90:10
genpmat
subroutine genpmat
Definition
genpmat.f90:7
init0
subroutine init0
Definition
init0.f90:10
init1
subroutine init1
Definition
init1.f90:10
linengy
subroutine linengy
Definition
linengy.f90:10
modmain
Definition
modmain.f90:6
modmpi
Definition
modmpi.f90:6
modmpi::mp_mpi
logical mp_mpi
Definition
modmpi.f90:17
modramdisk
Definition
modramdisk.f90:6
modramdisk::wrtdsk0
logical wrtdsk0
Definition
modramdisk.f90:15
modramdisk::wrtdsk
logical wrtdsk
Definition
modramdisk.f90:15
readstate
subroutine readstate
Definition
readstate.f90:10
writepmat
subroutine writepmat
Definition
writepmat.f90:10
writepmat.f90
Generated by
1.9.8