The Elk Code
Loading...
Searching...
No Matches
emdplot.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2014 D. Ernsting, S. Dugdale 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
emdplot
7
use
modmain
8
use
modpw
9
implicit none
10
! allocatable arrays
11
real
(4),
allocatable
:: emds(:,:)
12
if
(any(abs(
vkloff
(:)) >
epslat
))
then
13
write
(*,*)
14
write
(*,
'("Error(emdplot): use vkloff = 0 for the ground-state run")'
)
15
write
(*,*)
16
stop
17
end if
18
! initialise universal variables
19
call
init0
20
call
init1
21
call
init4
22
! read in the electron momentum density
23
allocate
(emds(
nhkmax
,
nkpt
))
24
call
reademd
(emds)
25
! write the density plot to file
26
select case
(
task
)
27
case
(171)
28
call
emdplot1d
(emds)
29
write
(*,*)
30
write
(*,
'("Info(emdplot): 1D electron momentum density written to &
31
&EMD1D.OUT")'
)
32
case
(172)
33
call
emdplot2d
(emds)
34
write
(*,*)
35
write
(*,
'("Info(emdplot): 2D electron momentum density written to &
36
&EMD2D.OUT")'
)
37
case
(173)
38
call
emdplot3d
(emds)
39
write
(*,*)
40
write
(*,
'("Info(emdplot): 3D electron momentum density written to &
41
&EMD3D.OUT")'
)
42
end select
43
deallocate
(emds)
44
end subroutine
45
emdplot1d
subroutine emdplot1d(emds)
Definition
emdplot1d.f90:7
emdplot2d
subroutine emdplot2d(emds)
Definition
emdplot2d.f90:7
emdplot3d
subroutine emdplot3d(emds)
Definition
emdplot3d.f90:7
emdplot
subroutine emdplot
Definition
emdplot.f90:7
init0
subroutine init0
Definition
init0.f90:10
init1
subroutine init1
Definition
init1.f90:10
init4
subroutine init4
Definition
init4.f90:7
modmain
Definition
modmain.f90:6
modmain::epslat
real(8) epslat
Definition
modmain.f90:24
modmain::nkpt
integer nkpt
Definition
modmain.f90:461
modmain::task
integer task
Definition
modmain.f90:1298
modmain::vkloff
real(8), dimension(3) vkloff
Definition
modmain.f90:450
modpw
Definition
modpw.f90:6
modpw::nhkmax
integer nhkmax
Definition
modpw.f90:42
reademd
subroutine reademd(emds)
Definition
reademd.f90:7
emdplot.f90
Generated by
1.9.8