File: eveqn.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 (72 lines) | stat: -rw-r--r-- 2,331 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
68
69
70
71
72

! Copyright (C) 2002-2007 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

!BOP
! !ROUTINE: eveqn
subroutine eveqn(ik,evalfv,evecfv,evecsv)
! !USES:
use modmain
use modulr
use modomp
! !INPUT/OUTPUT PARAMETERS:
!   ik     : k-point number (in,integer)
!   evalfv : first-variational eigenvalues (out,real(nstfv))
!   evecfv : first-variational eigenvectors (out,complex(nmatmax,nstfv))
!   evecsv : second-variational eigenvectors (out,complex(nstsv,nstsv))
! !DESCRIPTION:
!   Solves the first- and second-variational eigenvalue equations. See routines
!   {\tt match}, {\tt eveqnfv}, {\tt eveqnss} and {\tt eveqnsv}.
!
! !REVISION HISTORY:
!   Created March 2004 (JKD)
!EOP
!BOC
implicit none
! arguments
integer, intent(in) :: ik
real(8), intent(out) :: evalfv(nstfv,nspnfv)
complex(8), intent(out) :: evecfv(nmatmax,nstfv,nspnfv),evecsv(nstsv,nstsv)
! local variables
integer jspn,nthd
! allocatable arrays
complex(8), allocatable :: apwalm(:,:,:,:,:)
allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv))
! loop over first-variational spins (nspnfv=2 for spin-spirals only)
call omp_hold(nspnfv,nthd)
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP NUM_THREADS(nthd)
!$OMP DO
do jspn=1,nspnfv
! find the matching coefficients
  call match(ngk(jspn,ik),gkc(:,jspn,ik),tpgkc(:,:,jspn,ik), &
   sfacgk(:,:,jspn,ik),apwalm(:,:,:,:,jspn))
! solve the first-variational eigenvalue equation
  if (tefvit) then
! iteratively
    call eveqnit(nmat(jspn,ik),ngk(jspn,ik),igkig(:,jspn,ik),vkl(:,ik), &
     vgkl(:,:,jspn,ik),vgkc(:,:,jspn,ik),apwalm(:,:,:,:,jspn),evalfv(:,jspn), &
     evecfv(:,:,jspn))
  else
! directly
    call eveqnfv(nmat(jspn,ik),ngk(jspn,ik),igkig(:,jspn,ik),vkc(:,ik), &
     vgkc(:,:,jspn,ik),apwalm(:,:,:,:,jspn),evalfv(:,jspn),evecfv(:,:,jspn))
  end if
end do
!$OMP END DO
!$OMP END PARALLEL
call omp_free(nthd)
if (spinsprl) then
! solve the spin-spiral second-variational eigenvalue equation
  call eveqnss(ngk(:,ik),igkig(:,:,ik),apwalm,evalfv,evecfv,evalsv(:,ik),evecsv)
else
! solve the second-variational eigenvalue equation
  call eveqnsv(ngk(1,ik),igkig(:,1,ik),vgkc(:,:,1,ik),apwalm,evalfv,evecfv, &
   evalsv(:,ik),evecsv)
end if
deallocate(apwalm)
return
end subroutine
!EOC