The Elk Code
 
Loading...
Searching...
No Matches
findscq.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
6subroutine findscq(iq,nsc,vsc)
7use modmain
8use modphonon
9implicit none
10! arguments
11integer, intent(in) :: iq
12integer, intent(out) :: nsc
13real(8), intent(out) :: vsc(3,nqptnr)
14! local variables
15integer i1,i2,i3
16integer scl(3,3),i,n
17real(8) dmin,t1
18real(8) v1(3),v2(3)
19! check for Gamma-point phonon
20if (tphq0) then
21 scl(:,:)=0
22 scl(1,1)=1
23 scl(2,2)=1
24 scl(3,3)=1
25 nsc=1
26 goto 10
27end if
28! find the first lattice vector
29dmin=1.d8
30do i1=-ngridq(1),ngridq(1)
31 do i2=-ngridq(2),ngridq(2)
32 do i3=-ngridq(3),ngridq(3)
33 t1=dble(i1)*vql(1,iq)+dble(i2)*vql(2,iq)+dble(i3)*vql(3,iq)
34 if (abs(t1-nint(t1)) < epslat) then
35 v1(:)=dble(i1)*avec0(:,1)+dble(i2)*avec0(:,2)+dble(i3)*avec0(:,3)
36 t1=sqrt(v1(1)**2+v1(2)**2+v1(3)**2)
37 if ((t1 < dmin).and.(t1 > epslat)) then
38 scl(1,1)=i1
39 scl(2,1)=i2
40 scl(3,1)=i3
41 dmin=t1
42 end if
43 end if
44 end do
45 end do
46end do
47! find the second lattice vector
48dmin=1.d8
49do i1=-ngridq(1),ngridq(1)
50 do i2=-ngridq(2),ngridq(2)
51 do i3=-ngridq(3),ngridq(3)
52 t1=dble(i1)*vql(1,iq)+dble(i2)*vql(2,iq)+dble(i3)*vql(3,iq)
53 if (abs(t1-nint(t1)) < epslat) then
54! area defined by first two lattice vectors
55 n=(i2*scl(3,1)-i3*scl(2,1))**2 &
56 +(i3*scl(1,1)-i1*scl(3,1))**2 &
57 +(i1*scl(2,1)-i2*scl(1,1))**2
58 if (n /= 0) then
59 v1(:)=dble(i1)*avec0(:,1)+dble(i2)*avec0(:,2)+dble(i3)*avec0(:,3)
60 t1=v1(1)**2+v1(2)**2+v1(3)**2
61 if (t1 < dmin) then
62 scl(1,2)=i1
63 scl(2,2)=i2
64 scl(3,2)=i3
65 dmin=t1
66 end if
67 end if
68 end if
69 end do
70 end do
71end do
72! find the third lattice vector
73nsc=0
74dmin=1.d8
75do i1=-ngridq(1),ngridq(1)
76 do i2=-ngridq(2),ngridq(2)
77 do i3=-ngridq(3),ngridq(3)
78 t1=dble(i1)*vql(1,iq)+dble(i2)*vql(2,iq)+dble(i3)*vql(3,iq)
79 if (abs(t1-nint(t1)) < epslat) then
80! number of primitive unit cells in supercell
81 n=scl(1,2)*(i2*scl(3,1)-i3*scl(2,1)) &
82 +scl(2,2)*(i3*scl(1,1)-i1*scl(3,1)) &
83 +scl(3,2)*(i1*scl(2,1)-i2*scl(1,1))
84 if (n /= 0) then
85 v1(:)=dble(i1)*avec0(:,1)+dble(i2)*avec0(:,2)+dble(i3)*avec0(:,3)
86 t1=v1(1)**2+v1(2)**2+v1(3)**2
87 if (t1 < dmin) then
88 nsc=abs(n)
89 scl(1,3)=i1
90 scl(2,3)=i2
91 scl(3,3)=i3
92 dmin=t1
93 end if
94 end if
95 end if
96 end do
97 end do
98end do
99if (nsc == 0) goto 30
10010 continue
101! new lattice vectors
102do i=1,3
103 avec(:,i)=dble(scl(1,i))*avec0(:,1) &
104 +dble(scl(2,i))*avec0(:,2) &
105 +dble(scl(3,i))*avec0(:,3)
106end do
107! inverse of lattice vector matrix
108call r3minv(avec,ainv)
109! generate offset vectors for each primitive cell in the supercell
110n=1
111vsc(:,1)=0.d0
112do i1=-ngridq(1),ngridq(1)
113 do i2=-ngridq(2),ngridq(2)
114 do i3=-ngridq(3),ngridq(3)
115 if (n == nsc) return
116 v1(:)=dble(i1)*avec0(:,1)+dble(i2)*avec0(:,2)+dble(i3)*avec0(:,3)
117 call r3mv(ainv,v1,v2)
118 call r3frac(epslat,v2)
119 call r3mv(avec,v2,v1)
120 do i=1,n
121 t1=abs(v1(1)-vsc(1,i))+abs(v1(2)-vsc(2,i))+abs(v1(3)-vsc(3,i))
122 if (t1 < epslat) goto 20
123 end do
124 n=n+1
125 vsc(:,n)=v1(:)
12620 continue
127 end do
128 end do
129end do
13030 continue
131write(*,*)
132write(*,'("Error(findscq): unable to generate supercell")')
133write(*,*)
134stop
135end subroutine
136
subroutine findscq(iq, nsc, vsc)
Definition findscq.f90:7
real(8), dimension(3, 3) avec0
Definition modmain.f90:12
real(8), dimension(3, 3) avec
Definition modmain.f90:12
real(8) epslat
Definition modmain.f90:24
real(8), dimension(:,:), allocatable vql
Definition modmain.f90:545
real(8), dimension(3, 3) ainv
Definition modmain.f90:14
integer, dimension(3) ngridq
Definition modmain.f90:515
logical tphq0
Definition modphonon.f90:17
pure subroutine r3frac(eps, v)
Definition r3frac.f90:10
subroutine r3minv(a, b)
Definition r3minv.f90:10
pure subroutine r3mv(a, x, y)
Definition r3mv.f90:10