The Elk Code
 
Loading...
Searching...
No Matches
findsymlat.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2007 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3! This file is distributed under the terms of the GNU Lesser General Public
4! License. See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: findsymlat
8! !INTERFACE:
9subroutine findsymlat
10! !USES:
11use modmain
12use modtddft
13! !DESCRIPTION:
14! Finds the point group symmetries which leave the Bravais lattice invariant.
15! Let $A$ be the matrix consisting of the lattice vectors in columns, then
16! $$ g=A^{\rm T}A $$
17! is the metric tensor. Any $3\times 3$ matrix $S$ with elements $-1$, 0 or 1
18! is a point group symmetry of the lattice if $\det(S)$ is $-1$ or 1, and
19! $$ S^{\rm T}gS=g. $$
20! The first matrix in the set returned is the identity.
21!
22! !REVISION HISTORY:
23! Created January 2003 (JKD)
24! Removed arguments and simplified, April 2007 (JKD)
25!EOP
26!BOC
27implicit none
28! local variables
29integer md,sym(3,3),its,i,j
30integer i11,i12,i13,i21,i22,i23,i31,i32,i33
31real(8) s(3,3),g(3,3),sgs(3,3),sc(3,3)
32real(8) c(3,3),v(3),t1
33! determine metric tensor
34call r3mtm(avec,avec,g)
35! loop over all possible symmetry matrices
36nsymlat=0
37do i11=-1,1; do i12=-1,1; do i13=-1,1
38do i21=-1,1; do i22=-1,1; do i23=-1,1
39do i31=-1,1; do i32=-1,1; do i33=-1,1
40 sym(1,1)=i11; sym(1,2)=i12; sym(1,3)=i13
41 sym(2,1)=i21; sym(2,2)=i22; sym(2,3)=i23
42 sym(3,1)=i31; sym(3,2)=i32; sym(3,3)=i33
43! determinant of matrix
44 md=i3mdet(sym)
45! matrix should be orthogonal
46 if (abs(md) /= 1) goto 10
47! check invariance of metric tensor
48 s(:,:)=dble(sym(:,:))
49 call r3mtm(s,g,c)
50 call r3mm(c,s,sgs)
51 do j=1,3
52 do i=1,3
53 if (abs(sgs(i,j)-g(i,j)) > epslat) goto 10
54 end do
55 end do
56! check invariance of spin-spiral q-vector if required
57 if (spinsprl) then
58 call r3mtv(s,vqlss,v)
59 t1=abs(vqlss(1)-v(1))+abs(vqlss(2)-v(2))+abs(vqlss(3)-v(3))
60 if (t1 > epslat) goto 10
61 end if
62! check invariance of electric field if required
63 if (tefield) then
64 call r3mv(s,efieldl,v)
65 t1=abs(efieldl(1)-v(1))+abs(efieldl(2)-v(2))+abs(efieldl(3)-v(3))
66 if (t1 > epslat) goto 10
67 end if
68! check invariance of static A-field if required
69 if (tafield) then
70 call r3mv(s,afieldl,v)
71 t1=abs(afieldl(1)-v(1))+abs(afieldl(2)-v(2))+abs(afieldl(3)-v(3))
72 if (t1 > epslat) goto 10
73 end if
74! check invariance of time-dependent A-field if required
75 if (tafieldt) then
76! symmetry matrix in Cartesian coordinates
77 call r3mm(s,ainv,c)
78 call r3mm(avec,c,sc)
79 do its=1,ntimes
80 call r3mv(sc,afieldt(:,its),v)
81 t1=abs(afieldt(1,its)-v(1)) &
82 +abs(afieldt(2,its)-v(2)) &
83 +abs(afieldt(3,its)-v(3))
84 if (t1 > epslat) goto 10
85 end do
86 end if
88 if (nsymlat > 48) then
89 write(*,*)
90 write(*,'("Error(findsymlat): more than 48 symmetries found")')
91 write(*,'(" (lattice vectors may be linearly dependent)")')
92 write(*,*)
93 stop
94 end if
95! store the symmetry and determinant in global arrays
96 symlat(:,:,nsymlat)=sym(:,:)
9810 continue
99end do; end do; end do
100end do; end do; end do
101end do; end do; end do
102if (nsymlat == 0) then
103 write(*,*)
104 write(*,'("Error(findsymlat): no symmetries found")')
105 write(*,*)
106 stop
107end if
108! make the first symmetry the identity
109do i=1,nsymlat
110 if ((symlat(1,1,i) == 1).and.(symlat(1,2,i) == 0).and.(symlat(1,3,i) == 0) &
111 .and.(symlat(2,1,i) == 0).and.(symlat(2,2,i) == 1).and.(symlat(2,3,i) == 0) &
112 .and.(symlat(3,1,i) == 0).and.(symlat(3,2,i) == 0).and.(symlat(3,3,i) == 1)) &
113 then
114 sym(:,:)=symlat(:,:,1)
115 symlat(:,:,1)=symlat(:,:,i)
116 symlat(:,:,i)=sym(:,:)
117 md=symlatd(1)
118 symlatd(1)=symlatd(i)
119 symlatd(i)=md
120 exit
121 end if
122end do
123! index to the inverse of each operation
124do i=1,nsymlat
125 call i3minv(symlat(:,:,i),sym)
126 do j=1,nsymlat
127 if ((symlat(1,1,j) == sym(1,1)).and.(symlat(1,2,j) == sym(1,2)).and. &
128 (symlat(1,3,j) == sym(1,3)).and.(symlat(2,1,j) == sym(2,1)).and. &
129 (symlat(2,2,j) == sym(2,2)).and.(symlat(2,3,j) == sym(2,3)).and. &
130 (symlat(3,1,j) == sym(3,1)).and.(symlat(3,2,j) == sym(3,2)).and. &
131 (symlat(3,3,j) == sym(3,3))) then
132 isymlat(i)=j
133 goto 30
134 end if
135 end do
136 write(*,*)
137 write(*,'("Error(findsymlat): inverse operation not found")')
138 write(*,'(" for lattice symmetry ",I2)') i
139 write(*,*)
140 stop
14130 continue
142end do
143! determine the lattice symmetries in Cartesian coordinates
144do i=1,nsymlat
145 s(:,:)=dble(symlat(:,:,i))
146 call r3mm(s,ainv,c)
147 call r3mm(avec,c,symlatc(:,:,i))
148end do
149return
150
151contains
152
153pure integer function i3mdet(a)
154implicit none
155! arguments
156integer, intent(in) :: a(3,3)
157! determinant of integer matrix
158i3mdet=a(1,1)*(a(2,2)*a(3,3)-a(3,2)*a(2,3)) &
159 +a(2,1)*(a(3,2)*a(1,3)-a(1,2)*a(3,3)) &
160 +a(3,1)*(a(1,2)*a(2,3)-a(2,2)*a(1,3))
161end function
162
163end subroutine
164!EOC
165
pure integer function i3mdet(a)
subroutine findsymlat
subroutine i3minv(a, b)
Definition i3minv.f90:10
integer nsymlat
Definition modmain.f90:342
logical tafield
Definition modmain.f90:322
real(8), dimension(3, 3) avec
Definition modmain.f90:12
real(8), dimension(3) efieldl
Definition modmain.f90:314
logical tefield
Definition modmain.f90:310
real(8) epslat
Definition modmain.f90:24
integer, dimension(48) symlatd
Definition modmain.f90:346
logical spinsprl
Definition modmain.f90:283
real(8), dimension(3, 3) ainv
Definition modmain.f90:14
real(8), dimension(3) afieldl
Definition modmain.f90:327
real(8), dimension(3) vqlss
Definition modmain.f90:293
integer, dimension(3, 3, 48) symlat
Definition modmain.f90:344
logical tafieldt
Definition modtddft.f90:56
real(8), dimension(:,:), allocatable afieldt
Definition modtddft.f90:58
integer ntimes
Definition modtddft.f90:42
pure subroutine r3mm(a, b, c)
Definition r3mm.f90:10
pure subroutine r3mtm(a, b, c)
Definition r3mtm.f90:10
pure subroutine r3mtv(a, x, y)
Definition r3mtv.f90:10
pure subroutine r3mv(a, x, y)
Definition r3mv.f90:10