The Elk Code
genidxthc.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2024 J. K. Dewhurst and S. Sharma.
3
! This file is distributed under the terms of the GNU General Public License.
4
! See the file COPYING for license details.
5
6
subroutine
genidxthc
7
use
modmain
8
use
modtdhfc
9
implicit none
10
! local variables
11
integer
ik,jk,nst,ist
12
! determine the maximum number of second-variational states within the energy
13
! window over all k-points
14
nmaxthc
=0
15
do
ik=1,
nkpt
16
nst=0
17
do
ist=1,
nstsv
18
if
(abs(
evalsv
(ist,ik)-
efermi
) <
ecutthc
) nst=nst+1
19
end do
20
nmaxthc
=max(
nmaxthc
,nst)
21
end do
22
if
(
nmaxthc
== 0)
then
23
write
(*,*)
24
write
(*,
'("Error(genidxthc): no states within energy window ecutthc")'
)
25
write
(*,*)
26
stop
27
end if
28
! allocate global arrays
29
if
(
allocated
(
nstthc
))
deallocate
(
nstthc
)
30
allocate
(
nstthc
(
nkptnr
))
31
if
(
allocated
(
idxthc
))
deallocate
(
idxthc
)
32
allocate
(
idxthc
(
nmaxthc
,
nkptnr
))
33
if
(
allocated
(
ithc
))
deallocate
(
ithc
)
34
allocate
(
ithc
(
nmaxthc
,
nkptnr
))
35
! determine the number of second-variational states within the energy window for
36
! each k-point as well as the index to and total number of TDHFC states
37
nthc
=0
38
do
ik=1,
nkptnr
39
jk=
ivkik
(
ivk
(1,ik),
ivk
(2,ik),
ivk
(3,ik))
40
nst=0
41
do
ist=1,
nstsv
42
if
(abs(
evalsv
(ist,jk)-
efermi
) <
ecutthc
)
then
43
nst=nst+1
44
idxthc
(nst,ik)=ist
45
nthc
=
nthc
+1
46
ithc
(nst,ik)=
nthc
47
end if
48
end do
49
nstthc
(ik)=nst
50
end do
51
end subroutine
52
modmain::efermi
real(8) efermi
Definition:
modmain.f90:903
modmain::evalsv
real(8), dimension(:,:), allocatable evalsv
Definition:
modmain.f90:915
modtdhfc::ecutthc
real(8) ecutthc
Definition:
modtdhfc.f90:9
modmain::nkpt
integer nkpt
Definition:
modmain.f90:461
modmain::ivkik
integer, dimension(:,:,:), allocatable ivkik
Definition:
modmain.f90:467
modmain::nkptnr
integer nkptnr
Definition:
modmain.f90:463
modmain
Definition:
modmain.f90:6
modmain::nstsv
integer nstsv
Definition:
modmain.f90:885
modtdhfc::nthc
integer nthc
Definition:
modtdhfc.f90:14
modtdhfc::nmaxthc
integer nmaxthc
Definition:
modtdhfc.f90:12
modtdhfc::idxthc
integer, dimension(:,:), allocatable idxthc
Definition:
modtdhfc.f90:18
modtdhfc
Definition:
modtdhfc.f90:6
genidxthc
subroutine genidxthc
Definition:
genidxthc.f90:7
modtdhfc::ithc
integer, dimension(:,:), allocatable ithc
Definition:
modtdhfc.f90:20
modtdhfc::nstthc
integer, dimension(:), allocatable nstthc
Definition:
modtdhfc.f90:16
modmain::ivk
integer, dimension(:,:), allocatable ivk
Definition:
modmain.f90:465
genidxthc.f90
Generated by
1.8.14