File: writesym.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (111 lines) | stat: -rw-r--r-- 2,958 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111

! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

!BOP
! !ROUTINE: writesym
! !INTERFACE:
subroutine writesym
! !USES:
use modmain
! !DESCRIPTION:
!   Outputs the Bravais, crystal and site symmetry matrices to files
!   {\tt SYMLAT.OUT}, {\tt SYMCRYS.OUT} and {\tt SYMSITE.OUT}, respectively.
!   Also writes out equivalent atoms and related crystal symmetries to
!   {\tt EQATOMS.OUT}.
!
! !REVISION HISTORY:
!   Created October 2002 (JKD)
!EOP
!BOC
implicit none
! local variables
integer is,ia,ja,ias,i
integer isym,lspl,lspn
! output the Bravais lattice symmetries
open(50,file='SYMLAT'//trim(filext),form='FORMATTED')
write(50,'(I4," : nsymlat")') nsymlat
do isym=1,nsymlat
  write(50,*)
  write(50,'(I4)') isym
  do i=1,3
    write(50,'(3I4)') symlat(i,:,isym)
  end do
end do
close(50)
! output the crystal symmetries
open(50,file='SYMCRYS'//trim(filext),form='FORMATTED')
write(50,*)
write(50,'("(translation vectors and rotation matrices are in lattice &
 &coordinates)")')
write(50,*)
write(50,'(I4," : nsymcrys")') nsymcrys
do isym=1,nsymcrys
  write(50,*)
  write(50,'("Crystal symmetry : ",I4)') isym
  write(50,'(" spatial translation :")')
  write(50,'(3G18.10)') vtlsymc(:,isym)
  write(50,'(" spatial rotation :")')
  lspl=lsplsymc(isym)
  do i=1,3
    write(50,'(3I4)') symlat(i,:,lspl)
  end do
  write(50,'(" global spin rotation :")')
  lspn=lspnsymc(isym)
  do i=1,3
    write(50,'(3I4)') symlat(i,:,lspn)
  end do
end do
close(50)
! output the site symmetries
open(50,file='SYMSITE'//trim(filext),form='FORMATTED')
write(50,*)
write(50,'("(rotation matrices are in lattice coordinates)")')
do is=1,nspecies
  do ia=1,natoms(is)
    ias=idxas(ia,is)
    write(50,*)
    write(50,*)
    write(50,'("Species : ",I4," (",A,"), atom : ",I4)') is,trim(spsymb(is)),ia
    write(50,'(I4," : nsymsite")') nsymsite(ias)
    do isym=1,nsymsite(ias)
      write(50,*)
      write(50,'(" Site symmetry : ",I4)') isym
      write(50,'("  spatial rotation :")')
      lspl=lsplsyms(isym,ias)
      do i=1,3
        write(50,'(3I4)') symlat(i,:,lspl)
      end do
      write(50,'("  global spin rotation :")')
      lspn=lspnsyms(isym,ias)
      do i=1,3
        write(50,'(3I4)') symlat(i,:,lspn)
      end do
    end do
  end do
end do
close(50)
! output the equivalent atoms and related symmetries
open(50,file='EQATOMS'//trim(filext),form='FORMATTED')
do is=1,nspecies
  write(50,*)
  write(50,'("Species : ",I4," (",A,")")') is,trim(spsymb(is))
  do ia=1,natoms(is)
    write(50,'(" atom ",I4," is equivalent to atom(s)")') ia
    i=0
    do ja=1,natoms(is)
      if (eqatoms(ia,ja,is)) then
        if ((i.gt.0).and.(mod(i,20).eq.0)) write(50,*)
        write(50,'(I4)',advance='NO') ja
        i=i+1
      end if
    end do
    write(50,*)
  end do
end do
close(50)
return
end subroutine
!EOC