The Elk Code
eveqneph.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2019 Chung-Yu Wang, J. K. Dewhurst, S. Sharma and
3 ! E. K. U. Gross. This file is distributed under the terms of the GNU General
4 ! Public License. See the file COPYING for license details.
5 
6 subroutine eveqneph
7 use modmain
8 use modphonon
9 use modbog
10 use modmpi
11 use modomp
12 implicit none
13 ! local variables
14 integer iq,ik
15 integer n,lp,nthd
16 ! allocatable arrays
17 complex(8), allocatable :: dw(:,:,:),ex(:,:,:),fy(:,:)
18 complex(8), allocatable :: au(:,:,:),bv(:,:,:)
19 !------------------------------------!
20 ! phonon eigenvalue equation !
21 !------------------------------------!
22 allocate(dw(nbph,nbph,nqpt),ex(nbph,nbph,nqpt),fy(nbph,nqpt))
23 ! parallel loop over reduced q-point set
24 call holdthd(nqpt/np_mpi,nthd)
25 !$OMP PARALLEL DEFAULT(SHARED) &
26 !$OMP NUM_THREADS(nthd)
27 !$OMP DO SCHEDULE(DYNAMIC)
28 do iq=1,nqpt
29 ! distribute among MPI processes
30  if (mod(iq-1,np_mpi) /= lp_mpi) cycle
31 ! generate the matrices D and E
32  call hmlephde(iq,dw(:,:,iq),ex(:,:,iq))
33 ! zero the vector F
34  fy(:,iq)=0.d0
35 ! solve the phononic Bogoliubov equation
36  call eveqnwxy(nbph,pwxpsn,dw(:,:,iq),ex(:,:,iq),fy(:,iq),evalwx(:,iq))
37 end do
38 !$OMP END DO
39 !$OMP DO SCHEDULE(DYNAMIC)
40 do iq=1,nqpt
41 ! distribute among MPI processes
42  if (mod(iq-1,np_mpi) /= lp_mpi) cycle
43 ! compute the density matrices
44  call dmatwx(nbph,dw(:,:,iq),ex(:,:,iq),dxx(:,:,iq),dwx(:,:,iq),xnorm(:,iq))
45 ! write the eigenvalues, eigenvectors and X-norms to file
46  if (tlast) then
47  call putevalwx(iq,evalwx(:,iq))
48  call putevecwxy(iq,dw(:,:,iq),ex(:,:,iq),fy(:,iq))
49  end if
50 end do
51 !$OMP END DO
52 !$OMP END PARALLEL
53 call freethd(nthd)
54 deallocate(dw,ex,fy)
55 ! broadcast arrays to every MPI process
56 if (np_mpi > 1) then
57  n=nbph*nbph
58  do iq=1,nqpt
59  lp=mod(iq-1,np_mpi)
60  call mpi_bcast(evalwx(:,iq),nbph,mpi_double_precision,lp,mpicom,ierror)
61  call mpi_bcast(xnorm(:,iq),nbph,mpi_double_precision,lp,mpicom,ierror)
62  call mpi_bcast(dxx(:,:,iq),n,mpi_double_complex,lp,mpicom,ierror)
63  call mpi_bcast(dwx(:,:,iq),n,mpi_double_complex,lp,mpicom,ierror)
64  end do
65 end if
66 !--------------------------------------!
67 ! electron eigenvalue equation !
68 !--------------------------------------!
69 allocate(au(nstsv,nstsv,nkpt),bv(nstsv,nstsv,nkpt))
70 ! parallel loop over reduced k-point set
71 call holdthd(nkpt/np_mpi,nthd)
72 !$OMP PARALLEL DEFAULT(SHARED) &
73 !$OMP NUM_THREADS(nthd)
74 !$OMP DO SCHEDULE(DYNAMIC)
75 do ik=1,nkpt
76 ! distribute among MPI processes
77  if (mod(ik-1,np_mpi) /= lp_mpi) cycle
78 ! generate the matrix A
79  call hmlepha(ik,au(:,:,ik))
80 ! generate the matrix B
81  call hmlephb(ik,bv(:,:,ik))
82 ! solve the electronic Bogoliubov equation
83  call eveqnuv(nstsv,au(:,:,ik),bv(:,:,ik),evaluv(:,ik))
84 end do
85 !$OMP END DO
86 !$OMP DO SCHEDULE(DYNAMIC)
87 do ik=1,nkpt
88 ! distribute among MPI processes
89  if (mod(ik-1,np_mpi) /= lp_mpi) cycle
90 ! compute the density matrices
91  call dmatuv(nstsv,efermi,evalsv(:,ik),au(:,:,ik),bv(:,:,ik),dvv(:,:,ik), &
92  duv(:,:,ik),vnorm(:,ik))
93 ! write the eigenvalues and eigenvectors to file
94  if (tlast) then
95  call putevaluv(ik,evaluv(:,ik))
96  call putevecuv(ik,au(:,:,ik),bv(:,:,ik))
97  end if
98 end do
99 !$OMP END DO
100 !$OMP END PARALLEL
101 call freethd(nthd)
102 deallocate(au,bv)
103 ! broadcast arrays to every MPI process
104 if (np_mpi > 1) then
105  n=nstsv*nstsv
106  do ik=1,nkpt
107  lp=mod(ik-1,np_mpi)
108  call mpi_bcast(evaluv(:,ik),nstsv,mpi_double_precision,lp,mpicom,ierror)
109  call mpi_bcast(vnorm(:,ik),nstsv,mpi_double_precision,lp,mpicom,ierror)
110  call mpi_bcast(dvv(:,:,ik),n,mpi_double_complex,lp,mpicom,ierror)
111  call mpi_bcast(duv(:,:,ik),n,mpi_double_complex,lp,mpicom,ierror)
112  end do
113 end if
114 end subroutine
115 
subroutine hmlephb(ik, b)
Definition: hmlephb.f90:7
subroutine hmlepha(ik, a)
Definition: hmlepha.f90:7
real(8) efermi
Definition: modmain.f90:907
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:921
real(8), dimension(:,:), allocatable evalwx
Definition: modbog.f90:37
integer nqpt
Definition: modmain.f90:525
integer nkpt
Definition: modmain.f90:461
real(8), dimension(:,:), allocatable evaluv
Definition: modbog.f90:15
Definition: modomp.f90:6
subroutine putevalwx(iq, evalwxp)
Definition: putevalwx.f90:7
subroutine dmatwx(n, w, x, dxx, dwx, xn)
Definition: dmatwx.f90:7
subroutine eveqneph
Definition: eveqneph.f90:7
complex(8), dimension(:,:,:), pointer, contiguous dwx
Definition: modbog.f90:43
integer np_mpi
Definition: modmpi.f90:13
subroutine putevecwxy(iq, w, x, y)
Definition: putevecwxy.f90:7
integer nstsv
Definition: modmain.f90:889
logical tlast
Definition: modmain.f90:1053
subroutine eveqnwxy(n, p, dw, ex, fy, w)
Definition: eveqnwxy.f90:7
real(8), dimension(:,:), allocatable xnorm
Definition: modbog.f90:41
Definition: modbog.f90:6
integer nbph
Definition: modphonon.f90:13
subroutine dmatuv(n, ef, e, u, v, dvv, duv, vn)
Definition: dmatuv.f90:7
subroutine hmlephde(iq, d, e)
Definition: hmlephde.f90:7
complex(8), dimension(:,:,:), pointer, contiguous duv
Definition: modbog.f90:25
Definition: modmpi.f90:6
complex(8), dimension(:,:,:), pointer, contiguous dxx
Definition: modbog.f90:43
subroutine eveqnuv(n, au, bv, w)
Definition: eveqnuv.f90:7
integer lp_mpi
Definition: modmpi.f90:15
real(8), dimension(:,:), allocatable vnorm
Definition: modbog.f90:17
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer pwxpsn
Definition: modbog.f90:39
subroutine putevaluv(ik, evaluvp)
Definition: putevaluv.f90:7
complex(8), dimension(:,:,:), pointer, contiguous dvv
Definition: modbog.f90:25
subroutine putevecuv(ik, evecu, evecv)
Definition: putevecuv.f90:7
integer mpicom
Definition: modmpi.f90:11
integer ierror
Definition: modmpi.f90:19