The Elk Code
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:
9 subroutine rdmdkdc
10 ! !USES:
11 use modmain
12 use modrdm
13 use 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
22 implicit none
23 ! local variables
24 integer ik,nthd
25 ! allocatable arrays
26 complex(8), allocatable :: evecsv(:,:),kmat(:,:)
27 call holdthd(nkpt,nthd)
28 !$OMP PARALLEL DEFAULT(SHARED) &
29 !$OMP PRIVATE(evecsv,kmat) &
30 !$OMP NUM_THREADS(nthd)
31 allocate(evecsv(nstsv,nstsv),kmat(nstsv,nstsv))
32 !$OMP DO SCHEDULE(DYNAMIC)
33 do 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)
38 end do
39 !$OMP END DO
40 deallocate(evecsv,kmat)
41 !$OMP END PARALLEL
42 call freethd(nthd)
43 end subroutine
44 !EOC
45 
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
character(256) filext
Definition: modmain.f90:1301
integer nkpt
Definition: modmain.f90:461
Definition: modomp.f90:6
complex(8), parameter zone
Definition: modmain.f90:1240
subroutine getkmat(ik, kmat)
Definition: getkmat.f90:7
complex(8), dimension(:,:,:), allocatable dkdc
Definition: modrdm.f90:15
Definition: modrdm.f90:6
integer nstsv
Definition: modmain.f90:889
complex(8), parameter zzero
Definition: modmain.f90:1240
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
subroutine rdmdkdc
Definition: rdmdkdc.f90:10
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78