The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine genidxthc
7use modmain
8use modtdhfc
9implicit none
10! local variables
11integer ik,jk,nst,nmax,ist
12! determine the maximum number of states within energy window over all k-points
13nmax=0
14do ik=1,nkpt
15 nst=0
16 do ist=1,nstsv
17 if (abs(evalsv(ist,ik)-efermi) < ecutthc) nst=nst+1
18 end do
19 nmax=max(nmax,nst)
20end do
21if (nmax == 0) then
22 write(*,*)
23 write(*,'("Error(genidxthc): no states within energy window ecutthc")')
24 write(*,*)
25 stop
26end if
27! allocate global arrays
28if (allocated(nthck)) deallocate(nthck)
29allocate(nthck(nkpt))
30if (allocated(idxthc)) deallocate(idxthc)
31allocate(idxthc(nmax,nkpt))
32if (allocated(istthc)) deallocate(istthc)
33allocate(istthc(nmax,nkptnr))
34! determine the number of and index to used states
35do ik=1,nkpt
36 nst=0
37 do ist=1,nstsv
38 if (abs(evalsv(ist,ik)-efermi) < ecutthc) then
39 nst=nst+1
40 idxthc(nst,ik)=ist
41 end if
42 end do
43 nthck(ik)=nst
44end do
45! calculate the index to and total number of TDHFC states
46nthc=0
47do ik=1,nkptnr
48 jk=ivkik(ivk(1,ik),ivk(2,ik),ivk(3,ik))
49 do ist=1,nthck(jk)
50 nthc=nthc+1
51 istthc(ist,ik)=nthc
52 end do
53end do
54end subroutine
55
subroutine genidxthc
Definition genidxthc.f90:7
real(8) efermi
Definition modmain.f90:904
integer nkptnr
Definition modmain.f90:463
integer, dimension(:,:), allocatable ivk
Definition modmain.f90:465
integer nkpt
Definition modmain.f90:461
integer nstsv
Definition modmain.f90:886
integer, dimension(:,:,:), allocatable ivkik
Definition modmain.f90:467
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
real(8) ecutthc
Definition modtdhfc.f90:9
integer nthc
Definition modtdhfc.f90:11
integer, dimension(:,:), allocatable istthc
Definition modtdhfc.f90:17
integer, dimension(:), allocatable nthck
Definition modtdhfc.f90:13
integer, dimension(:,:), allocatable idxthc
Definition modtdhfc.f90:15