https://mooseframework.inl.gov
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
FieldSplitPreconditionerTempl< Base > Class Template Referenceabstract

Implements a preconditioner designed to map onto PETSc's PCFieldSplit. More...

#include <FieldSplitPreconditioner.h>

Inheritance diagram for FieldSplitPreconditionerTempl< Base >:
[legend]

Public Member Functions

 FieldSplitPreconditionerTempl (const InputParameters &parameters)
 
virtual void setupDM ()=0
 setup the data management data structure that manages the field split More...
 
virtual KSP getKSP ()=0
 

Static Public Member Functions

static InputParameters validParams ()
 Constructor. More...
 

Protected Member Functions

virtual const libMesh::DofMapBasedofMap () const =0
 
virtual const libMesh::Systemsystem () const =0
 
virtual std::string prefix () const =0
 
void createMooseDM (DM *dm)
 creates the MOOSE data management object More...
 

Protected Attributes

NonlinearSystemBase_nl
 The nonlinear system this FSP is associated with (convenience reference) More...
 
std::string _decomposition_split
 The decomposition split. More...
 

Detailed Description

template<typename Base>
class FieldSplitPreconditionerTempl< Base >

Implements a preconditioner designed to map onto PETSc's PCFieldSplit.

Definition at line 48 of file FieldSplitPreconditioner.h.

Constructor & Destructor Documentation

◆ FieldSplitPreconditionerTempl()

template<typename Base >
FieldSplitPreconditionerTempl< Base >::FieldSplitPreconditionerTempl ( const InputParameters parameters)

Definition at line 49 of file FieldSplitPreconditioner.C.

51  : Base(parameters),
52  _nl(this->_fe_problem.getNonlinearSystemBase(this->_nl_sys_num)),
53  _decomposition_split(this->template getParam<std::string>("topsplit"))
54 {
56  mooseError("The field split preconditioner can only be used with PETSc");
57 
58  // number of variables
59  unsigned int n_vars = _nl.nVariables();
60  // if we want to construct a full Jacobian?
61  // it is recommended to have a full Jacobian for using
62  // the fieldSplit preconditioner
63  bool full = this->template getParam<bool>("full");
64 
65  // how variables couple
66  std::unique_ptr<CouplingMatrix> cm = std::make_unique<CouplingMatrix>(n_vars);
67  if (!full)
68  {
69  if (this->isParamValid("off_diag_row") && this->isParamValid("off_diag_column"))
70  {
71 
72  const auto off_diag_rows =
73  this->template getParam<std::vector<NonlinearVariableName>>("off_diag_row");
74  const auto off_diag_columns =
75  this->template getParam<std::vector<NonlinearVariableName>>("off_diag_column");
76 
77  // put 1s on diagonal
78  for (unsigned int i = 0; i < n_vars; i++)
79  (*cm)(i, i) = 1;
80 
81  // off-diagonal entries
82  std::vector<std::vector<unsigned int>> off_diag(n_vars);
83  if (off_diag_rows.size() * off_diag_columns.size() != 0 &&
84  off_diag_rows.size() == off_diag_columns.size())
85  for (const auto i : index_range(off_diag_rows))
86  {
87  unsigned int row = _nl.getVariable(0, off_diag_rows[i]).number();
88  unsigned int column = _nl.getVariable(0, off_diag_columns[i]).number();
89  (*cm)(row, column) = 1;
90  }
91  }
92  }
93  else
94  {
95  for (unsigned int i = 0; i < n_vars; i++)
96  for (unsigned int j = 0; j < n_vars; j++)
97  (*cm)(i, j) = 1; // full coupling
98  }
99  this->setCouplingMatrix(std::move(cm));
100 
101  // turn on a flag
103 }
NonlinearSystemBase & _nl
The nonlinear system this FSP is associated with (convenience reference)
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:333
unsigned int number() const
Get variable number coming from libMesh.
virtual unsigned int nVariables() const
Get the number of variables in this system.
Definition: SystemBase.C:883
SolverPackage default_solver_package()
unsigned int n_vars
std::string _decomposition_split
The decomposition split.
MooseVariableFieldBase & getVariable(THREAD_ID tid, const std::string &var_name) const
Gets a reference to a variable of with specified name.
Definition: SystemBase.C:90
void useFieldSplitPreconditioner(FieldSplitPreconditionerBase *fsp)
If called with a non-null object true this system will use a field split preconditioner matrix...
auto index_range(const T &sizable)

Member Function Documentation

◆ createMooseDM()

template<typename Base >
void FieldSplitPreconditionerTempl< Base >::createMooseDM ( DM *  dm)
protected

creates the MOOSE data management object

Definition at line 107 of file FieldSplitPreconditioner.C.

108 {
109  LibmeshPetscCallA(
110  _nl.comm().get(),
112  LibmeshPetscCallA(
113  _nl.comm().get(),
115  LibmeshPetscCallA(_nl.comm().get(),
116  PetscObjectSetOptionsPrefix((PetscObject)*dm, prefix().c_str()));
117  LibmeshPetscCall(DMSetFromOptions(*dm));
118  LibmeshPetscCall(DMSetUp(*dm));
119 }
NonlinearSystemBase & _nl
The nonlinear system this FSP is associated with (convenience reference)
virtual const libMesh::DofMapBase & dofMap() const =0
PetscErrorCode PetscOptionItems *PetscErrorCode DM dm
const Parallel::Communicator & comm() const
PetscErrorCode DMCreateMoose(MPI_Comm comm, NonlinearSystemBase &nl, const libMesh::DofMapBase &dof_map, const libMesh::System &system, const std::string &dm_name, DM *dm)
Create a MOOSE DM.
std::string _decomposition_split
The decomposition split.
virtual std::string prefix() const =0
virtual const libMesh::System & system() const =0

◆ dofMap()

template<typename Base>
virtual const libMesh::DofMapBase& FieldSplitPreconditionerTempl< Base >::dofMap ( ) const
protectedpure virtual
Returns
The degree of freedom map to use for decomposition

Implemented in FieldSplitPreconditioner, and StaticCondensationFieldSplitPreconditioner.

◆ getKSP()

virtual KSP FieldSplitPreconditionerBase::getKSP ( )
pure virtualinherited
Returns
The KSP object associated with the field split preconditioner

Implemented in FieldSplitPreconditioner, and StaticCondensationFieldSplitPreconditioner.

◆ prefix()

template<typename Base>
virtual std::string FieldSplitPreconditionerTempl< Base >::prefix ( ) const
protectedpure virtual
Returns
The prefix to pass to PETSc for the DM

Implemented in FieldSplitPreconditioner, and StaticCondensationFieldSplitPreconditioner.

◆ setupDM()

virtual void FieldSplitPreconditionerBase::setupDM ( )
pure virtualinherited

setup the data management data structure that manages the field split

Implemented in FieldSplitPreconditioner, and StaticCondensationFieldSplitPreconditioner.

Referenced by NonlinearSystemBase::setupDM().

◆ system()

template<typename Base>
virtual const libMesh::System& FieldSplitPreconditionerTempl< Base >::system ( ) const
protectedpure virtual

◆ validParams()

template<typename Base >
InputParameters FieldSplitPreconditionerTempl< Base >::validParams ( )
static

Constructor.

Initializes SplitBasedPreconditioner data structures

Definition at line 32 of file FieldSplitPreconditioner.C.

Referenced by StaticCondensationFieldSplitPreconditioner::validParams(), and FieldSplitPreconditioner::validParams().

33 {
35  params.addClassDescription("Preconditioner designed to map onto PETSc's PCFieldSplit.");
36 
37  params.addRequiredParam<std::string>(
38  "topsplit", "Entrance to splits, the top split will specify how splits will go.");
39  // We should use full coupling Jacobian matrix by default
40  params.addParam<bool>("full",
41  true,
42  "Set to true if you want the full set of couplings between variables "
43  "simply for convenience so you don't have to set every off_diag_row "
44  "and off_diag_column combination.");
45  return params;
46 }
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
InputParameters validParams()
void addRequiredParam(const std::string &name, const std::string &doc_string)
This method adds a parameter and documentation string to the InputParameters object that will be extr...
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...

Member Data Documentation

◆ _decomposition_split

template<typename Base>
std::string FieldSplitPreconditionerTempl< Base >::_decomposition_split
protected

The decomposition split.

Definition at line 85 of file FieldSplitPreconditioner.h.

◆ _nl

template<typename Base>
NonlinearSystemBase& FieldSplitPreconditionerTempl< Base >::_nl
protected

The nonlinear system this FSP is associated with (convenience reference)

Definition at line 80 of file FieldSplitPreconditioner.h.

Referenced by FieldSplitPreconditionerTempl< MoosePreconditioner >::FieldSplitPreconditionerTempl().


The documentation for this class was generated from the following files: