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 logical tasv
19 integer ist,ispn,jspn
20 integer n,igp,j,k,nthd
21 real(8) t0
22 complex(8) z1
23 ! automatic arrays
24 complex(8) wfgp(ngkmax)
25 ! check if all the second-variational wavefunctions should be calculated
26 tasv=(idx(1) == 0)
27 t0=1.d0/sqrt(omega)
28 call holdthd(nst,nthd)
29 !$OMP PARALLEL DO DEFAULT(SHARED) &
30 !$OMP PRIVATE(wfgp,k,ispn,jspn) &
31 !$OMP PRIVATE(n,ist,z1,igp) &
32 !$OMP SCHEDULE(DYNAMIC) &
33 !$OMP NUM_THREADS(nthd)
34 do j=1,nst
35 ! index to state in evecsv
36  k=merge(j,idx(j),tasv)
37  if (tevecsv) then
38 ! generate spinor wavefunction from second-variational eigenvectors
39  do ispn=1,nspinor
40  jspn=jspnfv(ispn)
41  n=ngp(jspn)
42  if (tgp) then
43  wfir(1:n,ispn,j)=0.d0
44  else
45  wfgp(1:n)=0.d0
46  end if
47  do ist=1,nstfv
48  z1=evecsv((ispn-1)*nstfv+ist,k)
49  if (abs(z1%re)+abs(z1%im) < 1.d-12) cycle
50  if (tgp) then
51 ! wavefunction in G+p-space
52  wfir(1:n,ispn,j)=wfir(1:n,ispn,j)+z1*evecfv(1:n,ist,jspn)
53  else
54 ! wavefunction in real-space
55  z1=t0*z1
56  wfgp(1:n)=wfgp(1:n)+z1*evecfv(1:n,ist,jspn)
57  end if
58  end do
59 ! Fourier transform wavefunction to real-space if required
60  if (.not.tgp) then
61  wfir(1:ld,ispn,j)=0.d0
62  do igp=1,n
63  wfir(igfft_(igpig(igp,jspn)),ispn,j)=wfgp(igp)
64  end do
65  call zfftifc(3,ngridg_,1,wfir(:,ispn,j))
66  end if
67  end do
68  else
69 ! spin-unpolarised wavefunction
70  n=ngp(1)
71  if (tgp) then
72  wfir(1:n,1,j)=evecfv(1:n,k,1)
73  else
74  wfir(1:ld,1,j)=0.d0
75  do igp=1,n
76  wfir(igfft_(igpig(igp,1)),1,j)=t0*evecfv(igp,k,1)
77  end do
78  call zfftifc(3,ngridg_,1,wfir(:,1,j))
79  end if
80  end if
81 end do
82 !$OMP END PARALLEL DO
83 call freethd(nthd)
84 end subroutine
85 
logical tevecsv
Definition: modmain.f90:917
real(8) omega
Definition: modmain.f90:20
Definition: modomp.f90:6
subroutine zfftifc(nd, n, sgn, z)
Definition: zfftifc_fftw.f90:7
subroutine wfirsv(tgp, nst, idx, ngridg_, igfft_, ngp, igpig, evecfv, evecsv, ld, wfir)
Definition: wfirsv.f90:7
subroutine freethd(nthd)
Definition: modomp.f90:112
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer, dimension(2) jspnfv
Definition: modmain.f90:291