The Elk Code
gendmat.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2019 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
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 gendmat(tspndg,tlmdg,lmin,lmax,ld,dmat)
7 use modmain
8 use modmpi
9 use modomp
10 implicit none
11 ! arguments
12 logical, intent(in) :: tspndg,tlmdg
13 integer, intent(in) :: lmin,lmax,ld
14 complex(8), intent(out) :: dmat(ld,nspinor,ld,nspinor,natmtot)
15 ! local variables
16 integer ik,ispn
17 integer nst,ist,jst
18 integer ias,n,nthd
19 real(8) wo
20 ! automatic arrays
21 integer idx(nstsv)
22 ! allocatable arrays
23 complex(8), allocatable :: apwalm(:,:,:,:,:),evecfv(:,:,:),evecsv(:,:)
24 complex(8), allocatable :: dmatk(:,:,:,:,:)
25 ! zero the density matrix
26 dmat(:,:,:,:,:)=0.d0
27 ! begin parallel loop over k-points
28 call holdthd(nkpt/np_mpi,nthd)
29 !$OMP PARALLEL DEFAULT(SHARED) &
30 !$OMP PRIVATE(idx,apwalm,evecfv,evecsv,dmatk) &
31 !$OMP PRIVATE(ispn,nst,ist,jst,ias,wo) &
32 !$OMP NUM_THREADS(nthd)
33 allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv))
34 allocate(evecfv(nmatmax,nstfv,nspnfv),evecsv(nstsv,nstsv))
35 allocate(dmatk(ld,nspinor,ld,nspinor,nstsv))
36 !$OMP DO SCHEDULE(DYNAMIC)
37 do ik=1,nkpt
38 ! distribute among MPI processes
39  if (mod(ik-1,np_mpi) /= lp_mpi) cycle
40 ! find the matching coefficients
41  do ispn=1,nspnfv
42  call match(ngk(ispn,ik),vgkc(:,:,ispn,ik),gkc(:,ispn,ik), &
43  sfacgk(:,:,ispn,ik),apwalm(:,:,:,:,ispn))
44  end do
45 ! get the eigenvectors from file
46  call getevecfv(filext,ik,vkl(:,ik),vgkl(:,:,:,ik),evecfv)
47  call getevecsv(filext,ik,vkl(:,ik),evecsv)
48 ! index to occupied states
49  nst=0
50  do ist=1,nstsv
51  if (abs(occsv(ist,ik)) < epsocc) cycle
52  nst=nst+1
53  idx(nst)=ist
54  end do
55 ! begin loop over atoms and species
56  do ias=1,natmtot
57  call gendmatk(tspndg,tlmdg,lmin,lmax,ias,nst,idx,ngk(:,ik),apwalm,evecfv, &
58  evecsv,ld,dmatk)
59  do ist=1,nst
60  jst=idx(ist)
61  wo=occsv(jst,ik)*wkpt(ik)
62 !$OMP CRITICAL(gendmat_)
63  dmat(:,:,:,:,ias)=dmat(:,:,:,:,ias)+wo*dmatk(:,:,:,:,ist)
64 !$OMP END CRITICAL(gendmat_)
65  end do
66  end do
67 end do
68 !$OMP END DO
69 deallocate(apwalm,evecfv,evecsv,dmatk)
70 !$OMP END PARALLEL
71 call freethd(nthd)
72 ! add density matrices from each process and redistribute
73 if (np_mpi > 1) then
74  n=((ld*nspinor)**2)*natmtot
75  call mpi_allreduce(mpi_in_place,dmat,n,mpi_double_complex,mpi_sum,mpicom, &
76  ierror)
77 end if
78 ! symmetrise the density matrix
79 call symdmat(lmax,ld,dmat)
80 ! synchronise MPI processes
81 call mpi_barrier(mpicom,ierror)
82 end subroutine
83 
integer nmatmax
Definition: modmain.f90:853
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
character(256) filext
Definition: modmain.f90:1301
integer lmmaxapw
Definition: modmain.f90:199
integer nkpt
Definition: modmain.f90:461
integer ngkmax
Definition: modmain.f90:499
subroutine getevecfv(fext, ikp, vpl, vgpl, evecfv)
Definition: getevecfv.f90:10
subroutine match(ngp, vgpc, gpc, sfacgp, apwalm)
Definition: match.f90:10
Definition: modomp.f90:6
real(8) epsocc
Definition: modmain.f90:903
integer np_mpi
Definition: modmpi.f90:13
complex(8), dimension(:,:,:,:), allocatable sfacgk
Definition: modmain.f90:509
integer, dimension(:,:), allocatable ngk
Definition: modmain.f90:497
real(8), dimension(:), allocatable wkpt
Definition: modmain.f90:475
subroutine gendmat(tspndg, tlmdg, lmin, lmax, ld, dmat)
Definition: gendmat.f90:7
real(8), dimension(:,:,:,:), allocatable vgkl
Definition: modmain.f90:503
real(8), dimension(:,:), allocatable occsv
Definition: modmain.f90:905
subroutine symdmat(lmax, ld, dmat)
Definition: symdmat.f90:7
real(8), dimension(:,:,:,:), allocatable vgkc
Definition: modmain.f90:505
subroutine gendmatk(tspndg, tlmdg, lmin, lmax, ias, nst, idx, ngp, apwalm, evecfv, evecsv, ld, dmat)
Definition: gendmatk.f90:8
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
integer apwordmax
Definition: modmain.f90:760
Definition: modmpi.f90:6
real(8), dimension(:,:,:), allocatable gkc
Definition: modmain.f90:507
integer lp_mpi
Definition: modmpi.f90:15
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer nstfv
Definition: modmain.f90:887
integer mpicom
Definition: modmpi.f90:11
integer nspnfv
Definition: modmain.f90:289
integer ierror
Definition: modmpi.f90:19