The Elk Code
 
Loading...
Searching...
No Matches
eveqn.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2007 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: eveqn
8! !INTERFACE:
9subroutine eveqn(ik,evalfv,evecfv,evecsv)
10! !USES:
11use modmain
12! !INPUT/OUTPUT PARAMETERS:
13! ik : k-point number (in,integer)
14! evalfv : first-variational eigenvalues (out,real(nstfv))
15! evecfv : first-variational eigenvectors (out,complex(nmatmax,nstfv))
16! evecsv : second-variational eigenvectors (out,complex(nstsv,nstsv))
17! !DESCRIPTION:
18! Solves the first- and second-variational eigenvalue equations. See routines
19! {\tt match}, {\tt eveqnfv}, {\tt eveqnss} and {\tt eveqnsv}.
20!
21! !REVISION HISTORY:
22! Created March 2004 (JKD)
23!EOP
24!BOC
25implicit none
26! arguments
27integer, intent(in) :: ik
28real(8), intent(out) :: evalfv(nstfv,nspnfv)
29complex(8), intent(out) :: evecfv(nmatmax,nstfv,nspnfv),evecsv(nstsv,nstsv)
30! local variables
31integer jspn
32! allocatable arrays
33complex(8), allocatable :: apwalm(:,:,:,:,:)
34allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv))
35! loop over first-variational spins (nspnfv=2 for spin-spirals only)
36do jspn=1,nspnfv
37! find the matching coefficients
38 call match(ngk(jspn,ik),vgkc(:,:,jspn,ik),gkc(:,jspn,ik), &
39 sfacgk(:,:,jspn,ik),apwalm(:,:,:,:,jspn))
40! solve the first-variational eigenvalue equation
41 if (tefvit) then
42! iteratively
43 call eveqnit(nmat(jspn,ik),ngk(jspn,ik),igkig(:,jspn,ik),vkl(:,ik), &
44 vgkl(:,:,jspn,ik),vgkc(:,:,jspn,ik),apwalm(:,:,:,:,jspn),evalfv(:,jspn), &
45 evecfv(:,:,jspn))
46 else
47! directly
48 call eveqnfv(nmat(jspn,ik),ngk(jspn,ik),igkig(:,jspn,ik),vkc(:,ik), &
49 vgkc(:,:,jspn,ik),apwalm(:,:,:,:,jspn),evalfv(:,jspn),evecfv(:,:,jspn))
50 end if
51end do
52if (spinsprl) then
53! solve the spin-spiral second-variational eigenvalue equation
54 call eveqnss(ngk(:,ik),igkig(:,:,ik),apwalm,evalfv,evecfv,evalsv(:,ik),evecsv)
55else
56! solve the second-variational eigenvalue equation
57 call eveqnsv(ngk(1,ik),igkig(:,1,ik),vgkc(:,:,1,ik),apwalm,evalfv,evecfv, &
58 evalsv(:,ik),evecsv)
59end if
60deallocate(apwalm)
61end subroutine
62!EOC
63
subroutine eveqn(ik, evalfv, evecfv, evecsv)
Definition eveqn.f90:10
subroutine eveqnfv(nmatp, ngp, igpig, vpc, vgpc, apwalm, evalfv, evecfv)
Definition eveqnfv.f90:10
subroutine eveqnit(nmatp, ngp, igpig, vpl, vgpl, vgpc, apwalm, evalfv, evecfv)
Definition eveqnit.f90:7
subroutine eveqnss(ngp, igpig, apwalm, evalfv, evecfv, evalsv_, evecsv)
Definition eveqnss.f90:7
subroutine eveqnsv(ngp, igpig, vgpc, apwalm, evalfv, evecfv, evalsv_, evecsv)
Definition eveqnsv.f90:8
subroutine match(ngp, vgpc, gpc, sfacgp, apwalm)
Definition match.f90:10
real(8), dimension(:,:,:,:), allocatable vgkc
Definition modmain.f90:505
real(8), dimension(:,:,:), allocatable gkc
Definition modmain.f90:507
logical tefvit
Definition modmain.f90:868
integer, dimension(:,:), allocatable ngk
Definition modmain.f90:497
integer, dimension(:,:,:), allocatable igkig
Definition modmain.f90:501
integer lmmaxapw
Definition modmain.f90:199
integer apwordmax
Definition modmain.f90:760
integer natmtot
Definition modmain.f90:40
real(8), dimension(:,:,:,:), allocatable vgkl
Definition modmain.f90:503
logical spinsprl
Definition modmain.f90:283
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
integer ngkmax
Definition modmain.f90:499
integer, dimension(:,:), allocatable nmat
Definition modmain.f90:846
real(8), dimension(:,:), allocatable vkc
Definition modmain.f90:473
complex(8), dimension(:,:,:,:), allocatable sfacgk
Definition modmain.f90:509
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918