The Elk Code
 
Loading...
Searching...
No Matches
polar.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
6!BOP
7! !ROUTINE: polar
8! !INTERFACE:
9subroutine polar(pvl)
10! !USES:
11use modmain
12use modmpi
13use modomp
14! !INPUT/OUTPUT PARAMETERS:
15! pvl : polarisation vector modulo $2\pi$ (out,real(8))
16! !DESCRIPTION:
17! Calculates the polarisation vector modulo $2\pi$ in lattice coordinates
18! using the formula of R. D. King-Smith and David Vanderbilt [Phys. Rev. B
19! {\bf 47}, 1651(R) (1993)], namely
20! $$ P_l=\sum_{\bf k}{\rm Im}\ln\det\left(\langle
21! u_{i{\bf k}+\Delta{\bf k}_l}|u_{j{\bf k}}\rangle\right),$$
22! where $\Delta{\bf k}_l=(1/n_l){\bf B}_l$ and ${\bf B}_l$ is a reciprocal
23! lattice vector. The number of points $n_l$ is equal to that of the original
24! $k$-point grid in direction of ${\bf B}_l$, multiplied by {\tt nskpolar}.
25! See also the routines {\tt polark} and {\tt bornechg}.
26!
27! !REVISION HISTORY:
28! Created May 2020 (JKD)
29!EOP
30!BOC
31implicit none
32! arguments
33real(8), intent(out) :: pvl(3)
34! local variables
35integer ik,l,nthd
36real(8) vc(3),vgqc(3),gqc,pl
37! allocatable arrays
38real(8), allocatable :: jlgqr(:,:)
39complex(8), allocatable :: ylmgq(:),sfacgq(:),expqmt(:,:)
40! allocate local arrays
41allocate(jlgqr(njcmax,nspecies))
42allocate(ylmgq(lmmaxo),sfacgq(natmtot))
43allocate(expqmt(npcmtmax,natmtot))
44maxscl=1
45! loop over reciprocal lattice vectors
46do l=1,3
47! create fine k-point grid in direction l
48 ngridk(:)=ngridk0(:)
50! run one loop of the ground-state calculation
51 call gndstate
52! difference between adjacent k-vectors in this reciprocal lattice direction
53 vc(:)=bvec(:,l)/dble(ngridk(l))
54! calculate the phase factor function exp(iq.r)
55 call gengqf(1,vc,vgqc,gqc,jlgqr,ylmgq,sfacgq)
56 call genexpmt(1,jlgqr,ylmgq,1,sfacgq,expqmt)
57 pl=0.d0
58! parallel loop over non-reduced k-points
59 call holdthd(nkptnr/np_mpi,nthd)
60!$OMP PARALLEL DO DEFAULT(SHARED) &
61!$OMP REDUCTION(+:pl) &
62!$OMP SCHEDULE(DYNAMIC) &
63!$OMP NUM_THREADS(nthd)
64 do ik=1,nkptnr
65! distribute among MPI processes
66 if (mod(ik-1,np_mpi) /= lp_mpi) cycle
67 call polark(ik,l,expqmt,pl)
68 end do
69!$OMP END PARALLEL DO
70 call freethd(nthd)
71! add polarisation from each process and redistribute
72 if (np_mpi > 1) then
73 call mpi_allreduce(mpi_in_place,pl,1,mpi_double_precision,mpi_sum,mpicom, &
74 ierror)
75 end if
76 pvl(l)=pl
77end do
78! restore original input parameters
79ngridk(:)=ngridk0(:)
81deallocate(jlgqr,ylmgq,sfacgq,expqmt)
82end subroutine
83!EOC
84
subroutine genexpmt(ngp, jlgpr, ylmgp, ld, sfacgp, expmt)
Definition genexpmt.f90:7
subroutine gengqf(ng, vqpc, vgqc, gqc, jlgqr, ylmgq, sfacgq)
Definition gengqf.f90:7
subroutine gndstate
Definition gndstate.f90:10
integer njcmax
Definition modmain.f90:1170
integer nkptnr
Definition modmain.f90:463
integer, dimension(3) ngridk0
Definition modmain.f90:448
real(8), dimension(3, 3) bvec
Definition modmain.f90:16
integer nkspolar
Definition modmain.f90:485
integer natmtot
Definition modmain.f90:40
integer npcmtmax
Definition modmain.f90:216
integer, dimension(3) ngridk
Definition modmain.f90:448
integer maxscl0
Definition modmain.f90:1046
integer maxscl
Definition modmain.f90:1046
integer lmmaxo
Definition modmain.f90:203
integer nspecies
Definition modmain.f90:34
integer lp_mpi
Definition modmpi.f90:15
integer ierror
Definition modmpi.f90:19
integer mpicom
Definition modmpi.f90:11
integer np_mpi
Definition modmpi.f90:13
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine polar(pvl)
Definition polar.f90:10
subroutine polark(ik, l, expqmt, pl)
Definition polark.f90:7