37#ifndef __H__UG__LIB_DISC__LOCAL_SHAPE_FUNCTION_SET__NEDELEC__NEDELEC__
38#define __H__UG__LIB_DISC__LOCAL_SHAPE_FUNCTION_SET__NEDELEC__NEDELEC__
41#include "../local_dof_set.h"
51template <
typename TRefElem>
56 static const int refDim = TRefElem::dim;
72 nsh = rRefElem.num(1);
75 for(
size_t i = 0; i <
nsh; ++i)
112template <
typename TRefElement>
122 NedelecLSFS<ReferenceTriangle>,
123 ReferenceTriangle::dim,
124 MathVector<ReferenceTriangle::dim>,
125 MathMatrix<ReferenceTriangle::dim, ReferenceTriangle::dim>
133 static const size_t order = 1;
136 static const int dim = reference_element_type::dim;
151 static const size_t nsh = reference_element_type::numEdges;
161 inline size_t num_sh()
const {
return nsh;}
168 case 0: pos[0] = 0.5;
171 case 1: pos[0] = 0.5;
174 case 2: pos[0] = 0.0;
177 default:
UG_THROW(
"NedelecLSFS: shape function "<<i<<
178 " not found. Only "<<nsh<<
" shapes present.");
185 UG_THROW (
"NedelecLSFS: Nedelec shapes cannot be computed in the reference space.");
191 UG_THROW (
"NedelecLSFS: Gradients of the Nedelec shapes cannot be computed in the reference space.");
202 NedelecLSFS<ReferenceTetrahedron>,
203 ReferenceTetrahedron::dim,
204 MathVector<ReferenceTetrahedron::dim>,
205 MathMatrix<ReferenceTetrahedron::dim, ReferenceTetrahedron::dim>
213 static const size_t order = 1;
216 static const int dim = reference_element_type::dim;
231 static const size_t nsh = reference_element_type::numEdges;
241 inline size_t num_sh()
const {
return nsh;}
248 case 0: pos[0] = 0.5;
252 case 1: pos[0] = 0.5;
256 case 2: pos[0] = 0.0;
260 case 3: pos[0] = 0.0;
264 case 4: pos[0] = 0.5;
268 case 5: pos[0] = 0.0;
272 default:
UG_THROW(
"NedelecLSFS: shape function "<<i<<
273 " not found. Only "<<nsh<<
" shapes present.");
280 UG_THROW (
"NedelecLSFS: Nedelec shapes cannot be computed in the reference space.");
286 UG_THROW (
"NedelecLSFS: Gradients of the Nedelec shapes cannot be computed in the reference space.");
static interface for trial spaces
Definition local_shape_function_set.h:176
Definition local_dof_set.h:58
A class for fixed size, dense matrices.
Definition math_matrix.h:63
a mathematical Vector with N entries.
Definition math_vector.h:97
Nedelec, i.e. the edge local dof set.
Definition nedelec.h:53
size_t nsh
number of shapes (== number of edges)
Definition nedelec.h:100
size_t num_dof(ReferenceObjectID type) const
returns the number of DoFs on a sub-geometric object type
Definition nedelec.h:86
const LocalDoF & local_dof(size_t dof) const
returns the dof storage
Definition nedelec.h:93
bool exact_position_available() const
returns if the local dof position are exact
Definition nedelec.h:96
size_t num_dof() const
returns the total number of DoFs on the finite element
Definition nedelec.h:83
std::vector< LocalDoF > m_vLocalDoF
association to elements
Definition nedelec.h:103
static const int refDim
dimension of reference element
Definition nedelec.h:56
ReferenceObjectID roid() const
returns the type of reference element
Definition nedelec.h:80
NedelecLDS()
constructor
Definition nedelec.h:60
base_type::grad_type grad_type
Gradient type.
Definition nedelec.h:227
void grad(MathMatrix< dim, dim > &g, const size_t i, const MathVector< dim > &x) const
evaluates the gradient of the shape function
Definition nedelec.h:284
MathVector< dim > shape(const size_t i, const MathVector< dim > &x) const
evaluates the shape function
Definition nedelec.h:278
base_type::shape_type shape_type
Shape type.
Definition nedelec.h:224
NedelecLSFS()
Constructor.
Definition nedelec.h:235
ReferenceTetrahedron reference_element_type
Reference Element type.
Definition nedelec.h:210
bool position(size_t i, MathVector< dim > &pos) const
local position of DoF i
Definition nedelec.h:244
size_t num_sh() const
Definition nedelec.h:241
bool continuous() const
returns if space constructs continuous functions
Definition nedelec.h:238
BaseLSFS< NedelecLSFS< reference_element_type >, dim, MathVector< dim >, MathMatrix< dim, dim > > base_type
Base class.
Definition nedelec.h:220
NedelecLSFS()
Constructor.
Definition nedelec.h:155
MathVector< dim > shape(const size_t i, const MathVector< dim > &x) const
evaluates the shape function
Definition nedelec.h:183
size_t num_sh() const
Definition nedelec.h:161
ReferenceTriangle reference_element_type
Reference Element type.
Definition nedelec.h:130
base_type::grad_type grad_type
Gradient type.
Definition nedelec.h:147
bool continuous() const
returns if space constructs continuous functions
Definition nedelec.h:158
BaseLSFS< NedelecLSFS< reference_element_type >, dim, MathVector< dim >, MathMatrix< dim, dim > > base_type
Base class.
Definition nedelec.h:140
bool position(size_t i, MathVector< dim > &pos) const
local position of DoF i
Definition nedelec.h:164
base_type::shape_type shape_type
Shape type.
Definition nedelec.h:144
void grad(MathMatrix< dim, dim > &g, const size_t i, const MathVector< dim > &x) const
evaluates the gradient of the shape function
Definition nedelec.h:189
static TClass & get()
returns a singleton based on the identifier
Definition provider.h:54
Definition reference_element.h:494
Definition reference_element.h:394
#define UG_THROW(msg)
Definition error.h:57
ReferenceObjectID
these ids are used to identify the shape of a geometric object.
Definition grid_base_objects.h:74
int ReferenceElementDimension(ReferenceObjectID roid)
returns the reference element dimension at run-time
Definition reference_element_util.h:40