24 #include "smith/numerics/nonlinear_convergence.hpp"
26 #include "smith/numerics/petsc_solvers.hpp"
56 EquationSolver(std::unique_ptr<mfem::NewtonSolver> nonlinear_solver, std::unique_ptr<mfem::Solver> linear_solver,
69 MPI_Comm comm = MPI_COMM_WORLD);
85 void solve(mfem::Vector& x)
const;
133 void attachConvergenceManager()
const;
134 void initializeConvergenceManager(
double abs_tol,
double rel_tol, MPI_Comm comm)
const;
139 std::unique_ptr<mfem::Solver> preconditioner_;
144 std::unique_ptr<mfem::Solver> lin_solver_;
149 std::unique_ptr<mfem::NewtonSolver> nonlin_solver_;
156 bool nonlin_solver_set_solver_called_ =
false;
158 mutable std::shared_ptr<EquationSolverConvergenceManager> convergence_manager_ =
nullptr;
173 superlu_solver_.SetColumnPermutation(mfem::superlu::PARMETIS);
174 if (print_level == 0) {
175 superlu_solver_.SetPrintStatistics(
false);
185 void Mult(
const mfem::Vector& input, mfem::Vector& output)
const;
202 mutable std::unique_ptr<mfem::SuperLURowLocMatrix> superlu_mat_;
208 mutable std::unique_ptr<mfem::HypreParMatrix> monolithic_mat_;
215 mfem::SuperLUSolver superlu_solver_;
218 #ifdef MFEM_USE_STRUMPACK
222 class StrumpackSolver :
public mfem::Solver {
229 StrumpackSolver(
int print_level, MPI_Comm comm) : strumpack_solver_(comm)
231 strumpack_solver_.SetKrylovSolver(strumpack::KrylovSolver::DIRECT);
232 strumpack_solver_.SetReorderingStrategy(strumpack::ReorderingStrategy::METIS);
234 if (print_level == 1) {
235 strumpack_solver_.SetPrintFactorStatistics(
true);
236 strumpack_solver_.SetPrintSolveStatistics(
true);
246 void Mult(
const mfem::Vector& input, mfem::Vector& output)
const;
254 void SetOperator(
const mfem::Operator& op);
261 mutable std::unique_ptr<mfem::STRUMPACKRowLocMatrix> strumpack_mat_;
267 mutable std::unique_ptr<mfem::HypreParMatrix> monolithic_mat_;
274 mfem::STRUMPACKSolver strumpack_solver_;
295 const LinearSolverOptions& linear_opts,
296 mfem::Solver& preconditioner, MPI_Comm comm = MPI_COMM_WORLD);
306 LinearSolverOptions linear_opts = {}, MPI_Comm comm = MPI_COMM_WORLD);
321 [[maybe_unused]] MPI_Comm comm = MPI_COMM_WORLD);
331 std::unique_ptr<mfem::AmgXSolver> buildAMGX(
const AMGXOptions& options,
const MPI_Comm comm);
342 struct FromInlet<
smith::LinearSolverOptions> {
353 struct FromInlet<
smith::NonlinearSolverOptions> {
364 struct FromInlet<
smith::EquationSolver> {
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.
Accelerator functionality.
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.