Skip to content
Snippets Groups Projects
Commit b3d37f45 authored by Sander, Oliver's avatar Sander, Oliver
Browse files

Make the actual projection explicit in the LocalProjectedFEFunction

Previously, the actual projection happened by calling the constructor
of one of the classes that implements points on a manifold.  This patch
introduces an new method projectOnto, which makes the projection more
explicit.
parent 11a57898
No related branches found
No related tags found
No related merge requests found
...@@ -122,7 +122,7 @@ Dune::FieldMatrix< K, m, p > operator* ( const Dune::FieldMatrix< K, m, n > &A, ...@@ -122,7 +122,7 @@ Dune::FieldMatrix< K, m, p > operator* ( const Dune::FieldMatrix< K, m, n > &A,
for (size_t i=0; i<coefficients_.size(); i++) for (size_t i=0; i<coefficients_.size(); i++)
c.axpy(w[i][0], coefficients_[i].globalCoordinates()); c.axpy(w[i][0], coefficients_[i].globalCoordinates());
return TargetSpace(c); return TargetSpace::projectOnto(c);
} }
template <int dim, class ctype, class LocalFiniteElement, class TargetSpace> template <int dim, class ctype, class LocalFiniteElement, class TargetSpace>
......
...@@ -176,6 +176,15 @@ public: ...@@ -176,6 +176,15 @@ public:
return EmbeddedTangentVector(0); return EmbeddedTangentVector(0);
} }
/** \brief Projection from the embedding space onto the manifold
*
* For RealTuples this is simply the identity map
*/
static RealTuple<T,N> projectOnto(const CoordinateType& p)
{
return RealTuple<T,N>(p);
}
/** \brief Derivative of the projection from the embedding space onto the manifold /** \brief Derivative of the projection from the embedding space onto the manifold
* *
* For RealTuples this is simply the identity * For RealTuples this is simply the identity
......
...@@ -403,6 +403,13 @@ public: ...@@ -403,6 +403,13 @@ public:
return result; return result;
} }
static UnitVector<T,N> projectOnto(const CoordinateType& p)
{
UnitVector<T,N> result(p);
result.data_ /= result.data_.two_norm();
return result;
}
static DerivativeOfProjection derivativeOfProjection(const Dune::FieldVector<T,N>& p) static DerivativeOfProjection derivativeOfProjection(const Dune::FieldVector<T,N>& p)
{ {
Dune::FieldMatrix<T,N,N> result; Dune::FieldMatrix<T,N,N> result;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment