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:
9subroutine rdmvaryc
10! !USES:
11use modmain
12use modrdm
13use 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
22implicit none
23! local variables
24integer ik
25! allocatable arrays
26complex(8), allocatable :: dedc(:,:,:),evecsv(:,:)
27! compute the derivative w.r.t. evecsv
28allocate(dedc(nstsv,nstsv,nkpt))
29call rdmdedc(dedc)
30allocate(evecsv(nstsv,nstsv))
31do 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
43end do
44deallocate(dedc,evecsv)
45! synchronise MPI processes
46call mpi_barrier(mpicom,ierror)
47end subroutine
48!EOC
49
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition getevecsv.f90:7
character(256) filext
Definition modmain.f90:1300
integer nkpt
Definition modmain.f90:461
integer nstsv
Definition modmain.f90:886
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
integer lp_mpi
Definition modmpi.f90:15
integer ierror
Definition modmpi.f90:19
integer mpicom
Definition modmpi.f90:11
integer np_mpi
Definition modmpi.f90:13
real(8) taurdmc
Definition modrdm.f90:19
subroutine putevecsv(fext, ik, evecsv)
Definition putevecsv.f90:7
subroutine rdmdedc(dedc)
Definition rdmdedc.f90:10
subroutine rdmvaryc
Definition rdmvaryc.f90:10
subroutine unitary(n, a)
Definition unitary.f90:10