freegsnke.nk_solver_H module
Implements the core Newton Krylov nonlinear solver used by both static GS solver and evolutive solver.
Copyright 2025 UKAEA, UKRI-STFC, and The Authors, as per the COPYRIGHT and README files.
This file is part of FreeGSNKE.
FreeGSNKE is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
FreeGSNKE is free software: you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.
You should have received a copy of the GNU Lesser General Public License along with FreeGSNKE. If not, see <http://www.gnu.org/licenses/>.
- class freegsnke.nk_solver_H.nksolver(problem_dimension, verbose=False)[source]
Bases:
object
Implementation of Newton Krylow algorithm for solving a generic root problem of the type F(x, other args) = 0 in the variable x – F(x) should have the same dimensions as x. Problem must be formulated so that x is a 1d np.array.
In practice, given a guess x_0 and F(x_0) = R_0 it aims to find the best step dx such that F(x_0 + dx) is minimum.
- Arnoldi_iteration(x0, dx, R0, F_function, args, step_size, scaling_with_n, target_relative_unexplained_residual, max_n_directions, clip)[source]
Performs the iteration of the NK solution method: 1) explores direction dx 2) checks what fraction of the residual can be (linearly) canceled 3) restarts if not satisfied The best candidate step combining all explored directions is stored at self.dx
- Parameters:
x0 (1d np.array, np.shape(x0) = self.problem_dimension) – The expansion point x_0
dx (1d np.array, np.shape(dx) = self.problem_dimension) – The first direction to be explored. This will be sized appropriately.
R0 (1d np.array, np.shape(R0) = self.problem_dimension) – Residual of the root problem F_function at expansion point x_0
F_function (1d np.array, np.shape(x0) = self.problem_dimension) – Function representing the root problem at hand
args (list) – Additional arguments for using function F F = F(x, *args)
step_size (float) – l2 norm of proposed step in units of the residual norm
scaling_with_n (float) – allows to further scale dx candidate steps as a function of the iteration number n_it, by a factor (1 + self.n_it)**scaling_with_n
target_relative_explained_residual (float between 0 and 1) – terminates iteration when such a fraction of the initial residual R0 can be (linearly) cancelled
max_n_directions (int) – terminates iteration even though condition on explained residual is not met
clip (float) – maximum step size for each explored direction, in units of exploratory step dx_i
- Arnoldi_unit(x0, dx, R0, F_function, args)[source]
Explores direction dx and proposes new direction for next exploration.
- Parameters:
x0 (1d np.array, np.shape(x0) = self.problem_dimension) – The expansion point x_0
dx (1d np.array, np.shape(dx) = self.problem_dimension) – The first direction to be explored. This will be sized appropriately.
R0 (1d np.array, np.shape(R0) = self.problem_dimension) – Residual of the root problem F_function at expansion point x_0
F_function (1d np.array, np.shape(x0) = self.problem_dimension) – Function representing the root problem at hand
args (list) – Additional arguments for using function F F = F(x, *args)
- Returns:
new_candidate_step – The direction to be explored next
- Return type:
1d np.array, with same self.problem_dimension