The Elk Code
rdmminc.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: rdmminc
8 ! !INTERFACE:
9 subroutine rdmminc
10 ! !USES:
11 use modmain
12 use modrdm
13 use modmpi
14 ! !DESCRIPTION:
15 ! Minimizes the total energy with respect to the second-variational
16 ! coefficients {\tt evecsv}. The steepest-descent algorithm is used.
17 !
18 ! !REVISION HISTORY:
19 ! Created 2008 (Sharma)
20 !EOP
21 !BOC
22 implicit none
23 ! local variables
24 integer it
25 if (maxitc < 1) return
26 ! begin iteration loop
27 do it=1,maxitc
28  if (mp_mpi) then
29  write(*,'("Info(rdmminc): iteration ",I4," of ",I4)') it,maxitc
30  end if
31 ! generate the density and magnetisation
32  call rhomag
33 ! calculate the Coulomb potential
34  call potcoul
35 ! calculate Coulomb potential matrix elements
37 ! calculate derivative of kinetic energy w.r.t. evecsv
38  call rdmdkdc
39 ! write the Coulomb matrix elements to file
40  call writevcl1223
41 ! update evecsv, orthogonalise and write to file
42  call rdmvaryc
43 ! calculate the energy
44  call rdmenergy
45 ! write energy to file
46  if (mp_mpi) then
47  write(62,'(I6,G18.10)') it,engytot
48  flush(62)
49  end if
50 ! end iteration loop
51 end do
52 if (mp_mpi) then
53  write(60,*)
54  write(60,'("Natural orbital minimisation done")')
55  write(62,*)
56  if (spinpol) write(64,*)
57 end if
58 end subroutine
59 !EOC
60 
logical mp_mpi
Definition: modmpi.f90:17
logical spinpol
Definition: modmain.f90:228
real(8), dimension(:,:), allocatable vclmt
Definition: modmain.f90:624
subroutine writevcl1223
subroutine genvmat(vmt, vir, vmat)
Definition: genvmat.f90:7
subroutine rhomag
Definition: rhomag.f90:7
Definition: modrdm.f90:6
integer maxitc
Definition: modrdm.f90:27
real(8) engytot
Definition: modmain.f90:983
real(8), dimension(:), allocatable vclir
Definition: modmain.f90:624
subroutine rdmdkdc
Definition: rdmdkdc.f90:10
Definition: modmpi.f90:6
subroutine rdmminc
Definition: rdmminc.f90:10
complex(8), dimension(:,:,:), allocatable vclmat
Definition: modrdm.f90:13
subroutine rdmenergy
Definition: rdmenergy.f90:10
subroutine potcoul
Definition: potcoul.f90:10
subroutine rdmvaryc
Definition: rdmvaryc.f90:10