The Elk Code
genwfsvp_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 
6 subroutine genwfsvp_sp(tsh,tgp,nst,idx,ngridg_,igfft_,vpl,ngp,igpig,wfmt,ld, &
7  wfir)
8 use modmain
9 implicit none
10 ! arguments
11 logical, intent(in) :: tsh,tgp
12 integer, intent(in) :: nst,idx(*),ngridg_(3),igfft_(*)
13 real(8), intent(in) :: vpl(3)
14 integer, intent(out) :: ngp(nspnfv),igpig(ngkmax,nspnfv)
15 complex(4), intent(out) :: wfmt(npcmtmax,natmtot,nspinor,nst)
16 integer, intent(in) :: ld
17 complex(4), intent(out) :: wfir(ld,nspinor,nst)
18 ! local variables
19 integer ispn
20 real(8) vl(3),vc(3)
21 ! automatic arrays
22 real(8) vgpl(3,ngkmax,nspnfv),vgpc(3,ngkmax),gpc(ngkmax)
23 ! allocatable arrays
24 complex(8), allocatable :: sfacgp(:,:),apwalm(:,:,:,:,:)
25 complex(8), allocatable :: evecfv(:,:,:),evecsv(:,:)
26 allocate(sfacgp(ngkmax,natmtot))
27 allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv))
28 ! loop over first-variational spins
29 do ispn=1,nspnfv
30  vl(1:3)=vpl(1:3)
31  vc(1:3)=bvec(1:3,1)*vpl(1)+bvec(1:3,2)*vpl(2)+bvec(1:3,3)*vpl(3)
32 ! spin-spiral case
33  if (spinsprl) then
34  if (ispn == 1) then
35  vl(1:3)=vl(1:3)+0.5d0*vqlss(1:3)
36  vc(1:3)=vc(1:3)+0.5d0*vqcss(1:3)
37  else
38  vl(1:3)=vl(1:3)-0.5d0*vqlss(1:3)
39  vc(1:3)=vc(1:3)-0.5d0*vqcss(1:3)
40  end if
41  end if
42 ! generate the G+p-vectors
43  call gengkvec(ngvc,ivg,vgc,vl,vc,gkmax,ngkmax,ngp(ispn),igpig(:,ispn), &
44  vgpl(:,:,ispn),vgpc,gpc)
45 ! generate structure factors for G+p-vectors
46  call gensfacgp(ngp(ispn),vgpc,ngkmax,sfacgp)
47 ! find the matching coefficients
48  call match(ngp(ispn),vgpc,gpc,sfacgp,apwalm(:,:,:,:,ispn))
49 end do
50 deallocate(sfacgp)
51 ! get the first- and second-variational eigenvectors from file
52 allocate(evecfv(nmatmax,nstfv,nspnfv),evecsv(nstsv,nstsv))
53 call getevecfv(filext,0,vpl,vgpl,evecfv)
54 call getevecsv(filext,0,vpl,evecsv)
55 ! calculate the second-variational wavefunctions
56 call genwfsv_sp(tsh,tgp,nst,idx,ngridg_,igfft_,ngp,igpig,apwalm,evecfv,evecsv, &
57  wfmt,ld,wfir)
58 deallocate(apwalm,evecfv,evecsv)
59 end subroutine
60 
integer nmatmax
Definition: modmain.f90:853
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
character(256) filext
Definition: modmain.f90:1301
pure subroutine gensfacgp(ngp, vgpc, ld, sfacgp)
Definition: gensfacgp.f90:10
integer lmmaxapw
Definition: modmain.f90:199
subroutine genwfsvp_sp(tsh, tgp, nst, idx, ngridg_, igfft_, vpl, ngp, igpig, wfmt, ld, wfir)
Definition: genwfsvp_sp.f90:8
subroutine getevecfv(fext, ikp, vpl, vgpl, evecfv)
Definition: getevecfv.f90:10
subroutine match(ngp, vgpc, gpc, sfacgp, apwalm)
Definition: match.f90:10
logical spinsprl
Definition: modmain.f90:283
real(8), dimension(3) vqlss
Definition: modmain.f90:293
integer ngvc
Definition: modmain.f90:398
real(8), dimension(:,:), allocatable vgc
Definition: modmain.f90:420
integer nstsv
Definition: modmain.f90:889
subroutine genwfsv_sp(tsh, tgp, nst, idx, ngridg_, igfft_, ngp, igpig, apwalm, evecfv, evecsv, wfmt, ld, wfir)
Definition: genwfsv_sp.f90:8
real(8), dimension(3, 3) bvec
Definition: modmain.f90:16
integer, dimension(:,:), allocatable ivg
Definition: modmain.f90:400
real(8), dimension(3) vqcss
Definition: modmain.f90:295
pure subroutine gengkvec(ngv, ivg, vgc, vkl, vkc, gkmax, ngkmax, ngk, igkig, vgkl, vgkc, gkc)
Definition: gengkvec.f90:11
integer apwordmax
Definition: modmain.f90:760
real(8) gkmax
Definition: modmain.f90:495
integer nstfv
Definition: modmain.f90:887