Loading...
Searching...
No Matches
fgmres.hpp
Go to the documentation of this file.
1#pragma once
2
4#include "identity_solver.hpp"
9#include "util/table.hpp"
10
11namespace terra::linalg::solvers {
12
13template < typename ScalarType >
15{
16 /// @brief Number of inner iterations before restart (FGMRES(m)).
17 int restart = 30;
18
19 /// @brief Relative residual tolerance for convergence.
20 ScalarType relative_residual_tolerance = 1e-6;
21
22 /// @brief Absolute residual tolerance for convergence.
23 ScalarType absolute_residual_tolerance = 1e-6;
24
25 /// @brief Maximum number of total iterations (across all restarts).
26 int max_iterations = 100;
27};
28
29/// @brief Flexible GMRES (FGMRES) iterative solver for nonsymmetric linear systems.
30///
31/// FGMRES allows for varying (flexible) preconditioning at each iteration,
32/// making it suitable for use with inexact or variable preconditioners.
33///
34/// Reference:
35/// @code
36/// Saad, Y. (1993).
37/// A flexible inner-outer preconditioned GMRES algorithm.
38/// SIAM Journal on Scientific Computing, 14(2), 461-469.
39/// @endcode
40///
41/// Satisfies the SolverLike concept (see solver.hpp).
42/// Supports optional right preconditioning via the flexible framework.
43///
44/// @tparam OperatorT Operator type (must satisfy OperatorLike).
45/// @tparam PreconditionerT Preconditioner type (must satisfy SolverLike, defaults to IdentitySolver).
46template < OperatorLike OperatorT, SolverLike PreconditionerT = IdentitySolver< OperatorT > >
47class FGMRES
48{
49 public:
50 /// @brief Operator type to be solved.
51 using OperatorType = OperatorT;
52
53 /// @brief Solution vector type.
55
56 /// @brief Right-hand side vector type.
58
59 /// @brief Scalar type for computations.
60 using ScalarType = typename SolutionVectorType::ScalarType;
61
62 /// @brief Construct an FGMRES solver with a custom preconditioner.
63 ///
64 /// @param tmp Temporary vectors for workspace. Must contain at least 2*restart + 4 vectors:
65 /// - tmp[0]: residual vector r
66 /// - tmp[1 .. restart+1]: Arnoldi basis vectors V_0..V_restart
67 /// - tmp[restart+2 .. 2*restart+1]: preconditioned directions Z_0..Z_{restart-1}
68 /// - tmp[2*restart+2]: work vector w for A*z
69 /// - tmp[2*restart+3]: accumulator vector for solution update
70 /// @param options FGMRES solver parameters (restart, tolerances, max iterations).
71 /// @param statistics Shared pointer to statistics table for logging iteration progress (optional).
72 /// @param preconditioner Preconditioner solver (defaults to identity).
74 const std::vector< SolutionVectorType >& tmp,
75 const FGMRESOptions< ScalarType >& options = {},
76 const std::shared_ptr< util::Table >& statistics = nullptr,
77 const PreconditionerT preconditioner = IdentitySolver< OperatorT >() )
78 : tag_( "fgmres_solver" )
79 , tmp_( tmp )
80 , options_( options )
81 , statistics_( statistics )
82 , preconditioner_( preconditioner )
83 , skip_preconditioner_in_case_of_nan_or_infs_( true )
84 {}
85
86 /// @brief Set a tag string for statistics output identification.
87 /// @param tag Tag string to identify this solver instance in statistics.
88 void set_tag( const std::string& tag ) { tag_ = tag; }
89
90 /// @brief Set the number of inner iterations before restart (FGMRES(m)).
91 /// @param m Restart parameter (must be at least 1).
92 void set_restart( int m ) { options_.restart = std::max( 1, m ); }
93
94 /// @brief Solve the linear system \f$ Ax = b \f$ using flexible GMRES with restarts.
95 ///
96 /// Uses right preconditioning: solves \f$ A M^{-1} y = b \f$ where \f$ x = M^{-1} y \f$.
97 /// The preconditioner M can vary at each iteration (flexibility).
98 ///
99 /// The method builds an Arnoldi basis via modified Gram-Schmidt orthogonalization,
100 /// applies Givens rotations to minimize the least-squares residual in the Krylov subspace,
101 /// and updates the solution. If convergence is not achieved within the restart window,
102 /// the process restarts with the updated solution until max_iterations is reached or
103 /// convergence criteria are satisfied.
104 ///
105 /// @param A Operator (matrix) to solve with.
106 /// @param x Solution vector (input: initial guess, output: final solution).
107 /// @param b Right-hand side vector (input).
109 {
110 // Compute initial residual r = b - A*x
111 auto& r = tmp_[0];
112 apply( A, x, r );
113 lincomb( r, { 1.0, -1.0 }, { b, r } );
114
115 ScalarType beta0 = std::sqrt( dot( r, r ) );
116 const ScalarType initial_residual = beta0;
117
118 if ( statistics_ )
119 {
120 statistics_->add_row(
121 { { "tag", tag_ }, { "iteration", 0 }, { "relative_residual", 1.0 }, { "absolute_residual", beta0 } } );
122 }
123
124 if ( beta0 <= options_.absolute_residual_tolerance )
125 {
126 return;
127 }
128
129 // Validate workspace size.
130 // Required: 1 (r) + (m+1) (V) + m (Z) + 1 (w) + 1 (acc) = 2m + 4
131 const int nTmp = static_cast< int >( tmp_.size() );
132 const int m_req = options_.restart;
133 const int needed = 2 * m_req + 4;
134 if ( nTmp < needed )
135 {
136 std::cerr << "FGMRES: insufficient tmp vectors. Provided " << nTmp
137 << ", required at least (2*restart + 4) = " << needed << " for restart m = " << m_req
138 << std::endl;
139 Kokkos::abort( "FGMRES: insufficient tmp vectors" );
140 }
141
142 // Workspace layout offsets:
143 const int offV = 1; // V_0 .. V_m stored in tmp[1 .. m+1]
144 const int offZ = offV + ( m_req + 1 ); // Z_0 .. Z_{m-1} stored in tmp[offV + m+1 .. offV + 2m]
145 const int idxW = offZ + m_req; // w = A*z_j
146 const int idxAcc = idxW + 1; // accumulator for solution update
147
148 // Allocate small dense arrays for Hessenberg matrix and Givens rotations.
149 Eigen::Matrix< ScalarType, Eigen::Dynamic, Eigen::Dynamic > H( m_req + 1, m_req );
150 Eigen::Matrix< ScalarType, Eigen::Dynamic, 1 > g( m_req + 1 );
151 std::vector< ScalarType > cs( m_req + 1, 0 ), sn( m_req + 1, 0 );
152
153 int total_iters = 0;
154
155 // Outer restart loop
156 while ( total_iters < options_.max_iterations )
157 {
158 // Initialize Arnoldi: V_0 = r / ||r||
159 auto& V0 = tmp_[offV + 0];
160 lincomb( V0, { 1.0 / beta0 }, { r } );
161
162 // Reset dense workspace for this restart cycle
163 H.setZero();
164 g.setZero();
165 std::fill( cs.begin(), cs.end(), ScalarType( 0 ) );
166 std::fill( sn.begin(), sn.end(), ScalarType( 0 ) );
167 g( 0 ) = beta0;
168
169 int inner_its = 0;
170
171 // Inner Arnoldi iteration (build Krylov subspace up to dimension m_req)
172 for ( int j = 0; j < m_req && total_iters < options_.max_iterations; ++j, ++total_iters )
173 {
174 // Apply preconditioner: z_j = M^{-1} v_j
175 auto& vj = tmp_[offV + j];
176 auto& zj = tmp_[offZ + j];
177 assign( zj, 0 );
178 solve( preconditioner_, A, zj, vj );
179
180 const bool preconditioner_result_contains_nan_or_inf = has_nan_or_inf( zj );
181
182 if ( skip_preconditioner_in_case_of_nan_or_infs_ && preconditioner_result_contains_nan_or_inf )
183 {
185 << "FGMRES: The preconditioner appears to have produced a vector that has NaN or Inf entries.\n"
186 " This may be a result of the preconditioner or of the input provided by FGMRES.\n"
187 " To at least provide a fix for the first case, we keep the input to the preconditioner\n"
188 " and write it to the output (equivalent of skipping the preconditioner in this iteration).\n"
189 " (Details: total_iters = "
190 << total_iters << ", j = " << j << ")" << std::endl;
191
192 assign( zj, vj );
193 }
194
195 // Apply operator: w = A * z_j
196 auto& w = tmp_[idxW];
197 apply( A, zj, w );
198
199 // Modified Gram-Schmidt orthogonalization against V_0..V_j
200 for ( int i = 0; i <= j; ++i )
201 {
202 auto& vi = tmp_[offV + i];
203 const ScalarType hij = dot( w, vi );
204 H( i, j ) = hij;
205 lincomb( w, { 1.0, -hij }, { w, vi } );
206 }
207
208 // Compute h_{j+1,j} and normalize to get v_{j+1}
209 const ScalarType h_jp1j = std::sqrt( dot( w, w ) );
210
211 // SAFE BREAK 1: Arnoldi Breakdown
212 if ( h_jp1j < std::numeric_limits< ScalarType >::epsilon() * initial_residual )
213 {
214 H( j + 1, j ) = 0;
215 inner_its = j + 1;
216
217 util::logroot << "FGMRES: Arnoldi breakdown. Restarting.\n"
218 " (Details: total_iters = "
219 << total_iters << ", j = " << j << ")" << std::endl;
220
221 break; // Subspace is "full" or converged; exit to solve current least squares
222 }
223
224 H( j + 1, j ) = h_jp1j;
225
226 if ( h_jp1j > ScalarType( 0 ) )
227 {
228 auto& vjp1 = tmp_[offV + j + 1];
229 lincomb( vjp1, { 1.0 / h_jp1j }, { w } );
230 }
231
232 // Apply previous Givens rotations to column j of H
233 for ( int i = 0; i < j; ++i )
234 {
235 const ScalarType temp = cs[i] * H( i, j ) + sn[i] * H( i + 1, j );
236 H( i + 1, j ) = -sn[i] * H( i, j ) + cs[i] * H( i + 1, j );
237 H( i, j ) = temp;
238 }
239
240 // Compute new Givens rotation to zero H(j+1, j)
241 {
242 const ScalarType a = H( j, j );
243 const ScalarType b = H( j + 1, j );
244 const ScalarType r_ab = std::hypot( a, b );
245 cs[j] = ( r_ab == ScalarType( 0 ) ) ? ScalarType( 1 ) : a / r_ab;
246 sn[j] = ( r_ab == ScalarType( 0 ) ) ? ScalarType( 0 ) : b / r_ab;
247 H( j, j ) = r_ab;
248 H( j + 1, j ) = 0;
249
250 // SAFE BREAK 2: Stagnation / Singularity
251 if ( std::abs( H( j, j ) ) < std::numeric_limits< ScalarType >::epsilon() )
252 {
253 inner_its = j; // Don't include this column in the back-solve
254
255 util::logroot << "FGMRES: Stagnation after Givens rotation. Restarting.\n"
256 " (Details: total_iters = "
257 << total_iters << ", j = " << j << ")" << std::endl;
258
259 break;
260 }
261
262 // Apply new Givens rotation to g
263 const ScalarType gj = g( j );
264 const ScalarType gj1 = g( j + 1 );
265 g( j ) = cs[j] * gj + sn[j] * gj1;
266 g( j + 1 ) = -sn[j] * gj + cs[j] * gj1;
267 }
268
269 // Estimate residual norm from least-squares problem
270 const ScalarType abs_res = std::abs( g( j + 1 ) );
271 const ScalarType rel_res = abs_res / initial_residual;
272
273 if ( statistics_ )
274 {
275 statistics_->add_row(
276 { { "tag", tag_ },
277 { "iteration", total_iters + 1 },
278 { "relative_residual", rel_res },
279 { "absolute_residual", abs_res } } );
280 }
281
282 inner_its = j + 1;
283
284 // Check for convergence
285 if ( rel_res <= options_.relative_residual_tolerance || abs_res < options_.absolute_residual_tolerance )
286 {
287 break;
288 }
289 }
290
291 // Solve upper-triangular system R*y = g(0..inner_its-1) via back-substitution
292 Eigen::Matrix< ScalarType, Eigen::Dynamic, 1 > y( inner_its );
293 y.setZero();
294 for ( int row = inner_its - 1; row >= 0; --row )
295 {
296 ScalarType sum = g( row );
297 for ( int col = row + 1; col < inner_its; ++col )
298 sum -= H( row, col ) * y( col );
299 y( row ) = sum / H( row, row );
300 }
301
302 // Update solution: x += sum_{i=0}^{inner_its-1} y_i * z_i
303 auto& acc = tmp_[idxAcc];
304 assign( acc, 0 );
305 for ( int i = 0; i < inner_its; ++i )
306 {
307 auto& zi = tmp_[offZ + i];
308 lincomb( acc, { 1.0, y( i ) }, { acc, zi } );
309 }
310 lincomb( x, { 1.0, 1.0 }, { x, acc } );
311
312 // Recompute residual r = b - A*x for next restart
313 apply( A, x, r );
314 lincomb( r, { 1.0, -1.0 }, { b, r } );
315 beta0 = std::sqrt( dot( r, r ) );
316
317 // Check for final convergence
318 if ( beta0 <= options_.absolute_residual_tolerance ||
319 beta0 / initial_residual <= options_.relative_residual_tolerance )
320 {
321 return;
322 }
323 }
324 }
325
326 private:
327 std::string tag_; ///< Tag for statistics output identification.
328
329 std::vector< SolutionVectorType > tmp_; ///< Temporary workspace vectors.
330
331 FGMRESOptions< ScalarType > options_; ///< Solver parameters (restart, tolerances, max_iterations).
332
333 std::shared_ptr< util::Table > statistics_; ///< Statistics table for iteration logging (optional).
334
335 PreconditionerT preconditioner_; ///< Preconditioner solver.
336
337 bool skip_preconditioner_in_case_of_nan_or_infs_;
338};
339
340/// @brief Static assertion: FGMRES satisfies SolverLike concept.
341static_assert(
342 SolverLike<
343 FGMRES< linalg::detail::
344 DummyOperator< linalg::detail::DummyVector< double >, linalg::detail::DummyVector< double > > > > );
345
346} // namespace terra::linalg::solvers
Dummy vector class for concept checks and testing. Implements required vector operations as no-ops.
Definition vector.hpp:210
Flexible GMRES (FGMRES) iterative solver for nonsymmetric linear systems.
Definition fgmres.hpp:48
OperatorT OperatorType
Operator type to be solved.
Definition fgmres.hpp:51
SrcOf< OperatorType > SolutionVectorType
Solution vector type.
Definition fgmres.hpp:54
void set_restart(int m)
Set the number of inner iterations before restart (FGMRES(m)).
Definition fgmres.hpp:92
void set_tag(const std::string &tag)
Set a tag string for statistics output identification.
Definition fgmres.hpp:88
FGMRES(const std::vector< SolutionVectorType > &tmp, const FGMRESOptions< ScalarType > &options={}, const std::shared_ptr< util::Table > &statistics=nullptr, const PreconditionerT preconditioner=IdentitySolver< OperatorT >())
Construct an FGMRES solver with a custom preconditioner.
Definition fgmres.hpp:73
typename SolutionVectorType::ScalarType ScalarType
Scalar type for computations.
Definition fgmres.hpp:60
void solve_impl(OperatorType &A, SolutionVectorType &x, const RHSVectorType &b)
Solve the linear system using flexible GMRES with restarts.
Definition fgmres.hpp:108
DstOf< OperatorType > RHSVectorType
Right-hand side vector type.
Definition fgmres.hpp:57
"Identity solver" for linear systems.
Definition identity_solver.hpp:21
Definition block_preconditioner_2x2.hpp:7
void solve(Solver &solver, Operator &A, SolutionVector &x, const RHSVector &b)
Solve a linear system using the given solver and operator. Calls the solver's solve_impl method.
Definition solver.hpp:51
void lincomb(Vector &y, const std::vector< ScalarOf< Vector > > &c, const std::vector< Vector > &x, const ScalarOf< Vector > &c0)
Compute a linear combination of vectors. Implements: .
Definition vector.hpp:72
Operator::SrcVectorType SrcOf
Alias for the source vector type of an operator.
Definition operator.hpp:145
ScalarOf< Vector > dot(const Vector &y, const Vector &x)
Compute the dot product of two vectors. Implements: .
Definition vector.hpp:118
void apply(LinearForm &L, typename LinearForm::DstVectorType &dst)
Apply a linear form and write to a destination vector.
Definition linear_form.hpp:37
Operator::DstVectorType DstOf
Alias for the destination vector type of an operator.
Definition operator.hpp:149
bool has_nan_or_inf(const Vector &y)
Check if a vector contains NaN or inf entries. Returns true if any entry of is NaN or inf.
Definition vector.hpp:189
void assign(Vector &y, const ScalarOf< Vector > &c0)
Assign a scalar value to a vector. Implements: .
Definition vector.hpp:97
detail::PrefixCout logroot([]() { return detail::log_prefix();})
std::ostream subclass that just logs on root and adds a timestamp for each line.
ScalarType absolute_residual_tolerance
Absolute residual tolerance for convergence.
Definition fgmres.hpp:23
ScalarType relative_residual_tolerance
Relative residual tolerance for convergence.
Definition fgmres.hpp:20
int max_iterations
Maximum number of total iterations (across all restarts).
Definition fgmres.hpp:26
int restart
Number of inner iterations before restart (FGMRES(m)).
Definition fgmres.hpp:17