The Elk Code
 
Loading...
Searching...
No Matches
rdmdkdc.f90
Go to the documentation of this file.
1
2! Copyright (C) 2007 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: rdmdkdc
8! !INTERFACE:
9subroutine rdmdkdc
10! !USES:
11use modmain
12use modrdm
13use modomp
14! !DESCRIPTION:
15! Calculates the derivative of kinetic energy w.r.t. the second-variational
16! coefficients {\tt evecsv}.
17!
18! !REVISION HISTORY:
19! Created October 2008 (Sharma)
20!EOP
21!BOC
22implicit none
23! local variables
24integer ik,nthd
25! allocatable arrays
26complex(8), allocatable :: evecsv(:,:),kmat(:,:)
27call holdthd(nkpt,nthd)
28!$OMP PARALLEL DEFAULT(SHARED) &
29!$OMP PRIVATE(evecsv,kmat) &
30!$OMP NUM_THREADS(nthd)
31allocate(evecsv(nstsv,nstsv),kmat(nstsv,nstsv))
32!$OMP DO SCHEDULE(DYNAMIC)
33do ik=1,nkpt
34 call getevecsv(filext,ik,vkl(:,ik),evecsv)
35 call getkmat(ik,kmat)
36 call zgemm('N','N',nstsv,nstsv,nstsv,zone,kmat,nstsv,evecsv,nstsv,zzero, &
37 dkdc(:,:,ik),nstsv)
38end do
39!$OMP END DO
40deallocate(evecsv,kmat)
41!$OMP END PARALLEL
42call freethd(nthd)
43end subroutine
44!EOC
45
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition getevecsv.f90:7
subroutine getkmat(ik, kmat)
Definition getkmat.f90:7
complex(8), parameter zzero
Definition modmain.f90:1238
character(256) filext
Definition modmain.f90:1300
complex(8), parameter zone
Definition modmain.f90:1238
integer nkpt
Definition modmain.f90:461
integer nstsv
Definition modmain.f90:886
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
complex(8), dimension(:,:,:), allocatable dkdc
Definition modrdm.f90:15
subroutine rdmdkdc
Definition rdmdkdc.f90:10