The Elk Code
writegeom.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: writegeom
8 ! !INTERFACE:
9 subroutine writegeom(fnum)
10 ! !USES:
11 use modmain
12 ! !INPUT/OUTPUT PARAMETERS:
13 ! fnum : file number for writing output (in,integer)
14 ! !DESCRIPTION:
15 ! Outputs the lattice vectors and atomic positions to file, in a format
16 ! which may be then used directly in {\tt elk.in}.
17 !
18 ! !REVISION HISTORY:
19 ! Created January 2004 (JKD)
20 !EOP
21 !BOC
22 implicit none
23 ! arguments
24 integer, intent(in) :: fnum
25 ! local variables
26 integer is,ia,i
27 real(8) v1(3),v2(3)
28 write(fnum,*)
29 write(fnum,'("scale")')
30 write(fnum,'(" 1.0")')
31 write(fnum,*)
32 write(fnum,'("scale1")')
33 write(fnum,'(" 1.0")')
34 write(fnum,*)
35 write(fnum,'("scale2")')
36 write(fnum,'(" 1.0")')
37 write(fnum,*)
38 write(fnum,'("scale3")')
39 write(fnum,'(" 1.0")')
40 write(fnum,*)
41 write(fnum,'("avec")')
42 write(fnum,'(3G18.10)') avec(:,1)
43 write(fnum,'(3G18.10)') avec(:,2)
44 write(fnum,'(3G18.10)') avec(:,3)
45 if (molecule) then
46  write(fnum,*)
47  write(fnum,'("molecule")')
48  write(fnum,'(" ",L1)') molecule
49 end if
50 write(fnum,*)
51 write(fnum,'("atoms")')
52 write(fnum,'(I4,T40," : nspecies")') nspecies
53 do is=1,nspecies
54  write(fnum,'("''",A,"''",T40," : spfname")') trim(spfname(is))
55  write(fnum,'(I4,T40," : natoms; atpos, bfcmt below")') natoms(is)
56  do ia=1,natoms(is)
57  if (molecule) then
58 ! map lattice coordinates to [-0.5,0.5) and write as Cartesian coordinates
59  v1(:)=atposl(:,ia,is)
60  do i=1,3
61  if (v1(i) > 0.5d0) v1(i)=v1(i)-1.d0
62  end do
63  call r3mv(avec,v1,v2)
64  else
65 ! otherwise write lattice coordinates
66  v2(:)=atposl(:,ia,is)
67  end if
68  write(fnum,'(3F14.8," ",3F12.8)') v2(:),bfcmt(:,ia,is)
69  end do
70 end do
71 end subroutine
72 !EOC
73 
character(256), dimension(maxspecies) spfname
Definition: modmain.f90:74
real(8), dimension(3, maxatoms, maxspecies) atposl
Definition: modmain.f90:51
real(8), dimension(3, 3) avec
Definition: modmain.f90:12
subroutine writegeom(fnum)
Definition: writegeom.f90:10
real(8), dimension(3, maxatoms, maxspecies) bfcmt
Definition: modmain.f90:273
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
integer nspecies
Definition: modmain.f90:34
pure subroutine r3mv(a, x, y)
Definition: r3mv.f90:10
logical molecule
Definition: modmain.f90:47