The Elk Code
wfirsv.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2021 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 wfirsv(tgp,nst,idx,ngridg_,igfft_,ngp,igpig,evecfv,evecsv,ld,wfir)
7 use modmain
8 use modomp
9 implicit none
10 ! arguments
11 logical, intent(in) :: tgp
12 integer, intent(in) :: nst,idx(*),ngridg_(3),igfft_(*)
13 integer, intent(in) :: ngp(nspnfv),igpig(ngkmax,nspnfv)
14 complex(8), intent(in) :: evecfv(nmatmax,nstfv,nspnfv),evecsv(nstsv,nstsv)
15 integer, intent(in) :: ld
16 complex(8), intent(out) :: wfir(ld,nspinor,nst)
17 ! local variables
18 integer ist,ispn,jspn
19 integer n,igp,j,k,nthd
20 real(8) t0
21 complex(8) z1
22 ! automatic arrays
23 complex(8) wfgp(ngkmax)
24 t0=1.d0/sqrt(omega)
25 call holdthd(nst,nthd)
26 !$OMP PARALLEL DO DEFAULT(SHARED) &
27 !$OMP PRIVATE(wfgp,k,ispn,jspn) &
28 !$OMP PRIVATE(n,ist,z1,igp) &
29 !$OMP SCHEDULE(DYNAMIC) &
30 !$OMP NUM_THREADS(nthd)
31 do j=1,nst
32 ! index to state in evecsv
33  if (idx(1) == 0) then; k=j; else; k=idx(j); end if
34  if (tevecsv) then
35 ! generate spinor wavefunction from second-variational eigenvectors
36  do ispn=1,nspinor
37  jspn=jspnfv(ispn)
38  n=ngp(jspn)
39  if (tgp) then
40  wfir(1:n,ispn,j)=0.d0
41  else
42  wfgp(1:n)=0.d0
43  end if
44  do ist=1,nstfv
45  z1=evecsv((ispn-1)*nstfv+ist,k)
46  if (abs(z1%re)+abs(z1%im) < epswf) cycle
47  if (tgp) then
48 ! wavefunction in G+p-space
49  wfir(1:n,ispn,j)=wfir(1:n,ispn,j)+z1*evecfv(1:n,ist,jspn)
50  else
51 ! wavefunction in real-space
52  z1=t0*z1
53  wfgp(1:n)=wfgp(1:n)+z1*evecfv(1:n,ist,jspn)
54  end if
55  end do
56 ! Fourier transform wavefunction to real-space if required
57  if (.not.tgp) then
58  wfir(1:ld,ispn,j)=0.d0
59  do igp=1,n
60  wfir(igfft_(igpig(igp,jspn)),ispn,j)=wfgp(igp)
61  end do
62  call zfftifc(3,ngridg_,1,wfir(:,ispn,j))
63  end if
64  end do
65  else
66 ! spin-unpolarised wavefunction
67  n=ngp(1)
68  if (tgp) then
69  wfir(1:n,1,j)=evecfv(1:n,k,1)
70  else
71  wfir(1:ld,1,j)=0.d0
72  do igp=1,n
73  wfir(igfft_(igpig(igp,1)),1,j)=t0*evecfv(igp,k,1)
74  end do
75  call zfftifc(3,ngridg_,1,wfir(:,1,j))
76  end if
77  end if
78 end do
79 !$OMP END PARALLEL DO
80 call freethd(nthd)
81 end subroutine
82 
logical tevecsv
Definition: modmain.f90:923
real(8) omega
Definition: modmain.f90:20
Definition: modomp.f90:6
subroutine zfftifc(nd, n, sgn, z)
Definition: zfftifc_fftw.f90:7
real(8), parameter epswf
Definition: modmain.f90:845
subroutine wfirsv(tgp, nst, idx, ngridg_, igfft_, ngp, igpig, evecfv, evecsv, ld, wfir)
Definition: wfirsv.f90:7
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer, dimension(2) jspnfv
Definition: modmain.f90:291