The Elk Code
Loading...
Searching...
No Matches
sumrule.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
!BOP
7
! !ROUTINE: sumrule
8
! !INTERFACE:
9
subroutine
sumrule
10
! !DESCRIPTION:
11
! Applies the same correction to all the dynamical matrices such that the
12
! matrix for ${\bf q}=0$ satisfies the acoustic sum rule. In other words, the
13
! matrices are updated with
14
! $$ D_{ij}^{\bf q}\rightarrow D_{ij}^{\bf q}-\sum_{k=1}^3 \omega_k^0
15
! v_{k;i}^0 v_{k;j}^0 $$
16
! for all ${\bf q}$, where $\omega_k^0$ is the $k$th eigenvalue of the
17
! ${\bf q}=0$ dynamical matrix and $v_{k;i}^0$ the $i$th component of its
18
! eigenvector. The eigenvalues are assumed to be arranged in ascending order.
19
! This ensures that the ${\bf q}=0$ dynamical matrix has 3 zero eigenvalues,
20
! which the uncorrected matrix may not have due to the finite
21
! exchange-correlation grid.
22
!
23
! !REVISION HISTORY:
24
! Created May 2005 (JKD)
25
!EOP
26
!BOC
27
use
modmain
28
use
modphonon
29
implicit none
30
! local variables
31
integer
iq,i,j,k
32
! automatic arrays
33
real
(8) wq0(nbph)
34
complex(8)
ev(nbph,nbph)
35
! compute the eigenvalues and vectors of the q = 0 dynamical matrix
36
do
i=1,nbph
37
do
j=i,nbph
38
ev(i,j)=0.5d0*(
dynq
(i,j,1)+conjg(
dynq
(j,i,1)))
39
end do
40
end do
41
call
eveqnzh
(nbph,nbph,ev,wq0)
42
! subtract outer products of 3 lowest eigenvectors for q = 0 from all the
43
! dynamical matrices
44
do
iq=1,
nqpt
45
do
i=1,nbph
46
do
j=1,nbph
47
do
k=1,3
48
dynq
(i,j,iq)=
dynq
(i,j,iq)-wq0(k)*ev(i,k)*conjg(ev(j,k))
49
end do
50
end do
51
end do
52
end do
53
end subroutine
54
!EOC
55
eveqnzh
subroutine eveqnzh(n, ld, a, w)
Definition
eveqnzh.f90:7
modmain
Definition
modmain.f90:6
modmain::nqpt
integer nqpt
Definition
modmain.f90:525
modphonon
Definition
modphonon.f90:6
modphonon::dynq
complex(8), dimension(:,:,:), allocatable dynq
Definition
modphonon.f90:27
sumrule
subroutine sumrule
Definition
sumrule.f90:10
sumrule.f90
Generated by
1.9.8