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
|
! Copyright (C) 2007 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: rdmminc
! !INTERFACE:
subroutine rdmminc
! !USES:
use modmain
use modrdm
use modmpi
! !DESCRIPTION:
! Minimizes the total energy with respect to the second-variational
! coefficients {\tt evecsv}. The steepest-descent algorithm is used.
!
! !REVISION HISTORY:
! Created 2008 (Sharma)
!EOP
!BOC
implicit none
! local variables
integer it
if (maxitc < 1) return
! begin iteration loop
do it=1,maxitc
if (mp_mpi) then
write(*,'("Info(rdmminc): iteration ",I0," of ",I0)') it,maxitc
end if
! generate the density and magnetisation
call rhomag
! calculate the Coulomb potential
call potcoul
! calculate Coulomb potential matrix elements
call genvmat(vclmt,vclir,vclmat)
! calculate derivative of kinetic energy w.r.t. evecsv
call rdmdkdc
! write the Coulomb matrix elements to file
call writevcl1223
! update evecsv, orthogonalise and write to file
call rdmvaryc
! calculate the energy
call rdmenergy
! write energy to file
if (mp_mpi) then
write(62,'(I6,G18.10)') it,engytot
flush(62)
end if
! end iteration loop
end do
if (mp_mpi) then
write(60,*)
write(60,'("Natural orbital minimisation done")')
write(62,*)
if (spinpol) write(64,*)
end if
end subroutine
!EOC
|