38 void nedelec_one_nodal_soln(
const Elem * elem,
40 const std::vector<Number> & elem_soln,
42 std::vector<Number> & nodal_soln)
44 const unsigned int n_nodes = elem->n_nodes();
45 const ElemType elem_type = elem->type();
47 const Order totalorder =
static_cast<Order>(order+elem->p_level());
61 libmesh_assert_equal_to (elem_soln.size(), 3);
62 libmesh_assert_equal_to (nodal_soln.size(), 6*2);
68 libmesh_assert_equal_to (elem_soln.size(), 4);
70 if (elem_type ==
QUAD8)
71 libmesh_assert_equal_to (nodal_soln.size(), 8*2);
73 libmesh_assert_equal_to (nodal_soln.size(), 9*2);
78 libmesh_assert_equal_to (elem_soln.size(), 6);
79 libmesh_assert_equal_to (nodal_soln.size(), 10*3);
81 libmesh_not_implemented();
90 libmesh_assert_equal_to (elem_soln.size(), 12);
92 if (elem_type ==
HEX20)
93 libmesh_assert_equal_to (nodal_soln.size(), 20*3);
95 libmesh_assert_equal_to (nodal_soln.size(), 27*3);
101 libmesh_error_msg(
"ERROR: Invalid ElemType " <<
Utility::enum_to_string(elem_type) <<
" selected for NEDELEC_ONE FE family!");
105 const unsigned int n_sf =
108 std::vector<Point> refspace_nodes;
110 libmesh_assert_equal_to (refspace_nodes.size(),
n_nodes);
117 const std::vector<std::vector<RealGradient>> & vis_phi = vis_fe->get_phi();
119 vis_fe->reinit(elem,&refspace_nodes);
121 for (
unsigned int n = 0; n <
n_nodes; n++)
123 libmesh_assert_equal_to (elem_soln.size(), n_sf);
126 for (
int d = 0; d < dim; d++)
128 nodal_soln[dim*n+d] = 0;
132 for (
unsigned int i=0; i<n_sf; i++)
134 for (
int d = 0; d < dim; d++)
136 nodal_soln[dim*n+d] += elem_soln[i]*(vis_phi[i][n](d));
145 libmesh_error_msg(
"ERROR: Invalid total order " <<
Utility::enum_to_string(totalorder) <<
" selected for NEDELEC_ONE FE family!");
153 unsigned int nedelec_one_n_dofs(
const ElemType t,
const Order o)
184 libmesh_error_msg(
"ERROR: Invalid Order " <<
Utility::enum_to_string(o) <<
" selected for NEDELEC_ONE FE family!");
191 unsigned int nedelec_one_n_dofs_at_node(
const ElemType t,
193 const unsigned int n)
215 libmesh_error_msg(
"ERROR: Invalid node ID " << n);
234 libmesh_error_msg(
"ERROR: Invalid node ID " << n);
254 libmesh_error_msg(
"ERROR: Invalid node ID " << n);
275 libmesh_error_msg(
"ERROR: Invalid node ID " << n);
307 libmesh_error_msg(
"ERROR: Invalid node ID " << n);
345 libmesh_error_msg(
"ERROR: Invalid node ID " << n);
358 libmesh_error_msg(
"ERROR: Invalid Order " <<
Utility::enum_to_string(o) <<
" selected for NEDELEC_ONE FE family!");
364 #ifdef LIBMESH_ENABLE_AMR 365 void nedelec_one_compute_constraints (DofConstraints & ,
368 const Elem * libmesh_dbg_var(elem),
375 libmesh_assert(elem);
377 libmesh_not_implemented();
480 #endif // #ifdef LIBMESH_ENABLE_AMR 484 #define NEDELEC_LOW_D_ERROR_MESSAGE \ 485 libmesh_error_msg("ERROR: This method makes no sense for low-D elements!"); 493 const std::vector<Number> &,
494 std::vector<Number> &)
495 { NEDELEC_LOW_D_ERROR_MESSAGE }
500 const std::vector<Number> &,
501 std::vector<Number> &)
502 { NEDELEC_LOW_D_ERROR_MESSAGE }
507 const std::vector<Number> & elem_soln,
508 std::vector<Number> & nodal_soln)
509 { nedelec_one_nodal_soln(elem, order, elem_soln, 2 , nodal_soln); }
514 const std::vector<Number> & elem_soln,
515 std::vector<Number> & nodal_soln)
516 { nedelec_one_nodal_soln(elem, order, elem_soln, 3 , nodal_soln); }
561 #ifdef LIBMESH_ENABLE_AMR 567 { NEDELEC_LOW_D_ERROR_MESSAGE }
574 { NEDELEC_LOW_D_ERROR_MESSAGE }
579 const unsigned int variable_number,
581 { nedelec_one_compute_constraints(constraints, dof_map, variable_number, elem, 2); }
586 const unsigned int variable_number,
588 { nedelec_one_compute_constraints(constraints, dof_map, variable_number, elem, 3); }
589 #endif // LIBMESH_ENABLE_AMR 594 { NEDELEC_LOW_D_ERROR_MESSAGE }
597 { NEDELEC_LOW_D_ERROR_MESSAGE }
600 const unsigned int,
const Point &)
601 { NEDELEC_LOW_D_ERROR_MESSAGE }
604 const unsigned int,
const Point &)
605 { NEDELEC_LOW_D_ERROR_MESSAGE }
608 const unsigned int,
const Point &)
609 { NEDELEC_LOW_D_ERROR_MESSAGE }
612 const unsigned int,
const Point &)
613 { NEDELEC_LOW_D_ERROR_MESSAGE }
617 { NEDELEC_LOW_D_ERROR_MESSAGE }
620 { NEDELEC_LOW_D_ERROR_MESSAGE }
623 const unsigned int,
const Point &)
624 { NEDELEC_LOW_D_ERROR_MESSAGE }
627 const unsigned int,
const Point &)
628 { NEDELEC_LOW_D_ERROR_MESSAGE }
631 const unsigned int,
const Point &)
632 { NEDELEC_LOW_D_ERROR_MESSAGE }
635 const unsigned int,
const Point &)
636 { NEDELEC_LOW_D_ERROR_MESSAGE }
static unsigned int n_dofs(const ElemType t, const Order o)
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
The base class for all geometric element types.
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)
virtual bool shapes_need_reinit() const override
virtual bool is_hierarchic() const override
Manages the degrees of freedom (DOFs) in a simulation.
const dof_id_type n_nodes
static unsigned int n_shape_functions(const unsigned int dim, const FEType &fe_t, const ElemType t)
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)
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)
static void nodal_soln(const Elem *elem, const Order o, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln)
A geometric point in (x,y,z) space.
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)