From b3d37f45dc5d51fdc8913853d46343e7212d59e0 Mon Sep 17 00:00:00 2001 From: Oliver Sander <oliver.sander@tu-dresden.de> Date: Mon, 15 Jan 2018 13:49:32 +0100 Subject: [PATCH] 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. --- dune/gfe/localprojectedfefunction.hh | 2 +- dune/gfe/realtuple.hh | 9 +++++++++ dune/gfe/unitvector.hh | 7 +++++++ 3 files changed, 17 insertions(+), 1 deletion(-) diff --git a/dune/gfe/localprojectedfefunction.hh b/dune/gfe/localprojectedfefunction.hh index 31db2a19..583eee2e 100644 --- a/dune/gfe/localprojectedfefunction.hh +++ b/dune/gfe/localprojectedfefunction.hh @@ -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++) c.axpy(w[i][0], coefficients_[i].globalCoordinates()); - return TargetSpace(c); + return TargetSpace::projectOnto(c); } template <int dim, class ctype, class LocalFiniteElement, class TargetSpace> diff --git a/dune/gfe/realtuple.hh b/dune/gfe/realtuple.hh index 64bb0ee9..5d47b0d1 100644 --- a/dune/gfe/realtuple.hh +++ b/dune/gfe/realtuple.hh @@ -176,6 +176,15 @@ public: 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 * * For RealTuples this is simply the identity diff --git a/dune/gfe/unitvector.hh b/dune/gfe/unitvector.hh index 168d3153..51a6351e 100644 --- a/dune/gfe/unitvector.hh +++ b/dune/gfe/unitvector.hh @@ -403,6 +403,13 @@ public: 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) { Dune::FieldMatrix<T,N,N> result; -- GitLab