The Elk Code
 
Loading...
Searching...
No Matches
rdmminn.f90
Go to the documentation of this file.
1
2! Copyright (C) 2007-2008 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: rdmminn
8! !INTERFACE:
9subroutine rdmminn
10! !USES:
11use modmain
12use modrdm
13use modmpi
14! !DESCRIPTION:
15! Minimizes the total energy w.r.t. occupation numbers. The steepest-descent
16! algorithm is used.
17!
18! !REVISION HISTORY:
19! Created 2008 (Sharma)
20!EOP
21!BOC
22implicit none
23! local variables
24integer it,n
25if (maxitn < 1) return
26! write the Coulomb matrix elements to file
27call writevcl1221
28! calculate derivative of kinetic energy w.r.t. evecsv
29call rdmdkdc
30! begin iteration loop
31do it=1,maxitn
32 if (mp_mpi) then
33 if (mod(it,10) == 0) then
34 write(*,'("Info(rdmminn): iteration ",I4," of ",I4)') it,maxitn
35 end if
36 end if
37! generate the density and magnetisation
38 call rhomag
39! calculate the Coulomb potential
40 call potcoul
41! calculate Coulomb potential matrix elements
43! update occupation numbers and write to file (MPI master process only)
44 if (mp_mpi) call rdmvaryn
45! broadcast occupation numbers to all other processes
46 n=nstsv*nkpt
47 call mpi_bcast(occsv,n,mpi_double_precision,0,mpicom,ierror)
48! calculate the energy
49 call rdmenergy
50! write energy to file
51 if (mp_mpi) then
52 write(61,'(I6,G18.10)') it,engytot
53 flush(61)
54 end if
55! end iteration loop
56end do
57if (mp_mpi) then
58 write(60,*)
59 write(60,'("Occupation number minimisation done")')
60 write(61,*)
61 if (spinpol) write(63,*)
62end if
63end subroutine
64!EOC
65
subroutine genvmat(vmt, vir, vmat)
Definition genvmat.f90:7
logical spinpol
Definition modmain.f90:228
integer nkpt
Definition modmain.f90:461
real(8), dimension(:), allocatable vclir
Definition modmain.f90:624
integer nstsv
Definition modmain.f90:886
real(8), dimension(:,:), allocatable vclmt
Definition modmain.f90:624
real(8) engytot
Definition modmain.f90:980
real(8), dimension(:,:), allocatable occsv
Definition modmain.f90:902
integer ierror
Definition modmpi.f90:19
integer mpicom
Definition modmpi.f90:11
logical mp_mpi
Definition modmpi.f90:17
complex(8), dimension(:,:,:), allocatable vclmat
Definition modrdm.f90:13
integer maxitn
Definition modrdm.f90:25
subroutine potcoul
Definition potcoul.f90:10
subroutine rdmdkdc
Definition rdmdkdc.f90:10
subroutine rdmenergy
Definition rdmenergy.f90:10
subroutine rdmminn
Definition rdmminn.f90:10
subroutine rdmvaryn
Definition rdmvaryn.f90:10
subroutine rhomag
Definition rhomag.f90:7
subroutine writevcl1221