www.mooseframework.org
PhysicsBasedPreconditioner.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
11 
12 // MOOSE includes
14 #include "FEProblem.h"
15 #include "MooseEnum.h"
16 #include "MooseVariableFE.h"
17 #include "NonlinearSystem.h"
18 #include "PetscSupport.h"
19 
20 #include "libmesh/libmesh_common.h"
21 #include "libmesh/equation_systems.h"
22 #include "libmesh/nonlinear_implicit_system.h"
23 #include "libmesh/nonlinear_solver.h"
24 #include "libmesh/linear_implicit_system.h"
25 #include "libmesh/transient_system.h"
26 #include "libmesh/numeric_vector.h"
27 #include "libmesh/sparse_matrix.h"
28 #include "libmesh/string_to_enum.h"
29 #include "libmesh/coupling_matrix.h"
30 
32 
34 
37 {
39 
40  params.addRequiredParam<std::vector<std::string>>(
41  "solve_order",
42  "The order the block rows will be solved in. Put the name of variables here "
43  "to stand for solving that variable's block row. A variable may appear more "
44  "than once (to create cylces if you like).");
45  params.addRequiredParam<std::vector<std::string>>("preconditioner", "TODO: docstring");
46 
47  params.addParam<std::vector<std::string>>(
48  "off_diag_row",
49  "The off diagonal row you want to add into the matrix, it will be associated "
50  "with an off diagonal column from the same position in off_diag_colum.");
51  params.addParam<std::vector<std::string>>("off_diag_column",
52  "The off diagonal column you want to add into the "
53  "matrix, it will be associated with an off diagonal "
54  "row from the same position in off_diag_row.");
55 
56  return params;
57 }
58 
60  : MoosePreconditioner(params),
61  Preconditioner<Number>(MoosePreconditioner::_communicator),
62  _nl(_fe_problem.getNonlinearSystemBase()),
63  _init_timer(registerTimedSection("init", 2)),
64  _apply_timer(registerTimedSection("apply", 1))
65 {
66  unsigned int num_systems = _nl.system().n_vars();
67  _systems.resize(num_systems);
68  _preconditioners.resize(num_systems);
69  _off_diag.resize(num_systems);
70  _off_diag_mats.resize(num_systems);
71  _pre_type.resize(num_systems);
72 
73  { // Setup the Coupling Matrix so MOOSE knows what we're doing
75  unsigned int n_vars = nl.nVariables();
76 
77  // The coupling matrix is held and released by FEProblemBase, so it is not released in this
78  // object
79  std::unique_ptr<CouplingMatrix> cm = libmesh_make_unique<CouplingMatrix>(n_vars);
80 
81  bool full = false; // getParam<bool>("full"); // TODO: add a FULL option for PBP
82 
83  if (!full)
84  {
85  // put 1s on diagonal
86  for (unsigned int i = 0; i < n_vars; i++)
87  (*cm)(i, i) = 1;
88 
89  // off-diagonal entries
90  std::vector<std::vector<unsigned int>> off_diag(n_vars);
91  for (unsigned int i = 0; i < getParam<std::vector<std::string>>("off_diag_row").size(); i++)
92  {
93  unsigned int row =
94  nl.getVariable(0, getParam<std::vector<std::string>>("off_diag_row")[i]).number();
95  unsigned int column =
96  nl.getVariable(0, getParam<std::vector<std::string>>("off_diag_column")[i]).number();
97  (*cm)(row, column) = 1;
98  }
99 
100  // TODO: handle coupling entries between NL-vars and SCALAR-vars
101  }
102  else
103  {
104  for (unsigned int i = 0; i < n_vars; i++)
105  for (unsigned int j = 0; j < n_vars; j++)
106  (*cm)(i, j) = 1;
107  }
108 
109  _fe_problem.setCouplingMatrix(std::move(cm));
110  }
111 
112  // PC types
113  const std::vector<std::string> & pc_types = getParam<std::vector<std::string>>("preconditioner");
114  for (unsigned int i = 0; i < num_systems; i++)
115  _pre_type[i] = Utility::string_to_enum<PreconditionerType>(pc_types[i]);
116 
117  // solve order
118  const std::vector<std::string> & solve_order = getParam<std::vector<std::string>>("solve_order");
119  _solve_order.resize(solve_order.size());
120  for (unsigned int i = 0; i < solve_order.size(); i++)
121  _solve_order[i] = _nl.system().variable_number(solve_order[i]);
122 
123  // diag and off-diag systems
124  unsigned int n_vars = _nl.system().n_vars();
125 
126  // off-diagonal entries
127  const std::vector<std::string> & odr = getParam<std::vector<std::string>>("off_diag_row");
128  const std::vector<std::string> & odc = getParam<std::vector<std::string>>("off_diag_column");
129  std::vector<std::vector<unsigned int>> off_diag(n_vars);
130  for (unsigned int i = 0; i < odr.size(); i++)
131  {
132  unsigned int row = _nl.system().variable_number(odr[i]);
133  unsigned int column = _nl.system().variable_number(odc[i]);
134  off_diag[row].push_back(column);
135  }
136  // Add all of the preconditioning systems
137  for (unsigned int var = 0; var < n_vars; var++)
138  addSystem(var, off_diag[var], _pre_type[var]);
139 
140  _nl.nonlinearSolver()->attach_preconditioner(this);
141 
143  mooseError("PBP must be used with JFNK solve type");
144 }
145 
147 
148 void
150  std::vector<unsigned int> off_diag,
151  PreconditionerType type)
152 {
153  std::string var_name = _nl.system().variable_name(var);
154 
155  LinearImplicitSystem & precond_system =
156  _fe_problem.es().add_system<LinearImplicitSystem>(var_name + "_system");
157  precond_system.assemble_before_solve = false;
158 
159  const std::set<SubdomainID> * active_subdomains = _nl.getVariableBlocks(var);
160  precond_system.add_variable(
161  var_name + "_prec", _nl.system().variable(var).type(), active_subdomains);
162 
163  _systems[var] = &precond_system;
164  _pre_type[var] = type;
165 
166  _off_diag_mats[var].resize(off_diag.size());
167  for (unsigned int i = 0; i < off_diag.size(); i++)
168  {
169  // Add the matrix to hold the off-diagonal piece
170  _off_diag_mats[var][i] = &precond_system.add_matrix(_nl.system().variable_name(off_diag[i]));
171  }
172 
173  _off_diag[var] = off_diag;
174 }
175 
176 void
178 {
179  TIME_SECTION(_init_timer);
180 
181  // Tell libMesh that this is initialized!
182  _is_initialized = true;
183 
184  const unsigned int num_systems = _systems.size();
185 
186  // If no order was specified, just solve them in increasing order
187  if (_solve_order.size() == 0)
188  {
189  _solve_order.resize(num_systems);
190  for (unsigned int i = 0; i < num_systems; i++)
191  _solve_order[i] = i;
192  }
193 
194  // Loop over variables
195  for (unsigned int system_var = 0; system_var < num_systems; system_var++)
196  {
197  LinearImplicitSystem & u_system = *_systems[system_var];
198 
199  if (!_preconditioners[system_var])
200  _preconditioners[system_var] =
201  Preconditioner<Number>::build_preconditioner(MoosePreconditioner::_communicator);
202 
203  // we have to explicitly set the matrix in the preconditioner, because h-adaptivity could have
204  // changed it and we have to work with the current one
205  Preconditioner<Number> * preconditioner = _preconditioners[system_var].get();
206  preconditioner->set_matrix(*u_system.matrix);
207  preconditioner->set_type(_pre_type[system_var]);
208 
209  preconditioner->init();
210  }
211 }
212 
213 void
215 {
216  const unsigned int num_systems = _systems.size();
217 
218  std::vector<JacobianBlock *> blocks;
219 
220  // Loop over variables
221  for (unsigned int system_var = 0; system_var < num_systems; system_var++)
222  {
223  LinearImplicitSystem & u_system = *_systems[system_var];
224 
225  {
226  JacobianBlock * block = new JacobianBlock(u_system, *u_system.matrix, system_var, system_var);
227  blocks.push_back(block);
228  }
229 
230  for (unsigned int diag = 0; diag < _off_diag[system_var].size(); diag++)
231  {
232  unsigned int coupled_var = _off_diag[system_var][diag];
233  std::string coupled_name = _nl.system().variable_name(coupled_var);
234 
235  JacobianBlock * block =
236  new JacobianBlock(u_system, *_off_diag_mats[system_var][diag], system_var, coupled_var);
237  blocks.push_back(block);
238  }
239  }
240 
242 
243  // cleanup
244  for (auto & block : blocks)
245  delete block;
246 }
247 
248 void
249 PhysicsBasedPreconditioner::apply(const NumericVector<Number> & x, NumericVector<Number> & y)
250 {
251  TIME_SECTION(_apply_timer);
252 
253  const unsigned int num_systems = _systems.size();
254 
255  MooseMesh & mesh = _fe_problem.mesh();
256 
257  // Zero out the solution vectors
258  for (unsigned int sys = 0; sys < num_systems; sys++)
259  _systems[sys]->solution->zero();
260 
261  // Loop over solve order
262  for (unsigned int i = 0; i < _solve_order.size(); i++)
263  {
264  unsigned int system_var = _solve_order[i];
265 
266  LinearImplicitSystem & u_system = *_systems[system_var];
267 
268  // Copy rhs from the big system into the small one
270  mesh, _nl.system().number(), system_var, x, u_system.number(), 0, *u_system.rhs);
271 
272  // Modify the RHS by subtracting off the matvecs of the solutions for the other preconditioning
273  // systems with the off diagonal blocks in this system.
274  for (unsigned int diag = 0; diag < _off_diag[system_var].size(); diag++)
275  {
276  unsigned int coupled_var = _off_diag[system_var][diag];
277  LinearImplicitSystem & coupled_system = *_systems[coupled_var];
278  SparseMatrix<Number> & off_diag = *_off_diag_mats[system_var][diag];
279  NumericVector<Number> & rhs = *u_system.rhs;
280 
281  // This next bit computes rhs -= A*coupled_solution
282  // It does what it does because there is no vector_mult_sub()
283  rhs.close();
284  rhs.scale(-1.0);
285  rhs.close();
286  off_diag.vector_mult_add(rhs, *coupled_system.solution);
287  rhs.close();
288  rhs.scale(-1.0);
289  rhs.close();
290  }
291 
292  // Apply the preconditioner to the small system
293  _preconditioners[system_var]->apply(*u_system.rhs, *u_system.solution);
294 
295  // Copy solution from small system into the big one
296  // copyVarValues(mesh,system,0,*u_system.solution,0,system_var,y);
297  }
298 
299  // Copy the solutions out
300  for (unsigned int system_var = 0; system_var < num_systems; system_var++)
301  {
302  LinearImplicitSystem & u_system = *_systems[system_var];
303 
305  mesh, u_system.number(), 0, *u_system.solution, _nl.system().number(), system_var, y);
306  }
307 
308  y.close();
309 }
310 
311 void
313 {
314 }
type
MatType type
Definition: PetscDMMoose.C:1477
SystemBase::getVariable
MooseVariableFEBase & getVariable(THREAD_ID tid, const std::string &var_name)
Gets a reference to a variable of with specified name.
Definition: SystemBase.C:109
PhysicsBasedPreconditioner::setup
virtual void setup()
This is called every time the "operator might have changed".
Definition: PhysicsBasedPreconditioner.C:214
MoosePreconditioner::validParams
static InputParameters validParams()
Definition: MoosePreconditioner.C:21
FEProblem.h
PhysicsBasedPreconditioner::_systems
std::vector< LinearImplicitSystem * > _systems
List of linear system that build up the preconditioner.
Definition: PhysicsBasedPreconditioner.h:80
MooseObject::mooseError
void mooseError(Args &&... args) const
Definition: MooseObject.h:141
PhysicsBasedPreconditioner::_off_diag
std::vector< std::vector< unsigned int > > _off_diag
Holds which off diagonal blocks to compute.
Definition: PhysicsBasedPreconditioner.h:88
PhysicsBasedPreconditioner::_solve_order
std::vector< unsigned int > _solve_order
Holds the order the blocks are solved for.
Definition: PhysicsBasedPreconditioner.h:84
MoosePreconditioner::_fe_problem
FEProblemBase & _fe_problem
Subproblem this preconditioner is part of.
Definition: MoosePreconditioner.h:56
MoosePreconditioner::copyVarValues
static void copyVarValues(MeshBase &mesh, const unsigned int from_system, const unsigned int from_var, const NumericVector< Number > &from_vector, const unsigned int to_system, const unsigned int to_var, NumericVector< Number > &to_vector)
Helper function for copying values associated with variables in vectors from two different systems.
Definition: MoosePreconditioner.C:52
nl
NonlinearSystemBase * nl
Definition: PetscDMMoose.C:1267
MooseObject::type
const std::string & type() const
Get the type of this object.
Definition: MooseObject.h:63
PhysicsBasedPreconditioner::_pre_type
std::vector< PreconditionerType > _pre_type
Which preconditioner to use for each solve.
Definition: PhysicsBasedPreconditioner.h:86
PhysicsBasedPreconditioner::_nl
NonlinearSystemBase & _nl
The nonlinear system this PBP is associated with (convenience reference)
Definition: PhysicsBasedPreconditioner.h:78
ComputeJacobianBlocksThread.h
InputParameters::addParam
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an option parameter and a documentation string to the InputParameters object.
Definition: InputParameters.h:1198
MoosePreconditioner
Base class for MOOSE preconditioners.
Definition: MoosePreconditioner.h:34
SystemBase::nVariables
virtual unsigned int nVariables() const
Get the number of variables in this system.
Definition: SystemBase.C:734
defineLegacyParams
defineLegacyParams(PhysicsBasedPreconditioner)
PhysicsBasedPreconditioner::_init_timer
PerfID _init_timer
Timers.
Definition: PhysicsBasedPreconditioner.h:100
PhysicsBasedPreconditioner::~PhysicsBasedPreconditioner
virtual ~PhysicsBasedPreconditioner()
Definition: PhysicsBasedPreconditioner.C:146
PhysicsBasedPreconditioner::apply
virtual void apply(const NumericVector< Number > &x, NumericVector< Number > &y)
Computes the preconditioned vector "y" based on input "x".
Definition: PhysicsBasedPreconditioner.C:249
InputParameters
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system.
Definition: InputParameters.h:53
x
static PetscErrorCode Vec x
Definition: PetscDMMoose.C:1263
NonlinearSystem.h
registerMooseObjectAliased
registerMooseObjectAliased("MooseApp", PhysicsBasedPreconditioner, "PBP")
PhysicsBasedPreconditioner::addSystem
void addSystem(unsigned int var, std::vector< unsigned int > off_diag, PreconditionerType type=AMG_PRECOND)
Add a diagonal system + possibly off-diagonals ones as well, also specifying type of preconditioning.
Definition: PhysicsBasedPreconditioner.C:149
FEProblemBase::computeJacobianBlocks
virtual void computeJacobianBlocks(std::vector< JacobianBlock * > &blocks)
Computes several Jacobian blocks simultaneously, summing their contributions into smaller preconditio...
Definition: FEProblemBase.C:5276
PhysicsBasedPreconditioner::PhysicsBasedPreconditioner
PhysicsBasedPreconditioner(const InputParameters &params)
Definition: PhysicsBasedPreconditioner.C:59
PhysicsBasedPreconditioner::_off_diag_mats
std::vector< std::vector< SparseMatrix< Number > * > > _off_diag_mats
Holds pointers to the off-diagonal matrices.
Definition: PhysicsBasedPreconditioner.h:97
PhysicsBasedPreconditioner::_apply_timer
PerfID _apply_timer
Definition: PhysicsBasedPreconditioner.h:101
NonlinearSystemBase::system
virtual System & system() override
Get the reference to the libMesh system.
Definition: NonlinearSystemBase.h:618
MooseVariableFE.h
PhysicsBasedPreconditioner::clear
virtual void clear()
Release all memory and clear data structures.
Definition: PhysicsBasedPreconditioner.C:312
PhysicsBasedPreconditioner::validParams
static InputParameters validParams()
Constructor.
Definition: PhysicsBasedPreconditioner.C:36
FEProblemBase::solverParams
SolverParams & solverParams()
Get the solver parameters.
Definition: FEProblemBase.C:6313
PhysicsBasedPreconditioner::_preconditioners
std::vector< std::unique_ptr< Preconditioner< Number > > > _preconditioners
Holds one Preconditioner object per small system to solve.
Definition: PhysicsBasedPreconditioner.h:82
PhysicsBasedPreconditioner::init
virtual void init()
Initialize data structures if not done so already.
Definition: PhysicsBasedPreconditioner.C:177
PhysicsBasedPreconditioner
Implements a segregated solve preconditioner.
Definition: PhysicsBasedPreconditioner.h:33
MooseEnum.h
SystemBase::getVariableBlocks
virtual const std::set< SubdomainID > * getVariableBlocks(unsigned int var_number)
Get the block where a variable of this system is defined.
Definition: SystemBase.C:161
JacobianBlock
Helper class for holding the preconditioning blocks to fill.
Definition: ComputeJacobianBlocksThread.h:20
Moose::ST_JFNK
Jacobian-Free Newton Krylov.
Definition: MooseTypes.h:727
MooseMesh
MooseMesh wraps a libMesh::Mesh object and enhances its capabilities by caching additional data and s...
Definition: MooseMesh.h:74
FEProblemBase::getNonlinearSystemBase
NonlinearSystemBase & getNonlinearSystemBase()
Definition: FEProblemBase.h:560
FEProblemBase::es
virtual EquationSystems & es() override
Definition: FEProblemBase.h:147
FEProblemBase::setCouplingMatrix
void setCouplingMatrix(std::unique_ptr< CouplingMatrix > cm)
Set custom coupling matrix.
Definition: FEProblemBase.C:4462
NonlinearSystemBase::nonlinearSolver
virtual NonlinearSolver< Number > * nonlinearSolver()=0
NonlinearSystemBase
Nonlinear system to be solved.
Definition: NonlinearSystemBase.h:61
PhysicsBasedPreconditioner.h
InputParameters::addRequiredParam
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...
Definition: InputParameters.h:1176
FEProblemBase::mesh
virtual MooseMesh & mesh() override
Definition: FEProblemBase.h:148
PetscSupport.h
SolverParams::_type
Moose::SolveType _type
Definition: SolverParams.h:19
MooseVariableBase::number
unsigned int number() const
Get variable number coming from libMesh.
Definition: MooseVariableBase.h:48