File: genkmat.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (79 lines) | stat: -rw-r--r-- 2,251 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79

! Copyright (C) 2007-2010 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

!BOP
! !ROUTINE: genkmat
! !INTERFACE:
subroutine genkmat(tfv,tvclcr)
! !USES:
use modmain
use modmpi
use modomp
! !INPUT/OUTPUT PARAMETERS:
!   tfv    : .true. if the matrix elements are to be expressed in the
!            first-variational basis; second-variational otherwise (in,logical)
!   tvclvr : .true. if the non-local Coulomb potential from the core states is
!            to be included in the kinetic matrix elements (in,logical)
! !DESCRIPTION:
!   Computes the kinetic matrix elements in the first- or second-variational
!   basis and stores them in the file {\tt KMAT.OUT}. See routine {\tt putkmat}.
!
! !REVISION HISTORY:
!   Created January 2007 (JKD)
!EOP
!BOC
implicit none
! arguments
logical, intent(in) :: tfv,tvclcr
! local variables
integer ik,is,ias,nthd
! allocatable arrays
real(8), allocatable :: vmt(:,:),vir(:),rfmt(:)
allocate(vmt(npcmtmax,natmtot),vir(ngtot))
! convert muffin-tin Kohn-Sham potential to spherical coordinates on a coarse
! radial mesh
call omp_hold(natmtot,nthd)
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP PRIVATE(rfmt,is) &
!$OMP NUM_THREADS(nthd)
!$OMP DO
do ias=1,natmtot
  allocate(rfmt(npcmtmax))
  is=idxis(ias)
  call rfmtftoc(nrmt(is),nrmti(is),vsmt(:,ias),rfmt)
  call rbsht(nrcmt(is),nrcmti(is),rfmt,vmt(:,ias))
  deallocate(rfmt)
end do
!$OMP END DO
!$OMP END PARALLEL
call omp_free(nthd)
! multiply Kohn-Sham interstitial potential by characteristic function
vir(:)=vsir(:)*cfunir(:)
if (mp_mpi) write(*,*)
! synchronise MPI processes
call mpi_barrier(mpicom,ierror)
! loop over k-points
call omp_hold(nkpt/np_mpi,nthd)
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP NUM_THREADS(nthd)
!$OMP DO
do ik=1,nkpt
! distribute among MPI processes
  if (mod(ik-1,np_mpi).ne.lp_mpi) cycle
!$OMP CRITICAL(genkmat_)
  write(*,'("Info(genkmat): ",I6," of ",I6," k-points")') ik,nkpt
!$OMP END CRITICAL(genkmat_)
  call putkmat(tfv,tvclcr,ik,vmt,vir)
end do
!$OMP END DO
!$OMP END PARALLEL
call omp_free(nthd)
deallocate(vmt,vir)
! synchronise MPI processes
call mpi_barrier(mpicom,ierror)
return
end subroutine
!EOC