The Elk Code
eveqnfv.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 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: eveqnfv
8 ! !INTERFACE:
9 subroutine eveqnfv(nmatp,ngp,igpig,vpc,vgpc,apwalm,evalfv,evecfv)
10 ! !USES:
11 use modmain
12 use modomp
13 ! !INPUT/OUTPUT PARAMETERS:
14 ! nmatp : order of overlap and Hamiltonian matrices (in,integer)
15 ! ngp : number of G+p-vectors (in,integer)
16 ! igpig : index from G+p-vectors to G-vectors (in,integer(ngkmax))
17 ! vpc : p-vector in Cartesian coordinates (in,real(3))
18 ! vgpc : G+p-vectors in Cartesian coordinates (in,real(3,ngkmax))
19 ! apwalm : APW matching coefficients
20 ! (in,complex(ngkmax,apwordmax,lmmaxapw,natmtot))
21 ! evalfv : first-variational eigenvalues (out,real(nstfv))
22 ! evecfv : first-variational eigenvectors (out,complex(nmatmax,nstfv))
23 ! !DESCRIPTION:
24 ! Solves the eigenvalue equation,
25 ! $$ (H-\epsilon O)b=0, $$
26 ! for the all the first-variational states of the input $p$-point.
27 !
28 ! !REVISION HISTORY:
29 ! Created March 2004 (JKD)
30 !EOP
31 !BOC
32 implicit none
33 ! arguments
34 integer, intent(in) :: nmatp,ngp,igpig(ngkmax)
35 real(8), intent(in) :: vpc(3),vgpc(3,ngkmax)
36 complex(8), intent(in) :: apwalm(ngkmax,apwordmax,lmmaxapw,natmtot)
37 real(8), intent(out) :: evalfv(nstfv)
38 complex(8), intent(out) :: evecfv(nmatmax,nstfv)
39 ! local variables
40 integer nthd
41 real(8) ts0,ts1
42 ! allocatable arrays
43 complex(8), allocatable :: h(:,:),o(:,:)
44 !-----------------------------------------------!
45 ! Hamiltonian and overlap matrix set up !
46 !-----------------------------------------------!
47 call timesec(ts0)
48 allocate(h(nmatp,nmatp),o(nmatp,nmatp))
49 call holdthd(2,nthd)
50 !$OMP PARALLEL SECTIONS DEFAULT(SHARED) &
51 !$OMP NUM_THREADS(nthd)
52 !$OMP SECTION
53 ! Hamiltonian
54 call hmlfv(nmatp,ngp,igpig,vgpc,apwalm,h)
55 !$OMP SECTION
56 ! overlap
57 call olpfv(nmatp,ngp,igpig,apwalm,o)
58 !$OMP END PARALLEL SECTIONS
59 call freethd(nthd)
60 call timesec(ts1)
61 !$OMP ATOMIC
62 timemat=timemat+ts1-ts0
63 !---------------------------------------!
64 ! solve the eigenvalue equation !
65 !---------------------------------------!
66 if (tefvr) then
67 ! system has inversion symmetry: use real symmetric matrix eigen solver
68  call eveqnfvr(nmatp,ngp,vpc,h,o,evalfv,evecfv)
69 else
70 ! no inversion symmetry: use complex Hermitian matrix eigen solver
71  call eveqnfvz(nmatp,h,o,evalfv,evecfv)
72 end if
73 deallocate(h,o)
74 end subroutine
75 !EOC
76 
subroutine olpfv(nmatp, ngp, igpig, apwalm, o)
Definition: olpfv.f90:7
subroutine eveqnfvz(nmatp, h, o, evalfv, evecfv)
Definition: eveqnfvz.f90:7
Definition: modomp.f90:6
real(8) timemat
Definition: modmain.f90:1217
logical tefvr
Definition: modmain.f90:870
subroutine eveqnfvr(nmatp, ngp, vpc, h_, o_, evalfv, evecfv)
Definition: eveqnfvr.f90:10
subroutine timesec(ts)
Definition: timesec.f90:10
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
subroutine eveqnfv(nmatp, ngp, igpig, vpc, vgpc, apwalm, evalfv, evecfv)
Definition: eveqnfv.f90:10
subroutine hmlfv(nmatp, ngp, igpig, vgpc, apwalm, h)
Definition: hmlfv.f90:7