Serac  0.1
Serac is an implicit thermal strucural mechanics simulation code.
equation_solver.hpp
Go to the documentation of this file.
1 // Copyright (c) 2019-2024, Lawrence Livermore National Security, LLC and
2 // other Serac 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 
19 #include "mfem.hpp"
20 
23 
24 namespace serac {
25 
43 public:
44  // _equationsolver_constructor_start
52  EquationSolver(std::unique_ptr<mfem::NewtonSolver> nonlinear_solver, std::unique_ptr<mfem::Solver> linear_solver,
53  std::unique_ptr<mfem::Solver> preconditioner = nullptr);
54  // _equationsolver_constructor_end
55 
56  // _build_equationsolver_start
64  EquationSolver(NonlinearSolverOptions nonlinear_opts = {}, LinearSolverOptions lin_opts = {},
65  MPI_Comm comm = MPI_COMM_WORLD);
66  // _build_equationsolver_end
67 
74  void setOperator(const mfem::Operator& op);
75 
81  void solve(mfem::Vector& x) const;
82 
87  mfem::NewtonSolver& nonlinearSolver() { return *nonlin_solver_; }
88 
92  const mfem::NewtonSolver& nonlinearSolver() const { return *nonlin_solver_; }
93 
98  mfem::Solver& linearSolver() { return *lin_solver_; }
99 
103  const mfem::Solver& linearSolver() const { return *lin_solver_; }
104 
110  mfem::Solver* preconditioner() { return preconditioner_.get(); }
111 
115  const mfem::Solver* preconditioner() const { return preconditioner_.get(); }
116 
120  static void defineInputFileSchema(axom::inlet::Container& container);
121 
122 private:
126  std::unique_ptr<mfem::Solver> preconditioner_;
127 
131  std::unique_ptr<mfem::Solver> lin_solver_;
132 
136  std::unique_ptr<mfem::NewtonSolver> nonlin_solver_;
137 
143  bool nonlin_solver_set_solver_called_ = false;
144 };
145 
149 class SuperLUSolver : public mfem::Solver {
150 public:
156  SuperLUSolver(int print_level, MPI_Comm comm) : superlu_solver_(comm)
157  {
158  superlu_solver_.SetColumnPermutation(mfem::superlu::PARMETIS);
159  if (print_level == 0) {
160  superlu_solver_.SetPrintStatistics(false);
161  }
162  }
163 
170  void Mult(const mfem::Vector& input, mfem::Vector& output) const;
171 
180  void SetOperator(const mfem::Operator& op);
181 
182 private:
187  mutable std::unique_ptr<mfem::SuperLURowLocMatrix> superlu_mat_;
188 
194  mfem::SuperLUSolver superlu_solver_;
195 };
196 
197 #ifdef MFEM_USE_STRUMPACK
201 class StrumpackSolver : public mfem::Solver {
202 public:
208  StrumpackSolver(int print_level, MPI_Comm comm) : strumpack_solver_(comm)
209  {
210  strumpack_solver_.SetKrylovSolver(strumpack::KrylovSolver::DIRECT);
211  strumpack_solver_.SetReorderingStrategy(strumpack::ReorderingStrategy::METIS);
212 
213  if (print_level == 1) {
214  strumpack_solver_.SetPrintFactorStatistics(true);
215  strumpack_solver_.SetPrintSolveStatistics(true);
216  }
217  }
218 
225  void Mult(const mfem::Vector& input, mfem::Vector& output) const;
226 
233  void SetOperator(const mfem::Operator& op);
234 
235 private:
240  mutable std::unique_ptr<mfem::STRUMPACKRowLocMatrix> strumpack_mat_;
241 
247  mfem::STRUMPACKSolver strumpack_solver_;
248 };
249 
250 #endif
251 
259 std::unique_ptr<mfem::NewtonSolver> buildNonlinearSolver(NonlinearSolverOptions nonlinear_opts = {},
260  MPI_Comm comm = MPI_COMM_WORLD);
261 
269 std::pair<std::unique_ptr<mfem::Solver>, std::unique_ptr<mfem::Solver>> buildLinearSolverAndPreconditioner(
270  LinearSolverOptions linear_opts = {}, MPI_Comm comm = MPI_COMM_WORLD);
271 
280 std::unique_ptr<mfem::Solver> buildPreconditioner(Preconditioner preconditioner, int print_level = 0,
281  [[maybe_unused]] MPI_Comm comm = MPI_COMM_WORLD);
282 
283 #ifdef MFEM_USE_AMGX
291 std::unique_ptr<mfem::AmgXSolver> buildAMGX(const AMGXOptions& options, const MPI_Comm comm);
292 #endif
293 } // namespace serac
294 
300 template <>
301 struct FromInlet<serac::LinearSolverOptions> {
303  serac::LinearSolverOptions operator()(const axom::inlet::Container& base);
304 };
305 
311 template <>
312 struct FromInlet<serac::NonlinearSolverOptions> {
314  serac::NonlinearSolverOptions operator()(const axom::inlet::Container& base);
315 };
316 
322 template <>
323 struct FromInlet<serac::EquationSolver> {
325  serac::EquationSolver operator()(const axom::inlet::Container& base);
326 };
This class manages the objects typically required to solve a nonlinear set of equations arising from ...
mfem::NewtonSolver & nonlinearSolver()
const mfem::Solver & linearSolver() const
void setOperator(const mfem::Operator &op)
void solve(mfem::Vector &x) const
mfem::Solver * preconditioner()
static void defineInputFileSchema(axom::inlet::Container &container)
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
const mfem::NewtonSolver & nonlinearSolver() const
mfem::Solver & linearSolver()
A wrapper class for using the MFEM SuperLU solver with a HypreParMatrix.
void Mult(const mfem::Vector &input, mfem::Vector &output) const
Factor and solve the linear system y = Op^{-1} x using DSuperLU.
SuperLUSolver(int print_level, MPI_Comm comm)
Constructs a wrapper over an mfem::SuperLUSolver.
void SetOperator(const mfem::Operator &op)
Set the underlying matrix operator to use in the solution algorithm.
This file contains the all the necessary functions for reading input files.
Accelerator functionality.
Definition: serac.cpp:38
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.
Preconditioner
The type of preconditioner to be used.
std::unique_ptr< mfem::Solver > buildPreconditioner(Preconditioner preconditioner, int print_level, [[maybe_unused]] MPI_Comm comm)
Build a preconditioner from the available options.
std::unique_ptr< mfem::NewtonSolver > buildNonlinearSolver(NonlinearSolverOptions nonlinear_opts, MPI_Comm comm)
Build a nonlinear solver using the nonlinear option struct.
This file contains enumerations and record types for physics solver configuration.
Parameters for an iterative linear solution scheme.
Nonlinear solution scheme parameters.