33 #ifndef __H__UG__LIB_DISC__LOCAL_FINITE_ELEMENT__LOCAL_SHAPE_FUCNTION_SET__
34 #define __H__UG__LIB_DISC__LOCAL_FINITE_ELEMENT__LOCAL_SHAPE_FUCNTION_SET__
68 typename TGrad = MathVector<TDim> >
73 static const int dim = TDim;
127 virtual void shapes(std::vector<std::vector<shape_type> >& vvShape,
158 virtual void grads(std::vector<std::vector<grad_type> >& vvGrad,
172 template <
typename TImpl,
int TDim,
174 typename TGrad = MathVector<TDim> >
182 static const int dim = TDim;
204 getImpl().grad(g, i, x);
return;
222 for(
size_t sh = 0; sh <
getImpl().num_sh(); ++sh)
223 vShape[sh] =
shape(sh, x);
233 inline void shapes(std::vector<std::vector<shape_type> >& vvShape,
236 vvShape.resize(vLocPos.size());
237 for(
size_t ip = 0; ip < vLocPos.size(); ++ip)
238 shapes(vvShape[ip], vLocPos[ip]);
244 for(
size_t sh = 0; sh <
getImpl().num_sh(); ++sh)
245 grad(vGrad[sh], sh, x);
255 inline void grads(std::vector<std::vector<grad_type> >& vvGrad,
258 vvGrad.resize(vLocPos.size());
259 for(
size_t ip = 0; ip < vLocPos.size(); ++ip)
260 grads(vvGrad[ip], vLocPos[ip]);
static interface for trial spaces
Definition: local_shape_function_set.h:176
void grad(grad_type &g, size_t i, const MathVector< dim > &x) const
evaluates the gradient of the shape function
Definition: local_shape_function_set.h:202
void shapes(std::vector< shape_type > &vShape, const MathVector< dim > &x) const
returns all shape functions evaluated at a point
Definition: local_shape_function_set.h:227
shape_type shape(size_t i, const MathVector< dim > &x) const
evaluates the shape function
Definition: local_shape_function_set.h:196
TImpl ImplType
type of implementation
Definition: local_shape_function_set.h:179
ImplType & getImpl()
access to implementation
Definition: local_shape_function_set.h:265
void shape(shape_type &sh, size_t i, const MathVector< dim > &x) const
evaluates the shape function
Definition: local_shape_function_set.h:213
void grads(std::vector< std::vector< grad_type > > &vvGrad, const std::vector< MathVector< dim > > &vLocPos) const
returns all gradients evaluated at a point
Definition: local_shape_function_set.h:255
void shapes(shape_type *vShape, const MathVector< dim > &x) const
returns all shape functions evaluated at a point
Definition: local_shape_function_set.h:220
void grads(std::vector< grad_type > &vGrad, const MathVector< dim > &x) const
returns all gradients evaluated at a point
Definition: local_shape_function_set.h:249
TShape shape_type
Shape type.
Definition: local_shape_function_set.h:185
void shapes(std::vector< std::vector< shape_type > > &vvShape, const std::vector< MathVector< dim > > &vLocPos) const
returns all shape functions evaluated at a point
Definition: local_shape_function_set.h:233
static const int dim
dimension of reference element
Definition: local_shape_function_set.h:182
const ImplType & getImpl() const
const access to implementation
Definition: local_shape_function_set.h:268
void grads(grad_type *vGrad, const MathVector< dim > &x) const
returns all gradients evaluated at a point
Definition: local_shape_function_set.h:242
TGrad grad_type
Gradient type.
Definition: local_shape_function_set.h:188
Definition: local_dof_set.h:158
virtual base class for local shape function sets
Definition: local_shape_function_set.h:70
virtual void grads(std::vector< std::vector< grad_type > > &vvGrad, const std::vector< MathVector< dim > > &vLocPos) const =0
returns all gradients evaluated at a several points
TGrad grad_type
Gradient type.
Definition: local_shape_function_set.h:82
TShape shape_type
Shape type.
Definition: local_shape_function_set.h:79
MathVector< dim > position_type
Domain position type.
Definition: local_shape_function_set.h:76
virtual void grads(grad_type *vGrad, const MathVector< dim > &x) const =0
returns all gradients evaluated at a point
virtual shape_type shape(size_t i, const MathVector< dim > &x) const =0
evaluates the shape function
virtual void shapes(shape_type *vShape, const MathVector< dim > &x) const =0
returns all shape functions evaluated at a point
virtual void shapes(std::vector< shape_type > &vShape, const MathVector< dim > &x) const =0
returns all shape functions evaluated at a point
virtual void shapes(std::vector< std::vector< shape_type > > &vvShape, const std::vector< MathVector< dim > > &vLocPos) const =0
returns all shape functions evaluated at several point
virtual void grad(grad_type &g, size_t i, const MathVector< dim > &x) const =0
evaluates the gradient of the shape function
static const int dim
Dimension, where shape functions are defined.
Definition: local_shape_function_set.h:73
virtual bool continuous() const =0
returns if space constructs continuous functions
virtual void shape(shape_type &shape, size_t i, const MathVector< dim > &x) const =0
evaluates the shape function
virtual ~LocalShapeFunctionSet()
virtual destructor
Definition: local_shape_function_set.h:162
virtual void grads(std::vector< grad_type > &vGrad, const MathVector< dim > &x) const =0
returns all gradients evaluated at a point
double number
Definition: types.h:124