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:
9subroutine genkmat(tfv,tvclcr)
10! !USES:
11use modmain
12use modmpi
13use 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
27implicit none
28! arguments
29logical, intent(in) :: tfv,tvclcr
30! local variables
31integer ik,nthd
32! allocatable arrays
33real(8), allocatable :: vmt(:,:),vir(:),bmt(:,:,:),bir(:,:)
34allocate(vmt(npcmtmax,natmtot),vir(ngtc))
35if (spinpol) then
36 allocate(bmt(npcmtmax,natmtot,ndmag),bir(ngtc,ndmag))
37end 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
41call vblocal(vmt,vir,bmt,bir)
42! synchronise MPI processes
43call mpi_barrier(mpicom,ierror)
44if (mp_mpi) write(*,*)
45! loop over k-points
46call holdthd(nkpt/np_mpi,nthd)
47!$OMP PARALLEL DO DEFAULT(SHARED) &
48!$OMP SCHEDULE(DYNAMIC) &
49!$OMP NUM_THREADS(nthd)
50do 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)
57end do
58!$OMP END PARALLEL DO
59call freethd(nthd)
60deallocate(vmt,vir)
61if (spinpol) deallocate(bmt,bir)
62! synchronise MPI processes
63call mpi_barrier(mpicom,ierror)
64end subroutine
65!EOC
66
subroutine genkmat(tfv, tvclcr)
Definition genkmat.f90:10
logical spinpol
Definition modmain.f90:228
integer ngtc
Definition modmain.f90:392
integer nkpt
Definition modmain.f90:461
integer natmtot
Definition modmain.f90:40
integer npcmtmax
Definition modmain.f90:216
integer ndmag
Definition modmain.f90:238
integer lp_mpi
Definition modmpi.f90:15
integer ierror
Definition modmpi.f90:19
integer mpicom
Definition modmpi.f90:11
integer np_mpi
Definition modmpi.f90:13
logical mp_mpi
Definition modmpi.f90:17
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine putkmat(tfv, tvclcr, ik, vmt, vir, bmt, bir)
Definition putkmat.f90:7
subroutine vblocal(vmt, vir, bmt, bir)
Definition vblocal.f90:7