Skip to content
Snippets Groups Projects
Commit f7cb26cb authored by Elias Pipping's avatar Elias Pipping Committed by Elias Pipping
Browse files

Drop constness for state updates

parent f61840a5
No related branches found
No related tags found
No related merge requests found
......@@ -91,7 +91,7 @@ class FrictionPotential : public FrictionPotentialWrapper {
private:
double const coefficientProduct;
double const V_m;
double V_m;
};
class TrivialFunction : public FrictionPotentialWrapper {
......
......@@ -31,7 +31,7 @@ class GlobalNonlinearity {
/*
Return a restriction of the outer function to the i'th node.
*/
virtual shared_ptr<LocalFriction<dim> const> restriction(int i) const = 0;
virtual shared_ptr<LocalFriction<dim>> restriction(int i) const = 0;
void addHessian(VectorType const &v, MatrixType &hessian) const {
for (size_t i = 0; i < v.size(); ++i) {
......
......@@ -27,26 +27,26 @@ class GlobalRuinaNonlinearity
dataref nodalIntegrals, FrictionData const &fd,
dataref state)
: restrictions(nodalIntegrals.size()) {
auto trivialNonlinearity = make_shared<LocalFriction<dim> const>(
make_shared<TrivialFunction const>());
auto trivialNonlinearity =
make_shared<LocalFriction<dim>>(make_shared<TrivialFunction>());
for (size_t i = 0; i < restrictions.size(); ++i) {
restrictions[i] = frictionalNodes[i][0]
? trivialNonlinearity
: make_shared<LocalFriction<dim> const>(
make_shared<FrictionPotential const>(
nodalIntegrals[i], fd, state[i]));
restrictions[i] =
frictionalNodes[i][0]
? trivialNonlinearity
: make_shared<LocalFriction<dim>>(make_shared<FrictionPotential>(
nodalIntegrals[i], fd, state[i]));
}
}
/*
Return a restriction of the outer function to the i'th node.
*/
virtual shared_ptr<LocalFriction<dim> const> restriction(int i) const {
shared_ptr<LocalFriction<dim>> restriction(int i) const override {
return restrictions[i];
}
private:
std::vector<shared_ptr<LocalFriction<dim> const>> restrictions;
std::vector<shared_ptr<LocalFriction<dim>>> restrictions;
};
}
#endif
......@@ -18,7 +18,7 @@ template <int dimension> class LocalFriction {
using VectorType = FieldVector<double, dimension>;
using MatrixType = FieldMatrix<double, dimension, dimension>;
LocalFriction(shared_ptr<FrictionPotentialWrapper const> func)
LocalFriction(shared_ptr<FrictionPotentialWrapper> func)
: func(func), smp(func->smallestPositivePoint()) {}
double operator()(VectorType const &x) const {
......@@ -125,8 +125,8 @@ template <int dimension> class LocalFriction {
}
private:
shared_ptr<FrictionPotentialWrapper const> const func;
double const smp;
shared_ptr<FrictionPotentialWrapper> const func;
double smp;
};
}
#endif
......@@ -48,14 +48,14 @@ assemble_frictional(GridView const &gridView, AssemblerType const &assembler,
}
template <class MatrixType, class VectorType>
Dune::shared_ptr<Dune::GlobalNonlinearity<MatrixType, VectorType> const>
Dune::shared_ptr<Dune::GlobalNonlinearity<MatrixType, VectorType>>
assemble_nonlinearity(
Dune::BitSetVector<1> const &frictionalNodes,
Dune::BlockVector<Dune::FieldVector<double, 1>> const &nodalIntegrals,
FrictionData const &fd,
Dune::BlockVector<Dune::FieldVector<double, 1>> const &state) {
return Dune::make_shared<
Dune::GlobalRuinaNonlinearity<MatrixType, VectorType> const>(
Dune::GlobalRuinaNonlinearity<MatrixType, VectorType>>(
frictionalNodes, nodalIntegrals, fd, state);
}
......
......@@ -24,7 +24,7 @@ assemble_frictional(GridView const &gridView, AssemblerType const &assembler,
Dune::BitSetVector<1> const &frictionalNodes);
template <class MatrixType, class VectorType>
Dune::shared_ptr<Dune::GlobalNonlinearity<MatrixType, VectorType> const>
Dune::shared_ptr<Dune::GlobalNonlinearity<MatrixType, VectorType>>
assemble_nonlinearity(
Dune::BitSetVector<1> const &frictionalNodes,
Dune::BlockVector<Dune::FieldVector<double, 1>> const &nodalIntegrals,
......
......@@ -31,8 +31,7 @@ assemble_frictional<GridView, SmallVector, AssemblerType>(
GridView const &gridView, AssemblerType const &assembler,
Dune::BitSetVector<1> const &frictionalNodes);
template Dune::shared_ptr<
Dune::GlobalNonlinearity<MatrixType, VectorType> const>
template Dune::shared_ptr<Dune::GlobalNonlinearity<MatrixType, VectorType>>
assemble_nonlinearity<MatrixType, VectorType>(
Dune::BitSetVector<1> const &frictionalNodes,
Dune::BlockVector<Dune::FieldVector<double, 1>> const &nodalIntegrals,
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment