libMesh
fem_system.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 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 
69  virtual ~FEMSystem ();
70 
75 
80 
92  virtual void assembly (bool get_residual,
93  bool get_jacobian,
94  bool apply_heterogeneous_constraints = false,
95  bool apply_no_constraints = false) override;
96 
105  virtual void solve () override;
106 
111  void mesh_position_get();
112 
117  void mesh_position_set();
118 
127  virtual std::unique_ptr<DiffContext> build_context() override;
128 
129  /*
130  * Prepares the result of a build_context() call for use.
131  *
132  * Most FEMSystem-based problems will need to reimplement this in order to
133  * call FE::get_*() as their particular physics requires.
134  */
135  virtual void init_context(DiffContext &) override;
136 
141  virtual void postprocess () override;
142 
151  virtual void assemble_qoi (const QoISet & indices = QoISet()) override;
152 
160  virtual void assemble_qoi_derivative (const QoISet & qoi_indices = QoISet(),
161  bool include_liftfunc = true,
162  bool apply_constraints = true) override;
163 
171 
182 
192  Real numerical_jacobian_h_for_var(unsigned int var_num) const;
193 
194  void set_numerical_jacobian_h_for_var(unsigned int var_num, Real new_h);
195 
210 
214  typedef bool (TimeSolver::*TimeSolverResPtr)(bool, DiffContext &);
215 
221  FEMContext & context) const;
222 
228  void numerical_elem_jacobian (FEMContext & context) const;
229 
235  void numerical_side_jacobian (FEMContext & context) const;
236 
242  void numerical_nonlocal_jacobian (FEMContext & context) const;
243 
244 protected:
249  virtual void init_data () override;
250 
251 private:
252  std::vector<Real> _numerical_jacobian_h_for_var;
253 };
254 
255 // --------------------------------------------------------------
256 // FEMSystem inline methods
257 inline
258 Real
259 FEMSystem::numerical_jacobian_h_for_var(unsigned int var_num) const
260 {
261  if ((var_num >= _numerical_jacobian_h_for_var.size()) ||
262  _numerical_jacobian_h_for_var[var_num] == Real(0))
263  return numerical_jacobian_h;
264 
265  return _numerical_jacobian_h_for_var[var_num];
266 }
267 
268 inline
270  Real new_h)
271 {
272  if (_numerical_jacobian_h_for_var.size() <= var_num)
273  _numerical_jacobian_h_for_var.resize(var_num+1,Real(0));
274 
275  libmesh_assert_greater(new_h, 0);
276 
277  _numerical_jacobian_h_for_var[var_num] = new_h;
278 }
279 
280 } // namespace libMesh
281 
282 
283 #endif // LIBMESH_FEM_SYSTEM_H
libMesh::FEMSystem::solve
virtual void solve() override
Invokes the solver associated with the system.
Definition: fem_system.C:1046
libMesh::FEMSystem::init_context
virtual void init_context(DiffContext &) override
Definition: fem_system.C:1338
libMesh::FEMSystem::assemble_qoi_derivative
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:1147
libMesh::FEMSystem::set_numerical_jacobian_h_for_var
void set_numerical_jacobian_h_for_var(unsigned int var_num, Real new_h)
Definition: fem_system.h:269
libMesh::FEMSystem::init_data
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:843
libMesh::FEMSystem::sys_type
FEMSystem sys_type
The type of system.
Definition: fem_system.h:74
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::System::number
unsigned int number() const
Definition: system.h:2075
libMesh::DifferentiableSystem
This class provides a specific system class.
Definition: diff_system.h:53
libMesh::FEMPhysics
This class provides a specific system class.
Definition: fem_physics.h:44
libMesh::FEMSystem::FEMSystem
FEMSystem(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
Definition: fem_system.C:826
libMesh::FEMSystem
This class provides a specific system class.
Definition: fem_system.h:53
libMesh::FEMSystem::numerical_jacobian_h
Real numerical_jacobian_h
If calculating numeric jacobians is required, the FEMSystem will perturb each solution vector entry b...
Definition: fem_system.h:181
libMesh::FEMSystem::numerical_jacobian
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:1178
libMesh::FEMSystem::Parent
DifferentiableSystem Parent
The type of the parent.
Definition: fem_system.h:79
libMesh::QoISet
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
Definition: qoi_set.h:45
libMesh::FEMSystem::postprocess
virtual void postprocess() override
Runs a postprocessing loop over all elements, and if postprocess_sides is true over all sides.
Definition: fem_system.C:1097
libMesh::FEMSystem::numerical_side_jacobian
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:1298
libMesh::FEMSystem::~FEMSystem
virtual ~FEMSystem()
Destructor.
Definition: fem_system.C:837
libMesh::FEMSystem::mesh_position_get
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:1385
libMesh::FEMSystem::assembly
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:850
libMesh::DiffContext
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::FEMSystem::build_context
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:1314
libMesh::FEMSystem::_numerical_jacobian_h_for_var
std::vector< Real > _numerical_jacobian_h_for_var
Definition: fem_system.h:252
libMesh::FEMSystem::mesh_position_set
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:1057
libMesh::FEMSystem::verify_analytic_jacobians
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:209
libMesh::System::name
const std::string & name() const
Definition: system.h:2067
libMesh::FEMSystem::fe_reinit_during_postprocess
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:170
libMesh::TimeSolver
This is a generic class that defines a solver to handle time integration of DifferentiableSystems.
Definition: time_solver.h:58
libMesh::FEMSystem::assemble_qoi
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:1117
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::FEMSystem::numerical_jacobian_h_for_var
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),...
Definition: fem_system.h:259
libMesh::FEMSystem::numerical_nonlocal_jacobian
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:1306
libMesh::FEMSystem::numerical_elem_jacobian
void numerical_elem_jacobian(FEMContext &context) const
Uses the results of multiple element_residual() calls to numerically differentiate the corresponding ...
Definition: fem_system.C:1290
libMesh::FEMSystem::TimeSolverResPtr
bool(TimeSolver::* TimeSolverResPtr)(bool, DiffContext &)
Syntax sugar to make numerical_jacobian() declaration easier.
Definition: fem_system.h:214
libMesh::FEMContext
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62