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