The Elk Code
Loading...
Searching...
No Matches
rdmvaryc.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: rdmvaryc
8
! !INTERFACE:
9
subroutine
rdmvaryc
10
! !USES:
11
use
modmain
12
use
modrdm
13
use
modmpi
14
! !DESCRIPTION:
15
! Calculates new {\tt evecsv} from old by using the derivatives of the total
16
! energy w.r.t. {\tt evecsv}. A single step of steepest-descent is made.
17
!
18
! !REVISION HISTORY:
19
! Created 2009 (Sharma)
20
!EOP
21
!BOC
22
implicit none
23
! local variables
24
integer
ik
25
! allocatable arrays
26
complex(8)
,
allocatable
:: dedc(:,:,:),evecsv(:,:)
27
! compute the derivative w.r.t. evecsv
28
allocate
(dedc(
nstsv
,
nstsv
,
nkpt
))
29
call
rdmdedc
(dedc)
30
allocate
(evecsv(
nstsv
,
nstsv
))
31
do
ik=1,
nkpt
32
! distribute among MPI processes
33
if
(mod(ik-1,
np_mpi
) /=
lp_mpi
) cycle
34
! get the eigenvectors from file
35
call
getevecsv
(
filext
,ik,
vkl
(:,ik),evecsv)
36
! calculate new evecsv
37
evecsv(:,:)=evecsv(:,:)-
taurdmc
*dedc(:,:,ik)
38
! othogonalise evecsv
39
call
unitary
(
nstsv
,evecsv)
40
! write new evecsv to file
41
call
putevecsv
(
filext
,ik,evecsv)
42
! end loop over k-points
43
end do
44
deallocate
(dedc,evecsv)
45
! synchronise MPI processes
46
call
mpi_barrier(
mpicom
,
ierror
)
47
end subroutine
48
!EOC
49
getevecsv
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition
getevecsv.f90:7
modmain
Definition
modmain.f90:6
modmain::filext
character(256) filext
Definition
modmain.f90:1300
modmain::nkpt
integer nkpt
Definition
modmain.f90:461
modmain::nstsv
integer nstsv
Definition
modmain.f90:886
modmain::vkl
real(8), dimension(:,:), allocatable vkl
Definition
modmain.f90:471
modmpi
Definition
modmpi.f90:6
modmpi::lp_mpi
integer lp_mpi
Definition
modmpi.f90:15
modmpi::ierror
integer ierror
Definition
modmpi.f90:19
modmpi::mpicom
integer mpicom
Definition
modmpi.f90:11
modmpi::np_mpi
integer np_mpi
Definition
modmpi.f90:13
modrdm
Definition
modrdm.f90:6
modrdm::taurdmc
real(8) taurdmc
Definition
modrdm.f90:19
putevecsv
subroutine putevecsv(fext, ik, evecsv)
Definition
putevecsv.f90:7
rdmdedc
subroutine rdmdedc(dedc)
Definition
rdmdedc.f90:10
rdmvaryc
subroutine rdmvaryc
Definition
rdmvaryc.f90:10
unitary
subroutine unitary(n, a)
Definition
unitary.f90:10
rdmvaryc.f90
Generated by
1.9.8