23 #include "libmesh/string_to_enum.h" 26 #include "libmesh/enum_solver_package.h" 29 #include "libmesh/mesh.h" 30 #include "libmesh/mesh_generation.h" 31 #include "libmesh/mesh_modification.h" 34 #include "libmesh/dense_matrix.h" 35 #include "libmesh/sparse_matrix.h" 36 #include "libmesh/dense_vector.h" 37 #include "libmesh/numeric_vector.h" 40 #include "libmesh/fe.h" 41 #include "libmesh/fe_interface.h" 42 #include "libmesh/elem.h" 45 #include "libmesh/quadrature_gauss.h" 46 #include "libmesh/enum_quadrature_type.h" 49 #include "libmesh/dof_map.h" 52 #include "libmesh/equation_systems.h" 53 #include "libmesh/nonlinear_implicit_system.h" 54 #include "libmesh/nonlinear_solver.h" 55 #include "libmesh/static_condensation.h" 58 #include "libmesh/exact_solution.h" 59 #include "libmesh/enum_norm_type.h" 63 #include "libmesh/getpot.h" 64 #include "libmesh/exodusII_io.h" 69 #if defined(LIBMESH_HAVE_EIGEN_DENSE) && defined(LIBMESH_HAVE_PETSC) 74 main(
int argc,
char ** argv)
81 "--enable-petsc, --enable-trilinos, or --enable-eigen");
84 GetPot infile(
"vector_fe_ex9.in");
87 infile.parse_command_line(argc, argv);
90 const unsigned int dimension = 2;
91 const unsigned int grid_size = infile(
"grid_size", 2);
92 const bool mms = infile(
"mms",
true);
93 const Real nu = infile(
"nu", 1.);
94 const bool cavity = infile(
"cavity",
false);
97 libmesh_example_requires(dimension <= LIBMESH_DIM, dimension <<
"D support");
106 const std::string elem_str = infile(
"element_type", std::string(
"TRI6"));
108 libmesh_error_msg_if(elem_str !=
"TRI6" && elem_str !=
"TRI7",
111 <<
" but this example must be run with TRI6, TRI7, QUAD8, or QUAD9 in 2d" 112 <<
" or with TET14, or HEX27 in 3d.");
116 mesh, grid_size, grid_size, 0., 2., -1, 1., Utility::string_to_enum<ElemType>(elem_str));
119 mesh, grid_size, grid_size, -1, 1, -1, 1., Utility::string_to_enum<ElemType>(elem_str));
128 Utility::string_to_enum<ElemType>(elem_str));
155 if (system.has_static_condensation())
157 sc = &system.get_static_condensation();
164 hdg.
dof_map = &system.get_dof_map();
175 system.nonlinear_solver->residual_object = &hdg;
176 system.nonlinear_solver->jacobian_object = &hdg;
181 equation_systems.
init();
200 exact_sol.compute_error(
"HDG",
"vel_x");
201 exact_sol.compute_error(
"HDG",
"vel_y");
202 exact_sol.compute_error(
"HDG",
"pressure");
205 libMesh::out <<
"L2 error for u is: " << exact_sol.l2_error(
"HDG",
"vel_x") << std::endl;
206 libMesh::out <<
"L2 error for v is: " << exact_sol.l2_error(
"HDG",
"vel_y") << std::endl;
207 libMesh::out <<
"L2 error for pressure is: " << exact_sol.l2_error(
"HDG",
"pressure")
211 #ifdef LIBMESH_HAVE_EXODUS_API 216 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
This class handles the computation of the L2 and/or H1 error for the Systems in the EquationSystems o...
std::unique_ptr< QBase > qface
This is the EquationSystems class.
int main(int argc, char **argv)
std::unique_ptr< FEVectorBase > vector_fe
std::unique_ptr< FEBase > lm_fe_face
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Number compute_error(const Point &p, const Parameters ¶ms, const std::string &, const std::string &unknown_name)
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
std::unique_ptr< FEBase > scalar_fe
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
The libMesh namespace provides an interface to certain functionality in the library.
void dont_condense_vars(const std::unordered_set< unsigned int > &vars)
Add vars to the list of variables not to condense.
SolverPackage default_solver_package()
virtual void write_equation_systems(const std::string &fname, const EquationSystems &es, const std::set< std::string > *system_names=nullptr) override
Writes out the solution for no specific time or timestep.
void print_info(std::ostream &os=libMesh::out, const unsigned int verbosity=0, const bool global=true) const
Prints relevant information about the mesh.
std::unique_ptr< QBase > qrule
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and non-linear solv...
std::unique_ptr< FEBase > scalar_fe_face
unsigned int add_variable(std::string_view var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
T & set(const std::string &)
Parameters parameters
Data structure holding arbitrary parameters.
static std::unique_ptr< QBase > build(std::string_view name, const unsigned int dim, const Order order=INVALID_ORDER)
Builds a specific quadrature rule based on the name string.
void attach_exact_value(unsigned int sys_num, FunctionBase< Number > *f)
Clone and attach an arbitrary functor which computes the exact value of the system sys_num solution a...
virtual void init()
Initialize all the systems.
virtual System & add_system(std::string_view system_type, std::string_view name)
Add the system of type system_type named name to the systems array.
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
std::unique_ptr< FEVectorBase > vector_fe_face