The Elk Code
geomplot.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2007 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
subroutine
geomplot
7
use
modmain
8
implicit none
9
! local variables
10
integer
is,ia
11
real(8)
v1(3),v2(3),v3(3),v4(3),t1
12
real(8)
dxx,dyx,dyy,dzx,dzy,dzz
13
! initialise universal variables
14
call
init0
15
!------------------------------------------------!
16
! write the XCrysden file to crystal.xsf !
17
!------------------------------------------------!
18
open
(50,file=
'crystal.xsf'
,form=
'FORMATTED'
)
19
write
(50,*)
20
write
(50,
'("CRYSTAL")'
)
21
write
(50,*)
22
write
(50,
'("PRIMVEC")'
)
23
write
(50,
'(3G18.10)'
)
avec
(:,1)*
br_ang
24
write
(50,
'(3G18.10)'
)
avec
(:,2)*
br_ang
25
write
(50,
'(3G18.10)'
)
avec
(:,3)*
br_ang
26
write
(50,*)
27
write
(50,
'("PRIMCOORD")'
)
28
write
(50,
'(2I8)'
)
natmtot
,1
29
do
is=1,
nspecies
30
do
ia=1,
natoms
(is)
31
write
(50,
'(A,3G18.10)'
) trim(
spsymb
(is)),
atposc
(:,ia,is)*
br_ang
32
end do
33
end do
34
close
(50)
35
write
(*,*)
36
write
(*,
'("Info(geomplot):")'
)
37
write
(*,
'(" XCrysDen file written to crystal.xsf")'
)
38
!-----------------------------------------------!
39
! write the V_Sim file to crystal.ascii !
40
!-----------------------------------------------!
41
! determine coordinate system vectors
42
t1=sqrt(
avec
(1,1)**2+
avec
(2,1)**2+
avec
(3,1)**2)
43
v1(:)=
avec
(:,1)/t1
44
t1=sqrt(
avec
(1,2)**2+
avec
(2,2)**2+
avec
(3,2)**2)
45
v2(:)=
avec
(:,2)/t1
46
call
r3cross
(v1,v2,v3)
47
t1=sqrt(v3(1)**2+v3(2)**2+v3(3)**2)
48
v3(:)=v3(:)/t1
49
call
r3cross
(v3,v1,v2)
50
t1=sqrt(v2(1)**2+v2(2)**2+v2(3)**2)
51
v2(:)=v2(:)/t1
52
dxx=dot_product(
avec
(:,1),v1(:))
53
dyx=dot_product(
avec
(:,2),v1(:))
54
dyy=dot_product(
avec
(:,2),v2(:))
55
dzx=dot_product(
avec
(:,3),v1(:))
56
dzy=dot_product(
avec
(:,3),v2(:))
57
dzz=dot_product(
avec
(:,3),v3(:))
58
open
(50,file=
'crystal.ascii'
,form=
'FORMATTED'
)
59
write
(50,*)
60
write
(50,
'(3G18.10)'
) dxx,dyx,dyy
61
write
(50,
'(3G18.10)'
) dzx,dzy,dzz
62
write
(50,*)
63
do
is=1,
nspecies
64
do
ia=1,
natoms
(is)
65
v4(1)=dot_product(
atposc
(:,ia,is),v1(:))
66
v4(2)=dot_product(
atposc
(:,ia,is),v2(:))
67
v4(3)=dot_product(
atposc
(:,ia,is),v3(:))
68
write
(50,
'(3G18.10," ",A)'
) v4,trim(
spsymb
(is))
69
end do
70
end do
71
close
(50)
72
write
(*,*)
73
write
(*,
'("Info(geomplot):")'
)
74
write
(*,
'(" V_Sim file written to crystal.ascii")'
)
75
end subroutine
76
modmain
Definition:
modmain.f90:6
modmain::avec
real(8), dimension(3, 3) avec
Definition:
modmain.f90:12
modmain::natoms
integer, dimension(maxspecies) natoms
Definition:
modmain.f90:36
modmain::br_ang
real(8), parameter br_ang
Definition:
modmain.f90:1269
modmain::nspecies
integer nspecies
Definition:
modmain.f90:34
geomplot
subroutine geomplot
Definition:
geomplot.f90:7
init0
subroutine init0
Definition:
init0.f90:10
r3cross
pure subroutine r3cross(x, y, z)
Definition:
r3cross.f90:10
modmain::natmtot
integer natmtot
Definition:
modmain.f90:40
modmain::atposc
real(8), dimension(3, maxatoms, maxspecies) atposc
Definition:
modmain.f90:54
modmain::spsymb
character(64), dimension(maxspecies) spsymb
Definition:
modmain.f90:78
geomplot.f90
Generated by
1.8.14