The Elk Code
Loading...
Searching...
No Matches
genkmat.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2007-2010 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
!BOP
7
! !ROUTINE: genkmat
8
! !INTERFACE:
9
subroutine
genkmat
(tfv,tvclcr)
10
! !USES:
11
use
modmain
12
use
modmpi
13
use
modomp
14
! !INPUT/OUTPUT PARAMETERS:
15
! tfv : .true. if the matrix elements are to be expressed in the
16
! first-variational basis; second-variational otherwise (in,logical)
17
! tvclvr : .true. if the non-local Coulomb potential from the core states is
18
! to be included in the kinetic matrix elements (in,logical)
19
! !DESCRIPTION:
20
! Computes the kinetic matrix elements in the first- or second-variational
21
! basis and stores them in the file {\tt KMAT.OUT}. See routine {\tt putkmat}.
22
!
23
! !REVISION HISTORY:
24
! Created January 2007 (JKD)
25
!EOP
26
!BOC
27
implicit none
28
! arguments
29
logical
,
intent(in)
:: tfv,tvclcr
30
! local variables
31
integer
ik,nthd
32
! allocatable arrays
33
real
(8),
allocatable
:: vmt(:,:),vir(:),bmt(:,:,:),bir(:,:)
34
allocate
(vmt(
npcmtmax
,
natmtot
),vir(
ngtc
))
35
if
(
spinpol
)
then
36
allocate
(bmt(
npcmtmax
,
natmtot
,
ndmag
),bir(
ngtc
,
ndmag
))
37
end if
38
! generate the Kohn-Sham potential and magnetic field in spherical coordinates
39
! and multiply by the radial integration weights; also multiply the interstitial
40
! potential with the characteristic function
41
call
vblocal
(vmt,vir,bmt,bir)
42
! synchronise MPI processes
43
call
mpi_barrier(
mpicom
,
ierror
)
44
if
(
mp_mpi
)
write
(*,*)
45
! loop over k-points
46
call
holdthd
(
nkpt
/
np_mpi
,nthd)
47
!$OMP PARALLEL DO DEFAULT(SHARED) &
48
!$OMP SCHEDULE(DYNAMIC) &
49
!$OMP NUM_THREADS(nthd)
50
do
ik=1,
nkpt
51
! distribute among MPI processes
52
if
(mod(ik-1,
np_mpi
) /=
lp_mpi
) cycle
53
!$OMP CRITICAL(genkmat_)
54
write
(*,
'("Info(genkmat): ",I6," of ",I6," k-points")'
) ik,
nkpt
55
!$OMP END CRITICAL(genkmat_)
56
call
putkmat
(tfv,tvclcr,ik,vmt,vir,bmt,bir)
57
end do
58
!$OMP END PARALLEL DO
59
call
freethd
(nthd)
60
deallocate
(vmt,vir)
61
if
(
spinpol
)
deallocate
(bmt,bir)
62
! synchronise MPI processes
63
call
mpi_barrier(
mpicom
,
ierror
)
64
end subroutine
65
!EOC
66
genkmat
subroutine genkmat(tfv, tvclcr)
Definition
genkmat.f90:10
modmain
Definition
modmain.f90:6
modmain::spinpol
logical spinpol
Definition
modmain.f90:228
modmain::ngtc
integer ngtc
Definition
modmain.f90:392
modmain::nkpt
integer nkpt
Definition
modmain.f90:461
modmain::natmtot
integer natmtot
Definition
modmain.f90:40
modmain::npcmtmax
integer npcmtmax
Definition
modmain.f90:216
modmain::ndmag
integer ndmag
Definition
modmain.f90:238
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
modmpi::mp_mpi
logical mp_mpi
Definition
modmpi.f90:17
modomp
Definition
modomp.f90:6
modomp::holdthd
subroutine holdthd(nloop, nthd)
Definition
modomp.f90:78
modomp::freethd
subroutine freethd(nthd)
Definition
modomp.f90:106
putkmat
subroutine putkmat(tfv, tvclcr, ik, vmt, vir, bmt, bir)
Definition
putkmat.f90:7
vblocal
subroutine vblocal(vmt, vir, bmt, bir)
Definition
vblocal.f90:7
genkmat.f90
Generated by
1.9.8