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"    46 void raviart_thomas_nodal_soln(
const Elem * elem,
    48                                const std::vector<Number> & elem_soln,
    51                                std::vector<Number> & nodal_soln,
    52                                const bool add_p_level)
    57   const Order totalorder = order + add_p_level*elem->
p_level();
    59   nodal_soln.resize(
n_nodes*vdim);
    63   if (elem_type != 
TRI6  && elem_type != 
TRI7  &&
    66     libmesh_error_msg(
"ERROR: Invalid ElemType " << 
Utility::enum_to_string(elem_type) << 
" selected for RAVIART_THOMAS FE family!");
    70   std::vector<Point> refspace_nodes;
    72   libmesh_assert_equal_to (refspace_nodes.size(), 
n_nodes);
    73   libmesh_assert_equal_to (elem_soln.size(), n_sf);
    79   const std::vector<std::vector<RealGradient>> & vis_phi = vis_fe->get_phi();
    81   vis_fe->reinit(elem,&refspace_nodes);
    84   std::fill(nodal_soln.begin(), nodal_soln.end(), 0);
    86   for (
unsigned int n = 0; n < 
n_nodes; n++)
    88     for (
unsigned int i=0; i<n_sf; i++)
    89       for (
int d = 0; d < vdim; d++)
    90         nodal_soln[vdim*n+d] += elem_soln[i]*(vis_phi[i][n](d));
    97 unsigned int raviart_thomas_n_dofs(
const ElemType t, 
const Order o)
    99   libmesh_assert_greater (o, 0);
   109       return o*(o+1)*(o+3)/2;
   115       libmesh_error_msg(
"ERROR: Invalid ElemType " << 
Utility::enum_to_string(t) << 
" selected for RAVIART_THOMAS FE family!");
   121 unsigned int raviart_thomas_n_dofs(
const Elem * e, 
const Order o)
   124   return raviart_thomas_n_dofs(e->
type(), o);
   129 unsigned int raviart_thomas_n_dofs_at_node(
const ElemType t,
   131                                            const unsigned int n)
   133   libmesh_assert_greater (o, 0);
   150             libmesh_assert_equal_to(t, 
TRI7);
   153             libmesh_error_msg(
"ERROR: Invalid node ID " << n);
   172             libmesh_assert_equal_to(t, 
QUAD9);
   175             libmesh_error_msg(
"ERROR: Invalid node ID " << n);
   199             libmesh_error_msg(
"ERROR: Invalid node ID " << n);
   237             libmesh_error_msg(
"ERROR: Invalid node ID " << n);
   243       libmesh_error_msg(
"ERROR: Invalid ElemType " << 
Utility::enum_to_string(t) << 
" selected for RAVIART_THOMAS FE family!");
   249 unsigned int raviart_thomas_n_dofs_at_node(
const Elem & e,
   251                                            const unsigned int n)
   253   return raviart_thomas_n_dofs_at_node(e.
type(), o, n);
   258 unsigned int raviart_thomas_n_dofs_per_elem(
const ElemType t,
   261   libmesh_assert_greater (o, 0);
   271       return (o+1)*o*(o-1)/2;
   277       libmesh_error_msg(
"ERROR: Invalid ElemType " << 
Utility::enum_to_string(t) << 
" selected for RAVIART_THOMAS FE family!");
   283 unsigned int raviart_thomas_n_dofs_per_elem(
const Elem & e,
   286   return raviart_thomas_n_dofs_per_elem(e.
type(), o);
   290 #ifdef LIBMESH_ENABLE_AMR   294                                          const Elem * libmesh_dbg_var(elem),
   303   libmesh_not_implemented();
   305 #endif // #ifdef LIBMESH_ENABLE_AMR   309 #define RAVIART_LOW_D_ERROR_MESSAGE                                     \   310   libmesh_error_msg("ERROR: This method makes no sense for low-D elements!");   318                                       const std::vector<Number> &,
   319                                       std::vector<Number> &,
   322 { RAVIART_LOW_D_ERROR_MESSAGE }
   327                                       const std::vector<Number> &,
   328                                       std::vector<Number> &,
   331 { RAVIART_LOW_D_ERROR_MESSAGE }
   336                                       const std::vector<Number> & elem_soln,
   337                                       std::vector<Number> & nodal_soln,
   338                                       const bool add_p_level,
   340 { raviart_thomas_nodal_soln(elem, order, elem_soln, 2 , vdim, nodal_soln, add_p_level); }
   345                                       const std::vector<Number> & elem_soln,
   346                                       std::vector<Number> & nodal_soln,
   347                                       const bool add_p_level,
   349 { raviart_thomas_nodal_soln(elem, order, elem_soln, 3 , 3 , nodal_soln, add_p_level); }
   354                                          const std::vector<Number> &,
   355                                          std::vector<Number> &,
   358 { RAVIART_LOW_D_ERROR_MESSAGE }
   363                                          const std::vector<Number> &,
   364                                          std::vector<Number> &,
   367 { RAVIART_LOW_D_ERROR_MESSAGE }
   372                                          const std::vector<Number> & elem_soln,
   373                                          std::vector<Number> & nodal_soln,
   374                                          const bool add_p_level,
   376 { raviart_thomas_nodal_soln(elem, order, elem_soln, 2 , vdim, nodal_soln, add_p_level); }
   381                                          const std::vector<Number> & elem_soln,
   382                                          std::vector<Number> & nodal_soln,
   383                                          const bool add_p_level,
   385 { raviart_thomas_nodal_soln(elem, order, elem_soln, 3 , 3 , nodal_soln, add_p_level); }
   489 #ifdef LIBMESH_ENABLE_AMR   495 { RAVIART_LOW_D_ERROR_MESSAGE }
   502 { RAVIART_LOW_D_ERROR_MESSAGE }
   507                                                 const unsigned int variable_number,
   509 { raviart_thomas_compute_constraints(constraints, dof_map, variable_number, elem, 2); }
   514                                                 const unsigned int variable_number,
   516 { raviart_thomas_compute_constraints(constraints, dof_map, variable_number, elem, 3); }
   523 { RAVIART_LOW_D_ERROR_MESSAGE }
   530 { RAVIART_LOW_D_ERROR_MESSAGE }
   535                                                    const unsigned int variable_number,
   537 { raviart_thomas_compute_constraints(constraints, dof_map, variable_number, elem, 2); }
   542                                                    const unsigned int variable_number,
   544 { raviart_thomas_compute_constraints(constraints, dof_map, variable_number, elem, 3); }
   545 #endif // LIBMESH_ENABLE_AMR   550 { RAVIART_LOW_D_ERROR_MESSAGE }
   553 { RAVIART_LOW_D_ERROR_MESSAGE }
   556 { RAVIART_LOW_D_ERROR_MESSAGE }
   559 { RAVIART_LOW_D_ERROR_MESSAGE }
   562                                                const unsigned int,
const Point &)
   563 { RAVIART_LOW_D_ERROR_MESSAGE }
   566                                                const unsigned int,
const Point &,
const bool)
   567 { RAVIART_LOW_D_ERROR_MESSAGE }
   570                                                   const unsigned int,
const Point &)
   571 { RAVIART_LOW_D_ERROR_MESSAGE }
   574                                                   const unsigned int,
const Point &,
const bool)
   575 { RAVIART_LOW_D_ERROR_MESSAGE }
   577 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES   580                                                       const unsigned int,
const Point &)
   581 { RAVIART_LOW_D_ERROR_MESSAGE }
   584                                                       const unsigned int,
const Point &,
const bool)
   585 { RAVIART_LOW_D_ERROR_MESSAGE }
   588                                                          const unsigned int,
const Point &)
   589 { RAVIART_LOW_D_ERROR_MESSAGE }
   592                                                          const unsigned int,
const Point &,
const bool)
   593 { RAVIART_LOW_D_ERROR_MESSAGE }
   599 { RAVIART_LOW_D_ERROR_MESSAGE }
   602 { RAVIART_LOW_D_ERROR_MESSAGE }
   605 { RAVIART_LOW_D_ERROR_MESSAGE }
   608 { RAVIART_LOW_D_ERROR_MESSAGE }
   611                                                const unsigned int,
const Point &)
   612 { RAVIART_LOW_D_ERROR_MESSAGE }
   615                                                const unsigned int,
const Point &,
const bool)
   616 { RAVIART_LOW_D_ERROR_MESSAGE }
   619                                                   const unsigned int,
const Point &)
   620 { RAVIART_LOW_D_ERROR_MESSAGE }
   623                                                   const unsigned int,
const Point &,
const bool)
   624 { RAVIART_LOW_D_ERROR_MESSAGE }
   626 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES   629                                                       const unsigned int,
const Point &)
   630 { RAVIART_LOW_D_ERROR_MESSAGE }
   633                                                       const unsigned int,
const Point &,
const bool)
   634 { RAVIART_LOW_D_ERROR_MESSAGE }
   637                                                          const unsigned int,
const Point &)
   638 { RAVIART_LOW_D_ERROR_MESSAGE }
   641                                                          const unsigned int,
const Point &,
const bool)
   642 { RAVIART_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)