The Elk Code
 
Loading...
Searching...
No Matches
gengclq.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2006 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: gengclq
8! !INTERFACE:
9subroutine gengclq
10! !USES:
11use modmain
12use modtest
13! !DESCRIPTION:
14! The Fock matrix elements
15! $$ V_{ij{\bf k}}\equiv\sum_{l{\bf k'}}\int
16! \frac{\Psi^{\dag}_{i{\bf k}}({\bf r})\cdot\Psi_{l{\bf k}'}({\bf r})
17! \Psi^{\dag}_{l{\bf k}'}({\bf r}')\cdot\Psi_{j{\bf k}}({\bf r}')}
18! {|{\bf r}-{\bf r'}|}\,d^3r\,d^3r' $$
19! contain a divergent term in the sum over ${\bf k}'$ which behaves as
20! $1/q^2$, where ${\bf q}\equiv{\bf k}-{\bf k}'$ is in the first Brillouin
21! zone. The resulting convergence with respect to the number of discrete
22! $q$-points, $N_q$, is very slow. This routine computes the regularised
23! Coulomb Green's function
24! \begin{align}
25! g({\bf q}_i)=\frac{4\pi}{V}\int_{V_i}\frac{1}{q^2}\,d^3q,
26! \end{align}
27! where the integral is over the small parallelepiped with volume
28! $V=\Omega_{\rm BZ}/N_q$ and centered on the discrete point ${\bf q}_i$.
29! This dramatically increases the rate of convergence of methods which involve
30! a summation over the $1/q^2$ part of the Coulomb interaction. The above
31! integral is evaluated numerically on increasingly finer grids and then
32! extrapolated to the continuum.
33!
34! !REVISION HISTORY:
35! Created August 2004 (JKD,SS)
36! Changed from genwiq2, July 2017 (JKD)
37!EOP
38!BOC
39implicit none
40! local variables
41integer, parameter :: np=5,ns0=10,nss=20
42integer ns,iq,i1,i2,i3,ip
43real(8) d(3),sm,t1,t2
44real(8) v1(3),v2(3),v3(3)
45real(8) xa(np),ya(np)
46! external functions
47real(8), external :: polynm
48! allocate global gclq array
49if (allocated(gclq)) deallocate(gclq)
50allocate(gclq(nqpt))
51! begin loop over q-points, note that the vectors vqc are assumed to be in the
52! first Brillouin zone
53do iq=1,nqpt
54! loop over different subdivisions
55 ns=ns0
56 do ip=1,np
57! subdivision vectors in lattice coordinates
58 d(1:3)=1.d0/dble(ngridq(1:3)*2*ns)
59! compute the integral of 1/q²
60 sm=0.d0
61 do i1=-ns,ns-1
62 t1=dble(i1)*d(1)
63 v1(1:3)=vqc(1:3,iq)+t1*bvec(1:3,1)
64 do i2=-ns,ns-1
65 t1=dble(i2)*d(2)
66 v2(1:3)=v1(1:3)+t1*bvec(1:3,2)
67 do i3=-ns,ns-1
68 t1=dble(i3)*d(3)
69 v3(1:3)=v2(1:3)+t1*bvec(1:3,3)
70 t2=v3(1)**2+v3(2)**2+v3(3)**2
71 if (t2 > 1.d-14) sm=sm+1.d0/t2
72 end do
73 end do
74 end do
75 t1=1.d0/dble(2*ns)
76 xa(ip)=t1
77 ya(ip)=fourpi*sm*t1**3
78! increment number of subdivisions
79 ns=ns+nss
80 end do
81! extrapolate the volume element to zero with a polynomial
82 gclq(iq)=polynm(0,np,xa,ya,0.d0)
83end do
84! zero the Green's function at q = 0 if required
85if (t0gclq0) gclq(1)=0.d0
86! write gclq to test file
87call writetest(700,"regularised Coulomb Green''s function (gclq)",nv=nqpt, &
88 tol=1.d-8,rva=gclq)
89end subroutine
90!EOC
91
subroutine gengclq
Definition gengclq.f90:10
real(8), dimension(:), allocatable gclq
Definition modmain.f90:553
real(8), dimension(3, 3) bvec
Definition modmain.f90:16
real(8), parameter fourpi
Definition modmain.f90:1231
real(8), dimension(:,:), allocatable vqc
Definition modmain.f90:547
integer nqpt
Definition modmain.f90:525
integer, dimension(3) ngridq
Definition modmain.f90:515
logical t0gclq0
Definition modmain.f90:555
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition modtest.f90:16
pure real(8) function polynm(m, np, xa, ya, x)
Definition polynm.f90:10