21 #include "libmesh/dof_map.h" 22 #include "libmesh/elem.h" 23 #include "libmesh/enum_to_string.h" 24 #include "libmesh/fe.h" 25 #include "libmesh/fe_interface.h" 26 #include "libmesh/fe_macro.h" 27 #include "libmesh/tensor_value.h" 42 void nedelec_one_nodal_soln(
const Elem * elem,
44 const std::vector<Number> & elem_soln,
47 std::vector<Number> & nodal_soln,
48 const bool add_p_level)
53 const Order totalorder = order + add_p_level*elem->
p_level();
55 nodal_soln.resize(
n_nodes*vdim);
59 if (elem_type !=
TRI6 && elem_type !=
TRI7 &&
63 libmesh_error_msg(
"ERROR: Invalid ElemType " <<
Utility::enum_to_string(elem_type) <<
" selected for NEDELEC_ONE FE family!");
67 std::vector<Point> refspace_nodes;
69 libmesh_assert_equal_to (refspace_nodes.size(),
n_nodes);
70 libmesh_assert_equal_to (elem_soln.size(), n_sf);
76 const std::vector<std::vector<RealGradient>> & vis_phi = vis_fe->get_phi();
78 vis_fe->reinit(elem,&refspace_nodes);
81 std::fill(nodal_soln.begin(), nodal_soln.end(), 0);
83 for (
unsigned int n = 0; n <
n_nodes; n++)
85 for (
unsigned int i=0; i<n_sf; i++)
86 for (
int d = 0; d < vdim; d++)
87 nodal_soln[vdim*n+d] += elem_soln[i]*(vis_phi[i][n](d));
94 unsigned int nedelec_one_n_dofs(
const ElemType t,
const Order o)
96 libmesh_assert_greater (o, 0);
106 libmesh_assert_less (o, 2);
107 libmesh_fallthrough();
109 return o*(o+2)*(o+3)/2;
111 libmesh_assert_less (o, 2);
112 libmesh_fallthrough();
114 return 3*o*(o+1)*(o+1);
118 libmesh_error_msg(
"ERROR: Invalid ElemType " <<
Utility::enum_to_string(t) <<
" selected for NEDELEC_ONE FE family!");
124 unsigned int nedelec_one_n_dofs(
const Elem * e,
const Order o)
127 return nedelec_one_n_dofs(e->
type(), o);
132 unsigned int nedelec_one_n_dofs_at_node(
const ElemType t,
134 const unsigned int n)
136 libmesh_assert_greater (o, 0);
153 libmesh_assert_equal_to(t,
TRI7);
156 libmesh_error_msg(
"ERROR: Invalid node ID " << n);
175 libmesh_assert_equal_to(t,
QUAD9);
178 libmesh_error_msg(
"ERROR: Invalid node ID " << n);
182 libmesh_assert_less (o, 2);
183 libmesh_fallthrough();
204 libmesh_assert_equal_to(t,
TET14);
207 libmesh_error_msg(
"ERROR: Invalid node ID " << n);
211 libmesh_assert_less (o, 2);
212 libmesh_fallthrough();
245 libmesh_assert_equal_to(t,
HEX27);
248 libmesh_assert_equal_to(t,
HEX27);
251 libmesh_error_msg(
"ERROR: Invalid node ID " << n);
257 libmesh_error_msg(
"ERROR: Invalid ElemType " <<
Utility::enum_to_string(t) <<
" selected for NEDELEC_ONE FE family!");
263 unsigned int nedelec_one_n_dofs_at_node(
const Elem & e,
265 const unsigned int n)
267 return nedelec_one_n_dofs_at_node(e.
type(), o, n);
272 unsigned int nedelec_one_n_dofs_per_elem(
const ElemType t,
275 libmesh_assert_greater (o, 0);
285 libmesh_assert_less (o, 2);
286 libmesh_fallthrough();
288 return o*(o-1)*(o-2)/2;
290 libmesh_assert_less (o, 2);
291 libmesh_fallthrough();
293 return 3*o*(o-1)*(o-1);
297 libmesh_error_msg(
"ERROR: Invalid ElemType " <<
Utility::enum_to_string(t) <<
" selected for NEDELEC_ONE FE family!");
303 unsigned int nedelec_one_n_dofs_per_elem(
const Elem & e,
306 return nedelec_one_n_dofs_per_elem(e.
type(), o);
311 #ifdef LIBMESH_ENABLE_AMR 315 const Elem * libmesh_dbg_var(elem),
324 libmesh_not_implemented();
326 #endif // #ifdef LIBMESH_ENABLE_AMR 330 #define NEDELEC_LOW_D_ERROR_MESSAGE \ 331 libmesh_error_msg("ERROR: This method makes no sense for low-D elements!"); 339 const std::vector<Number> &,
340 std::vector<Number> &,
343 { NEDELEC_LOW_D_ERROR_MESSAGE }
348 const std::vector<Number> &,
349 std::vector<Number> &,
352 { NEDELEC_LOW_D_ERROR_MESSAGE }
357 const std::vector<Number> & elem_soln,
358 std::vector<Number> & nodal_soln,
359 const bool add_p_level,
361 { nedelec_one_nodal_soln(elem, order, elem_soln, 2 , vdim, nodal_soln, add_p_level); }
366 const std::vector<Number> & elem_soln,
367 std::vector<Number> & nodal_soln,
368 const bool add_p_level,
370 { nedelec_one_nodal_soln(elem, order, elem_soln, 3 , 3 , nodal_soln, add_p_level); }
425 #ifdef LIBMESH_ENABLE_AMR 431 { NEDELEC_LOW_D_ERROR_MESSAGE }
438 { NEDELEC_LOW_D_ERROR_MESSAGE }
443 const unsigned int variable_number,
445 { nedelec_one_compute_constraints(constraints, dof_map, variable_number, elem, 2); }
450 const unsigned int variable_number,
452 { nedelec_one_compute_constraints(constraints, dof_map, variable_number, elem, 3); }
453 #endif // LIBMESH_ENABLE_AMR 458 { NEDELEC_LOW_D_ERROR_MESSAGE }
461 { NEDELEC_LOW_D_ERROR_MESSAGE }
464 const unsigned int,
const Point &)
465 { NEDELEC_LOW_D_ERROR_MESSAGE }
468 const unsigned int,
const Point &,
const bool)
469 { NEDELEC_LOW_D_ERROR_MESSAGE }
471 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 474 const unsigned int,
const Point &)
475 { NEDELEC_LOW_D_ERROR_MESSAGE }
478 const unsigned int,
const Point &,
const bool)
479 { NEDELEC_LOW_D_ERROR_MESSAGE }
485 { NEDELEC_LOW_D_ERROR_MESSAGE }
488 { NEDELEC_LOW_D_ERROR_MESSAGE }
491 const unsigned int,
const Point &)
492 { NEDELEC_LOW_D_ERROR_MESSAGE }
495 const unsigned int,
const Point &,
const bool)
496 { NEDELEC_LOW_D_ERROR_MESSAGE }
498 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 501 const unsigned int,
const Point &)
502 { NEDELEC_LOW_D_ERROR_MESSAGE }
505 const unsigned int,
const Point &,
const bool)
506 { NEDELEC_LOW_D_ERROR_MESSAGE }
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
static unsigned int n_dofs(const ElemType t, const Order o)
ElemType
Defines an enum for geometric element types.
Order
defines an enum for polynomial orders.
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
This is the base class from which all geometric element types are derived.
static unsigned int n_dofs_at_node(const ElemType t, const Order o, const unsigned int n)
static OutputShape shape_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
unsigned int p_level() const
virtual bool shapes_need_reinit() const override
The libMesh namespace provides an interface to certain functionality in the library.
virtual bool is_hierarchic() const override
This class handles the numbering of degrees of freedom on a mesh.
LIBMESH_DEFAULT_VECTORIZED_FE(template<>Real FE< 0, BERNSTEIN)
LIBMESH_FE_SIDE_NODAL_SOLN(HIERARCHIC_VEC)
const dof_id_type n_nodes
static unsigned int n_shape_functions(const unsigned int dim, const FEType &fe_t, const ElemType t)
virtual unsigned int n_nodes() const =0
static void get_refspace_nodes(const ElemType t, std::vector< Point > &nodes)
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
static unsigned int n_dofs_per_elem(const ElemType t, const Order o)
virtual FEContinuity get_continuity() const override
std::string enum_to_string(const T e)
static void compute_constraints(DofConstraints &constraints, DofMap &dof_map, const unsigned int variable_number, const Elem *elem)
Computes the constraint matrix contributions (for non-conforming adapted meshes) corresponding to var...
static void nodal_soln(const Elem *elem, const Order o, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln, bool add_p_level=true, const unsigned vdim=1)
Build the nodal soln from the element soln.
FEContinuity
defines an enum for finite element types to libmesh_assert a certain level (or type? Hcurl?) of continuity.
virtual ElemType type() const =0
A Point defines a location in LIBMESH_DIM dimensional Real space.
The constraint matrix storage format.
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)