Octopus
root_solver.F90
Go to the documentation of this file.
1!! Copyright (C) 2005-2006 Heiko Appel
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
22 use debug_oct_m
23 use global_oct_m
24 use, intrinsic :: iso_fortran_env
28 use parser_oct_m
31
32 implicit none
33
34 private
35 public :: &
40
41 integer, public, parameter :: &
42 ROOT_NEWTON = 3
43
44 type root_solver_t
45 private
46 integer :: solver_type
47 integer :: dim
48 integer :: maxiter
49 integer :: usediter
50 real(real64) :: abs_tolerance
51 real(real64) :: rel_tolerance
52 end type root_solver_t
53
54contains
55
56 ! ---------------------------------------------------------
57 subroutine root_solver_init(rs, namespace, dimensionality, solver_type, maxiter, &
58 rel_tolerance, abs_tolerance)
59 type(root_solver_t), intent(out) :: rs
60 type(namespace_t), intent(in) :: namespace
61 integer, intent(in) :: dimensionality
62 integer, optional, intent(in) :: solver_type, maxiter
63 real(real64), optional, intent(in) :: rel_tolerance, abs_tolerance
64
65 ! no push_sub, called too often
66
67 ! Set dimension
68 rs%dim = dimensionality
69
70 ! Read config parameters for the solver
71 call root_solver_read(rs, namespace)
72
73 if (present(solver_type)) rs%solver_type = solver_type
74 if (present(maxiter)) rs%maxiter = maxiter
75 if (present(rel_tolerance)) rs%rel_tolerance = rel_tolerance
76 if (present(abs_tolerance)) rs%abs_tolerance = abs_tolerance
77
78 end subroutine root_solver_init
79
80
81 ! ---------------------------------------------------------
82 subroutine root_solver_read(rs, namespace)
83 type(root_solver_t), intent(inout) :: rs
84 type(namespace_t), intent(in) :: namespace
85
86 push_sub(root_solver_read)
87
88 !%Variable RootSolver
89 !%Type integer
90 !%Default root_newton
91 !%Section Math::RootSolver
92 !%Description
93 !% Specifies what kind of root solver will be used.
94 !%Option root_newton 3
95 !% Newton method.
96 !%End
97 call parse_variable(namespace, 'RootSolver', root_newton, rs%solver_type)
98 if (.not. varinfo_valid_option('RootSolver', rs%solver_type)) then
99 call messages_input_error(namespace, 'RootSolver')
100 end if
101
102 !%Variable RootSolverMaxIter
103 !%Type integer
104 !%Default 500
105 !%Section Math::RootSolver
106 !%Description
107 !% In case of an iterative root solver, this variable determines the maximum number
108 !% of iteration steps.
109 !%End
110 call parse_variable(namespace, 'RootSolverMaxIter', 500, rs%maxiter)
111
112 !%Variable RootSolverRelTolerance
113 !%Type float
114 !%Default 1e-10
115 !%Section Math::RootSolver
116 !%Description
117 !% Relative tolerance for the root-finding process.
118 !%End
119 call parse_variable(namespace, 'RootSolverRelTolerance', 1e-10_real64, rs%rel_tolerance)
120
121 !%Variable RootSolverAbsTolerance
122 !%Type float
123 !%Default 1e-10
124 !%Section Math::RootSolver
125 !%Description
126 !% Relative tolerance for the root-finding process.
127 !%End
128 call parse_variable(namespace, 'RootSolverAbsTolerance', 1e-10_real64, rs%abs_tolerance)
129
130 pop_sub(root_solver_read)
131 end subroutine root_solver_read
132
133 ! ---------------------------------------------------------
134 subroutine droot_solver_run(rs, func, root, success, startval)
135 type(root_solver_t), intent(in) :: rs
136 real(real64), intent(out) :: root(:)
137 logical, intent(out) :: success
138 real(real64), optional, intent(in) :: startval(:)
139 interface
140 subroutine func(z, f, jf)
141 import real64
142 implicit none
143 real(real64), intent(in) :: z(:)
144 real(real64), intent(out) :: f(:), jf(:, :)
145 end subroutine func
146 end interface
147
148 ! no push_sub, called too often
149
150 ! Initializations
151 root = m_zero
152 success = .false.
153
154 select case (rs%solver_type)
155 case (root_newton)
156 call droot_newton(rs, func, root, startval, success)
157 case default
158 write(message(1), '(a,i4,a)') "Error in droot_solver_run: '", rs%solver_type, "' is not a valid root solver"
159 call messages_fatal(1)
160 end select
161
162 end subroutine
163
164 ! ---------------------------------------------------------
168 subroutine droot_newton(rs, func, root, startval, success)
169 type(root_solver_t), intent(in) :: rs
170 real(real64), intent(out) :: root(:)
171 real(real64), optional, intent(in) :: startval(:)
172 logical, intent(out) :: success
173 interface
174 subroutine func(z, f, jf)
175 import real64
176 implicit none
177 real(real64), intent(in) :: z(:)
178 real(real64), intent(out) :: f(:), jf(:, :)
179 end subroutine func
180 end interface
181
182 integer :: iter
183 real(real64) :: err, prev_err
184 real(real64), allocatable :: f(:), jf(:, :), delta(:, :), rhs(:, :), root_trial(:)
185
186 real(real64) :: mm, prev_m, alpha, err_trial
187
188 ! no push_sub, called too often
189
190 safe_allocate(f(1:rs%dim))
191 safe_allocate(root_trial(1:rs%dim))
192 safe_allocate(jf(1:rs%dim, 1:rs%dim))
193 safe_allocate(delta(1:rs%dim, 1))
194 safe_allocate(rhs(1:rs%dim, 1))
195
196 if (present(startval)) then
197 root = startval
198 else
199 root = m_zero
200 end if
201
202 call func(root, f, jf)
203 err = norm2(f)
204 prev_err = err
205 mm = m_one
206
207 success = .true.
208 iter = 0
209 do while(err > rs%abs_tolerance)
210 rhs(1:rs%dim, 1) = -mm * f(1:rs%dim)
211
212 call lalg_linsyssolve(rs%dim, 1, jf, rhs, delta)
213
214 ! Backtracking damping
215 alpha = 1.0_real64
216 do
217 root_trial = root + alpha * delta(:,1)
218 call func(root_trial, f, jf)
219 err_trial = norm2(f)
220
221 if (err_trial < err) exit ! accept
222 alpha = 0.5_real64 * alpha
223 if (alpha < 1e-12_real64) then
224 success = .false.
225 exit
226 end if
227 end do
228 if (.not. success) exit ! Damping failed
229
230 root = root_trial
231 prev_err = err
232 err = err_trial
233
234 iter = iter + 1
235 if (iter > rs%maxiter) then
236 success = .false.
237 exit
238 end if
239
240 ! Let us assume we have a linear convergence due to a multiplicity m of the solution,
241 ! Then we have |err_i+1| = (m-1)/m |err_i|
242 ! Or (m-1)/m = |err_i+1|/|err_i| = A, where A is the asymptotic error constant
243 ! for i sufficiently large.
244 !
245 ! Then we can estimate m from the first two iterations and use it to accelerate the convergence
246 !
247 ! See https://web.archive.org/web/20190524083302/http://mathfaculty.fullerton.edu/mathews/n2003/NewtonAccelerateMod.html
248 ! as well as
249 ! https://en.wikipedia.org/wiki/Newton%27s_method#Slow_convergence_for_roots_of_multiplicity_greater_than_1
250 ! for some more details
251 if (iter == 1) then
252 prev_m = floor(m_one / (m_one - err/prev_err))
253 end if
254 if (iter == 2) then
255 mm = floor(m_one / (m_one - err/prev_err))
256 ! If we get the same multiplicity, we use it for the next iterations
257 if (abs(mm - prev_m) > 0.5_real64) mm = m_one
258 ! If the multiplicity is going to infinity, we do not want to use this value
259 if ( mm > 10._real64 .or. mm < 1.0_real64) mm = m_one
260 end if
262 end do
263
264 safe_deallocate_a(f)
265 safe_deallocate_a(jf)
266 safe_deallocate_a(delta)
267 safe_deallocate_a(rhs)
268
269 end subroutine droot_newton
270
271end module root_solver_oct_m
272
273!! Local Variables:
274!! mode: f90
275!! coding: utf-8
276!! End:
double floor(double __x) __attribute__((__nothrow__
real(real64) function func(r1, rn, n, a)
Definition: logrid.F90:356
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_one
Definition: global.F90:189
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:714
subroutine droot_newton(rs, func, root, startval, success)
Newton-Raphson scheme can only be used in the real case.
subroutine, public root_solver_init(rs, namespace, dimensionality, solver_type, maxiter, rel_tolerance, abs_tolerance)
subroutine, public droot_solver_run(rs, func, root, success, startval)
subroutine, public root_solver_read(rs, namespace)
static double f(double w, void *p)
int true(void)