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