The Elk Code
Loading...
Searching...
No Matches
genevfsv.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2012 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
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
genevfsv
7
use
modmain
8
use
modmpi
9
use
modomp
10
implicit none
11
! local variables
12
integer
ik,lp,nthd
13
! automatic arrays
14
real
(8) evalfv(nstfv,nspnfv)
15
! allocatable arrays
16
complex(8)
,
allocatable
:: evecfv(:,:,:),evecsv(:,:)
17
! begin parallel loop over k-points
18
call
holdthd
(
nkpt
/
np_mpi
,nthd)
19
!$OMP PARALLEL DEFAULT(SHARED) &
20
!$OMP PRIVATE(evalfv,evecfv,evecsv) &
21
!$OMP NUM_THREADS(nthd)
22
allocate
(evecfv(
nmatmax
,nstfv,nspnfv),evecsv(
nstsv
,
nstsv
))
23
!$OMP DO SCHEDULE(DYNAMIC)
24
do
ik=1,
nkpt
25
! distribute among MPI processes
26
if
(mod(ik-1,
np_mpi
) /=
lp_mpi
) cycle
27
! solve the first- and second-variational eigenvalue equations
28
call
eveqn
(ik,evalfv,evecfv,evecsv)
29
! write the eigenvalues/vectors to file
30
call
putevalfv
(
filext
,ik,evalfv)
31
call
putevalsv
(
filext
,ik,
evalsv
(:,ik))
32
call
putevecfv
(
filext
,ik,evecfv)
33
call
putevecsv
(
filext
,ik,evecsv)
34
end do
35
!$OMP END DO
36
deallocate
(evecfv,evecsv)
37
!$OMP END PARALLEL
38
call
freethd
(nthd)
39
! broadcast eigenvalue array to every MPI process
40
if
(
np_mpi
> 1)
then
41
do
ik=1,
nkpt
42
lp=mod(ik-1,
np_mpi
)
43
call
mpi_bcast(
evalsv
(:,ik),
nstsv
,mpi_double_precision,lp,
mpicom
,
ierror
)
44
end do
45
end if
46
end subroutine
47
eveqn
subroutine eveqn(ik, evalfv, evecfv, evecsv)
Definition
eveqn.f90:10
genevfsv
subroutine genevfsv
Definition
genevfsv.f90:7
modmain
Definition
modmain.f90:6
modmain::filext
character(256) filext
Definition
modmain.f90:1300
modmain::nkpt
integer nkpt
Definition
modmain.f90:461
modmain::nmatmax
integer nmatmax
Definition
modmain.f90:848
modmain::nstsv
integer nstsv
Definition
modmain.f90:886
modmain::evalsv
real(8), dimension(:,:), allocatable evalsv
Definition
modmain.f90:918
modmpi
Definition
modmpi.f90:6
modmpi::lp_mpi
integer lp_mpi
Definition
modmpi.f90:15
modmpi::ierror
integer ierror
Definition
modmpi.f90:19
modmpi::mpicom
integer mpicom
Definition
modmpi.f90:11
modmpi::np_mpi
integer np_mpi
Definition
modmpi.f90:13
modomp
Definition
modomp.f90:6
modomp::holdthd
subroutine holdthd(nloop, nthd)
Definition
modomp.f90:78
modomp::freethd
subroutine freethd(nthd)
Definition
modomp.f90:106
putevalfv
subroutine putevalfv(fext, ik, evalfv)
Definition
putevalfv.f90:7
putevalsv
subroutine putevalsv(fext, ik, evalsv_)
Definition
putevalsv.f90:7
putevecfv
subroutine putevecfv(fext, ik, evecfv)
Definition
putevecfv.f90:7
putevecsv
subroutine putevecsv(fext, ik, evecsv)
Definition
putevecsv.f90:7
genevfsv.f90
Generated by
1.9.8