From bc617adfdcfe07fa50425fc55eb64969adcb1b99 Mon Sep 17 00:00:00 2001 From: Oliver Sander <sander@igpm.rwth-aachen.de> Date: Thu, 12 Dec 2013 21:32:21 +0000 Subject: [PATCH] New hacky method to assemble the embedded Hessian into a SymmetricMatrix [[Imported from SVN: r9584]] --- dune/gfe/averagedistanceassembler.hh | 14 ++++++++++++++ 1 file changed, 14 insertions(+) diff --git a/dune/gfe/averagedistanceassembler.hh b/dune/gfe/averagedistanceassembler.hh index 3b0e7ccd..c426a6ba 100644 --- a/dune/gfe/averagedistanceassembler.hh +++ b/dune/gfe/averagedistanceassembler.hh @@ -3,6 +3,8 @@ #include <vector> +#include <dune/gfe/symmetricmatrix.hh> + /** \tparam TargetSpace The manifold that we are mapping to */ template <class TargetSpace, class WeightType=double> class AverageDistanceAssembler @@ -75,6 +77,18 @@ public: matrix.axpy(weights_[i], TargetSpace::secondDerivativeOfDistanceSquaredWRTSecondArgument(coefficients_[i], x)); } + void assembleEmbeddedHessian(const TargetSpace& x, + Dune::SymmetricMatrix<ctype,embeddedSize>& matrix) const + { + matrix = 0; + for (size_t i=0; i<coefficients_.size(); i++) { + Dune::FieldMatrix<ctype,embeddedSize,embeddedSize> tmp = TargetSpace::secondDerivativeOfDistanceSquaredWRTSecondArgument(coefficients_[i], x); + for (size_t j=0; j<embeddedSize; j++) + for (size_t k=0; k<=j; k++) + matrix(j,k) += weights_[i] * tmp[j][k]; + } + } + void assembleHessian(const TargetSpace& x, Dune::FieldMatrix<ctype,size,size>& matrix) const { -- GitLab