libMesh
fem_system.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 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 #ifndef LIBMESH_FEM_SYSTEM_H
21 #define LIBMESH_FEM_SYSTEM_H
22 
23 // Local Includes
24 #include "libmesh/diff_system.h"
25 #include "libmesh/fem_physics.h"
26 
27 // C++ includes
28 #include <cstddef>
29 
30 namespace libMesh
31 {
32 
33 // Forward Declarations
34 class DiffContext;
35 class FEMContext;
36 
37 
54  public FEMPhysics
55 {
56 public:
57 
63  const std::string & name,
64  const unsigned int number);
65 
71  FEMSystem (const FEMSystem &) = delete;
72  FEMSystem & operator= (const FEMSystem &) = delete;
73  FEMSystem (FEMSystem &&) = delete;
74  FEMSystem & operator= (FEMSystem &&) = delete;
75  virtual ~FEMSystem ();
76 
81 
86 
98  virtual void assembly (bool get_residual,
99  bool get_jacobian,
100  bool apply_heterogeneous_constraints = false,
101  bool apply_no_constraints = false) override;
102 
111  virtual void solve () override;
112 
117  void mesh_position_get();
118 
123  void mesh_position_set();
124 
133  virtual std::unique_ptr<DiffContext> build_context() override;
134 
135  /*
136  * Prepares the result of a build_context() call for use.
137  *
138  * Most FEMSystem-based problems will need to reimplement this in order to
139  * call FE::get_*() as their particular physics requires.
140  */
141  virtual void init_context(DiffContext &) override;
142 
147  virtual void postprocess () override;
148 
157  virtual void assemble_qoi (const QoISet & indices = QoISet()) override;
158 
166  virtual void assemble_qoi_derivative (const QoISet & qoi_indices = QoISet(),
167  bool include_liftfunc = true,
168  bool apply_constraints = true) override;
169 
177 
188 
198  Real numerical_jacobian_h_for_var(unsigned int var_num) const;
199 
200  void set_numerical_jacobian_h_for_var(unsigned int var_num, Real new_h);
201 
216 
220  typedef bool (TimeSolver::*TimeSolverResPtr)(bool, DiffContext &);
221 
227  FEMContext & context) const;
228 
234  void numerical_elem_jacobian (FEMContext & context) const;
235 
241  void numerical_side_jacobian (FEMContext & context) const;
242 
248  void numerical_nonlocal_jacobian (FEMContext & context) const;
249 
250 protected:
255  virtual void init_data () override;
256 
257 private:
258  std::vector<Real> _numerical_jacobian_h_for_var;
259 };
260 
261 // --------------------------------------------------------------
262 // FEMSystem inline methods
263 inline
264 Real
265 FEMSystem::numerical_jacobian_h_for_var(unsigned int var_num) const
266 {
267  if ((var_num >= _numerical_jacobian_h_for_var.size()) ||
268  _numerical_jacobian_h_for_var[var_num] == Real(0))
269  return numerical_jacobian_h;
270 
271  return _numerical_jacobian_h_for_var[var_num];
272 }
273 
274 inline
276  Real new_h)
277 {
278  if (_numerical_jacobian_h_for_var.size() <= var_num)
279  _numerical_jacobian_h_for_var.resize(var_num+1,Real(0));
280 
281  libmesh_assert_greater(new_h, Real(0));
282 
283  _numerical_jacobian_h_for_var[var_num] = new_h;
284 }
285 
286 } // namespace libMesh
287 
288 
289 #endif // LIBMESH_FEM_SYSTEM_H
virtual std::unique_ptr< DiffContext > build_context() override
Builds a FEMContext object with enough information to do evaluations on each element.
Definition: fem_system.C:1346
This is the EquationSystems class.
Real verify_analytic_jacobians
If verify_analytic_jacobian is equal to zero (as it is by default), no numeric jacobians will be calc...
Definition: fem_system.h:215
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
virtual void assemble_qoi_derivative(const QoISet &qoi_indices=QoISet(), bool include_liftfunc=true, bool apply_constraints=true) override
Runs a qoi derivative assembly loop over all elements, and if assemble_qoi_sides is true over all sid...
Definition: fem_system.C:1179
FEMSystem sys_type
The type of system.
Definition: fem_system.h:80
This is a generic class that defines a solver to handle time integration of DifferentiableSystems.
Definition: time_solver.h:63
This class provides a specific system class.
Definition: fem_physics.h:44
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
Definition: qoi_set.h:45
virtual void init_context(DiffContext &) override
Definition: fem_system.C:1370
Real numerical_jacobian_h
If calculating numeric jacobians is required, the FEMSystem will perturb each solution vector entry b...
Definition: fem_system.h:187
void numerical_elem_jacobian(FEMContext &context) const
Uses the results of multiple element_residual() calls to numerically differentiate the corresponding ...
Definition: fem_system.C:1322
FEMSystem & operator=(const FEMSystem &)=delete
bool(TimeSolver::* TimeSolverResPtr)(bool, DiffContext &)
Syntax sugar to make numerical_jacobian() declaration easier.
Definition: fem_system.h:220
The libMesh namespace provides an interface to certain functionality in the library.
This class provides a specific system class.
Definition: fem_system.h:53
void set_numerical_jacobian_h_for_var(unsigned int var_num, Real new_h)
Definition: fem_system.h:275
This class provides a specific system class.
Definition: diff_system.h:54
unsigned int number() const
Definition: system.h:2350
bool fe_reinit_during_postprocess
If fe_reinit_during_postprocess is true (it is true by default), FE objects will be reinit()ed with t...
Definition: fem_system.h:176
void numerical_nonlocal_jacobian(FEMContext &context) const
Uses the results of multiple side_residual() calls to numerically differentiate the corresponding jac...
Definition: fem_system.C:1338
virtual void init_data() override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used...
Definition: fem_system.C:873
DifferentiableSystem Parent
The type of the parent.
Definition: fem_system.h:85
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62
FEMSystem(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
Definition: fem_system.C:858
virtual void postprocess() override
Runs a postprocessing loop over all elements, and if postprocess_sides is true over all sides...
Definition: fem_system.C:1127
void numerical_side_jacobian(FEMContext &context) const
Uses the results of multiple side_residual() calls to numerically differentiate the corresponding jac...
Definition: fem_system.C:1330
std::vector< Real > _numerical_jacobian_h_for_var
Definition: fem_system.h:258
virtual void assemble_qoi(const QoISet &indices=QoISet()) override
Runs a qoi assembly loop over all elements, and if assemble_qoi_sides is true over all sides...
Definition: fem_system.C:1147
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void mesh_position_set()
Tells the FEMSystem to set the mesh nodal coordinates which should correspond to degree of freedom co...
Definition: fem_system.C:1087
virtual void solve() override
Invokes the solver associated with the system.
Definition: fem_system.C:1076
virtual ~FEMSystem()
Real numerical_jacobian_h_for_var(unsigned int var_num) const
If numerical_jacobian_h_for_var(var_num) is changed from its default value (numerical_jacobian_h), the FEMSystem will perturb solution vector entries for variable var_num by that amount when calculating finite differences with respect to that variable.
Definition: fem_system.h:265
const std::string & name() const
Definition: system.h:2342
void numerical_jacobian(TimeSolverResPtr res, FEMContext &context) const
Uses the results of multiple res calls to numerically differentiate the corresponding jacobian...
Definition: fem_system.C:1210
virtual void assembly(bool get_residual, bool get_jacobian, bool apply_heterogeneous_constraints=false, bool apply_no_constraints=false) override
Prepares matrix or rhs for matrix assembly.
Definition: fem_system.C:880
void mesh_position_get()
Tells the FEMSystem to set the degree of freedom coefficients which should correspond to mesh nodal c...
Definition: fem_system.C:1423