18 #ifndef LIBMESH_VARIATIONAL_SMOOTHER_SYSTEM_H 19 #define LIBMESH_VARIATIONAL_SMOOTHER_SYSTEM_H 22 #include "libmesh/enum_fe_family.h" 23 #include "libmesh/fem_function_base.h" 24 #include "libmesh/fem_system.h" 25 #include "libmesh/libmesh_common.h" 49 const std::string &
name,
63 virtual void assembly (
bool get_residual,
65 bool apply_heterogeneous_constraints =
false,
66 bool apply_no_constraints =
false)
override;
100 #endif // LIBMESH_VARIATIONAL_SMOOTHER_SYSTEM_H virtual void init_data() override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used...
This is the EquationSystems class.
This class provides all data required for a physics package (e.g.
virtual void assembly(bool get_residual, bool get_jacobian, bool apply_heterogeneous_constraints=false, bool apply_no_constraints=false) override
Assembly method to update the mesh based on the smoother solve.
The libMesh namespace provides an interface to certain functionality in the library.
This class provides a specific system class.
Real _dilation_weight
The relative weight to give the dilation metric. The distortion metric is given weight 1 - _dilation_...
unsigned int number() const
Real & get_dilation_weight()
virtual bool element_time_derivative(bool request_jacobian, libMesh::DiffContext &context) override
Adds the time derivative contribution on elem to elem_residual.
Real _ref_vol
The reference volume for each element.
void compute_element_reference_volume()
const Real _epsilon_squared
The small nonzero constant to prevent zero denominators (degenerate elements only) ...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
VariationalSmootherSystem(libMesh::EquationSystems &es, const std::string &name, const unsigned int number)
This is an FEMSystem to solve the optimization probelem posed by the VariationalMeshSmoother class...
virtual void init_context(libMesh::DiffContext &context) override
const std::string & name() const
~VariationalSmootherSystem() override