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:
9subroutine writew90win
10! !USES:
11use modmain
12use 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
23implicit none
24! local variables
25integer ik,is,ia,i
26character(256) fname
27fname=trim(seedname)//'.win'
28open(50,file=trim(fname),action='WRITE',form='FORMATTED')
29! write the number of Wannier functions and bands
30write(50,'("length_unit = bohr")')
31write(50,'("num_wann = ",I8)') num_wann
32write(50,'("num_bands = ",I8)') num_bands
33write(50,'("num_iter = ",I8)') num_iter
34write(50,'("dis_num_iter = ",I8)') dis_num_iter
35write(50,*)
36write(50,'("trial_step = ",G18.10)') trial_step
37! write spinors if required
38if (spinpol) then
39 write(50,*)
40 write(50,'("spinors = true")')
41 write(50,'("spn_formatted = true")')
42end if
43! write lattice vectors
44write(50,*)
45write(50,'("begin unit_cell_cart")')
46write(50,'("bohr")')
47write(50,'(3G18.10)') avec(:,1)
48write(50,'(3G18.10)') avec(:,2)
49write(50,'(3G18.10)') avec(:,3)
50write(50,'("end unit_cell_cart")')
51! writes atomic positions
52write(50,*)
53write(50,'("begin atoms_frac")')
54do is=1,nspecies
55 do ia=1,natoms(is)
56 write(50,'(A5,3G18.10)') trim(spsymb(is)),atposl(:,ia,is)
57 end do
58end do
59write(50,'("end atoms_frac")')
60! write the list of k-points
61write(50,*)
62write(50,'("mp_grid = ",3I6)') ngridk
63write(50,*)
64write(50,'("begin kpoints")')
65do ik=1,nkptnr
66 write(50,'(3G18.10)') vkl(:,ik)
67end do
68write(50,'("end kpoints")')
69write(50,*)
70! write the extra lines
71do i=1,nxlwin
72 write(50,'(A)') trim(xlwin(i))
73end do
74close(50)
75write(*,*)
76write(*,'("Info(writew90win): created file ",A)') trim(fname)
77end subroutine
78!EOC
79
integer nkptnr
Definition modmain.f90:463
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
logical spinpol
Definition modmain.f90:228
real(8), dimension(3, 3) avec
Definition modmain.f90:12
integer, dimension(3) ngridk
Definition modmain.f90:448
character(64), dimension(maxspecies) spsymb
Definition modmain.f90:78
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
integer nspecies
Definition modmain.f90:34
real(8), dimension(3, maxatoms, maxspecies) atposl
Definition modmain.f90:51
integer num_wann
Definition modw90.f90:18
integer dis_num_iter
Definition modw90.f90:26
character(256), dimension(:), allocatable xlwin
Definition modw90.f90:16
real(8) trial_step
Definition modw90.f90:28
integer nxlwin
Definition modw90.f90:14
integer num_bands
Definition modw90.f90:20
integer num_iter
Definition modw90.f90:24
character(256) seedname
Definition modw90.f90:12
subroutine writew90win