libMesh
linear_implicit_system.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2026 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 // C++ includes
21 
22 // Local includes
23 #include "libmesh/linear_implicit_system.h"
24 #include "libmesh/linear_solver.h"
25 #include "libmesh/equation_systems.h"
26 #include "libmesh/numeric_vector.h" // for parameter sensitivity calcs
27 //#include "libmesh/parameter_vector.h"
28 #include "libmesh/sparse_matrix.h" // for get_transpose
29 #include "libmesh/system_subset.h"
30 #include "libmesh/static_condensation.h"
31 #include "libmesh/static_condensation_preconditioner.h"
32 
33 namespace libMesh
34 {
35 
36 
38  const std::string & name_in,
39  const unsigned int number_in) :
40 
41  Parent (es, name_in, number_in),
42  _n_linear_iterations (0),
43  _final_linear_residual (1.e20),
44  _shell_matrix(nullptr),
45  _subset(nullptr),
46  _subset_solve_mode(SUBSET_ZERO)
47 {
48  // linear_solver is now in the ImplicitSystem base class, but we are
49  // going to keep using it basically the way we did before it was
50  // moved.
52 
53  if (this->has_static_condensation())
55 }
56 
57 
58 
60 
61 
62 
64 {
67 }
68 
69 
70 
72 {
73  // clear the linear solver
74  linear_solver->clear();
75 
76  this->restrict_solve_to(nullptr);
77 
78  // clear the parent data
79  Parent::clear();
80 
81  // And restore any StaticCondensation to defaults
82  if (this->has_static_condensation())
84 }
85 
86 
87 
89 {
90  // initialize parent data
92 
93  // re-initialize the linear solver interface
94  linear_solver->clear();
95 }
96 
97 
98 
100 {
101  // re-initialize the linear solver interface
102  linear_solver->clear();
103 
104  // initialize parent data
105  Parent::reinit();
106 }
107 
108 
109 
111  const SubsetSolveMode subset_solve_mode)
112 {
113  _subset = subset;
114  _subset_solve_mode = subset_solve_mode;
115 
116  if (subset != nullptr)
117  libmesh_assert_equal_to (&subset->get_system(), this);
118 }
119 
120 
121 
123 {
124  if (this->assemble_before_solve)
125  // Assemble the linear system
126  this->assemble ();
127 
128  // If the linear solver hasn't been initialized, we do so here.
129  if (this->prefix_with_name())
130  linear_solver->init(this->prefix().c_str());
131  else
132  linear_solver->init();
133 
134  linear_solver->init_systems(*this);
135 
136  // Get the user-specified linear solver tolerance
137  const auto [maxits, tol] = this->get_linear_solve_parameters();
138 
139  if (_subset != nullptr)
140  linear_solver->restrict_solve_to(&_subset->dof_ids(),_subset_solve_mode);
141 
142  // Solve the linear system. Several cases:
143  std::pair<unsigned int, Real> rval = std::make_pair(0,0.0);
144  if (_shell_matrix)
145  // 1.) Shell matrix with or without user-supplied preconditioner.
146  rval = linear_solver->solve(*_shell_matrix, this->request_matrix("Preconditioner"), *solution, *rhs, tol, maxits);
147  else
148  // 2.) No shell matrix, with or without user-supplied preconditioner
149  rval = linear_solver->solve (*matrix, this->request_matrix("Preconditioner"), *solution, *rhs, tol, maxits);
150 
151  if (_subset != nullptr)
152  linear_solver->restrict_solve_to(nullptr);
153 
154  // Store the number of linear iterations required to
155  // solve and the final residual.
156  _n_linear_iterations = rval.first;
157  _final_linear_residual = rval.second;
158 
159  // Update the system after the solve
160  this->update();
161 }
162 
163 
164 
166 {
167  _shell_matrix = shell_matrix;
168 }
169 
170 
171 /*
172  void LinearImplicitSystem::sensitivity_solve (const ParameterVector & parameters)
173  {
174  if (this->assemble_before_solve)
175  {
176  // Assemble the linear system
177  this->assemble ();
178 
179  // But now assemble right hand sides with the residual's
180  // parameter derivatives
181  this->assemble_residual_derivatives(parameters);
182  }
183 
184  // Get a reference to the EquationSystems
185  const EquationSystems & es =
186  this->get_equation_systems();
187 
188  // Get the user-specified linear solver tolerance
189  const Real tol =
190  es.parameters.get<Real>("sensitivity solver tolerance");
191 
192  // Get the user-specified maximum # of linear solver iterations
193  const unsigned int maxits =
194  es.parameters.get<unsigned int>("sensitivity solver maximum iterations");
195 
196  // Our iteration counts and residuals will be sums of the individual
197  // results
198  _n_linear_iterations = 0;
199  _final_linear_residual = 0.0;
200  std::pair<unsigned int, Real> rval = std::make_pair(0,0.0);
201 
202  // Solve the linear system.
203  SparseMatrix<Number> * pc = this->request_matrix("Preconditioner");
204  for (std::size_t p=0; p != parameters.size(); ++p)
205  {
206  rval = linear_solver->solve (*matrix, pc,
207  this->get_sensitivity_solution(p),
208  this->get_sensitivity_rhs(p), tol, maxits);
209  _n_linear_iterations += rval.first;
210  _final_linear_residual += rval.second;
211  }
212 
213  // Our matrix is the *negative* of the Jacobian for b-A*u, so our
214  // solutions are all inverted
215  for (std::size_t p=0; p != parameters.size(); ++p)
216  {
217  this->get_sensitivity_solution(p) *= -1.0;
218  }
219  }
220 
221 
222 
223  void LinearImplicitSystem::adjoint_solve (const QoISet & qoi_indices)
224  {
225  const unsigned int Nq = this->n_qois();
226 
227  // We currently don't support adjoint solves of shell matrices
228  // FIXME - we should let shell matrices support
229  // vector_transpose_mult so that we can use them here.
230  if (_shell_matrix!=nullptr)
231  libmesh_not_implemented();
232 
233  if (this->assemble_before_solve)
234  {
235  // Assemble the linear system
236  this->assemble ();
237 
238  // And take the adjoint
239  matrix->get_transpose(*matrix);
240 
241  // Including of any separate preconditioner
242  SparseMatrix<Number> * pc = this->request_matrix("Preconditioner");
243  if (pc)
244  pc->get_transpose(*pc);
245 
246  // But now replace the right hand sides with the quantity of
247  // interest functional's derivative
248  this->assemble_qoi_derivative(qoi_indices);
249  }
250 
251  // Get a reference to the EquationSystems
252  const EquationSystems & es =
253  this->get_equation_systems();
254 
255  // Get the user-specified linear solver tolerance
256  const Real tol =
257  es.parameters.get<Real>("adjoint solver tolerance");
258 
259  // Get the user-specified maximum # of linear solver iterations
260  const unsigned int maxits =
261  es.parameters.get<unsigned int>("adjoint solver maximum iterations");
262 
263  // Our iteration counts and residuals will be sums of the individual
264  // results
265  _n_linear_iterations = 0;
266  _final_linear_residual = 0.0;
267  std::pair<unsigned int, Real> rval = std::make_pair(0,0.0);
268 
269  // Solve the linear system.
270  SparseMatrix<Number> * pc = this->request_matrix("Preconditioner");
271  for (unsigned int i=0; i != Nq; ++i)
272  if (qoi_indices.has_index(i))
273  {
274  rval = linear_solver->solve (*matrix, pc,
275  this->add_adjoint_solution(i),
276  this->get_adjoint_rhs(i), tol, maxits);
277  _n_linear_iterations += rval.first;
278  _final_linear_residual += rval.second;
279  }
280  }
281 
282 
283 
284  void LinearImplicitSystem::forward_qoi_parameter_sensitivity
285  (const QoISet & qoi_indices,
286  const ParameterVector & parameters,
287  SensitivityData & sensitivities)
288  {
289  const unsigned int Np = parameters.size();
290  const unsigned int Nq = this->n_qois();
291 
292  // An introduction to the problem:
293  //
294  // A(p)*u(p) = b(p), where x is determined implicitly.
295  // Residual R(u(p),p) := b(p) - A(p)*u(p)
296  // partial R / partial u = -A
297  //
298  // This implies that:
299  // d/dp(R) = 0
300  // (partial b / partial p) -
301  // (partial A / partial p) * u -
302  // A * (partial u / partial p) = 0
303  // A * (partial u / partial p) = (partial R / partial p)
304  // = (partial b / partial p) - (partial A / partial p) * u
305 
306  // We first solve for (partial u / partial p) for each parameter:
307  // -A * (partial u / partial p) = - (partial R / partial p)
308 
309  this->sensitivity_solve(parameters);
310 
311  // Get ready to fill in sensitivities:
312  sensitivities.allocate_data(qoi_indices, *this, parameters);
313 
314  // We use the identity:
315  // dq/dp = (partial q / partial p) + (partial q / partial u) *
316  // (partial u / partial p)
317 
318  // We get (partial q / partial u) from the user
319  this->assemble_qoi_derivative(qoi_indices);
320 
321  for (unsigned int j=0; j != Np; ++j)
322  {
323  // We currently get partial derivatives via central differencing
324  Number delta_p = 1e-6;
325 
326  // (partial q / partial p) ~= (q(p+dp)-q(p-dp))/(2*dp)
327 
328  Number old_parameter = *parameters[j];
329 
330  *parameters[j] = old_parameter - delta_p;
331  this->assemble_qoi(qoi_indices);
332  std::vector<Number> qoi_minus = this->qoi;
333 
334  *parameters[j] = old_parameter + delta_p;
335  this->assemble_qoi(qoi_indices);
336  std::vector<Number> & qoi_plus = this->qoi;
337  std::vector<Number> partialq_partialp(Nq, 0);
338  for (unsigned int i=0; i != Nq; ++i)
339  if (qoi_indices.has_index(i))
340  partialq_partialp[i] = (qoi_plus[i] - qoi_minus[i]) / (2.*delta_p);
341 
342  for (unsigned int i=0; i != Nq; ++i)
343  if (qoi_indices.has_index(i))
344  sensitivities[i][j] = partialq_partialp[i] +
345  this->get_adjoint_rhs(i).dot(this->get_sensitivity_solution(i));
346  }
347 
348  // All parameters have been reset.
349  // Don't leave the qoi or system changed - principle of least
350  // surprise.
351  this->assemble();
352  this->rhs->close();
353  this->matrix->close();
354  this->assemble_qoi(qoi_indices);
355  }
356 */
357 
358 
359 
361 {
362  return linear_solver.get();
363 }
364 
365 
366 
368  bool,
369  bool,
370  bool)
371 {
372  // Residual R(u(p),p) := A(p)*u(p) - b(p)
373  // partial R / partial u = A
374 
375  this->assemble();
376  this->rhs->close();
377  this->matrix->close();
378 
379  *(this->rhs) *= -1.0;
380  this->rhs->add_vector(*this->solution, *this->matrix);
381 }
382 
383 } // namespace libMesh
virtual void restrict_solve_to(const SystemSubset *subset, const SubsetSolveMode subset_solve_mode=SUBSET_ZERO) override
After calling this method, any solve will be limited to the given subset.
void attach_shell_matrix(ShellMatrix< Number > *shell_matrix)
This function enables the user to provide a shell matrix, i.e.
virtual void assemble() override
Prepares matrix and _dof_map for matrix assembly.
virtual void create_static_condensation() override
Request that static condensation be performed for this system.
This is the EquationSystems class.
ShellMatrix< Number > * _shell_matrix
User supplies shell matrix or nullptr if no shell matrix is used.
static std::unique_ptr< LinearSolver< T > > build(const libMesh::Parallel::Communicator &comm_in, const SolverPackage solver_package=libMesh::default_solver_package())
Builds a LinearSolver using the linear solver package specified by solver_package.
Definition: linear_solver.C:59
unsigned int _n_linear_iterations
The number of linear iterations required to solve the linear system Ax=b.
virtual void create_static_condensation() override
Request that static condensation be performed for this system.
virtual std::pair< unsigned int, Real > get_linear_solve_parameters() const
virtual void reinit()
Reinitializes degrees of freedom and other required data on the current mesh.
Definition: system.C:442
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
Computes , where v is a pointer and each dof_indices[i] specifies where to add value v[i]...
virtual void clear() override
Clear all the data structures associated with the system.
virtual void init_data()
Initializes the data for the system.
Definition: system.C:207
LinearImplicitSystem(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
NumericVector< Number > * rhs
The system matrix.
const Parallel::Communicator & comm() const
virtual void reinit() override
Reinitializes the member data fields associated with the system, so that, e.g., assemble() may be use...
SubsetSolveMode _subset_solve_mode
If restrict-solve-to-subset mode is active, this member decides what happens with the dofs outside th...
virtual LinearSolver< Number > * get_linear_solver() const override
bool has_static_condensation() const
Definition: system.C:2664
The libMesh namespace provides an interface to certain functionality in the library.
virtual void solve() override
Assembles & solves the linear system A*x=b.
const SparseMatrix< Number > * request_matrix(std::string_view mat_name) const
Definition: system.C:1087
void setup_static_condensation_preconditioner(T &solver)
Sets up the static condensation preconditioner for the supplied solver.
virtual void init_data() override
Initializes new data members of the system.
virtual void assembly(bool get_residual, bool get_jacobian, bool apply_heterogeneous_constraints=false, bool apply_no_constraints=false) override
Assembles a residual in rhs and/or a jacobian in matrix, as requested.
std::string prefix() const
Definition: system.h:1980
This is a base class for classes which represent subsets of the dofs of a System. ...
Definition: system_subset.h:42
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:98
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1655
virtual const std::vector< unsigned int > & dof_ids() const =0
virtual void close()=0
Calls the NumericVector&#39;s internal assembly routines, ensuring that the values are consistent across ...
Real _final_linear_residual
The final residual for the linear system Ax=b.
virtual void update()
Update the local values to reflect the solution on neighboring processors.
Definition: system.C:498
virtual void close()=0
Calls the SparseMatrix&#39;s internal assembly routines, ensuring that the values are consistent across p...
SparseMatrix< Number > * matrix
The system matrix.
const System & get_system() const
Definition: system_subset.C:38
virtual void clear() override
Clear all the data structures associated with the system.
const SystemSubset * _subset
The current subset on which to solve (or nullptr if none).
bool assemble_before_solve
Flag which tells the system to whether or not to call the user assembly function during each call to ...
Definition: system.h:1609
std::unique_ptr< LinearSolver< Number > > linear_solver
This class handles all the details of interfacing with various linear algebra packages like PETSc or ...
SubsetSolveMode
defines an enum for the question what happens to the dofs outside the given subset when a system is s...
bool prefix_with_name() const
Definition: system.h:1974