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_ = 200
85 if (present(def_maximumiter)) def_maximumiter_ = def_maximumiter
86
87 call parse_variable(namespace, 'LRMaximumIter', def_maximumiter_, this%max_iter)
88
89 !%Variable LRConvAbsDens
90 !%Type float
91 !%Default 1e-5
92 !%Section Linear Response::SCF in LR calculations
93 !%Description
94 !% The tolerance in the absolute variation of the density response, to determine if
95 !% the SCF for linear response is converged.
96 !% <math>\varepsilon = \int {\rm d}^3r \left| \rho^{out}(\bf r) -\rho^{inp}(\bf r) \right|</math>.
97 !% A zero value means do not use this criterion.
98 !%End
99
100 call parse_variable(namespace, 'LRConvAbsDens', 1e-5_real64, this%conv_abs_dens)
101
102 !%Variable LRConvRelDens
103 !%Type float
104 !%Default 0.0
105 !%Section Linear Response::SCF in LR calculations
106 !%Description
107 !% The tolerance in the relative variation of the density response, to determine if
108 !% the SCF for linear response is converged.
109 !% <math>\varepsilon = \frac{1}{N_{\rm elec}}</math> <tt>LRConvAbsDens</tt>.
110 !% A zero value means do not use this criterion.
111 !%End
112
113 call parse_variable(namespace, 'LRConvRelDens', m_zero, this%conv_rel_dens)
115 ! value to use for adaptive tol scheme
116 if (this%conv_abs_dens <= m_zero) then
117 this%conv_threshold_use = this%conv_rel_dens * qtot
118 else if (this%conv_rel_dens <= m_zero) then
119 this%conv_threshold_use = this%conv_abs_dens
120 else
121 this%conv_threshold_use = min(this%conv_abs_dens, this%conv_rel_dens * qtot)
122 end if
123
124 if (this%conv_abs_dens <= m_zero .and. this%conv_rel_dens <= m_zero) then
125 message(1) = "Input: Not all convergence criteria can be <= 0"
126 message(2) = "Please set one of the following to a nonzero value:"
127 message(3) = "LRConvAbsDens | LRConvRelDens"
128 call messages_fatal(3, namespace=namespace)
129 end if
130
131 !%Variable LRTolScheme
132 !%Type integer
133 !%Default tol_adaptive
134 !%Section Linear Response::SCF in LR calculations
135 !%Description
136 !% The scheme used to adjust the tolerance of the solver during
137 !% the SCF iteration. For <tt>kdotp</tt> and magnetic <tt>em_resp</tt> modes, or
138 !% whenever <tt>HamiltonianVariation = V_ext_only</tt>, the
139 !% scheme is set to <tt>tol_fixed</tt>, and this variable is ignored.
140 !%Option tol_fixed 0
141 !% The solver tolerance is fixed for all the iterations; this
142 !% improves convergence but increases the computational cost
143 !%Option tol_adaptive 1
144 !% The tolerance is increased according to the level of
145 !% convergence of the SCF.
146 !%Option tol_linear 2
147 !% The tolerance decreases linearly for the first <tt>LRTolIterWindow</tt> iterations.
148 !%Option tol_exp 3
149 !% The tolerance decreases exponentially for the first <tt>LRTolIterWindow</tt> iterations.
150 !%End
151 if (present(tol_scheme)) then
152 this%scheme = tol_scheme
153 else
154 call parse_variable(namespace, 'LRTolScheme', scf_tol_adaptive, this%scheme)
155 end if
156 if (.not. varinfo_valid_option('LRTolScheme', this%scheme)) then
157 call messages_input_error(namespace, 'LRTolScheme')
158 end if
160 !%Variable LRTolInitTol
161 !%Type float
162 !%Default 1e-2
163 !%Section Linear Response::Solver
164 !%Description
165 !% This is the tolerance to determine that the linear solver has converged,
166 !% for the first SCF iteration. Ignored if <tt>LRTolScheme = fixed</tt>.
167 !%End
168
169 call parse_variable(namespace, 'LRTolInitTol', 1e-2_real64, this%initial_tol)
170 this%current_tol = this%initial_tol
171
172 !%Variable LRTolFinalTol
173 !%Type float
174 !%Default 1e-6
175 !%Section Linear Response::Solver
176 !%Description
177 !% This is the tolerance to determine that the linear solver has converged.
178 !%End
179
180 call parse_variable(namespace, 'LRTolFinalTol', 1e-6_real64, this%final_tol)
181
182 if (this%scheme == scf_tol_adaptive) then
183 !%Variable LRTolAdaptiveFactor
184 !%Type float
185 !%Default 0.1
186 !%Section Linear Response::SCF in LR calculations
187 !%Description
188 !% This factor controls how much the tolerance is decreased
189 !% during the self-consistency process. Smaller values mean that
190 !% tolerance is decreased faster.
191 !%End
192
193 call parse_variable(namespace, 'LRTolAdaptiveFactor', 0.1_real64, this%dynamic_tol_factor)
194 end if
195
196 if (this%scheme == scf_tol_linear .or. this%scheme == scf_tol_exp) then
197 !%Variable LRTolIterWindow
198 !%Type float
199 !%Default 10
200 !%Section Linear Response::SCF in LR calculations
201 !%Description
202 !% Number of iterations necessary to reach the final tolerance.
203 !%End
204
205 call parse_variable(namespace, 'LRTolIterWindow', 10, this%iter_window)
206 end if
207
208 pop_sub(scf_tol_init)
209
210 end subroutine scf_tol_init
211
212
213 !-----------------------------------------------------------------
214 real(real64) function scf_tol_step(this, iter, scf_res) result(r)
215 type(scf_tol_t), intent(inout) :: this
216 integer, intent(in) :: iter
217 real(real64), intent(in) :: scf_res
218
219 real(real64) :: logi, logf
220
221 push_sub(scf_tol_step)
222
223 if (iter == 0) this%current_tol = m_huge
224
225 select case (this%scheme)
226 case (scf_tol_fixed)
227 r = this%final_tol
228
229 case (scf_tol_adaptive)
230 if (iter == 0) then
231 r = this%initial_tol
232 else
233 r = this%dynamic_tol_factor * (this%final_tol/this%conv_threshold_use)*scf_res
234 end if
235
236 case (scf_tol_linear)
237 r = this%initial_tol + (this%final_tol - this%initial_tol) * &
238 real(iter, real64) / real(this%iter_window, real64)
239
240 case (scf_tol_exp)
241 logi = log(this%initial_tol)
242 logf = log(this%final_tol)
243 r = logi + (logf - logi) * &
244 real(iter, real64) / real(this%iter_window, real64)
245 r = exp(r)
246 end select
247
248 ! tolerance can never be larger than final tolerance
249 r = max(r, this%final_tol)
250 ! tolerance always has to decrease
251 r = min(r, this%current_tol)
252
253 this%current_tol = r
254
255 pop_sub(scf_tol_step)
256 end function scf_tol_step
257
258
259 !-----------------------------------------------------------------
260 subroutine scf_tol_stop(this)
261 type(scf_tol_t), intent(inout) :: this
262
263 push_sub(scf_tol_stop)
264 this%current_tol = m_zero
265
266 pop_sub(scf_tol_stop)
267 end subroutine scf_tol_stop
268
269
270 !-----------------------------------------------------------------
271 subroutine scf_tol_end(this)
272 type(scf_tol_t), intent(inout) :: this
273
274 push_sub(scf_tol_end)
275 this%current_tol = m_zero
276
277 pop_sub(scf_tol_end)
278 end subroutine scf_tol_end
279
280
281 !-----------------------------------------------------------------
282 real(real64) function scf_tol_final(this)
283 type(scf_tol_t), intent(in) :: this
284
285 scf_tol_final = this%final_tol
286
287 end function scf_tol_final
288
289 !-----------------------------------------------------------------
290
291 subroutine scf_tol_obsolete_variables(namespace, old_prefix, new_prefix)
292 type(namespace_t), intent(in) :: namespace
293 character(len=*), intent(in) :: old_prefix
294 character(len=*), intent(in) :: new_prefix
295
297
298 call messages_obsolete_variable(namespace, trim(old_prefix)//'LRMaximumIter', trim(new_prefix)//'LRMaximumIter')
299 call messages_obsolete_variable(namespace, trim(old_prefix)//'LRConvAbsDens', trim(new_prefix)//'LRConvAbsDens')
300 call messages_obsolete_variable(namespace, trim(old_prefix)//'LRTolScheme', trim(new_prefix)//'LRTolScheme')
301 call messages_obsolete_variable(namespace, trim(old_prefix)//'LRTolInitTol', trim(new_prefix)//'LRTolInitTol')
302 call messages_obsolete_variable(namespace, trim(old_prefix)//'LRTolFinalTol', trim(new_prefix)//'LRTolFinalTol')
303 call messages_obsolete_variable(namespace, trim(old_prefix)//'LRTolAdaptiveFactor', trim(new_prefix)//'LRTolAdaptiveFactor')
304 call messages_obsolete_variable(namespace, trim(old_prefix)//'LRTolIterWindow', trim(new_prefix)//'LRTolIterWindow')
305
308
309end module scf_tol_oct_m
310
311!! Local Variables:
312!! mode: f90
313!! coding: utf-8
314!! End:
double log(double __x) __attribute__((__nothrow__
double exp(double __x) __attribute__((__nothrow__
real(real64), parameter, public m_huge
Definition: global.F90:205
real(real64), parameter, public m_zero
Definition: global.F90:187
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:723
real(real64) function, public scf_tol_step(this, iter, scf_res)
Definition: scf_tol.F90:308
integer, parameter, public scf_tol_adaptive
Definition: scf_tol.F90:127
subroutine, public scf_tol_init(this, namespace, qtot, def_maximumiter, tol_scheme)
Definition: scf_tol.F90:160
subroutine, public scf_tol_end(this)
Definition: scf_tol.F90:365
real(real64) function, public scf_tol_final(this)
Definition: scf_tol.F90:376
integer, parameter, public scf_tol_exp
Definition: scf_tol.F90:127
integer, parameter, public scf_tol_linear
Definition: scf_tol.F90:127
subroutine, public scf_tol_obsolete_variables(namespace, old_prefix, new_prefix)
Definition: scf_tol.F90:385
subroutine, public scf_tol_stop(this)
Definition: scf_tol.F90:354