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