Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

FE API deprecation/changes #4077

Open
wants to merge 18 commits into
base: devel
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
4 changes: 4 additions & 0 deletions include/base/dof_object.h
Original file line number Diff line number Diff line change
Expand Up @@ -181,12 +181,14 @@ class DofObject : public ReferenceCountedObject<DofObject>
*/
unique_id_type unique_id () const;

#ifdef LIBMESH_ENABLE_DEPRECATED
/**
* \returns The globally \p unique_id for this \p DofObject as a
* writable reference. Deprecated; use the API taking an input
* instead.
*/
unique_id_type & set_unique_id ();
#endif // LIBMESH_ENABLE_DEPRECATED

/**
* Sets the \p unique_id for this \p DofObject
Expand Down Expand Up @@ -851,6 +853,7 @@ unique_id_type DofObject::unique_id () const



#ifdef LIBMESH_ENABLE_DEPRECATED
inline
unique_id_type & DofObject::set_unique_id ()
{
Expand All @@ -861,6 +864,7 @@ unique_id_type & DofObject::set_unique_id ()
libmesh_not_implemented();
#endif
}
#endif // LIBMESH_ENABLE_DEPRECATED



Expand Down
112 changes: 112 additions & 0 deletions include/fe/fe.h
Original file line number Diff line number Diff line change
Expand Up @@ -414,6 +414,9 @@ class FE : public FEGenericBase<typename FEOutputType<T>::type>
* a finite element of type \p t and approximation order \p o.
*
* On a p-refined element, \p o should be the total order of the element.
*
* This method does not support all finite element types; e.g. the
* Polygon1 type may differ from element to element.
*/
static unsigned int n_shape_functions (const ElemType t,
const Order o)
Expand All @@ -424,29 +427,69 @@ class FE : public FEGenericBase<typename FEOutputType<T>::type>
* finite element.
*
* On a p-refined element, \p o should be the total order of the element.
*
* This method does not support all finite element types; e.g. the
* Polygon1 type may differ from element to element.
*/
static unsigned int n_dofs(const ElemType t,
const Order o);

/**
* \returns The number of shape functions associated with this
* finite element.
*
* On a p-refined element, \p o should be the total order of the element.
*
* \p e should only be a null pointer if using a FE family like
* SCALAR that has degrees of freedom independent of any element.
*/
static unsigned int n_dofs(const Elem * e,
const Order o);

/**
* \returns The number of dofs at node \p n for a finite element
* of type \p t and order \p o.
*
* On a p-refined element, \p o should be the total order of the element.
*
* This method does not support all finite element types; e.g. the
* Polygon1 type may differ from element to element.
*/
static unsigned int n_dofs_at_node(const ElemType t,
const Order o,
const unsigned int n);

/**
* \returns The number of dofs at node \p n for a finite element
* of type \p t and order \p o.
*
* On a p-refined element, \p o should be the total order of the element.
*/
static unsigned int n_dofs_at_node(const Elem & e,
const Order o,
const unsigned int n);

/**
* \returns The number of dofs interior to the element,
* not associated with any interior nodes.
*
* On a p-refined element, \p o should be the total order of the element.
*
* This method may not support all finite element types, e.g. higher
* order polygons or polyhedra may differ from element to element.
*/
static unsigned int n_dofs_per_elem(const ElemType t,
const Order o);

/**
* \returns The number of dofs interior to the element,
* not associated with any interior nodes.
*
* On a p-refined element, \p o should be the total order of the element.
*/
static unsigned int n_dofs_per_elem(const Elem & e,
const Order o);

/**
* \returns The continuity level of the finite element.
*/
Expand Down Expand Up @@ -1528,6 +1571,8 @@ void lagrange_nodal_soln(const Elem * elem,
*/
unsigned int monomial_n_dofs(const ElemType t, const Order o);

unsigned int monomial_n_dofs(const Elem * e, const Order o);

/**
* Helper functions for rational basis functions.
*/
Expand Down Expand Up @@ -1582,6 +1627,73 @@ void rational_all_shape_derivs (const Elem & elem,

} // namespace libMesh


// Full specialization of all n_dofs type functions, for every
// dimension, with both original ElemType and new Elem signatures
#define LIBMESH_DEFAULT_NDOFS(MyType) \
template <> unsigned int FE<0,MyType>::n_dofs(const ElemType t, const Order o) { return MyType##_n_dofs(t, o); } \
template <> unsigned int FE<1,MyType>::n_dofs(const ElemType t, const Order o) { return MyType##_n_dofs(t, o); } \
template <> unsigned int FE<2,MyType>::n_dofs(const ElemType t, const Order o) { return MyType##_n_dofs(t, o); } \
template <> unsigned int FE<3,MyType>::n_dofs(const ElemType t, const Order o) { return MyType##_n_dofs(t, o); } \
\
template <> unsigned int FE<0,MyType>::n_dofs(const Elem * e, const Order o) { return MyType##_n_dofs(e, o); } \
template <> unsigned int FE<1,MyType>::n_dofs(const Elem * e, const Order o) { return MyType##_n_dofs(e, o); } \
template <> unsigned int FE<2,MyType>::n_dofs(const Elem * e, const Order o) { return MyType##_n_dofs(e, o); } \
template <> unsigned int FE<3,MyType>::n_dofs(const Elem * e, const Order o) { return MyType##_n_dofs(e, o); } \
\
template <> unsigned int FE<0,MyType>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(t, o, n); } \
template <> unsigned int FE<1,MyType>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(t, o, n); } \
template <> unsigned int FE<2,MyType>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(t, o, n); } \
template <> unsigned int FE<3,MyType>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(t, o, n); } \
\
template <> unsigned int FE<0,MyType>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(e, o, n); } \
template <> unsigned int FE<1,MyType>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(e, o, n); } \
template <> unsigned int FE<2,MyType>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(e, o, n); } \
template <> unsigned int FE<3,MyType>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(e, o, n); } \
\
template <> unsigned int FE<0,MyType>::n_dofs_per_elem(const ElemType t, const Order o) { return MyType##_n_dofs_per_elem(t, o); } \
template <> unsigned int FE<1,MyType>::n_dofs_per_elem(const ElemType t, const Order o) { return MyType##_n_dofs_per_elem(t, o); } \
template <> unsigned int FE<2,MyType>::n_dofs_per_elem(const ElemType t, const Order o) { return MyType##_n_dofs_per_elem(t, o); } \
template <> unsigned int FE<3,MyType>::n_dofs_per_elem(const ElemType t, const Order o) { return MyType##_n_dofs_per_elem(t, o); } \
\
template <> unsigned int FE<0,MyType>::n_dofs_per_elem(const Elem & e, const Order o) { return MyType##_n_dofs_per_elem(e, o); } \
template <> unsigned int FE<1,MyType>::n_dofs_per_elem(const Elem & e, const Order o) { return MyType##_n_dofs_per_elem(e, o); } \
template <> unsigned int FE<2,MyType>::n_dofs_per_elem(const Elem & e, const Order o) { return MyType##_n_dofs_per_elem(e, o); } \
template <> unsigned int FE<3,MyType>::n_dofs_per_elem(const Elem & e, const Order o) { return MyType##_n_dofs_per_elem(e, o); }


#define LIBMESH_DEFAULT_VEC_NDOFS(MyType) \
template <> unsigned int FE<0,MyType##_VEC>::n_dofs(const ElemType t, const Order o) { return FE<0,MyType>::n_dofs(t, o); } \
template <> unsigned int FE<1,MyType##_VEC>::n_dofs(const ElemType t, const Order o) { return FE<1,MyType>::n_dofs(t, o); } \
template <> unsigned int FE<2,MyType##_VEC>::n_dofs(const ElemType t, const Order o) { return 2*FE<2,MyType>::n_dofs(t, o); } \
template <> unsigned int FE<3,MyType##_VEC>::n_dofs(const ElemType t, const Order o) { return 3*FE<3,MyType>::n_dofs(t, o); } \
\
template <> unsigned int FE<0,MyType##_VEC>::n_dofs(const Elem * e, const Order o) { return FE<0,MyType>::n_dofs(e, o); } \
template <> unsigned int FE<1,MyType##_VEC>::n_dofs(const Elem * e, const Order o) { return FE<1,MyType>::n_dofs(e, o); } \
template <> unsigned int FE<2,MyType##_VEC>::n_dofs(const Elem * e, const Order o) { return 2*FE<2,MyType>::n_dofs(e, o); } \
template <> unsigned int FE<3,MyType##_VEC>::n_dofs(const Elem * e, const Order o) { return 3*FE<3,MyType>::n_dofs(e, o); } \
\
template <> unsigned int FE<0,MyType##_VEC>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return FE<0,MyType>::n_dofs_at_node(t, o, n); } \
template <> unsigned int FE<1,MyType##_VEC>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return FE<1,MyType>::n_dofs_at_node(t, o, n); } \
template <> unsigned int FE<2,MyType##_VEC>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return 2*FE<2,MyType>::n_dofs_at_node(t, o, n); } \
template <> unsigned int FE<3,MyType##_VEC>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return 3*FE<3,MyType>::n_dofs_at_node(t, o, n); } \
\
template <> unsigned int FE<0,MyType##_VEC>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return FE<0,MyType>::n_dofs_at_node(e.type(), o, n); } \
template <> unsigned int FE<1,MyType##_VEC>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return FE<1,MyType>::n_dofs_at_node(e.type(), o, n); } \
template <> unsigned int FE<2,MyType##_VEC>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return 2*FE<2,MyType>::n_dofs_at_node(e.type(), o, n); } \
template <> unsigned int FE<3,MyType##_VEC>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return 3*FE<3,MyType>::n_dofs_at_node(e.type(), o, n); } \
\
template <> unsigned int FE<0,MyType##_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return FE<0,MyType>::n_dofs_per_elem(t, o); } \
template <> unsigned int FE<1,MyType##_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return FE<1,MyType>::n_dofs_per_elem(t, o); } \
template <> unsigned int FE<2,MyType##_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return 2*FE<2,MyType>::n_dofs_per_elem(t, o); } \
template <> unsigned int FE<3,MyType##_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return 3*FE<3,MyType>::n_dofs_per_elem(t, o); } \
\
template <> unsigned int FE<0,MyType##_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return FE<0,MyType>::n_dofs_per_elem(e.type(), o); } \
template <> unsigned int FE<1,MyType##_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return FE<1,MyType>::n_dofs_per_elem(e.type(), o); } \
template <> unsigned int FE<2,MyType##_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return 2*FE<2,MyType>::n_dofs_per_elem(e.type(), o); } \
template <> unsigned int FE<3,MyType##_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return 3*FE<3,MyType>::n_dofs_per_elem(e.type(), o); }


#define LIBMESH_DEFAULT_VECTORIZED_FE(MyDim, MyType) \
template<> \
void FE<MyDim,MyType>::all_shapes \
Expand Down
7 changes: 7 additions & 0 deletions include/fe/fe_abstract.h
Original file line number Diff line number Diff line change
Expand Up @@ -201,16 +201,23 @@ class FEAbstract : public ReferenceCountedObject<FEAbstract>
const std::vector<Point> & reference_side_points,
std::vector<Point> & reference_points) = 0;

#ifdef LIBMESH_ENABLE_DEPRECATED
/**
* \returns \p true if the point p is located on the reference element
* for element type t, false otherwise. Since we are doing floating
* point comparisons here the parameter \p eps can be specified to
* indicate a tolerance. For example, \f$ x \le 1 \f$ becomes
* \f$ x \le 1 + \epsilon \f$.
*
* This method overload is now deprecated, since it cannot support
* all finite element types; e.g. the Polygon1 type may differ from
* element to element. Use \p Elem::on_reference_element() instead.
*/
static bool on_reference_element(const Point & p,
const ElemType t,
const Real eps = TOLERANCE);
#endif // LIBMESH_ENABLE_DEPRECATED

/**
* \returns The reference space coordinates of \p nodes based on the
* element type.
Expand Down
Loading