16 #include "libmesh/implicit_system.h" 17 #include "libmesh/static_condensation.h" 26 params.
set<
bool>(
"full") =
true;
28 params.
addParam<std::vector<NonlinearVariableName>>(
31 "A list of variables for whom to not statically condense their degrees of freedom out of the " 32 "system. By default all degrees of freedom on element interiors are condensed out.");
42 auto check_param = [
this](
const auto & param_name)
46 "This class prefixes every PETSc option so that it applies to the condensed " 47 "system. Given that, there are multiple issues with setting '",
49 "': 1) it applies to the nonlinear solver algorithm whereas the prefix this class " 50 "applies makes PETSc options conceptually applicable only to the linear solve of " 51 "the statically condensed system 2) these are singleton MOOSE-wrapped PETSc " 52 "options. Consequently even if having multiple prefixes for a system's nonlinear " 53 "solver options made sense, we don't support it. E.g. if you specify '",
55 "' in both this object's block and the Executioner block, then there will be a " 58 check_param(
"mffd_type");
59 check_param(
"solve_type");
64 "Static condensation preconditioning should use a PJFNK solve type. This is because it " 65 "provides a means to compute the action of the Jacobian on vectors, which we otherwise " 66 "would not have because when using static condensation, the Jacobian is never formed. Note " 67 "that actions of the Jacobian on a vector are necessary for things like: GMRES, printing " 68 "of linear residuals, or application of line searches like cubic backtracking. These " 69 "particular operations can be avoided by using '-ksp_type preonly', disabling printing of " 70 "linear residuals, and using the 'cp' or 'basic' line search.");
74 mooseError(
"Static condensation can only be used with implicit systems");
75 implicit_sys->create_static_condensation();
76 auto & sc = implicit_sys->get_static_condensation();
77 std::unordered_set<unsigned int> uncondensed_vars;
78 for (
auto & nl_var_name :
getParam<std::vector<NonlinearVariableName>>(
"dont_condense_vars"))
80 sc.dont_condense_vars(uncondensed_vars);
NonlinearSystemBase & _nl
The nonlinear system whose linearization this preconditioner should be applied to.
virtual void initialSetup() override
Perform some setup tasks such as storing the PETSc options.
unsigned int number() const
Get variable number coming from libMesh.
static InputParameters validParams()
FEProblemBase & _fe_problem
Subproblem this preconditioner is part of.
static InputParameters validParams()
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
virtual const std::string & name() const
const std::string _type
The type of this class.
const T & getParam(const std::string &name) const
Retrieve a parameter for the object.
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
void paramError(const std::string ¶m, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
unsigned int number() const
Gets the number of this system.
Static condensation preconditioner.
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
SolverParams & solverParams(unsigned int solver_sys_num=0)
Get the solver parameters.
Single matrix preconditioner.
Preconditioned Jacobian-Free Newton Krylov.
registerMooseObjectAliased("MooseApp", MooseStaticCondensationPreconditioner, "StaticCondensation")
MooseVariableFieldBase & getVariable(THREAD_ID tid, const std::string &var_name) const
Gets a reference to a variable of with specified name.
void storePetscOptions(FEProblemBase &fe_problem, const std::string &prefix, const ParallelParamObject ¶m_object)
Stores the PETSc options supplied from the parameter object on the problem.
MooseStaticCondensationPreconditioner(const InputParameters ¶ms)
virtual libMesh::System & system() override
Get the reference to the libMesh system.