42 SolverBase() : l2_0(1.0), l2_n(1.0), n(0), rho(0.0), n10(0)
44 timer_tree = std::make_shared<TimerTree>();
55 update_preconditioner(
bool const update_preconditioner)
const = 0;
57 template<
typename Control>
59 compute_performance_metrics(Control
const & solver_control)
const
62 this->l2_0 = solver_control.initial_value();
63 this->l2_n = solver_control.last_value();
64 this->n = solver_control.last_step();
69 this->rho = std::pow(l2_n / l2_0, 1.0 / n);
70 this->n10 = -10.0 * std::log(10.0) / std::log(rho);
74 virtual std::shared_ptr<TimerTree>
83 mutable unsigned int n;
88 std::shared_ptr<TimerTree> timer_tree;
95 solver_tolerance_abs(1.e-20),
96 solver_tolerance_rel(1.e-6),
97 use_preconditioner(
false),
98 compute_performance_metrics(
false)
102 unsigned int max_iter;
103 double solver_tolerance_abs;
104 double solver_tolerance_rel;
105 bool use_preconditioner;
106 bool compute_performance_metrics;
110class SolverCG :
public SolverBase<VectorType>
113 SolverCG(Operator
const & underlying_operator_in,
114 Preconditioner & preconditioner_in,
116 : underlying_operator(underlying_operator_in),
117 preconditioner(preconditioner_in),
118 solver_data(solver_data_in)
123 update_preconditioner(
bool const update_preconditioner)
const override
125 if(solver_data.use_preconditioner)
127 if(preconditioner.needs_update() or update_preconditioner)
129 preconditioner.update();
139 dealii::ReductionControl solver_control(solver_data.max_iter,
140 solver_data.solver_tolerance_abs,
141 solver_data.solver_tolerance_rel);
143 dealii::SolverCG<VectorType> solver(solver_control);
145 if(solver_data.use_preconditioner ==
false)
147 solver.solve(underlying_operator, dst, rhs, dealii::PreconditionIdentity());
151 solver.solve(underlying_operator, dst, rhs, preconditioner);
154 AssertThrow(std::isfinite(solver_control.last_value()),
155 dealii::ExcMessage(
"Last iteration step contained NaN or Inf values."));
157 if(solver_data.compute_performance_metrics)
158 this->compute_performance_metrics(solver_control);
160 this->timer_tree->insert({
"SolverCG"}, timer.wall_time());
162 return solver_control.last_step();
165 std::shared_ptr<TimerTree>
166 get_timings()
const override
168 if(solver_data.use_preconditioner)
170 this->timer_tree->insert({
"SolverCG"}, preconditioner.get_timings());
173 return this->timer_tree;
177 Operator
const & underlying_operator;
178 Preconditioner & preconditioner;
199struct SolverDataGMRES
203 solver_tolerance_abs(1.e-20),
204 solver_tolerance_rel(1.e-6),
205 use_preconditioner(
false),
206 max_n_tmp_vectors(30),
207 compute_eigenvalues(
false),
208 compute_performance_metrics(
false)
212 unsigned int max_iter;
213 double solver_tolerance_abs;
214 double solver_tolerance_rel;
215 bool use_preconditioner;
216 unsigned int max_n_tmp_vectors;
217 bool compute_eigenvalues;
218 bool compute_performance_metrics;
222class SolverGMRES :
public SolverBase<VectorType>
225 SolverGMRES(Operator
const & underlying_operator_in,
226 Preconditioner & preconditioner_in,
228 MPI_Comm
const & mpi_comm_in)
229 : underlying_operator(underlying_operator_in),
230 preconditioner(preconditioner_in),
231 solver_data(solver_data_in),
232 mpi_comm(mpi_comm_in)
236 virtual ~SolverGMRES()
241 update_preconditioner(
bool const update_preconditioner)
const override
243 if(solver_data.use_preconditioner)
245 if(preconditioner.needs_update() or update_preconditioner)
247 preconditioner.update();
257 dealii::ReductionControl solver_control(solver_data.max_iter,
258 solver_data.solver_tolerance_abs,
259 solver_data.solver_tolerance_rel);
261 typename dealii::SolverGMRES<VectorType>::AdditionalData additional_data;
262 additional_data.max_n_tmp_vectors = solver_data.max_n_tmp_vectors;
263 additional_data.right_preconditioning =
true;
264 dealii::SolverGMRES<VectorType> solver(solver_control, additional_data);
266 if(solver_data.compute_eigenvalues ==
true)
268 solver.connect_eigenvalues_slot(std::bind(output_eigenvalues<std::complex<double>>,
269 std::placeholders::_1,
275 if(solver_data.use_preconditioner ==
false)
277 solver.solve(underlying_operator, dst, rhs, dealii::PreconditionIdentity());
281 solver.solve(this->underlying_operator, dst, rhs, this->preconditioner);
284 AssertThrow(std::isfinite(solver_control.last_value()),
285 dealii::ExcMessage(
"Last iteration step contained NaN or Inf values."));
287 if(solver_data.compute_performance_metrics)
288 this->compute_performance_metrics(solver_control);
290 this->timer_tree->insert({
"SolverGMRES"}, timer.wall_time());
292 return solver_control.last_step();
295 std::shared_ptr<TimerTree>
296 get_timings()
const override
298 if(solver_data.use_preconditioner)
300 this->timer_tree->insert({
"SolverGMRES"}, preconditioner.get_timings());
303 return this->timer_tree;
307 Operator
const & underlying_operator;
308 Preconditioner & preconditioner;
311 MPI_Comm
const mpi_comm;
314struct SolverDataFGMRES
318 solver_tolerance_abs(1.e-20),
319 solver_tolerance_rel(1.e-6),
320 use_preconditioner(
false),
321 max_n_tmp_vectors(30),
322 compute_performance_metrics(
false)
326 unsigned int max_iter;
327 double solver_tolerance_abs;
328 double solver_tolerance_rel;
329 bool use_preconditioner;
330 unsigned int max_n_tmp_vectors;
331 bool compute_performance_metrics;
335class SolverFGMRES :
public SolverBase<VectorType>
338 SolverFGMRES(Operator
const & underlying_operator_in,
339 Preconditioner & preconditioner_in,
341 : underlying_operator(underlying_operator_in),
342 preconditioner(preconditioner_in),
343 solver_data(solver_data_in)
347 virtual ~SolverFGMRES()
352 update_preconditioner(
bool const update_preconditioner)
const override
354 if(solver_data.use_preconditioner)
356 if(preconditioner.needs_update() or update_preconditioner)
358 preconditioner.update();
368 dealii::ReductionControl solver_control(solver_data.max_iter,
369 solver_data.solver_tolerance_abs,
370 solver_data.solver_tolerance_rel);
372 typename dealii::SolverFGMRES<VectorType>::AdditionalData additional_data;
373 additional_data.max_basis_size = solver_data.max_n_tmp_vectors;
376 dealii::SolverFGMRES<VectorType> solver(solver_control, additional_data);
378 if(solver_data.use_preconditioner ==
false)
380 solver.solve(underlying_operator, dst, rhs, dealii::PreconditionIdentity());
384 solver.solve(underlying_operator, dst, rhs, preconditioner);
387 AssertThrow(std::isfinite(solver_control.last_value()),
388 dealii::ExcMessage(
"Last iteration step contained NaN or Inf values."));
390 if(solver_data.compute_performance_metrics)
391 this->compute_performance_metrics(solver_control);
393 this->timer_tree->insert({
"SolverFGMRES"}, timer.wall_time());
395 return solver_control.last_step();
398 std::shared_ptr<TimerTree>
399 get_timings()
const override
401 if(solver_data.use_preconditioner)
403 this->timer_tree->insert({
"SolverFGMRES"}, preconditioner.get_timings());
406 return this->timer_tree;
410 Operator
const & underlying_operator;
411 Preconditioner & preconditioner;