File: eveqnfv.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 (78 lines) | stat: -rw-r--r-- 2,506 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
73
74
75
76
77
78

! Copyright (C) 2002-2005 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: eveqnfv
! !INTERFACE:
subroutine eveqnfv(nmatp,ngp,igpig,vpc,vgpc,apwalm,evalfv,evecfv)
! !USES:
use modmain
use modomp
! !INPUT/OUTPUT PARAMETERS:
!   nmatp  : order of overlap and Hamiltonian matrices (in,integer)
!   ngp    : number of G+p-vectors (in,integer)
!   igpig  : index from G+p-vectors to G-vectors (in,integer(ngkmax))
!   vpc    : p-vector in Cartesian coordinates (in,real(3))
!   vgpc   : G+p-vectors in Cartesian coordinates (in,real(3,ngkmax))
!   apwalm : APW matching coefficients
!            (in,complex(ngkmax,apwordmax,lmmaxapw,natmtot))
!   evalfv : first-variational eigenvalues (out,real(nstfv))
!   evecfv : first-variational eigenvectors (out,complex(nmatmax,nstfv))
! !DESCRIPTION:
!   Solves the eigenvalue equation,
!   $$ (H-\epsilon O)b=0, $$
!   for the all the first-variational states of the input $p$-point.
!
! !REVISION HISTORY:
!   Created March 2004 (JKD)
!EOP
!BOC
implicit none
! arguments
integer, intent(in) :: nmatp,ngp,igpig(ngkmax)
real(8), intent(in) :: vpc(3),vgpc(3,ngkmax)
complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw,natmtot)
real(8), intent(out) :: evalfv(nstfv)
complex(8), intent(out) :: evecfv(nmatmax,nstfv)
! local variables
integer nthd
real(8) ts0,ts1
! allocatable arrays
complex(8), allocatable :: h(:,:),o(:,:)
!-----------------------------------------------!
!     Hamiltonian and overlap matrix set up     !
!-----------------------------------------------!
call timesec(ts0)
allocate(h(nmatp,nmatp),o(nmatp,nmatp))
call omp_hold(2,nthd)
!$OMP PARALLEL SECTIONS DEFAULT(SHARED) &
!$OMP NUM_THREADS(nthd)
!$OMP SECTION
! Hamiltonian
call hmlfv(nmatp,ngp,igpig,vgpc,apwalm,h)
!$OMP SECTION
! overlap
call olpfv(nmatp,ngp,igpig,apwalm,o)
!$OMP END PARALLEL SECTIONS
call omp_free(nthd)
call timesec(ts1)
!$OMP CRITICAL(eveqnfv_)
timemat=timemat+ts1-ts0
!$OMP END CRITICAL(eveqnfv_)
!---------------------------------------!
!     solve the eigenvalue equation     !
!---------------------------------------!
if (tefvr) then
! system has inversion symmetry: use real symmetric matrix eigen solver
  call eveqnfvr(nmatp,ngp,vpc,h,o,evalfv,evecfv)
else
! no inversion symmetry: use complex Hermitian matrix eigen solver
  call eveqnfvz(nmatp,h,o,evalfv,evecfv)
end if
deallocate(h,o)
return
end subroutine
!EOC