From 361dfae4055ed047bf2cffd89f72ca3742d4fad2 Mon Sep 17 00:00:00 2001 From: Oliver Sander <sander@igpm.rwth-aachen.de> Date: Thu, 22 Sep 2011 14:59:43 +0000 Subject: [PATCH] remove some more now-obsolete template parameters [[Imported from SVN: r7838]] --- dune/gfe/averagedistanceassembler.hh | 12 ++++-------- dune/gfe/localgeodesicfefunction.hh | 20 +++++++++++++------- dune/gfe/targetspacertrsolver.cc | 10 +++++----- dune/gfe/targetspacertrsolver.hh | 7 +++---- 4 files changed, 25 insertions(+), 24 deletions(-) diff --git a/dune/gfe/averagedistanceassembler.hh b/dune/gfe/averagedistanceassembler.hh index 33319bc7..bd675e9e 100644 --- a/dune/gfe/averagedistanceassembler.hh +++ b/dune/gfe/averagedistanceassembler.hh @@ -3,8 +3,8 @@ #include <vector> -/** \tparam N Number of coefficients (i.e., simplex corners) */ -template <class TargetSpace, int N> +/** \tparam TargetSpace The manifold that we are mapping to */ +template <class TargetSpace> class AverageDistanceAssembler { static const int size = TargetSpace::TangentVector::dimension; @@ -12,12 +12,8 @@ class AverageDistanceAssembler public: - AverageDistanceAssembler(const Dune::array<TargetSpace,N>& coefficients, - const Dune::array<double,N>& weights) DUNE_DEPRECATED - : coefficients_(coefficients.begin(), coefficients.end()), - weights_(weights.begin(), weights.end()) - {} - + /** \brief Constructor with given coefficients \f$ v_i \f$ and weights \f$ w_i \f$ + */ AverageDistanceAssembler(const std::vector<TargetSpace>& coefficients, const std::vector<double>& weights) : coefficients_(coefficients), diff --git a/dune/gfe/localgeodesicfefunction.hh b/dune/gfe/localgeodesicfefunction.hh index 43ee577b..5f4f7299 100644 --- a/dune/gfe/localgeodesicfefunction.hh +++ b/dune/gfe/localgeodesicfefunction.hh @@ -156,10 +156,10 @@ evaluate(const Dune::FieldVector<ctype, dim>& local) const for (size_t i=0; i<w.size(); i++) w[i] = wNested[i][0]; - /** \todo The 'dim+1' parameter is a dummy here */ - AverageDistanceAssembler<TargetSpace,dim+1> assembler(coefficients_, w); + // The energy functional whose mimimizer is the value of the geodesic interpolation + AverageDistanceAssembler<TargetSpace> assembler(coefficients_, w); - TargetSpaceRiemannianTRSolver<TargetSpace,dim+1> solver; + TargetSpaceRiemannianTRSolver<TargetSpace> solver; solver.setup(&assembler, coefficients_[0], // initial iterate @@ -223,7 +223,7 @@ evaluateDerivative(const Dune::FieldVector<ctype, dim>& local) const for (size_t i=0; i<w.size(); i++) w[i] = wNested[i][0]; - AverageDistanceAssembler<TargetSpace,1> assembler(coefficients_, w); + AverageDistanceAssembler<TargetSpace> assembler(coefficients_, w); Dune::FieldMatrix<ctype,embeddedDim,embeddedDim> dFdq(0); assembler.assembleEmbeddedHessian(q,dFdq); @@ -314,8 +314,14 @@ evaluateDerivativeOfValueWRTCoefficient(const Dune::FieldVector<ctype, dim>& loc TargetSpace q = evaluate(local); // dFdq - Dune::array<ctype,dim+1> w = barycentricCoordinates(local); - AverageDistanceAssembler<TargetSpace,dim+1> assembler(coefficients_, w); + std::vector<Dune::FieldVector<ctype,1> > wNested; + localFiniteElement_.localBasis().evaluateFunction(local,wNested); + + std::vector<ctype> w(wNested.size()); + for (size_t i=0; i<w.size(); i++) + w[i] = wNested[i][0]; + + AverageDistanceAssembler<TargetSpace> assembler(coefficients_, w); Dune::FieldMatrix<ctype,embeddedDim,embeddedDim> dFdq(0); assembler.assembleEmbeddedHessian(q,dFdq); @@ -429,7 +435,7 @@ evaluateDerivativeOfGradientWRTCoefficient(const Dune::FieldVector<ctype, dim>& for (size_t i=0; i<w.size(); i++) w[i] = wNested[i][0]; - AverageDistanceAssembler<TargetSpace,dim+1> assembler(coefficients_, w); + AverageDistanceAssembler<TargetSpace> assembler(coefficients_, w); Dune::FieldMatrix<ctype,embeddedDim,embeddedDim> dFdq(0); assembler.assembleEmbeddedHessian(q,dFdq); diff --git a/dune/gfe/targetspacertrsolver.cc b/dune/gfe/targetspacertrsolver.cc index aa2a8ade..d59f2900 100644 --- a/dune/gfe/targetspacertrsolver.cc +++ b/dune/gfe/targetspacertrsolver.cc @@ -6,9 +6,9 @@ #include "maxnormtrustregion.hh" -template <class TargetSpace, int N> -void TargetSpaceRiemannianTRSolver<TargetSpace,N>:: -setup(const AverageDistanceAssembler<TargetSpace,N>* assembler, +template <class TargetSpace> +void TargetSpaceRiemannianTRSolver<TargetSpace>:: +setup(const AverageDistanceAssembler<TargetSpace>* assembler, const TargetSpace& x, double tolerance, int maxTrustRegionSteps, @@ -45,8 +45,8 @@ setup(const AverageDistanceAssembler<TargetSpace,N>* assembler, } -template <class TargetSpace, int N> -void TargetSpaceRiemannianTRSolver<TargetSpace,N>::solve() +template <class TargetSpace> +void TargetSpaceRiemannianTRSolver<TargetSpace>::solve() { MaxNormTrustRegion<blocksize> trustRegion(1, // we have only one block initialTrustRegionRadius_); diff --git a/dune/gfe/targetspacertrsolver.hh b/dune/gfe/targetspacertrsolver.hh index b721ae74..16e776a8 100644 --- a/dune/gfe/targetspacertrsolver.hh +++ b/dune/gfe/targetspacertrsolver.hh @@ -10,9 +10,8 @@ /** \brief Riemannian trust-region solver for geodesic finite-element problems \tparam TargetSpace The manifold that our functions take values in - \tparam N Number of simplex vertices */ -template <class TargetSpace, int N> +template <class TargetSpace> class TargetSpaceRiemannianTRSolver : public NumProc { @@ -34,7 +33,7 @@ public: {} /** \brief Set up the solver using a monotone multigrid method as the inner solver */ - void setup(const AverageDistanceAssembler<TargetSpace,N>* assembler, + void setup(const AverageDistanceAssembler<TargetSpace>* assembler, const TargetSpace& x, double tolerance, int maxTrustRegionSteps, @@ -71,7 +70,7 @@ protected: double innerTolerance_; /** \brief The assembler for the average-distance functional */ - const AverageDistanceAssembler<TargetSpace,N>* assembler_; + const AverageDistanceAssembler<TargetSpace>* assembler_; /** \brief The solver for the quadratic inner problems */ std::auto_ptr< ::LoopSolver<CorrectionType> > innerSolver_; -- GitLab