The Elk Code
dynsymapp.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2005-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 
6 subroutine dynsymapp(isym,vpl,dq,dqs)
7 use modmain
8 use modphonon
9 implicit none
10 ! arguments
11 integer, intent(in) :: isym
12 real(8), intent(in) :: vpl(3)
13 complex(8), intent(in) :: dq(nbph,nbph)
14 complex(8), intent(inout) :: dqs(nbph,nbph)
15 ! local variables
16 integer is,ia,ja,ias,jas
17 integer lspl,ilspl,i,j,k,l,m,n
18 real(8) sl(3,3),sic(3,3),v(3),t1
19 real(8) a(3,3),b(3,3),c(3,3)
20 complex(8) z1
21 ! automatic arrays
22 integer map(natmtot)
23 complex(8) zph(natmtot)
24 ! index to spatial rotation in lattice point group
25 lspl=lsplsymc(isym)
26 ! the inverse of the spatial symmetry
27 ilspl=isymlat(lspl)
28 ! symmetry matrix in lattice coordinates
29 sl(:,:)=dble(symlat(:,:,lspl))
30 ! inverse symmetry matrix in Cartesian coordinates
31 sic(:,:)=symlatc(:,:,ilspl)
32 ! operate with symmetry matrix on vpl
33 call r3mtv(sl,vpl,v)
34 do is=1,nspecies
35  do ia=1,natoms(is)
36  ias=idxas(ia,is)
37 ! equivalent atom with this symmetry
38  ja=ieqatom(ia,is,isym)
39  jas=idxas(ja,is)
40  map(ias)=jas
41 ! phase factor
42  t1=twopi*dot_product(vpl(:),atposl(:,ia,is))
43  z1=cmplx(cos(t1),sin(t1),8)
44  zph(ias)=z1
45  t1=-twopi*dot_product(v(:),atposl(:,ja,is))
46  zph(ias)=z1*cmplx(cos(t1),sin(t1),8)
47  end do
48 end do
49 ! rotate and phase-shift dynamical matrix with symmetry
50 do ias=1,natmtot
51  i=3*(ias-1)
52  k=3*(map(ias)-1)
53  do jas=1,natmtot
54  j=3*(jas-1)
55  l=3*(map(jas)-1)
56  do m=1,3
57  do n=1,3
58  a(m,n)=dble(dq(i+m,j+n))
59  b(m,n)=aimag(dq(i+m,j+n))
60  end do
61  end do
62  call r3mm(sic,a,c)
63  call r3mmt(c,sic,a)
64  call r3mm(sic,b,c)
65  call r3mmt(c,sic,b)
66  z1=zph(ias)*conjg(zph(jas))
67  do m=1,3
68  do n=1,3
69  dqs(k+m,l+n)=dqs(k+m,l+n)+z1*cmplx(a(m,n),b(m,n),8)
70  end do
71  end do
72  end do
73 end do
74 end subroutine
75 
subroutine dynsymapp(isym, vpl, dq, dqs)
Definition: dynsymapp.f90:7
pure subroutine r3mmt(a, b, c)
Definition: r3mmt.f90:10
pure subroutine r3mtv(a, x, y)
Definition: r3mtv.f90:10
real(8), parameter twopi
Definition: modmain.f90:1233
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
integer, dimension(48) isymlat
Definition: modmain.f90:348
integer, dimension(3, 3, 48) symlat
Definition: modmain.f90:344
integer, dimension(:,:,:), allocatable ieqatom
Definition: modmain.f90:368
real(8), dimension(3, maxatoms, maxspecies) atposl
Definition: modmain.f90:51
integer, dimension(maxsymcrys) lsplsymc
Definition: modmain.f90:364
real(8), dimension(3, 3, 48) symlatc
Definition: modmain.f90:350
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
integer nspecies
Definition: modmain.f90:34
pure subroutine r3mm(a, b, c)
Definition: r3mm.f90:10