The Elk Code
rdmft.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: rdmft
8 ! !INTERFACE:
9 subroutine rdmft
10 ! !USES:
11 use modmain
12 use modrdm
13 use modmpi
14 ! !DESCRIPTION:
15 ! Main routine for one-body reduced density matrix functional theory (RDMFT).
16 !
17 ! !REVISION HISTORY:
18 ! Created 2008 (Sharma)
19 !EOP
20 !BOC
21 implicit none
22 ! local variables
23 character(64) str
24 ! initialise global variables
25 call init0
26 call init1
27 ! generate q-point set and gclq array
28 call init2
29 ! read density and potentials from file
30 call readstate
31 ! Fourier transform Kohn-Sham potential to G-space
32 call genvsig
33 ! generate the core wavefunctions and densities
34 call gencore
35 ! find the new linearisation energies
36 call linengy
37 ! generate the APW and local-orbital radial functions and integrals
38 call genapwlofr
39 ! generate the spin-orbit coupling radial functions
40 call gensocfr
41 ! compute the kinetic energy of the core
42 call energykncr
43 ! generate the first- and second-variational eigenvectors and eigenvalues
44 call genevfsv
45 ! find the occupation numbers
46 call occupy
47 ! generate the kinetic matrix elements in the first-variational basis
48 call genkmat(.true.,.false.)
49 ! open information files (MPI master process only)
50 if (mp_mpi) then
51  open(60,file='RDM_INFO.OUT',form='FORMATTED')
52  open(61,file='RDMN_ENERGY.OUT',form='FORMATTED')
53  open(62,file='RDMC_ENERGY.OUT',form='FORMATTED')
54  if (spinpol) then
55  open(63,file='RDMN_MOMENT.OUT',form='FORMATTED')
56  open(64,file='RDMC_MOMENT.OUT',form='FORMATTED')
57  end if
58  open(65,file='RDM_ENERGY.OUT',form='FORMATTED')
59 ! write out general information to RDM_INFO.OUT
60  call writeinfo(60)
61  call writebox(60,"Self-consistent loop started")
62 end if
63 ! begin main RDMFT self-consistent loop
64 do iscl=1,rdmmaxscl
65  if (mp_mpi) then
66  write(str,'("Loop number : ",I0)') iscl
67  call writebox(60,trim(str))
68  flush(60)
69  write(*,'("Info(rdmft): self-consistent loop number : ",I4)') iscl
70  end if
71 ! synchronise MPI processes
72  call mpi_barrier(mpicom,ierror)
73 ! minimisation over natural orbitals
74  call rdmminc
75 ! minimisation over occupation number
76  call rdmminn
77 ! compute the RDMFT 'eigenvalues'
78  call rdmeval
79  if (mp_mpi) then
80  call rdmwriteengy(60)
81  call writechg(60)
82  if (spinpol) call writemom(60)
83  call writeeval
84 ! write out the total energy
85  write(65,'(G18.10)') engytot
86  flush(65)
87  end if
88 end do
89 if (mp_mpi) then
90  call writebox(60,"Self-consistent loop stopped")
91 ! write density to STATE.OUT
92  call writestate
93  write(60,*)
94  write(60,'("Wrote STATE.OUT")')
95 ! close information files
96  close(60)
97  close(61)
98  close(62)
99  if (spinpol) then
100  close(63)
101  close(64)
102  end if
103  close(65)
104 end if
105 end subroutine
106 !EOC
107 
subroutine rdmminn
Definition: rdmminn.f90:10
subroutine genevfsv
Definition: genevfsv.f90:7
logical mp_mpi
Definition: modmpi.f90:17
subroutine occupy
Definition: occupy.f90:10
logical spinpol
Definition: modmain.f90:228
subroutine writeinfo(fnum)
Definition: writeinfo.f90:10
subroutine rdmwriteengy(fnum)
integer iscl
Definition: modmain.f90:1051
subroutine gencore
Definition: gencore.f90:10
subroutine writestate
Definition: writestate.f90:10
subroutine genvsig
Definition: genvsig.f90:10
subroutine genapwlofr
Definition: genapwlofr.f90:7
subroutine gensocfr
Definition: gensocfr.f90:7
subroutine linengy
Definition: linengy.f90:10
Definition: modrdm.f90:6
real(8) engytot
Definition: modmain.f90:983
subroutine energykncr
Definition: energykncr.f90:7
subroutine init2
Definition: init2.f90:7
subroutine init1
Definition: init1.f90:10
Definition: modmpi.f90:6
subroutine writechg(fnum)
Definition: writechg.f90:7
subroutine writebox(fnum, str)
Definition: writebox.f90:7
subroutine readstate
Definition: readstate.f90:10
subroutine rdmminc
Definition: rdmminc.f90:10
subroutine rdmft
Definition: rdmft.f90:10
subroutine genkmat(tfv, tvclcr)
Definition: genkmat.f90:10
subroutine init0
Definition: init0.f90:10
integer rdmmaxscl
Definition: modrdm.f90:23
subroutine rdmeval
Definition: rdmeval.f90:10
subroutine writeeval
Definition: writeeval.f90:10
subroutine writemom(fnum)
Definition: writemom.f90:7
integer mpicom
Definition: modmpi.f90:11
integer ierror
Definition: modmpi.f90:19