The Elk Code
Loading...
Searching...
No Matches
writew90win.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2015 Manh Duc Le, 2017-18 Arsenii Gerasimov, Yaroslav Kvashnin
3
! and Lars Nordstrom. This file is distributed under the terms of the GNU
4
! General Public License. See the file COPYING for license details.
5
6
!BOP
7
! !ROUTINE: writew90win
8
! !INTERFACE:
9
subroutine
writew90win
10
! !USES:
11
use
modmain
12
use
modw90
13
! !DESCRIPTION:
14
! Writes out a template {\tt seedname.win} file with the Wannier90 input
15
! parameters. Uses {\tt wannier} and {\tt wannierExtra} blocks.
16
!
17
! !REVISION HISTORY:
18
! Created January 2015 (Manh Duc Le)
19
! Modified, August 2018 (Arsenii Gerasimov)
20
! Modified, February 2019 (JKD)
21
!EOP
22
!BOC
23
implicit none
24
! local variables
25
integer
ik,is,ia,i
26
character(256)
fname
27
fname=trim(
seedname
)//
'.win'
28
open
(50,file=trim(fname),action=
'WRITE'
,form=
'FORMATTED'
)
29
! write the number of Wannier functions and bands
30
write
(50,
'("length_unit = bohr")'
)
31
write
(50,
'("num_wann = ",I8)'
)
num_wann
32
write
(50,
'("num_bands = ",I8)'
)
num_bands
33
write
(50,
'("num_iter = ",I8)'
)
num_iter
34
write
(50,
'("dis_num_iter = ",I8)'
)
dis_num_iter
35
write
(50,*)
36
write
(50,
'("trial_step = ",G18.10)'
)
trial_step
37
! write spinors if required
38
if
(
spinpol
)
then
39
write
(50,*)
40
write
(50,
'("spinors = true")'
)
41
write
(50,
'("spn_formatted = true")'
)
42
end if
43
! write lattice vectors
44
write
(50,*)
45
write
(50,
'("begin unit_cell_cart")'
)
46
write
(50,
'("bohr")'
)
47
write
(50,
'(3G18.10)'
)
avec
(:,1)
48
write
(50,
'(3G18.10)'
)
avec
(:,2)
49
write
(50,
'(3G18.10)'
)
avec
(:,3)
50
write
(50,
'("end unit_cell_cart")'
)
51
! writes atomic positions
52
write
(50,*)
53
write
(50,
'("begin atoms_frac")'
)
54
do
is=1,
nspecies
55
do
ia=1,
natoms
(is)
56
write
(50,
'(A5,3G18.10)'
) trim(
spsymb
(is)),
atposl
(:,ia,is)
57
end do
58
end do
59
write
(50,
'("end atoms_frac")'
)
60
! write the list of k-points
61
write
(50,*)
62
write
(50,
'("mp_grid = ",3I6)'
)
ngridk
63
write
(50,*)
64
write
(50,
'("begin kpoints")'
)
65
do
ik=1,
nkptnr
66
write
(50,
'(3G18.10)'
)
vkl
(:,ik)
67
end do
68
write
(50,
'("end kpoints")'
)
69
write
(50,*)
70
! write the extra lines
71
do
i=1,
nxlwin
72
write
(50,
'(A)'
) trim(
xlwin
(i))
73
end do
74
close
(50)
75
write
(*,*)
76
write
(*,
'("Info(writew90win): created file ",A)'
) trim(fname)
77
end subroutine
78
!EOC
79
modmain
Definition
modmain.f90:6
modmain::nkptnr
integer nkptnr
Definition
modmain.f90:463
modmain::natoms
integer, dimension(maxspecies) natoms
Definition
modmain.f90:36
modmain::spinpol
logical spinpol
Definition
modmain.f90:228
modmain::avec
real(8), dimension(3, 3) avec
Definition
modmain.f90:12
modmain::ngridk
integer, dimension(3) ngridk
Definition
modmain.f90:448
modmain::spsymb
character(64), dimension(maxspecies) spsymb
Definition
modmain.f90:78
modmain::vkl
real(8), dimension(:,:), allocatable vkl
Definition
modmain.f90:471
modmain::nspecies
integer nspecies
Definition
modmain.f90:34
modmain::atposl
real(8), dimension(3, maxatoms, maxspecies) atposl
Definition
modmain.f90:51
modw90
Definition
modw90.f90:6
modw90::num_wann
integer num_wann
Definition
modw90.f90:18
modw90::dis_num_iter
integer dis_num_iter
Definition
modw90.f90:26
modw90::xlwin
character(256), dimension(:), allocatable xlwin
Definition
modw90.f90:16
modw90::trial_step
real(8) trial_step
Definition
modw90.f90:28
modw90::nxlwin
integer nxlwin
Definition
modw90.f90:14
modw90::num_bands
integer num_bands
Definition
modw90.f90:20
modw90::num_iter
integer num_iter
Definition
modw90.f90:24
modw90::seedname
character(256) seedname
Definition
modw90.f90:12
writew90win
subroutine writew90win
Definition
writew90win.f90:10
writew90win.f90
Generated by
1.9.8