The Elk Code
dynev.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 
6 subroutine dynev(dq,wq,ev)
7 use modmain
8 use modphonon
9 implicit none
10 ! arguments
11 complex(8), intent(in) :: dq(nbph,nbph)
12 real(8), intent(out) :: wq(nbph)
13 complex(8), intent(out) :: ev(nbph,nbph)
14 ! local variables
15 integer is,ia,js,ja
16 integer ip,jp,i,j
17 real(8) t1
18 ev(:,:)=0.d0
19 i=0
20 do is=1,nspecies
21  do ia=1,natoms(is)
22  do ip=1,3
23  i=i+1
24  j=0
25  do js=1,nspecies
26 ! mass factor
27  if ((spmass(is) <= 0.d0).or.(spmass(js) <= 0.d0)) then
28 ! infinite mass
29  t1=0.d0
30  else
31  t1=1.d0/sqrt(spmass(is)*spmass(js))
32  end if
33  do ja=1,natoms(js)
34  do jp=1,3
35  j=j+1
36  if (i <= j) then
37 ! use Hermitian average of dynamical matrix
38  ev(i,j)=0.5d0*t1*(dq(i,j)+conjg(dq(j,i)))
39  end if
40  end do
41  end do
42  end do
43  end do
44  end do
45 end do
46 ! find the eigenvalues and eigenvectors of the dynamical matrix
47 call eveqnzh(nbph,nbph,ev,wq)
48 do i=1,nbph
49  if (wq(i) >= 0.d0) then
50  wq(i)=sqrt(wq(i))
51  else
52  wq(i)=-sqrt(abs(wq(i)))
53  end if
54 end do
55 end subroutine
56 
subroutine eveqnzh(n, ld, a, w)
Definition: eveqnzh.f90:7
subroutine dynev(dq, wq, ev)
Definition: dynev.f90:7
real(8), dimension(maxspecies) spmass
Definition: modmain.f90:101
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
integer nspecies
Definition: modmain.f90:34