Octopus
scf_tol.F90
Go to the documentation of this file.
1!! Copyright (C) 2004 E.S. Kadantsev, M. Marques
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
21module scf_tol_oct_m
22 use debug_oct_m
23 use global_oct_m
24 use, intrinsic :: iso_fortran_env
27 use parser_oct_m
29
30 implicit none
31
32 private
33
34 integer, public, parameter :: &
35 SCF_TOL_FIXED = 0, &
36 scf_tol_adaptive = 1, &
37 scf_tol_linear = 2, &
38 scf_tol_exp = 3
39
40 public :: &
41 scf_tol_t, &
48
49 type scf_tol_t
50 private
51 integer, public :: max_iter
52 integer :: scheme
53 real(real64), public :: conv_rel_dens
54 real(real64), public :: conv_abs_dens
55 real(real64) :: conv_threshold_use
56 real(real64) :: dynamic_tol_factor
57 real(real64) :: current_tol
58 real(real64) :: initial_tol
59 real(real64), public :: final_tol
60 integer :: iter_window
61 end type scf_tol_t
62
63contains
64
65 !-----------------------------------------------------------------
66 subroutine scf_tol_init(this, namespace, qtot, def_maximumiter, tol_scheme)
67 type(scf_tol_t), intent(out) :: this
68 type(namespace_t), intent(in) :: namespace
69 real(real64), intent(in) :: qtot
70 integer, optional, intent(in) :: def_maximumiter
71 integer, optional, intent(in) :: tol_scheme
72
73 integer :: def_maximumiter_
74
75 push_sub(scf_tol_init)
76
77 !%Variable LRMaximumIter
78 !%Type integer
79 !%Default 200
80 !%Section Linear Response::SCF in LR calculations
81 !%Description
82 !% The maximum number of SCF iterations to calculate response.
83 !%End
84 def_maximumiter_ = optional_default(def_maximumiter, 200)
85 call parse_variable(namespace, 'LRMaximumIter', def_maximumiter_, this%max_iter)
86
87 !%Variable LRConvAbsDens
88 !%Type float
89 !%Default 1e-5
90 !%Section Linear Response::SCF in LR calculations
91 !%Description
92 !% The tolerance in the absolute variation of the density response, to determine if
93 !% the SCF for linear response is converged.
94 !% <math>\varepsilon = \int {\rm d}^3r \left| \rho^{out}(\bf r) -\rho^{inp}(\bf r) \right|</math>.
95 !% A zero value means do not use this criterion.
96 !%End
97 call parse_variable(namespace, 'LRConvAbsDens', 1e-5_real64, this%conv_abs_dens)
98
99 !%Variable LRConvRelDens
100 !%Type float
101 !%Default 0.0
102 !%Section Linear Response::SCF in LR calculations
103 !%Description
104 !% The tolerance in the relative variation of the density response, to determine if
105 !% the SCF for linear response is converged.
106 !% <math>\varepsilon = \frac{1}{N_{\rm elec}}</math> <tt>LRConvAbsDens</tt>.
107 !% A zero value means do not use this criterion.
108 !%End
109 call parse_variable(namespace, 'LRConvRelDens', m_zero, this%conv_rel_dens)
110
111 ! value to use for adaptive tol scheme
112 if (this%conv_abs_dens <= m_zero) then
113 this%conv_threshold_use = this%conv_rel_dens * qtot
114 else if (this%conv_rel_dens <= m_zero) then
115 this%conv_threshold_use = this%conv_abs_dens
116 else
117 this%conv_threshold_use = min(this%conv_abs_dens, this%conv_rel_dens * qtot)
118 end if
119
120 if (this%conv_abs_dens <= m_zero .and. this%conv_rel_dens <= m_zero) then
121 message(1) = "Input: Not all convergence criteria can be <= 0"
122 message(2) = "Please set one of the following to a nonzero value:"
123 message(3) = "LRConvAbsDens | LRConvRelDens"
124 call messages_fatal(3, namespace=namespace)
125 end if
126
127 !%Variable LRTolScheme
128 !%Type integer
129 !%Default tol_adaptive
130 !%Section Linear Response::SCF in LR calculations
131 !%Description
132 !% The scheme used to adjust the tolerance of the solver during
133 !% the SCF iteration. For <tt>kdotp</tt> and magnetic <tt>em_resp</tt> modes, or
134 !% whenever <tt>HamiltonianVariation = V_ext_only</tt>, the
135 !% scheme is set to <tt>tol_fixed</tt>, and this variable is ignored.
136 !%Option tol_fixed 0
137 !% The solver tolerance is fixed for all the iterations; this
138 !% improves convergence but increases the computational cost
139 !%Option tol_adaptive 1
140 !% The tolerance is increased according to the level of
141 !% convergence of the SCF.
142 !%Option tol_linear 2
143 !% The tolerance decreases linearly for the first <tt>LRTolIterWindow</tt> iterations.
144 !%Option tol_exp 3
145 !% The tolerance decreases exponentially for the first <tt>LRTolIterWindow</tt> iterations.
146 !%End
147 if (present(tol_scheme)) then
148 this%scheme = tol_scheme
149 else
150 call parse_variable(namespace, 'LRTolScheme', scf_tol_adaptive, this%scheme)
151 end if
152 if (.not. varinfo_valid_option('LRTolScheme', this%scheme)) then
153 call messages_input_error(namespace, 'LRTolScheme')
154 end if
156 !%Variable LRTolInitTol
157 !%Type float
158 !%Default 1e-2
159 !%Section Linear Response::Solver
160 !%Description
161 !% This is the tolerance to determine that the linear solver has converged,
162 !% for the first SCF iteration. Ignored if <tt>LRTolScheme = fixed</tt>.
163 !%End
164
165 call parse_variable(namespace, 'LRTolInitTol', 1e-2_real64, this%initial_tol)
166 this%current_tol = this%initial_tol
167
168 !%Variable LRTolFinalTol
169 !%Type float
170 !%Default 1e-6
171 !%Section Linear Response::Solver
172 !%Description
173 !% This is the tolerance to determine that the linear solver has converged.
174 !%End
175
176 call parse_variable(namespace, 'LRTolFinalTol', 1e-6_real64, this%final_tol)
177
178 if (this%scheme == scf_tol_adaptive) then
179 !%Variable LRTolAdaptiveFactor
180 !%Type float
181 !%Default 0.1
182 !%Section Linear Response::SCF in LR calculations
183 !%Description
184 !% This factor controls how much the tolerance is decreased
185 !% during the self-consistency process. Smaller values mean that
186 !% tolerance is decreased faster.
187 !%End
188
189 call parse_variable(namespace, 'LRTolAdaptiveFactor', 0.1_real64, this%dynamic_tol_factor)
190 end if
191
192 if (this%scheme == scf_tol_linear .or. this%scheme == scf_tol_exp) then
193 !%Variable LRTolIterWindow
194 !%Type float
195 !%Default 10
196 !%Section Linear Response::SCF in LR calculations
197 !%Description
198 !% Number of iterations necessary to reach the final tolerance.
199 !%End
200
201 call parse_variable(namespace, 'LRTolIterWindow', 10, this%iter_window)
202 end if
203
204 pop_sub(scf_tol_init)
205
206 end subroutine scf_tol_init
207
208
209 !-----------------------------------------------------------------
210 real(real64) function scf_tol_step(this, iter, scf_res) result(r)
211 type(scf_tol_t), intent(inout) :: this
212 integer, intent(in) :: iter
213 real(real64), intent(in) :: scf_res
214
215 real(real64) :: logi, logf
216
217 push_sub(scf_tol_step)
218
219 if (iter == 0) this%current_tol = m_huge
220
221 select case (this%scheme)
222 case (scf_tol_fixed)
223 r = this%final_tol
224
225 case (scf_tol_adaptive)
226 if (iter == 0) then
227 r = this%initial_tol
228 else
229 r = this%dynamic_tol_factor * (this%final_tol/this%conv_threshold_use)*scf_res
230 end if
231
232 case (scf_tol_linear)
233 r = this%initial_tol + (this%final_tol - this%initial_tol) * &
234 real(iter, real64) / real(this%iter_window, real64)
235
236 case (scf_tol_exp)
237 logi = log(this%initial_tol)
238 logf = log(this%final_tol)
239 r = logi + (logf - logi) * &
240 real(iter, real64) / real(this%iter_window, real64)
241 r = exp(r)
242 end select
243
244 ! tolerance can never be larger than final tolerance
245 r = max(r, this%final_tol)
246 ! tolerance always has to decrease
247 r = min(r, this%current_tol)
248
249 this%current_tol = r
250
251 pop_sub(scf_tol_step)
252 end function scf_tol_step
253
254
255 !-----------------------------------------------------------------
256 subroutine scf_tol_stop(this)
257 type(scf_tol_t), intent(inout) :: this
258
259 push_sub(scf_tol_stop)
260 this%current_tol = m_zero
261
262 pop_sub(scf_tol_stop)
263 end subroutine scf_tol_stop
264
265
266 !-----------------------------------------------------------------
267 subroutine scf_tol_end(this)
268 type(scf_tol_t), intent(inout) :: this
269
270 push_sub(scf_tol_end)
271 this%current_tol = m_zero
272
273 pop_sub(scf_tol_end)
274 end subroutine scf_tol_end
275
276
277 !-----------------------------------------------------------------
278 real(real64) function scf_tol_final(this)
279 type(scf_tol_t), intent(in) :: this
280
281 scf_tol_final = this%final_tol
282
283 end function scf_tol_final
284
285 !-----------------------------------------------------------------
286
287 subroutine scf_tol_obsolete_variables(namespace, old_prefix, new_prefix)
288 type(namespace_t), intent(in) :: namespace
289 character(len=*), intent(in) :: old_prefix
290 character(len=*), intent(in) :: new_prefix
291
293
294 call messages_obsolete_variable(namespace, trim(old_prefix)//'LRMaximumIter', trim(new_prefix)//'LRMaximumIter')
295 call messages_obsolete_variable(namespace, trim(old_prefix)//'LRConvAbsDens', trim(new_prefix)//'LRConvAbsDens')
296 call messages_obsolete_variable(namespace, trim(old_prefix)//'LRTolScheme', trim(new_prefix)//'LRTolScheme')
297 call messages_obsolete_variable(namespace, trim(old_prefix)//'LRTolInitTol', trim(new_prefix)//'LRTolInitTol')
298 call messages_obsolete_variable(namespace, trim(old_prefix)//'LRTolFinalTol', trim(new_prefix)//'LRTolFinalTol')
299 call messages_obsolete_variable(namespace, trim(old_prefix)//'LRTolAdaptiveFactor', trim(new_prefix)//'LRTolAdaptiveFactor')
300 call messages_obsolete_variable(namespace, trim(old_prefix)//'LRTolIterWindow', trim(new_prefix)//'LRTolIterWindow')
301
303 end subroutine scf_tol_obsolete_variables
304
306
307!! Local Variables:
308!! mode: f90
309!! coding: utf-8
310!! End:
double log(double __x) __attribute__((__nothrow__
double exp(double __x) __attribute__((__nothrow__
real(real64), parameter, public m_huge
Definition: global.F90:208
real(real64), parameter, public m_zero
Definition: global.F90:190
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:416
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:697
real(real64) function, public scf_tol_step(this, iter, scf_res)
Definition: scf_tol.F90:306
integer, parameter, public scf_tol_adaptive
Definition: scf_tol.F90:129
subroutine, public scf_tol_init(this, namespace, qtot, def_maximumiter, tol_scheme)
Definition: scf_tol.F90:162
subroutine, public scf_tol_end(this)
Definition: scf_tol.F90:363
real(real64) function, public scf_tol_final(this)
Definition: scf_tol.F90:374
integer, parameter, public scf_tol_exp
Definition: scf_tol.F90:129
integer, parameter, public scf_tol_linear
Definition: scf_tol.F90:129
subroutine, public scf_tol_obsolete_variables(namespace, old_prefix, new_prefix)
Definition: scf_tol.F90:383
subroutine, public scf_tol_stop(this)
Definition: scf_tol.F90:352