Newer
Older
#ifndef HARMONIC_ENERGY_LOCAL_STIFFNESS_HH
#define HARMONIC_ENERGY_LOCAL_STIFFNESS_HH
#include <dune/common/fmatrix.hh>
#include <dune/geometry/quadraturerules.hh>
#include "localgeodesicfestiffness.hh"
#include "localgeodesicfefunction.hh"
#ifdef HARMONIC_ENERGY_FD_GRADIENT
#warning Finite-difference approximation of the energy gradient
#endif
template<class GridView, class LocalFiniteElement, class TargetSpace>
: public LocalGeodesicFEStiffness<GridView,LocalFiniteElement,TargetSpace>
{
// grid types
typedef typename GridView::Grid::ctype DT;
typedef typename TargetSpace::ctype RT;
typedef typename GridView::template Codim<0>::Entity Entity;
// some other sizes
enum { blocksize = TargetSpace::TangentVector::dimension };
/** \brief Assemble the energy for a single element */
const LocalFiniteElement& localFiniteElement,
const std::vector<TargetSpace>& localSolution) const;
#ifndef HARMONIC_ENERGY_FD_GRADIENT
// The finite difference gradient method is in the base class.
// If the cpp macro is not set we overload it here.

Oliver Sander
committed
/** \brief Assemble the gradient of the energy functional on one element */
virtual void assembleEmbeddedGradient(const Entity& element,
const LocalFiniteElement& localFiniteElement,
const std::vector<TargetSpace>& solution,

Oliver Sander
committed
std::vector<typename TargetSpace::EmbeddedTangentVector>& gradient) const;

Oliver Sander
committed
virtual void assembleGradient(const Entity& element,
const LocalFiniteElement& localFiniteElement,
const std::vector<TargetSpace>& localSolution,

Oliver Sander
committed
std::vector<typename TargetSpace::TangentVector>& localGradient) const;
#endif
template <class GridView, class LocalFiniteElement, class TargetSpace>
typename HarmonicEnergyLocalStiffness<GridView, LocalFiniteElement, TargetSpace>::RT
HarmonicEnergyLocalStiffness<GridView, LocalFiniteElement, TargetSpace>::
const LocalFiniteElement& localFiniteElement,
const std::vector<TargetSpace>& localSolution) const
assert(element.type() == localFiniteElement.type());
RT energy = 0;
LocalGeodesicFEFunction<gridDim, double, LocalFiniteElement, TargetSpace> localGeodesicFEFunction(localFiniteElement,
localSolution);
int quadOrder = (element.type().isSimplex()) ? (localFiniteElement.localBasis().order()-1) * 2
: (localFiniteElement.localBasis().order()-1) * 2 * gridDim;
const Dune::QuadratureRule<double, gridDim>& quad
= Dune::QuadratureRules<double, gridDim>::rule(element.type(), quadOrder);
for (size_t pt=0; pt<quad.size(); pt++) {
const Dune::FieldVector<double,gridDim>& quadPos = quad[pt].position();
const double integrationElement = element.geometry().integrationElement(quadPos);
const Dune::FieldMatrix<double,gridDim,gridDim>& jacobianInverseTransposed = element.geometry().jacobianInverseTransposed(quadPos);
double weight = quad[pt].weight() * integrationElement;
// The derivative of the local function defined on the reference element
Dune::FieldMatrix<double, TargetSpace::EmbeddedTangentVector::dimension, gridDim> referenceDerivative = localGeodesicFEFunction.evaluateDerivative(quadPos);
// The derivative of the function defined on the actual element
Dune::FieldMatrix<double, TargetSpace::EmbeddedTangentVector::dimension, gridDim> derivative(0);
for (size_t comp=0; comp<referenceDerivative.N(); comp++)
jacobianInverseTransposed.umv(referenceDerivative[comp], derivative[comp]);
// Add the local energy density
// The Frobenius norm is the correct norm here if the metric of TargetSpace is the identity.
// (And if the metric of the domain space is the identity, which it always is here.)
energy += weight * derivative.frobenius_norm2();

Oliver Sander
committed
#ifndef HARMONIC_ENERGY_FD_GRADIENT
template <class GridView, class LocalFiniteElement, class TargetSpace>
void HarmonicEnergyLocalStiffness<GridView, LocalFiniteElement, TargetSpace>::

Oliver Sander
committed
assembleEmbeddedGradient(const Entity& element,
const LocalFiniteElement& localFiniteElement,
const std::vector<TargetSpace>& localSolution,
std::vector<typename TargetSpace::EmbeddedTangentVector>& localGradient) const

Oliver Sander
committed
{
// initialize gradient
localGradient.resize(localSolution.size());

Oliver Sander
committed
std::fill(localGradient.begin(), localGradient.end(), typename TargetSpace::EmbeddedTangentVector(0));

Oliver Sander
committed
// Set up local gfe function from the local coefficients
LocalGeodesicFEFunction<gridDim, double, LocalFiniteElement, TargetSpace> localGeodesicFEFunction(localFiniteElement,
localSolution);

Oliver Sander
committed
// I am not sure about the correct quadrature order
int quadOrder = (element.type().isSimplex()) ? (localFiniteElement.localBasis().order()-1) * 2
: (localFiniteElement.localBasis().order()-1) * 2 * gridDim;

Oliver Sander
committed
// numerical quadrature loop
const Dune::QuadratureRule<double, gridDim>& quad
= Dune::QuadratureRules<double, gridDim>::rule(element.type(), quadOrder);
for (size_t pt=0; pt<quad.size(); pt++) {
// Local position of the quadrature point
const Dune::FieldVector<double,gridDim>& quadPos = quad[pt].position();
const double integrationElement = element.geometry().integrationElement(quadPos);
const Dune::FieldMatrix<double,gridDim,gridDim>& jacobianInverseTransposed = element.geometry().jacobianInverseTransposed(quadPos);
double weight = quad[pt].weight() * integrationElement;
// The derivative of the local function defined on the reference element
Dune::FieldMatrix<double, TargetSpace::EmbeddedTangentVector::dimension, gridDim> referenceDerivative = localGeodesicFEFunction.evaluateDerivative(quadPos);

Oliver Sander
committed
// The derivative of the function defined on the actual element
Dune::FieldMatrix<double, TargetSpace::EmbeddedTangentVector::dimension, gridDim> derivative;

Oliver Sander
committed
for (size_t comp=0; comp<referenceDerivative.N(); comp++)
jacobianInverseTransposed.mv(referenceDerivative[comp], derivative[comp]);
// loop over all the element's degrees of freedom and compute the gradient wrt it
for (size_t i=0; i<localSolution.size(); i++) {
Tensor3<double, TargetSpace::EmbeddedTangentVector::dimension,TargetSpace::EmbeddedTangentVector::dimension,gridDim> referenceDerivativeDerivative;
#ifdef HARMONIC_ENERGY_FD_INNER_GRADIENT
#warning Using finite differences to compute the inner gradients!
localGeodesicFEFunction.evaluateFDDerivativeOfGradientWRTCoefficient(quadPos, i, referenceDerivativeDerivative);
#else

Oliver Sander
committed
localGeodesicFEFunction.evaluateDerivativeOfGradientWRTCoefficient(quadPos, i, referenceDerivativeDerivative);

Oliver Sander
committed
// multiply the transformation from the reference element to the actual element
Tensor3<double, TargetSpace::EmbeddedTangentVector::dimension,TargetSpace::EmbeddedTangentVector::dimension,gridDim> derivativeDerivative;
for (int ii=0; ii<TargetSpace::EmbeddedTangentVector::dimension; ii++)
for (int jj=0; jj<TargetSpace::EmbeddedTangentVector::dimension; jj++)

Oliver Sander
committed
for (int kk=0; kk<gridDim; kk++) {
derivativeDerivative[ii][jj][kk] = 0;
for (int ll=0; ll<gridDim; ll++)
derivativeDerivative[ii][jj][kk] += referenceDerivativeDerivative[ii][jj][ll] * jacobianInverseTransposed[kk][ll];
}

Oliver Sander
committed
for (int j=0; j<derivative.rows; j++) {
for (int k=0; k<derivative.cols; k++) {
for (int l=0; l<TargetSpace::EmbeddedTangentVector::dimension; l++)

Oliver Sander
committed
localGradient[i][l] += weight*derivative[j][k] * derivativeDerivative[l][j][k];

Oliver Sander
committed
}
}
}

Oliver Sander
committed
}
template <class GridView, class LocalFiniteElement, class TargetSpace>
void HarmonicEnergyLocalStiffness<GridView, LocalFiniteElement, TargetSpace>::

Oliver Sander
committed
assembleGradient(const Entity& element,
const LocalFiniteElement& localFiniteElement,
const std::vector<TargetSpace>& localSolution,

Oliver Sander
committed
std::vector<typename TargetSpace::TangentVector>& localGradient) const
{
std::vector<typename TargetSpace::EmbeddedTangentVector> embeddedLocalGradient;
// first compute the gradient in embedded coordinates
assembleEmbeddedGradient(element, localFiniteElement, localSolution, embeddedLocalGradient);

Oliver Sander
committed
// transform to coordinates on the tangent space
localGradient.resize(embeddedLocalGradient.size());
for (size_t i=0; i<localGradient.size(); i++)
localSolution[i].orthonormalFrame().mv(embeddedLocalGradient[i], localGradient[i]);
}
#endif