File: genwfsvp.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (67 lines) | stat: -rw-r--r-- 2,296 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67

! Copyright (C) 2011 S. Sharma, J. K. Dewhurst and E. K. U. Gross.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

subroutine genwfsvp(tsh,tgp,nst,idx,ngdg,igf,vpl,ngp,igpig,wfmt,ld,wfir)
use modmain
implicit none
! arguments
logical, intent(in) :: tsh,tgp
integer, intent(in) :: nst,idx(nst),ngdg(3),igf(*)
real(8), intent(in) :: vpl(3)
integer, intent(out) :: ngp(nspnfv),igpig(ngkmax,nspnfv)
complex(8), intent(out) :: wfmt(npcmtmax,natmtot,nspinor,nst)
integer, intent(in) :: ld
complex(8), intent(out) :: wfir(ld,nspinor,nst)
! local variables
integer ispn,igp
real(8) vl(3),vc(3)
! allocatable arrays
real(8), allocatable :: vgpl(:,:,:),vgpc(:,:),gpc(:),tpgpc(:,:)
complex(8), allocatable :: sfacgp(:,:),apwalm(:,:,:,:,:)
complex(8), allocatable :: evecfv(:,:,:),evecsv(:,:)
allocate(vgpl(3,ngkmax,nspnfv),vgpc(3,ngkmax))
allocate(gpc(ngkmax),tpgpc(2,ngkmax))
allocate(sfacgp(ngkmax,natmtot))
allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv))
! loop over first-variational spins
do ispn=1,nspnfv
  vl(:)=vpl(:)
  vc(:)=bvec(:,1)*vpl(1)+bvec(:,2)*vpl(2)+bvec(:,3)*vpl(3)
! spin-spiral case
  if (spinsprl) then
    if (ispn.eq.1) then
      vl(:)=vl(:)+0.5d0*vqlss(:)
      vc(:)=vc(:)+0.5d0*vqcss(:)
    else
      vl(:)=vl(:)-0.5d0*vqlss(:)
      vc(:)=vc(:)-0.5d0*vqcss(:)
    end if
  end if
! generate the G+p-vectors
  call gengkvec(ngvec,ivg,vgc,vl,vc,gkmax,ngkmax,ngp(ispn),igpig(:,ispn), &
   vgpl(:,:,ispn),vgpc)
! generate the spherical coordinates of the G+p-vectors
  do igp=1,ngp(ispn)
    call sphcrd(vgpc(:,igp),gpc(igp),tpgpc(:,igp))
  end do
! generate structure factors for G+p-vectors
  call gensfacgp(ngp(ispn),vgpc,ngkmax,sfacgp)
! find the matching coefficients
  call match(ngp(ispn),gpc,tpgpc,sfacgp,apwalm(:,:,:,:,ispn))
end do
deallocate(vgpc,gpc,tpgpc,sfacgp)
! get the first- and second-variational eigenvectors from file
allocate(evecfv(nmatmax,nstfv,nspnfv))
call getevecfv(filext,0,vpl,vgpl,evecfv)
deallocate(vgpl)
allocate(evecsv(nstsv,nstsv))
call getevecsv(filext,0,vpl,evecsv)
! calculate the second-variational wavefunctions
call genwfsv(tsh,tgp,nst,idx,ngdg,igf,ngp,igpig,apwalm,evecfv,evecsv,wfmt,ld, &
 wfir)
deallocate(apwalm,evecfv,evecsv)
return
end subroutine