The Elk Code
dynsym.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2006-2008 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 dynsym(vpl,dq)
7 use modmain
8 use modphonon
9 implicit none
10 ! arguments
11 real(8), intent(in) :: vpl(3)
12 complex(8), intent(inout) :: dq(nbph,nbph)
13 ! local variables
14 integer isym,lspl,i,j,n
15 real(8) v1(3),v2(3),s(3,3),t1
16 ! automatic arrays
17 complex(8) dqs(nbph,nbph)
18 ! map input vector to first Brillouin zone
19 v1(:)=vpl(:)
20 call vecfbz(epslat,bvec,v1)
21 n=0
22 dqs(:,:)=0.d0
23 ! use the symmetries which leave vpl invariant
24 do isym=1,nsymcrys
25  if (.not.tv0symc(isym)) cycle
26  lspl=lsplsymc(isym)
27  s(:,:)=dble(symlat(:,:,lspl))
28  call r3mtv(s,v1,v2)
29  t1=abs(v1(1)-v2(1))+abs(v1(2)-v2(2))+abs(v1(3)-v2(3))
30  if (t1 < epslat) then
31  call dynsymapp(isym,v1,dq,dqs)
32  n=n+1
33  end if
34 end do
35 if (n == 0) then
36  write(*,*)
37  write(*,'("Error(dynsym): no symmetries leave vpl invariant")')
38  write(*,*)
39  stop
40 end if
41 t1=1.d0/dble(n)
42 dq(:,:)=t1*dqs(:,:)
43 ! make the dynamical matrix Hermitian
44 do j=1,nbph
45  do i=1,j-1
46  dq(i,j)=0.5d0*(dq(i,j)+conjg(dq(j,i)))
47  dq(j,i)=conjg(dq(i,j))
48  end do
49  dq(j,j)=dble(dq(j,j))
50 end do
51 end subroutine
52 
subroutine dynsymapp(isym, vpl, dq, dqs)
Definition: dynsymapp.f90:7
pure subroutine r3mtv(a, x, y)
Definition: r3mtv.f90:10
subroutine vecfbz(eps, bvec, vpl)
Definition: vecfbz.f90:10
integer nsymcrys
Definition: modmain.f90:358
integer, dimension(3, 3, 48) symlat
Definition: modmain.f90:344
logical, dimension(maxsymcrys) tv0symc
Definition: modmain.f90:362
integer, dimension(maxsymcrys) lsplsymc
Definition: modmain.f90:364
real(8), dimension(3, 3) bvec
Definition: modmain.f90:16
real(8) epslat
Definition: modmain.f90:24
subroutine dynsym(vpl, dq)
Definition: dynsym.f90:7