Smith  0.1
Smith is an implicit thermal structural mechanics simulation code.
equation_solver.hpp
Go to the documentation of this file.
1 // Copyright (c) Lawrence Livermore National Security, LLC and
2 // other Smith Project Developers. See the top-level LICENSE file for
3 // details.
4 //
5 // SPDX-License-Identifier: (BSD-3-Clause)
6 
13 #pragma once
14 
15 #include <memory>
16 #include <optional>
17 #include <variant>
18 #include <utility>
19 
20 #include "mpi.h"
21 #include "mfem.hpp"
22 
24 #include "smith/numerics/nonlinear_convergence.hpp"
26 #include "smith/numerics/petsc_solvers.hpp"
27 
28 namespace smith {
29 
47  public:
48  // _equationsolver_constructor_start
56  EquationSolver(std::unique_ptr<mfem::NewtonSolver> nonlinear_solver, std::unique_ptr<mfem::Solver> linear_solver,
57  std::unique_ptr<mfem::Solver> preconditioner = nullptr);
58  // _equationsolver_constructor_end
59 
60  // _build_equationsolver_start
68  EquationSolver(NonlinearSolverOptions nonlinear_opts = {}, LinearSolverOptions lin_opts = {},
69  MPI_Comm comm = MPI_COMM_WORLD);
70  // _build_equationsolver_end
71 
78  void setOperator(const mfem::Operator& op);
79 
85  void solve(mfem::Vector& x) const;
86 
88  void setConvergenceTolerances(double abs_tol, double rel_tol, MPI_Comm comm) const;
89 
91  void resetConvergenceState() const;
92 
97  mfem::NewtonSolver& nonlinearSolver() { return *nonlin_solver_; }
98 
102  const mfem::NewtonSolver& nonlinearSolver() const { return *nonlin_solver_; }
103 
108  mfem::Solver& linearSolver() { return *lin_solver_; }
109 
113  const mfem::Solver& linearSolver() const { return *lin_solver_; }
114 
120  mfem::Solver& preconditioner() { return *preconditioner_; }
121 
125  const mfem::Solver& preconditioner() const { return *preconditioner_; }
126 
130  static void defineInputFileSchema(axom::inlet::Container& container);
131 
132  private:
133  void attachConvergenceManager() const;
134  void initializeConvergenceManager(double abs_tol, double rel_tol, MPI_Comm comm) const;
135 
139  std::unique_ptr<mfem::Solver> preconditioner_;
140 
144  std::unique_ptr<mfem::Solver> lin_solver_;
145 
149  std::unique_ptr<mfem::NewtonSolver> nonlin_solver_;
150 
156  bool nonlin_solver_set_solver_called_ = false;
157 
158  mutable std::shared_ptr<EquationSolverConvergenceManager> convergence_manager_ = nullptr;
159 };
160 
164 class SuperLUSolver : public mfem::Solver {
165  public:
171  SuperLUSolver(int print_level, MPI_Comm comm) : superlu_solver_(comm)
172  {
173  superlu_solver_.SetColumnPermutation(mfem::superlu::PARMETIS);
174  if (print_level == 0) {
175  superlu_solver_.SetPrintStatistics(false);
176  }
177  }
178 
185  void Mult(const mfem::Vector& input, mfem::Vector& output) const;
186 
195  void SetOperator(const mfem::Operator& op);
196 
197  private:
202  mutable std::unique_ptr<mfem::SuperLURowLocMatrix> superlu_mat_;
203 
208  mutable std::unique_ptr<mfem::HypreParMatrix> monolithic_mat_;
209 
215  mfem::SuperLUSolver superlu_solver_;
216 };
217 
218 #ifdef MFEM_USE_STRUMPACK
222 class StrumpackSolver : public mfem::Solver {
223  public:
229  StrumpackSolver(int print_level, MPI_Comm comm) : strumpack_solver_(comm)
230  {
231  strumpack_solver_.SetKrylovSolver(strumpack::KrylovSolver::DIRECT);
232  strumpack_solver_.SetReorderingStrategy(strumpack::ReorderingStrategy::METIS);
233 
234  if (print_level == 1) {
235  strumpack_solver_.SetPrintFactorStatistics(true);
236  strumpack_solver_.SetPrintSolveStatistics(true);
237  }
238  }
239 
246  void Mult(const mfem::Vector& input, mfem::Vector& output) const;
247 
254  void SetOperator(const mfem::Operator& op);
255 
256  private:
261  mutable std::unique_ptr<mfem::STRUMPACKRowLocMatrix> strumpack_mat_;
262 
267  mutable std::unique_ptr<mfem::HypreParMatrix> monolithic_mat_;
268 
274  mfem::STRUMPACKSolver strumpack_solver_;
275 };
276 
277 #endif
278 
283 std::unique_ptr<mfem::HypreParMatrix> buildMonolithicMatrix(const mfem::BlockOperator&);
284 
294 std::unique_ptr<mfem::NewtonSolver> buildNonlinearSolver(NonlinearSolverOptions nonlinear_opts,
295  const LinearSolverOptions& linear_opts,
296  mfem::Solver& preconditioner, MPI_Comm comm = MPI_COMM_WORLD);
297 
305 std::pair<std::unique_ptr<mfem::Solver>, std::unique_ptr<mfem::Solver>> buildLinearSolverAndPreconditioner(
306  LinearSolverOptions linear_opts = {}, MPI_Comm comm = MPI_COMM_WORLD);
307 
312 bool requiresMonolithicOperator(const LinearSolverOptions& linear_opts);
313 
320 std::unique_ptr<mfem::Solver> buildPreconditioner(LinearSolverOptions linear_opts,
321  [[maybe_unused]] MPI_Comm comm = MPI_COMM_WORLD);
322 
323 #ifdef MFEM_USE_AMGX
331 std::unique_ptr<mfem::AmgXSolver> buildAMGX(const AMGXOptions& options, const MPI_Comm comm);
332 #endif
333 
334 } // namespace smith
335 
341 template <>
342 struct FromInlet<smith::LinearSolverOptions> {
344  smith::LinearSolverOptions operator()(const axom::inlet::Container& base);
345 };
346 
352 template <>
353 struct FromInlet<smith::NonlinearSolverOptions> {
355  smith::NonlinearSolverOptions operator()(const axom::inlet::Container& base);
356 };
357 
363 template <>
364 struct FromInlet<smith::EquationSolver> {
366  smith::EquationSolver operator()(const axom::inlet::Container& base);
367 };
This class manages the objects typically required to solve a nonlinear set of equations arising from ...
void resetConvergenceState() const
Reset nonlinear convergence state stored across iterations of a single solve.
mfem::NewtonSolver & nonlinearSolver()
const mfem::Solver & linearSolver() const
mfem::Solver & preconditioner()
void setConvergenceTolerances(double abs_tol, double rel_tol, MPI_Comm comm) const
Update scalar convergence tolerances for the managed nonlinear convergence path.
EquationSolver(std::unique_ptr< mfem::NewtonSolver > nonlinear_solver, std::unique_ptr< mfem::Solver > linear_solver, std::unique_ptr< mfem::Solver > preconditioner=nullptr)
const mfem::Solver & preconditioner() const
void solve(mfem::Vector &x) const
static void defineInputFileSchema(axom::inlet::Container &container)
const mfem::NewtonSolver & nonlinearSolver() const
mfem::Solver & linearSolver()
void setOperator(const mfem::Operator &op)
A wrapper class for using the MFEM SuperLU solver with a HypreParMatrix.
void SetOperator(const mfem::Operator &op)
Set the underlying matrix operator to use in the solution algorithm.
SuperLUSolver(int print_level, MPI_Comm comm)
Constructs a wrapper over an mfem::SuperLUSolver.
void Mult(const mfem::Vector &input, mfem::Vector &output) const
Factor and solve the linear system y = Op^{-1} x using DSuperLU.
This file contains the all the necessary functions for reading input files.
Accelerator functionality.
Definition: smith.cpp:36
std::unique_ptr< mfem::NewtonSolver > buildNonlinearSolver(NonlinearSolverOptions nonlinear_opts, const LinearSolverOptions &linear_opts, mfem::Solver &prec, MPI_Comm comm)
Build a nonlinear solver using the nonlinear option struct.
bool requiresMonolithicOperator(const LinearSolverOptions &linear_opts)
Return true if the configured linear solve stack requires block operators to be merged.
std::unique_ptr< mfem::HypreParMatrix > buildMonolithicMatrix(const mfem::BlockOperator &block_operator)
Build a monolithic HypreParMatrix from a BlockOperator.
std::pair< std::unique_ptr< mfem::Solver >, std::unique_ptr< mfem::Solver > > buildLinearSolverAndPreconditioner(LinearSolverOptions linear_opts, MPI_Comm comm)
Build the linear solver and its associated preconditioner given a linear options struct.
std::unique_ptr< mfem::Solver > buildPreconditioner(LinearSolverOptions linear_opts, [[maybe_unused]] MPI_Comm comm)
Build a preconditioner from the available options.
This file contains enumerations and record types for physics solver configuration.
Parameters for an iterative linear solution scheme.
Nonlinear solution scheme parameters.