The Elk Code
 
Loading...
Searching...
No Matches
hartfock.f90
Go to the documentation of this file.
1
2! Copyright (C) 2007 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine hartfock
7use modmain
8use modmpi
9use modomp
10implicit none
11! local variables
12integer ik,lp,nthd
13real(8) etp,de
14character(64) str
15! allocatable arrays
16real(8), allocatable :: vmt(:,:),vir(:)
17real(8), allocatable :: bmt(:,:,:),bir(:,:)
18complex(8), allocatable :: evecsv(:,:)
19! initialise universal variables
20call init0
21call init1
22call init2
23! read the charge density from file
24call readstate
25! Fourier transform Kohn-Sham potential to G-space
26call genvsig
27! generate the core wavefunctions and densities
28call gencore
29! find the new linearisation energies
30call linengy
31! generate the APW and local-orbital radial functions and integrals
32call genapwlofr
33! generate the spin-orbit coupling radial functions
34call gensocfr
35! generate the first- and second-variational eigenvectors and eigenvalues
36call genevfsv
37! find the occupation numbers and Fermi energy
38call occupy
39! generate the kinetic matrix elements in the first-variational basis
40call genkmat(.true.,.true.)
41! allocate local arrays
42allocate(vmt(npcmtmax,natmtot),vir(ngtc))
43if (hybrid.and.spinpol) then
44 allocate(bmt(npcmtmax,natmtot,ndmag),bir(ngtc,ndmag))
45end if
46! only the MPI master process should write files
47if (mp_mpi) then
48! open HF_INFO.OUT file
49 open(60,file='HF_INFO'//trim(filext),form='FORMATTED')
50! open TOTENERGY.OUT
51 open(61,file='TOTENERGY'//trim(filext),form='FORMATTED')
52! open FERMIDOS.OUT
53 open(62,file='FERMIDOS'//trim(filext),form='FORMATTED')
54! open MOMENT.OUT if required
55 if (spinpol) open(63,file='MOMENT'//trim(filext),form='FORMATTED')
56! open GAP.OUT
57 open(64,file='GAP'//trim(filext),form='FORMATTED')
58! open DTOTENERGY.OUT
59 open(66,file='DTOTENERGY'//trim(filext),form='FORMATTED')
60! open MOMENTM.OUT
61 if (spinpol) open(68,file='MOMENTM'//trim(filext),form='FORMATTED')
62! write out general information to HF_INFO.OUT
63 call writeinfo(60)
64end if
65! set last self-consistent loop flag
66tlast=.false.
67etp=0.d0
68! begin the self-consistent loop
69if (mp_mpi) then
70 call writebox(60,"Self-consistent loop started")
71end if
72do iscl=1,maxscl
73 if (mp_mpi) then
74 write(str,'("Loop number : ",I0)') iscl
75 call writebox(60,trim(str))
76 flush(60)
77 write(*,'("Info(hartfock): self-consistent loop number : ",I4)') iscl
78 end if
79 if (iscl >= maxscl) then
80 if (mp_mpi) then
81 write(60,*)
82 write(60,'("Reached self-consistent loops maximum")')
83 end if
84 tlast=.true.
85 end if
86! reset the OpenMP thread variables
87 call omp_reset
88! compute the Hartree-Fock local potentials
89 call hflocal(vmt,vir,bmt,bir)
90! synchronise MPI processes
91 call mpi_barrier(mpicom,ierror)
92 call holdthd(nkpt/np_mpi,nthd)
93!$OMP PARALLEL DEFAULT(SHARED) &
94!$OMP PRIVATE(evecsv) &
95!$OMP NUM_THREADS(nthd)
96 allocate(evecsv(nstsv,nstsv))
97!$OMP DO SCHEDULE(DYNAMIC)
98 do ik=1,nkpt
99! distribute among MPI processes
100 if (mod(ik-1,np_mpi) /= lp_mpi) cycle
101 call getevecsv(filext,ik,vkl(:,ik),evecsv)
102! solve the Hartree-Fock eigenvalue equation
103 call eveqnhf(ik,vmt,vir,bmt,bir,evecsv)
104! write the eigenvalues/vectors to file
105 call putevalsv(filext,ik,evalsv(:,ik))
106 call putevecsv(filext,ik,evecsv)
107 end do
108!$OMP END DO
109 deallocate(evecsv)
110!$OMP END PARALLEL
111 call freethd(nthd)
112! synchronise MPI processes
113 call mpi_barrier(mpicom,ierror)
114! broadcast eigenvalue array to every process
115 do ik=1,nkpt
116 lp=mod(ik-1,np_mpi)
117 call mpi_bcast(evalsv(:,ik),nstsv,mpi_double_precision,lp,mpicom,ierror)
118 end do
119! find the occupation numbers and Fermi energy
120 call occupy
121 if (mp_mpi) then
122! write the occupation numbers to file
123 do ik=1,nkpt
124 call putoccsv(filext,ik,occsv(:,ik))
125 end do
126! write out the eigenvalues and occupation numbers
127 call writeeval
128! write the Fermi energy to file
129 call writefermi
130 end if
131! generate the density and magnetisation
132 call rhomag
133! compute the energy components
134 call energy
135 if (mp_mpi) then
136! output energy components
137 call writeengy(60)
138 write(60,*)
139 write(60,'("Density of states at Fermi energy : ",G18.10)') fermidos
140 write(60,'(" (states/Hartree/unit cell)")')
141 write(60,*)
142 write(60,'("Estimated indirect band gap : ",G18.10)') bandgap(1)
143 write(60,'(" from k-point ",I6," to k-point ",I6)') ikgap(1),ikgap(2)
144 write(60,'("Estimated direct band gap : ",G18.10)') bandgap(2)
145 write(60,'(" at k-point ",I6)') ikgap(3)
146! write total energy to TOTENERGY.OUT
147 write(61,'(G24.14)') engytot
148 flush(61)
149! write DOS at Fermi energy to FERMIDOS.OUT
150 write(62,'(G18.10)') fermidos
151 flush(62)
152! output charges and moments
153 call writechg(60)
154 if (spinpol) then
155 call writemom(60)
156! write total moment to MOMENT.OUT
157 write(63,'(3G18.10)') momtot(1:ndmag)
158 flush(63)
159! write total moment magnitude to MOMENTM.OUT
160 write(68,'(G18.10)') momtotm
161 flush(68)
162 end if
163! write estimated Hartree-Fock indirect band gap
164 write(64,'(G24.14)') bandgap(1)
165 flush(64)
166 end if
167 if (tlast) goto 10
168! compute the change in total energy and check for convergence
169 if (iscl >= 2) then
170 de=abs(engytot-etp)
171 if (mp_mpi) then
172 write(60,*)
173 write(60,'("Absolute change in total energy (target) : ",G18.10," (",&
174 &G18.10,")")') de,epsengy
175 end if
176 if (de < epsengy) then
177 if (mp_mpi) then
178 write(60,*)
179 write(60,'("Energy convergence target achieved")')
180 end if
181 tlast=.true.
182 end if
183 if (mp_mpi) then
184 write(66,'(G18.10)') de
185 flush(66)
186 end if
187 end if
188 etp=engytot
189! check for STOP file
190 call checkstop
191 if (tstop) tlast=.true.
192! broadcast tlast from master process to all other processes
193 call mpi_bcast(tlast,1,mpi_logical,0,mpicom,ierror)
194end do
19510 continue
196if (mp_mpi) then
197 call writebox(60,"Self-consistent loop stopped")
198 if (maxscl > 1) then
199 call writestate
200 write(60,*)
201 write(60,'("Wrote STATE.OUT")')
202 end if
203end if
204!-----------------------!
205! compute forces !
206!-----------------------!
207if (tforce) then
208 call force
209! output forces to HF_INFO.OUT
210 if (mp_mpi) call writeforces(60)
211end if
212if (mp_mpi) then
213! close the HF_INFO.OUT file
214 close(60)
215! close the TOTENERGY.OUT file
216 close(61)
217! close the FERMIDOS.OUT file
218 close(62)
219! close the MOMENT.OUT and MOMENTM.OUT files
220 if (spinpol) then
221 close(63); close(68)
222 end if
223! close the GAP.OUT file
224 close(64)
225! close the DTOTENERGY.OUT file
226 close(66)
227end if
228deallocate(vmt,vir)
229if (hybrid.and.spinpol) deallocate(bmt,bir)
230end subroutine
231
subroutine checkstop
Definition checkstop.f90:7
subroutine energy
Definition energy.f90:10
subroutine eveqnhf(ikp, vmt, vir, bmt, bir, evecsvp)
Definition eveqnhf.f90:7
subroutine force
Definition force.f90:10
subroutine genapwlofr
Definition genapwlofr.f90:7
subroutine gencore
Definition gencore.f90:10
subroutine genevfsv
Definition genevfsv.f90:7
subroutine genkmat(tfv, tvclcr)
Definition genkmat.f90:10
subroutine gensocfr
Definition gensocfr.f90:7
subroutine genvsig
Definition genvsig.f90:10
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition getevecsv.f90:7
subroutine hartfock
Definition hartfock.f90:7
subroutine hflocal(vmt, vir, bmt, bir)
Definition hflocal.f90:7
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
subroutine init2
Definition init2.f90:7
subroutine linengy
Definition linengy.f90:10
real(8) epsengy
Definition modmain.f90:1060
character(256) filext
Definition modmain.f90:1300
logical spinpol
Definition modmain.f90:228
integer, dimension(3) ikgap
Definition modmain.f90:914
real(8), dimension(2) bandgap
Definition modmain.f90:912
logical hybrid
Definition modmain.f90:1149
integer ngtc
Definition modmain.f90:392
logical tstop
Definition modmain.f90:1052
integer nkpt
Definition modmain.f90:461
integer natmtot
Definition modmain.f90:40
real(8), dimension(3) momtot
Definition modmain.f90:738
integer npcmtmax
Definition modmain.f90:216
real(8) fermidos
Definition modmain.f90:910
logical tlast
Definition modmain.f90:1050
integer nstsv
Definition modmain.f90:886
real(8) momtotm
Definition modmain.f90:740
integer maxscl
Definition modmain.f90:1046
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
real(8) engytot
Definition modmain.f90:980
logical tforce
Definition modmain.f90:986
real(8), dimension(:,:), allocatable occsv
Definition modmain.f90:902
integer iscl
Definition modmain.f90:1048
integer ndmag
Definition modmain.f90:238
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
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
logical mp_mpi
Definition modmpi.f90:17
subroutine omp_reset
Definition modomp.f90:71
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine occupy
Definition occupy.f90:10
subroutine putevalsv(fext, ik, evalsv_)
Definition putevalsv.f90:7
subroutine putevecsv(fext, ik, evecsv)
Definition putevecsv.f90:7
subroutine putoccsv(fext, ik, occsvp)
Definition putoccsv.f90:7
subroutine readstate
Definition readstate.f90:10
subroutine rhomag
Definition rhomag.f90:7
subroutine writebox(fnum, str)
Definition writebox.f90:7
subroutine writechg(fnum)
Definition writechg.f90:7
subroutine writeengy(fnum)
Definition writeengy.f90:7
subroutine writeeval
Definition writeeval.f90:10
subroutine writefermi
subroutine writeforces(fnum)
subroutine writeinfo(fnum)
Definition writeinfo.f90:10
subroutine writemom(fnum)
Definition writemom.f90:7
subroutine writestate