Skip to content
Snippets Groups Projects
Commit 2591b977 authored by Oliver Sander's avatar Oliver Sander Committed by sander@FU-BERLIN.DE
Browse files

Rewrite the interface. Test functions form a vector space.

Hence we want an interface more like a local basis, rather
than individual functions.

[[Imported from SVN: r7901]]
parent f28dd691
Branches
No related tags found
No related merge requests found
...@@ -31,17 +31,23 @@ public: ...@@ -31,17 +31,23 @@ public:
/** \brief Constructor /** \brief Constructor
*/ */
LocalGFETestFunction(const LocalFiniteElement& localFiniteElement, LocalGFETestFunction(const LocalFiniteElement& localFiniteElement,
const std::vector<TargetSpace>& baseCoefficients, const std::vector<TargetSpace>& baseCoefficients)
const std::vector<typename TargetSpace::TangentVector>& tangentCoefficients) : localGFEFunction_(localFiniteElement, baseCoefficients)
: localGFEFunction_(localFiniteElement, baseCoefficients),
tangentCoefficients_(tangentCoefficients)
{} {}
/** \brief The number of Lagrange points, NOT the number of basis functions */
unsigned int size() const
{
return localGFEFunction_.size();
}
/** \brief Evaluate the function */ /** \brief Evaluate all shape functions at the given point */
typename TargetSpace::EmbeddedTangentVector evaluate(const Dune::FieldVector<ctype, dim>& local) const; void evaluateFunction(const Dune::FieldVector<ctype, dim>& local,
std::vector<typename TargetSpace::EmbeddedTangentVector>& out) const;
/** \brief Evaluate the derivative of the function */ /** \brief Evaluate the derivatives of all shape functions function */
Dune::FieldMatrix<ctype, EmbeddedTangentVector::dimension, dim> evaluateDerivative(const Dune::FieldVector<ctype, dim>& local) const; void evaluateJacobian(const Dune::FieldVector<ctype, dim>& local,
std::vector<Dune::FieldMatrix<ctype, EmbeddedTangentVector::dimension, dim> >& out) const;
private: private:
...@@ -49,55 +55,60 @@ private: ...@@ -49,55 +55,60 @@ private:
*/ */
const LocalGeodesicFEFunction<dim,ctype,LocalFiniteElement,TargetSpace> localGFEFunction_; const LocalGeodesicFEFunction<dim,ctype,LocalFiniteElement,TargetSpace> localGFEFunction_;
std::vector<typename TargetSpace::TangentVector> tangentCoefficients_;
}; };
template <int dim, class ctype, class LocalFiniteElement, class TargetSpace> template <int dim, class ctype, class LocalFiniteElement, class TargetSpace>
typename TargetSpace::EmbeddedTangentVector LocalGFETestFunction<dim,ctype,LocalFiniteElement,TargetSpace>:: void LocalGFETestFunction<dim,ctype,LocalFiniteElement,TargetSpace>::evaluateFunction(const Dune::FieldVector<ctype, dim>& local,
evaluate(const Dune::FieldVector<ctype, dim>& local) const std::vector<typename TargetSpace::EmbeddedTangentVector>& out) const
{ {
typename TargetSpace::EmbeddedTangentVector result(0); out.resize(size() * spaceDim);
for (size_t i=0; i<tangentCoefficients_.size(); i++) { for (size_t i=0; i<size(); i++) {
Dune::FieldMatrix< double, embeddedDim, embeddedDim > derivative; Dune::FieldMatrix< double, embeddedDim, embeddedDim > derivative;
localGFEFunction_.evaluateDerivativeOfValueWRTCoefficient (local, localGFEFunction_.evaluateDerivativeOfValueWRTCoefficient (local,
i, i,
derivative); derivative);
Dune::FieldMatrix<ctype,spaceDim,embeddedDim> basisVectors = localGFEFunction_.coefficients_[i];
derivative.umv(tangentCoefficients_, result); for (int j=0; j<spaceDim; j++)
derivative.mv(basisVectors[j], out[i*spaceDim + j]);
} }
return result;
} }
template <int dim, class ctype, class LocalFiniteElement, class TargetSpace> template <int dim, class ctype, class LocalFiniteElement, class TargetSpace>
Dune::FieldMatrix<ctype, TargetSpace::EmbeddedTangentVector::dimension, dim> LocalGFETestFunction<dim,ctype,LocalFiniteElement,TargetSpace>:: void LocalGFETestFunction<dim,ctype,LocalFiniteElement,TargetSpace>::evaluateJacobian(const Dune::FieldVector<ctype, dim>& local,
evaluateDerivative(const Dune::FieldVector<ctype, dim>& local) const std::vector<Dune::FieldMatrix<ctype, EmbeddedTangentVector::dimension, dim> >& out) const
{ {
typename Dune::FieldMatrix<ctype, TargetSpace::EmbeddedTangentVector::dimension, dim> result(0); out.resize(size() * spaceDim);
for (size_t i=0; i<tangentCoefficients_.size(); i++) { for (size_t i=0; i<size(); i++) {
Tensor3< double, embeddedDim, embeddedDim, dim > derivative; Tensor3< double, embeddedDim, embeddedDim, dim > derivative;
localGFEFunction_.evaluateDerivativeOfGradientWRTCoefficient (local, localGFEFunction_.evaluateDerivativeOfGradientWRTCoefficient (local,
i, i,
derivative); derivative);
// Contract the second index of the derivative with the tangent vector at the i-th Lagrange point. Dune::FieldMatrix<ctype,spaceDim,embeddedDim> basisVectors = localGFEFunction_.coefficients_[i];
// Add that to the result.
for (size_t j=0; j<embeddedDim; j++) for (int j=0; j<spaceDim; j++) {
out[i*spaceDim + j] = 0;
// Contract the second index of the derivative with the tangent vector at the i-th Lagrange point.
// Add that to the result.
for (size_t k=0; k<embeddedDim; k++) for (size_t k=0; k<embeddedDim; k++)
for (size_t l=0; l<dim; l++) for (size_t l=0; l<embeddedDim; l++)
result[j][l] += derivative[j][k][l] * tangentCoefficients_[i][k]; for (size_t m=0; m<dim; m++)
out[i*spaceDim+j][k][m] += derivative[k][l][m] * basisVectors[j][l];
}
} }
return result;
} }
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment