Go to the documentation of this file.
20 #ifndef LIBMESH_FEM_SYSTEM_H
21 #define LIBMESH_FEM_SYSTEM_H
24 #include "libmesh/diff_system.h"
25 #include "libmesh/fem_physics.h"
63 const std::string &
name,
64 const unsigned int number);
92 virtual void assembly (
bool get_residual,
94 bool apply_heterogeneous_constraints =
false,
95 bool apply_no_constraints =
false)
override;
105 virtual void solve ()
override;
127 virtual std::unique_ptr<DiffContext>
build_context()
override;
161 bool include_liftfunc =
true,
162 bool apply_constraints =
true)
override;
275 libmesh_assert_greater(new_h, 0);
283 #endif // LIBMESH_FEM_SYSTEM_H
virtual void solve() override
Invokes the solver associated with the system.
virtual void init_context(DiffContext &) override
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...
void set_numerical_jacobian_h_for_var(unsigned int var_num, Real new_h)
virtual void init_data() override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used.
FEMSystem sys_type
The type of system.
The libMesh namespace provides an interface to certain functionality in the library.
unsigned int number() const
This class provides a specific system class.
This class provides a specific system class.
FEMSystem(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
This class provides a specific system class.
Real numerical_jacobian_h
If calculating numeric jacobians is required, the FEMSystem will perturb each solution vector entry b...
void numerical_jacobian(TimeSolverResPtr res, FEMContext &context) const
Uses the results of multiple res calls to numerically differentiate the corresponding jacobian.
DifferentiableSystem Parent
The type of the parent.
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
virtual void postprocess() override
Runs a postprocessing loop over all elements, and if postprocess_sides is true over all sides.
void numerical_side_jacobian(FEMContext &context) const
Uses the results of multiple side_residual() calls to numerically differentiate the corresponding jac...
virtual ~FEMSystem()
Destructor.
void mesh_position_get()
Tells the FEMSystem to set the degree of freedom coefficients which should correspond to mesh nodal c...
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.
This class provides all data required for a physics package (e.g.
This is the EquationSystems class.
virtual std::unique_ptr< DiffContext > build_context() override
Builds a FEMContext object with enough information to do evaluations on each element.
std::vector< Real > _numerical_jacobian_h_for_var
void mesh_position_set()
Tells the FEMSystem to set the mesh nodal coordinates which should correspond to degree of freedom co...
Real verify_analytic_jacobians
If verify_analytic_jacobian is equal to zero (as it is by default), no numeric jacobians will be calc...
const std::string & name() const
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...
This is a generic class that defines a solver to handle time integration of DifferentiableSystems.
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.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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),...
void numerical_nonlocal_jacobian(FEMContext &context) const
Uses the results of multiple side_residual() calls to numerically differentiate the corresponding jac...
void numerical_elem_jacobian(FEMContext &context) const
Uses the results of multiple element_residual() calls to numerically differentiate the corresponding ...
bool(TimeSolver::* TimeSolverResPtr)(bool, DiffContext &)
Syntax sugar to make numerical_jacobian() declaration easier.
This class provides all data required for a physics package (e.g.