From b14c4cf547d198dd26052a6e1e7c22af8b94cda7 Mon Sep 17 00:00:00 2001
From: Oliver Sander <sander@igpm.rwth-aachen.de>
Date: Wed, 12 Oct 2011 15:45:42 +0000
Subject: [PATCH] simplify the code by making the AverageDistanceAssembler
 accept the interpolation weights in a std::vector of FieldVectors

[[Imported from SVN: r7885]]
---
 dune/gfe/averagedistanceassembler.hh | 18 +++++++++++++++-
 dune/gfe/localgeodesicfefunction.hh  | 32 +++++++---------------------
 2 files changed, 25 insertions(+), 25 deletions(-)

diff --git a/dune/gfe/averagedistanceassembler.hh b/dune/gfe/averagedistanceassembler.hh
index bd675e9e..17646d9a 100644
--- a/dune/gfe/averagedistanceassembler.hh
+++ b/dune/gfe/averagedistanceassembler.hh
@@ -20,6 +20,22 @@ public:
           weights_(weights)
     {}
 
+    /** \brief Constructor with given coefficients \f$ v_i \f$ and weights \f$ w_i \f$
+     * 
+     * The weights are given as a vector of length-1 FieldVectors instead of as a
+     * vector of doubles.  The reason is that these weights are actually the values
+     * of Lagrange shape functions, and the dune-localfunction interface returns
+     * shape function values this way.
+     */
+    AverageDistanceAssembler(const std::vector<TargetSpace>& coefficients,
+                             const std::vector<Dune::FieldVector<double,1> >& weights)
+        : coefficients_(coefficients),
+          weights_(weights.size())
+    {
+        for (size_t i=0; i<weights.size(); i++)
+            weights_[i] = weights[i][0];
+    }
+
     double value(const TargetSpace& x) const {
 
         double result = 0;
@@ -79,7 +95,7 @@ public:
 
     const std::vector<TargetSpace> coefficients_;
 
-    const std::vector<double> weights_;
+    std::vector<double> weights_;
 
 };
 
diff --git a/dune/gfe/localgeodesicfefunction.hh b/dune/gfe/localgeodesicfefunction.hh
index fd6b2765..d82c7772 100644
--- a/dune/gfe/localgeodesicfefunction.hh
+++ b/dune/gfe/localgeodesicfefunction.hh
@@ -139,12 +139,8 @@ TargetSpace LocalGeodesicFEFunction<dim,ctype,LocalFiniteElement,TargetSpace>::
 evaluate(const Dune::FieldVector<ctype, dim>& local) const
 {
     // Evaluate the weighting factors---these are the Lagrangian shape function values at 'local'
-    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];
+    std::vector<Dune::FieldVector<ctype,1> > w;
+    localFiniteElement_.localBasis().evaluateFunction(local,w);
 
     // The energy functional whose mimimizer is the value of the geodesic interpolation
     AverageDistanceAssembler<TargetSpace> assembler(coefficients_, w);
@@ -206,12 +202,8 @@ evaluateDerivative(const Dune::FieldVector<ctype, dim>& local) const
     Dune::Matrix<Dune::FieldMatrix<ctype,1,1> > RHS = dFdw * B;
 
     // the actual system matrix
-    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];
+    std::vector<Dune::FieldVector<ctype,1> > w;
+    localFiniteElement_.localBasis().evaluateFunction(local, w);
     
     AverageDistanceAssembler<TargetSpace> assembler(coefficients_, w);
     
@@ -304,12 +296,8 @@ evaluateDerivativeOfValueWRTCoefficient(const Dune::FieldVector<ctype, dim>& loc
     TargetSpace q = evaluate(local);
 
     // dFdq
-    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];
+    std::vector<Dune::FieldVector<ctype,1> > w;
+    localFiniteElement_.localBasis().evaluateFunction(local,w);
 
     AverageDistanceAssembler<TargetSpace> assembler(coefficients_, w);
     
@@ -418,12 +406,8 @@ evaluateDerivativeOfGradientWRTCoefficient(const Dune::FieldVector<ctype, dim>&
     Dune::Matrix<Dune::FieldMatrix<ctype,1,1> > dFdw = computeDFdw(q);
     
     // the actual system matrix
-    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];
+    std::vector<Dune::FieldVector<ctype,1> > w;
+    localFiniteElement_.localBasis().evaluateFunction(local,w);
 
     AverageDistanceAssembler<TargetSpace> assembler(coefficients_, w);
     
-- 
GitLab