The Elk Code
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 else
38  allocate(bmt(1,1,1),bir(1,1))
39 end if
40 ! generate the Kohn-Sham potential and magnetic field in spherical coordinates
41 ! and multiply by the radial integration weights; also multiply the interstitial
42 ! potential with the characteristic function
43 call vblocal(vmt,vir,bmt,bir)
44 ! synchronise MPI processes
45 call mpi_barrier(mpicom,ierror)
46 if (mp_mpi) write(*,*)
47 ! loop over k-points
48 call holdthd(nkpt/np_mpi,nthd)
49 !$OMP PARALLEL DO DEFAULT(SHARED) &
50 !$OMP SCHEDULE(DYNAMIC) &
51 !$OMP NUM_THREADS(nthd)
52 do ik=1,nkpt
53 ! distribute among MPI processes
54  if (mod(ik-1,np_mpi) /= lp_mpi) cycle
55 !$OMP CRITICAL(genkmat_)
56  write(*,'("Info(genkmat): ",I6," of ",I6," k-points")') ik,nkpt
57 !$OMP END CRITICAL(genkmat_)
58  call putkmat(tfv,tvclcr,ik,vmt,vir,bmt,bir)
59 end do
60 !$OMP END PARALLEL DO
61 call freethd(nthd)
62 deallocate(vmt,vir,bmt,bir)
63 ! synchronise MPI processes
64 call mpi_barrier(mpicom,ierror)
65 end subroutine
66 !EOC
67 
integer npcmtmax
Definition: modmain.f90:216
logical mp_mpi
Definition: modmpi.f90:17
integer ngtc
Definition: modmain.f90:392
logical spinpol
Definition: modmain.f90:228
integer nkpt
Definition: modmain.f90:461
integer ndmag
Definition: modmain.f90:238
Definition: modomp.f90:6
integer np_mpi
Definition: modmpi.f90:13
Definition: modmpi.f90:6
integer lp_mpi
Definition: modmpi.f90:15
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
subroutine genkmat(tfv, tvclcr)
Definition: genkmat.f90:10
integer natmtot
Definition: modmain.f90:40
subroutine vblocal(vmt, vir, bmt, bir)
Definition: vblocal.f90:7
subroutine putkmat(tfv, tvclcr, ik, vmt, vir, bmt, bir)
Definition: putkmat.f90:7
integer mpicom
Definition: modmpi.f90:11
integer ierror
Definition: modmpi.f90:19