The Elk Code
 
Loading...
Searching...
No Matches
polark.f90
Go to the documentation of this file.
1
2! Copyright (C) 2020 J. K. Dewhurst and S. Sharma.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine polark(ik,l,expqmt,pl)
7use modmain
8implicit none
9! arguments
10integer, intent(in) :: ik,l
11complex(8), intent(in) :: expqmt(npcmtmax,natmtot)
12real(8), intent(inout) :: pl
13! local variables
14integer jk,nst,ist
15real(8) vkql(3)
16complex(8) z1
17! automatic arrays
18integer idx(nstsv),ngp(nspnfv),igpig(ngkmax,nspnfv)
19! allocatable arrays
20complex(8), allocatable :: wfmt(:,:,:,:),wfir(:,:,:)
21complex(8), allocatable :: wfmtq(:,:,:,:),wfgkq(:,:,:)
22complex(8), allocatable :: oq(:,:)
23! external functions
24complex(8), external :: zmdet
25! equivalent reduced k-point
26jk=ivkik(ivk(1,ik),ivk(2,ik),ivk(3,ik))
27! find the adjacent k-point in lattice coordinates
28vkql(1:3)=vkl(1:3,ik)
29vkql(l)=vkql(l)+1.d0/dble(ngridk(l))
30! count and index the occupied states
31nst=0
32do ist=1,nstsv
33 if (evalsv(ist,jk) > efermi) cycle
34 nst=nst+1
35 idx(nst)=ist
36end do
37! generate the wavefunctions for occupied states at k
38allocate(wfmt(npcmtmax,natmtot,nspinor,nst),wfir(ngtc,nspinor,nst))
39call genwfsvp(.false.,.false.,nst,idx,ngdgc,igfc,vkl(:,ik),ngp,igpig,wfmt,ngtc,&
40 wfir)
41! generate the wavefunctions for occupied states at k+q
42allocate(wfmtq(npcmtmax,natmtot,nspinor,nst),wfgkq(ngkmax,nspinor,nst))
43call genwfsvp(.false.,.true.,nst,idx,ngdgc,igfc,vkql,ngp,igpig,wfmtq,ngkmax, &
44 wfgkq)
45! determine the overlap matrix for all occupied states
46allocate(oq(nst,nst))
47call genolpq(nst,expqmt,ngp,igpig,wfmt,wfir,wfmtq,wfgkq,oq)
48! compute the determinant of the matrix
49z1=zmdet(nst,oq)
50! determine the phase of the determinant and add to total polarisation
51pl=pl+atan2(z1%im,z1%re)
52deallocate(wfmt,wfir,wfmtq,wfgkq,oq)
53end subroutine
54
subroutine genolpq(nst, expqmt, ngpq, igpqig, wfmt, wfir, wfmtq, wfgpq, oq)
Definition genolpq.f90:7
subroutine genwfsvp(tsh, tgp, nst, idx, ngridg_, igfft_, vpl, ngp, igpig, wfmt, ld, wfir)
Definition genwfsvp.f90:7
real(8) efermi
Definition modmain.f90:904
integer nspinor
Definition modmain.f90:267
integer, dimension(3) ngdgc
Definition modmain.f90:388
integer, dimension(:,:), allocatable ivk
Definition modmain.f90:465
integer, dimension(:), allocatable igfc
Definition modmain.f90:410
integer ngtc
Definition modmain.f90:392
integer, dimension(3) ngridk
Definition modmain.f90:448
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
integer, dimension(:,:,:), allocatable ivkik
Definition modmain.f90:467
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
subroutine polark(ik, l, expqmt, pl)
Definition polark.f90:7