KASKADE 7 development version
base.hh
Go to the documentation of this file.
1/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
2/* */
3/* This file is part of the library KASKADE 7 */
4/* https://www.zib.de/research/projects/kaskade7-finite-element-toolbox */
5/* */
6/* Copyright (C) 2002-2012 Zuse Institute Berlin */
7/* */
8/* KASKADE 7 is distributed under the terms of the ZIB Academic License. */
9/* see $KASKADE/academic.txt */
10/* */
11/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
12
13#ifndef HERMITE_INTERPOLATION_BASE
14#define HERMITE_INTERPOLATION_BASE
15
16#include <vector>
17#include <dune/common/fvector.hh>
18
19namespace Kaskade
20{
29 template<class Scalar_, int dimension>
30 class HIPBase
31 {
32 public:
33 typedef Scalar_ Scalar;
34 static int const dim = dimension;
36
38
40 /*
41 * \param entity entity containing x
42 * \param x evaluation point in local coordinates
43 * \param shapeFunctionSet shape function set of the reference entity
44 * \result interpolation polynomial evaluated at x
45 */
46 template <class ShapeFunctionSet>
47 Vector evaluate(Vector const &x, ShapeFunctionSet const& shapeFunctionSet) const
48 {
49 Vector result(0);
50 for(int sfId=0; sfId<shapeFunctionSet.size(); ++sfId) result += values[sfId] * shapeFunctionSet[sfId].evaluateFunction( x )[0] * directions[sfId];
51 return result;
52 }
53
54 protected:
56 explicit HIPBase(int numberOfShapeFunctions)
57 {
58 values = std::vector<Scalar>(numberOfShapeFunctions);
59 directions = std::vector<Vector>(numberOfShapeFunctions,Vector(0));
60 }
61
63 void insertEntry(Scalar const value, Vector const& dir, int const id){
64 values[id] = value;
65 directions[id] = dir;
66 }
67
68 std::vector<Scalar> values;
69 std::vector<Vector> directions;
70 };
71
72
73 // class declaration
75 /*
76 * Only 2D and 3D implemented.
77 * In 2D and 3D for every vertex a "normal" must be given.
78 * In 3D also "normals" for the vertices are necessary.
79 * The given directions are interpreted as desired normals of the surface.
80 * The polynomials calculate the necessary deformations of the vertices. In 2D this is done exactly.
81 * In 3D the displacement is exact at the edges. On the faces it is approximated as deformation in direction
82 * of the actual normal of the face.
83 */
84 template<class Scalar, class GridView, class InterpolationPolicy, int dimension> class HIPImpl;
85} /* end of namespace Kaskade */
86#endif
Base class for hermite interpolation.
Definition: base.hh:31
std::vector< Vector > directions
Definition: base.hh:69
void insertEntry(Scalar const value, Vector const &dir, int const id)
Store value and direction associated with the id's shape function.
Definition: base.hh:63
Dune::FieldVector< Scalar, dim > Vector
Definition: base.hh:35
std::vector< Scalar > values
Definition: base.hh:68
HIPBase(int numberOfShapeFunctions)
Constructor.
Definition: base.hh:56
static int const dim
Definition: base.hh:34
Scalar_ Scalar
Definition: base.hh:33
Vector evaluate(Vector const &x, ShapeFunctionSet const &shapeFunctionSet) const
Evaluate interpolation polynomial at position x (in local coordinates).
Definition: base.hh:47
Implementation of interpolation polynomials from given vertex and/or edge directions ("normals").
Definition: base.hh:84
A set of shape functions.
virtual int size() const
Number of shape functions in the set.