Forked from
agnumpde / dune-tectonic
26 commits behind the upstream repository.
-
Elias Pipping authoredElias Pipping authored
Code owners
Assign users and groups as approvers for specific file changes. Learn more.
localfriction.hh 4.15 KiB
#ifndef DUNE_TECTONIC_LOCALFRICTION_HH
#define DUNE_TECTONIC_LOCALFRICTION_HH
#include <cmath>
#include <limits>
#include <dune/common/fvector.hh>
#include <dune/common/fmatrix.hh>
#include <dune/fufem/arithmetic.hh>
#include <dune/solvers/common/interval.hh>
#include <dune/tectonic/frictionpotential.hh>
template <size_t dimension> class LocalFriction {
public:
virtual ~LocalFriction() {}
using VectorType = Dune::FieldVector<double, dimension>;
using MatrixType = Dune::FieldMatrix<double, dimension, dimension>;
void virtual updateAlpha(double alpha) = 0;
double virtual regularity(VectorType const &x) const = 0;
double virtual coefficientOfFriction(VectorType const &x) const = 0;
void virtual directionalSubDiff(VectorType const &x, VectorType const &v,
Dune::Solvers::Interval<double> &D) const = 0;
void virtual addHessian(VectorType const &x, MatrixType &A) const = 0;
void virtual addGradient(VectorType const &x, VectorType &y) const = 0;
void virtual directionalDomain(
VectorType const &, VectorType const &,
Dune::Solvers::Interval<double> &dom) const = 0;
};
template <size_t dimension, class ScalarFriction>
class WrappedScalarFriction : public LocalFriction<dimension> {
using VectorType = typename LocalFriction<dimension>::VectorType;
using MatrixType = typename LocalFriction<dimension>::MatrixType;
public:
template <typename... Args>
WrappedScalarFriction(Args... args)
: func_(args...) {}
void updateAlpha(double alpha) override { func_.updateAlpha(alpha); }
double regularity(VectorType const &x) const override {
double const xnorm = x.two_norm();
if (xnorm <= 0.0)
return std::numeric_limits<double>::infinity();
return func_.regularity(xnorm);
}
double coefficientOfFriction(VectorType const &x) const override {
return func_.coefficientOfFriction(x.two_norm());
}
// directional subdifferential: at u on the line u + t*v
// u and v are assumed to be non-zero
void directionalSubDiff(VectorType const &x, VectorType const &v,
Dune::Solvers::Interval<double> &D) const override {
double const xnorm = x.two_norm();
if (xnorm <= 0.0)
D[0] = D[1] = func_.differential(0.0) * v.two_norm();
else
D[0] = D[1] = func_.differential(xnorm) * (x * v) / xnorm;
}
/** Formula for the derivative:
\f{align*}{
\frac {d^2}{dz^2} H(|z|)
&= \frac d{dz} \left( H'(|z|) \otimes \frac z{|z|} \right)\\
&= H''(|z|) \frac z{|z|} \otimes \frac z{|z|}
+ H'(|z|) \left( \frac {|z| \operatorname{id} - z \otimes z/|z|}{|z|^2}
\right)\\
&= \frac {H''(|z|)}{|z|^2} z \otimes z
+ \frac {H'(|z|)}{|z|} \operatorname{id}
- \frac {H'(|z|)}{|z|^3} z \otimes z\\
&= \left( \frac {H''(|z|)}{|z|^2} - \frac {H'(|z|)}{|z|^3} \right) z
\otimes z
+ \frac {H'(|z|)}{|z|} \operatorname{id}
\f}
*/
void addHessian(VectorType const &x, MatrixType &A) const override {
double const xnorm2 = x.two_norm2();
double const xnorm = std::sqrt(xnorm2);
if (xnorm2 <= 0.0)
return;
double const H1 = func_.differential(xnorm);
double const H2 = func_.second_deriv(xnorm);
double const tensorweight = (H2 - H1 / xnorm) / xnorm2;
double const idweight = H1 / xnorm;
for (size_t i = 0; i < dimension; ++i)
for (size_t j = 0; j < i; ++j) {
double const entry = tensorweight * x[i] * x[j];
A[i][j] += entry;
A[j][i] += entry;
}
for (size_t k = 0; k < dimension; ++k) {
double const entry = tensorweight * x[k] * x[k];
A[k][k] += entry + idweight;
}
}
void addGradient(VectorType const &x, VectorType &y) const override {
double const xnorm = x.two_norm();
if (std::isinf(func_.regularity(xnorm)))
return;
if (xnorm > 0.0)
Arithmetic::addProduct(y, func_.differential(xnorm) / xnorm, x);
}
void directionalDomain(VectorType const &, VectorType const &,
Dune::Solvers::Interval<double> &dom) const override {
dom[0] = -std::numeric_limits<double>::max();
dom[1] = std::numeric_limits<double>::max();
}
private:
ScalarFriction func_;
};
#endif