The Elk Code
 
Loading...
Searching...
No Matches
dynqtor.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 dynqtor(dq,dr)
7use modmain
8use modphonon
9implicit none
10! arguments
11complex(8), intent(in) :: dq(nbph,nbph,nqpt)
12real(8), intent(out) :: dr(nbph,nbph,nqptnr)
13! local variables
14integer intr(2,3),iq,n
15integer isym,lspl,ir
16integer i1,i2,i3,j1,j2,j3
17real(8) v1(3),v2(3),s(3,3),t1
18complex(8) z1
19! automatic arrays
20complex(8) dqs(nbph,nbph)
21intr(2,:)=ngridq(:)/2
22intr(1,:)=intr(2,:)-ngridq(:)+1
23dr(:,:,:)=0.d0
24! loop over q-vectors
25do j1=0,ngridq(1)-1
26 v1(1)=dble(j1)/dble(ngridq(1))
27 do j2=0,ngridq(2)-1
28 v1(2)=dble(j2)/dble(ngridq(2))
29 do j3=0,ngridq(3)-1
30 v1(3)=dble(j3)/dble(ngridq(3))
31 iq=ivqiq(j1,j2,j3)
32! map v1 to the first Brillouin zone
33 call vecfbz(epslat,bvec,v1)
34! rotate and add the dynamical matrix of the reduced q-point with all symmetries
35 n=0
36 dqs(:,:)=0.d0
37 do isym=1,nsymcrys
38 lspl=lsplsymc(isym)
39 s(:,:)=dble(symlat(:,:,lspl))
40 call r3mtv(s,vql(:,iq),v2)
41 call vecfbz(epslat,bvec,v2)
42 t1=abs(v1(1)-v2(1))+abs(v1(2)-v2(2))+abs(v1(3)-v2(3))
43 if (t1 < epslat) then
44 call dynsymapp(isym,vql(:,iq),dq(:,:,iq),dqs)
45 n=n+1
46 end if
47 end do
48 if (n == 0) then
49 write(*,*)
50 write(*,'("Error(dynqtor): vector ",3G18.10)') v1
51 write(*,'(" cannot be mapped to reduced q-point set")')
52 write(*,*)
53 stop
54 end if
55 t1=1.d0/dble(n)
56 dqs(:,:)=t1*dqs(:,:)
57! subtract the non-analytic term if required
58 if (tphnat) call dynqnat(-1,v1,dqs)
59! loop over R-vectors
60 ir=0
61 do i3=intr(1,3),intr(2,3)
62 do i2=intr(1,2),intr(2,2)
63 do i1=intr(1,1),intr(2,1)
64 ir=ir+1
65 t1=twopi*(v1(1)*dble(i1)+v1(2)*dble(i2)+v1(3)*dble(i3))
66 z1=cmplx(cos(t1),sin(t1),8)
67 dr(:,:,ir)=dr(:,:,ir)+dble(z1*dqs(:,:))
68 end do
69 end do
70 end do
71 end do
72 end do
73end do
74t1=1.d0/dble(nqptnr)
75dr(:,:,:)=t1*dr(:,:,:)
76end subroutine
77
subroutine dynqnat(sgn, vpl, dq)
Definition dynqnat.f90:7
subroutine dynqtor(dq, dr)
Definition dynqtor.f90:7
subroutine dynsymapp(isym, vpl, dq, dqs)
Definition dynsymapp.f90:7
real(8), dimension(3, 3) bvec
Definition modmain.f90:16
real(8), parameter twopi
Definition modmain.f90:1230
integer, dimension(:,:,:), allocatable ivqiq
Definition modmain.f90:531
real(8) epslat
Definition modmain.f90:24
real(8), dimension(:,:), allocatable vql
Definition modmain.f90:545
integer nsymcrys
Definition modmain.f90:358
integer, dimension(3, 3, 48) symlat
Definition modmain.f90:344
integer, dimension(3) ngridq
Definition modmain.f90:515
integer, dimension(maxsymcrys) lsplsymc
Definition modmain.f90:364
logical tphnat
Definition modphonon.f90:34
pure subroutine r3mtv(a, x, y)
Definition r3mtv.f90:10
subroutine vecfbz(eps, bvec, vpl)
Definition vecfbz.f90:10