113 lincomb( r, { 1.0, -1.0 }, { b, r } );
120 statistics_->add_row(
121 { {
"tag", tag_ }, {
"iteration", 0 }, {
"relative_residual", 1.0 }, {
"absolute_residual", beta0 } } );
124 if ( beta0 <= options_.absolute_residual_tolerance )
131 const int nTmp =
static_cast< int >( tmp_.size() );
132 const int m_req = options_.restart;
133 const int needed = 2 * m_req + 4;
136 std::cerr <<
"FGMRES: insufficient tmp vectors. Provided " << nTmp
137 <<
", required at least (2*restart + 4) = " << needed <<
" for restart m = " << m_req
139 Kokkos::abort(
"FGMRES: insufficient tmp vectors" );
144 const int offZ = offV + ( m_req + 1 );
145 const int idxW = offZ + m_req;
146 const int idxAcc = idxW + 1;
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 );
156 while ( total_iters < options_.max_iterations )
159 auto& V0 = tmp_[offV + 0];
160 lincomb( V0, { 1.0 / beta0 }, { r } );
165 std::fill( cs.begin(), cs.end(),
ScalarType( 0 ) );
166 std::fill( sn.begin(), sn.end(),
ScalarType( 0 ) );
172 for (
int j = 0; j < m_req && total_iters < options_.max_iterations; ++j, ++total_iters )
175 auto& vj = tmp_[offV + j];
176 auto& zj = tmp_[offZ + j];
178 solve( preconditioner_, A, zj, vj );
180 const bool preconditioner_result_contains_nan_or_inf =
has_nan_or_inf( zj );
182 if ( skip_preconditioner_in_case_of_nan_or_infs_ && preconditioner_result_contains_nan_or_inf )
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;
196 auto& w = tmp_[idxW];
200 for (
int i = 0; i <= j; ++i )
202 auto& vi = tmp_[offV + i];
205 lincomb( w, { 1.0, -hij }, { w, vi } );
212 if ( h_jp1j < std::numeric_limits< ScalarType >::epsilon() * initial_residual )
218 " (Details: total_iters = "
219 << total_iters <<
", j = " << j <<
")" << std::endl;
224 H( j + 1, j ) = h_jp1j;
228 auto& vjp1 = tmp_[offV + j + 1];
229 lincomb( vjp1, { 1.0 / h_jp1j }, { w } );
233 for (
int i = 0; i < j; ++i )
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 );
251 if ( std::abs( H( j, j ) ) < std::numeric_limits< ScalarType >::epsilon() )
255 util::logroot <<
"FGMRES: Stagnation after Givens rotation. Restarting.\n"
256 " (Details: total_iters = "
257 << total_iters <<
", j = " << j <<
")" << std::endl;
265 g( j ) = cs[j] * gj + sn[j] * gj1;
266 g( j + 1 ) = -sn[j] * gj + cs[j] * gj1;
270 const ScalarType abs_res = std::abs( g( j + 1 ) );
271 const ScalarType rel_res = abs_res / initial_residual;
275 statistics_->add_row(
277 {
"iteration", total_iters + 1 },
278 {
"relative_residual", rel_res },
279 {
"absolute_residual", abs_res } } );
285 if ( rel_res <= options_.relative_residual_tolerance || abs_res < options_.absolute_residual_tolerance )
292 Eigen::Matrix< ScalarType, Eigen::Dynamic, 1 > y( inner_its );
294 for (
int row = inner_its - 1; row >= 0; --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 );
303 auto& acc = tmp_[idxAcc];
305 for (
int i = 0; i < inner_its; ++i )
307 auto& zi = tmp_[offZ + i];
308 lincomb( acc, { 1.0, y( i ) }, { acc, zi } );
310 lincomb( x, { 1.0, 1.0 }, { x, acc } );
314 lincomb( r, { 1.0, -1.0 }, { b, r } );
315 beta0 = std::sqrt(
dot( r, r ) );
318 if ( beta0 <= options_.absolute_residual_tolerance ||
319 beta0 / initial_residual <= options_.relative_residual_tolerance )