The Elk Code
 
Loading...
Searching...
No Matches
geomopt.f90
Go to the documentation of this file.
1
2! Copyright (C) 2011 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
6subroutine geomopt
7use modmain
8use modmpi
9use moddelf
10use modvars
11implicit none
12! local variables
13integer istp,jstp,is,i
14real(8) ds
15! initialise global variables
16call init0
17call init1
18! store orginal volume
20! atomic forces are required
21tforce=.true.
22trdstate=(task == 3)
23! initial atomic step sizes
24if (allocated(tauatp)) deallocate(tauatp)
25allocate(tauatp(natmtot))
27! initialise the previous total force on each atom
28if (allocated(forcetotp)) deallocate(forcetotp)
29allocate(forcetotp(3,natmtot))
30forcetotp(:,:)=0.d0
31! initial lattice optimisation step size
33! initialise previous stress tensor
34stressp(:)=0.d0
35! open formatted files for writing
36if (mp_mpi) then
37 open(71,file='TOTENERGY_OPT.OUT',form='FORMATTED')
38 open(72,file='FORCEMAX.OUT',form='FORMATTED')
39 open(73,file='GEOMETRY_OPT.OUT',form='FORMATTED')
40 open(74,file='IADIST_OPT.OUT',form='FORMATTED')
41 open(75,file='FORCES_OPT.OUT',form='FORMATTED')
42 if (spinpol) then
43 open(78,file='MOMENTM_OPT.OUT',form='FORMATTED')
44 end if
45 if (latvopt /= 0) then
46 open(86,file='STRESSMAX.OUT',form='FORMATTED')
47 open(87,file='STRESS_OPT.OUT',form='FORMATTED')
48 open(88,file='OMEGA_OPT.OUT',form='FORMATTED')
49 end if
50end if
51! synchronise MPI processes
52call mpi_barrier(mpicom,ierror)
53do istp=1,maxlatvstp
54 do jstp=1,maxatpstp
55 if (atpopt == 0) exit
56 if (mp_mpi) then
57 write(*,'("Info(geomopt): atomic position optimisation step : ",I6)') jstp
58 end if
59! ground-state and forces calculation
60 call gndstate
61! check for stop signal
62 if (tstop) goto 10
63! subsequent calculations will read in the potential from STATE.OUT
64 trdstate=.true.
65! update the atomic positions
66 call atpstep
67! write total energy, forces, atomic positions, interatomic distances to file
68 if (mp_mpi) then
69 write(71,'(G24.14)') engytot
70 flush(71)
71 write(72,'(G18.10)') forcemax
72 flush(72)
73 write(73,*); write(73,*)
74 write(73,'("! Atomic position optimisation step : ",I6)') jstp
75 call writegeom(73)
76 flush(73)
77 write(74,*); write(74,*)
78 write(74,'("Atomic position optimisation step : ",I6)') jstp
79 call writeiad(74)
80 flush(74)
81 write(75,*); write(75,*)
82 write(75,'("Atomic position optimisation step : ",I6)') jstp
83 call writeforces(75)
84 write(75,*)
85 write(75,'("Maximum force magnitude over all atoms (target) : ",G18.10,&
86 &" (",G18.10,")")') forcemax,epsforce
87 flush(75)
88 if (spinpol) then
89 write(78,'(G24.14)') momtotm
90 flush(78)
91 end if
92 end if
93! broadcast forcemax from master process to all other processes
94 call mpi_bcast(forcemax,1,mpi_double_precision,0,mpicom,ierror)
95! check force convergence
96 if (forcemax <= epsforce) then
97 if (mp_mpi) then
98 write(75,*)
99 write(75,'("Force convergence target achieved")')
100 end if
101 exit
102 end if
103 if (mp_mpi.and.(jstp == maxatpstp)) then
104 write(*,*)
105 write(*,'("Warning(geomopt): atomic position optimisation failed to &
106 &converge in ",I6," steps")') maxatpstp
107 end if
108! store the current forces array
109 forcetotp(:,:)=forcetot(:,:)
110! end loop over atomic position optimisation
111 end do
112! exit lattice optimisation loop if required
113 if (latvopt == 0) exit
114 if (mp_mpi) then
115 write(*,'("Info(geomopt): lattice optimisation step : ",I6)') istp
116 end if
117! generate the stress tensor
118 call genstress
119! take average of current and previous stress tensors
120 stress(:)=0.5d0*(stress(:)+stressp(:))
121! check for stop signal
122 if (tstop) goto 10
123! update the lattice vectors
124 call latvstep
125! write stress tensor components and maximum magnitude to file
126 if (mp_mpi) then
127 write(71,'(G24.14)') engytot
128 flush(71)
129 write(73,*); write(73,*)
130 write(73,'("! Lattice optimisation step : ",I6)') istp
131 call writegeom(73)
132 flush(73)
133 write(74,*); write(74,*)
134 write(74,'("Lattice optimisation step : ",I6)') istp
135 call writeiad(74)
136 flush(74)
137 if (spinpol) then
138 write(78,'(G24.14)') momtotm
139 flush(78)
140 end if
141 write(86,'(G18.10)') stressmax
142 flush(86)
143 write(87,*)
144 write(87,'("Lattice optimisation step : ",I6)') istp
145 write(87,'("Derivative of total energy w.r.t. strain tensors :")')
146 do i=1,nstrain
147 write(87,'(G18.10)') stress(i)
148 end do
149 flush(87)
150 write(88,'(G18.10)') omega
151 flush(88)
152 end if
153! check for stress convergence
154 if (latvopt == 1) then
155 ds=sum(abs(stress(:)))
156 else
157! stress may be non-zero because of volume constraint; check change in stress
158! tensor instead
159 ds=sum(abs(stress(:)-stressp(:)))
160 end if
161! broadcast ds from master process to all other processes
162 call mpi_bcast(ds,1,mpi_double_precision,0,mpicom,ierror)
163 if ((istp >= 3).and.(ds <= epsstress*tau0latv)) then
164 if (mp_mpi) then
165 write(87,*)
166 write(87,'("Stress convergence target achieved")')
167 end if
168 exit
169 end if
170 if (mp_mpi.and.(istp == maxlatvstp)) then
171 write(*,*)
172 write(*,'("Warning(geomopt): lattice optimisation failed to converge in ",&
173 &I6," steps")') maxlatvstp
174 end if
176! delete the eigenvector files
177 call delfiles(evec=.true.)
178! end loop over lattice optimisation
179end do
18010 continue
181if (mp_mpi) then
182 close(71); close(72); close(73); close(74); close(75)
183 if (spinpol) close(78)
184 if (latvopt /= 0) then
185 close(86); close(87); close(88)
186 end if
187end if
188! ground-state should be run again after lattice optimisation
189if (latvopt /= 0) call gndstate
190! write optimised variables
191if (wrtvars) then
192 call writevars('avec_opt',nv=9,rva=avec)
193 call writevars('omega_opt',rv=omega)
194 do is=1,nspecies
195 call writevars('atposl_opt',n1=is,nv=3*natoms(is),rva=atposl(:,:,is))
196 end do
197 do is=1,nspecies
198 call writevars('atposc_opt',n1=is,nv=3*natoms(is),rva=atposc(:,:,is))
199 end do
200 call writevars('engytot_opt',rv=engytot)
201end if
202end subroutine
203
subroutine atpstep
Definition atpstep.f90:10
subroutine genstress
Definition genstress.f90:7
subroutine geomopt
Definition geomopt.f90:7
subroutine gndstate
Definition gndstate.f90:10
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
subroutine latvstep
Definition latvstep.f90:7
subroutine delfiles(evec, devec, eval, occ, pmat, epsi)
Definition moddelf.f90:25
integer nstrain
Definition modmain.f90:1012
real(8) epsforce
Definition modmain.f90:1062
real(8) forcemax
Definition modmain.f90:994
real(8), dimension(3, maxatoms, maxspecies) atposc
Definition modmain.f90:54
real(8) epsstress
Definition modmain.f90:1064
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
logical spinpol
Definition modmain.f90:228
real(8) omega
Definition modmain.f90:20
integer maxlatvstp
Definition modmain.f90:1036
integer atpopt
Definition modmain.f90:1004
real(8), dimension(9) stressp
Definition modmain.f90:1023
real(8), dimension(3, 3) avec
Definition modmain.f90:12
integer maxatpstp
Definition modmain.f90:1006
real(8), dimension(:), allocatable tauatp
Definition modmain.f90:1010
integer latvopt
Definition modmain.f90:1034
logical tstop
Definition modmain.f90:1052
real(8) tau0latv
Definition modmain.f90:1038
integer natmtot
Definition modmain.f90:40
real(8), dimension(9) stress
Definition modmain.f90:1021
logical trdstate
Definition modmain.f90:682
real(8) omega0
Definition modmain.f90:20
integer task
Definition modmain.f90:1298
real(8), dimension(:,:), allocatable forcetot
Definition modmain.f90:990
real(8) momtotm
Definition modmain.f90:740
real(8) tau0atp
Definition modmain.f90:1008
real(8) stressmax
Definition modmain.f90:1025
real(8), dimension(9) taulatv
Definition modmain.f90:1040
real(8) engytot
Definition modmain.f90:980
logical tforce
Definition modmain.f90:986
real(8), dimension(:,:), allocatable forcetotp
Definition modmain.f90:992
integer nspecies
Definition modmain.f90:34
real(8), dimension(3, maxatoms, maxspecies) atposl
Definition modmain.f90:51
integer ierror
Definition modmpi.f90:19
integer mpicom
Definition modmpi.f90:11
logical mp_mpi
Definition modmpi.f90:17
subroutine writevars(vname, n1, n2, n3, n4, n5, n6, nv, iv, iva, rv, rva, zv, zva, sv, sva)
Definition modvars.f90:16
logical wrtvars
Definition modvars.f90:9
subroutine writeforces(fnum)
subroutine writegeom(fnum)
Definition writegeom.f90:10
subroutine writeiad(fnum)
Definition writeiad.f90:10