Newer
Older
#ifndef SRC_PROGRAM_STATE_HH
#define SRC_PROGRAM_STATE_HH
#include <dune/common/parametertree.hh>
#include <dune/fufem/boundarypatch.hh>
#include <dune/tnnmg/nonlinearities/zerononlinearity.hh>
#include <dune/tnnmg/problem-classes/blocknonlineartnnmgproblem.hh>
#include <dune/tectonic/bodydata.hh>
#include "levelcontactnetwork.hh"
#include "../spatial-solving/solverfactory.hh"
#include "../utils/debugutils.hh"
template <class VectorTEMPLATE, class ScalarVectorTEMPLATE> class BodyState {
using Vector = VectorTEMPLATE;
using ScalarVector = ScalarVectorTEMPLATE;
BodyState(Vector * _u, Vector * _v, Vector * _a, ScalarVector * _alpha, ScalarVector * _weightedNormalStress)
: u(_u),
v(_v),
a(_a),
alpha(_alpha),
weightedNormalStress(_weightedNormalStress) {}
public:
Vector * const u;
Vector * const v;
Vector * const a;
ScalarVector * const alpha;
ScalarVector * const weightedNormalStress;
};
template <class VectorTEMPLATE, class ScalarVectorTEMPLATE> class ProgramState {
public:
using Vector = VectorTEMPLATE;
using ScalarVector = ScalarVectorTEMPLATE;
using BodyState = BodyState<Vector, ScalarVector>;
private:
using LocalVector = typename Vector::block_type;
//using LocalMatrix = typename Matrix::block_type;
const static int dims = LocalVector::dimension;
public:
ProgramState(const std::vector<size_t>& leafVertexCounts)
: bodyCount_(leafVertexCounts.size()),
u(bodyCount_),
v(bodyCount_),
a(bodyCount_),
alpha(bodyCount_),
weightedNormalStress(bodyCount_) {
for (size_t i=0; i<bodyCount_; i++) {
size_t leafVertexCount = leafVertexCounts[i];
u[i].resize(leafVertexCount),
v[i].resize(leafVertexCount),
a[i].resize(leafVertexCount),
alpha[i].resize(leafVertexCount),
weightedNormalStress[i].resize(leafVertexCount),
bodyStates[i] = new BodyState(&u[i], &v[i], &a[i], &alpha[i], &weightedNormalStress[i]);
// Set up initial conditions
template <class GridType>
void setupInitialConditions(const Dune::ParameterTree& parset, const LevelContactNetwork<GridType, dims>& levelContactNetwork) {
using Matrix = typename LevelContactNetwork<GridType, dims>::Matrix;
const auto& nBodyAssembler = levelContactNetwork.nBodyAssembler();
// Solving a linear problem with a multigrid solver
const Dune::BitSetVector<dims>& _dirichletNodes, const std::vector<std::shared_ptr<Matrix>>& _matrices,
std::vector<const Matrix*> matrices_ptr(_matrices.size());
for (size_t i=0; i<matrices_ptr.size(); i++) {
matrices_ptr[i] = _matrices[i].get();
}
/*std::vector<Matrix> matrices(velocityMatrices.size());
std::vector<Vector> rhs(velocityRHSs.size());
for (size_t i=0; i<globalFriction_.size(); i++) {
matrices[i] = velocityMatrices[i];
rhs[i] = velocityRHSs[i];
globalFriction_[i]->addHessian(v_rel[i], matrices[i]);
globalFriction_[i]->addGradient(v_rel[i], rhs[i]);
matrices_ptr[i] = &matrices[i];
}*/
// assemble full global contact problem
Matrix bilinearForm;
nBodyAssembler.assembleJacobian(matrices_ptr, bilinearForm);
Vector totalRhs;
nBodyAssembler.assembleRightHandSide(_rhs, totalRhs);
using DeformedGridType = typename LevelContactNetwork<GridType, dims>::DeformedGridType;
using LinearFactory = SolverFactory<DeformedGridType, GlobalFriction<Matrix, Vector>, Matrix, Vector>;
LinearFactory factory(parset.sub("solver.tnnmg"), nBodyAssembler, _dirichletNodes);
/* std::vector<Dune::BitSetVector<dims>> bodyDirichletNodes;
nBodyAssembler.postprocess(_dirichletNodes, bodyDirichletNodes);
for (size_t i=0; i<bodyDirichletNodes.size(); i++) {
print(bodyDirichletNodes[i], "bodyDirichletNodes_" + std::to_string(i) + ": ");
}*/
/* print(bilinearForm, "matrix: ");
auto multigridStep = factory.getStep();
multigridStep->setProblem(bilinearForm, totalX, totalRhs);
const EnergyNorm<Matrix, Vector> norm(bilinearForm);
LoopSolver<Vector> solver(
multigridStep.get(), _localParset.get<size_t>("maximumIterations"),
_localParset.get<double>("tolerance"), &norm,
_localParset.get<Solver::VerbosityMode>("verbosity"),
false); // absolute error
solver.preprocess();
solver.solve();
relativeTime = 0.0;
relativeTau = 1e-6;
std::vector<Vector> ell0(bodyCount_);
for (size_t i=0; i<bodyCount_; i++) {
// Initial displacement: Start from a situation of minimal stress,
// which is automatically attained in the case [v = 0 = a].
// Assuming dPhi(v = 0) = 0, we thus only have to solve Au = ell0
Dune::BitSetVector<dims> dirichletNodes;
levelContactNetwork.totalNodes("dirichlet", dirichletNodes);
for (size_t i=0; i<dirichletNodes.size(); i++) {
bool val = false;
for (size_t d=0; d<dims; d++) {
val = val || dirichletNodes[i][d];
}
dirichletNodes[i] = val;
for (size_t d=0; d<dims; d++) {
dirichletNodes[i][d] = val;
}
}
solveLinearProblem(dirichletNodes, levelContactNetwork.matrices().elasticity, ell0, u,
parset.sub("u0.solver"));
// Initial acceleration: Computed in agreement with Ma = ell0 - Au
// (without Dirichlet constraints), again assuming dPhi(v = 0) = 0
// Initial state
alpha[i] = parset.get<double>("boundary.friction.initialAlpha");
// Initial normal stress
const auto& body = levelContactNetwork.body(i);
std::vector<std::shared_ptr<typename LevelContactNetwork<GridType, dims>::Body::LeafBoundaryCondition>> frictionBoundaryConditions;
body->leafBoundaryConditions("friction", frictionBoundaryConditions);
for (size_t j=0; j<frictionBoundaryConditions.size(); j++) {
ScalarVector frictionBoundaryStress(weightedNormalStress[i].size());
body->assembler()->assembleWeightedNormalStress(
*frictionBoundaryConditions[j]->boundaryPatch(), frictionBoundaryStress, body->data()->getYoungModulus(),
body->data()->getPoissonRatio(), u[i]);
weightedNormalStress[i] += frictionBoundaryStress;
}
Dune::MatrixVector::subtractProduct(accelerationRHS[i], *body->matrices().elasticity, u[i]);
solveLinearProblem(noNodes, levelContactNetwork.matrices().mass, accelerationRHS, a,
std::vector<BodyState* > bodyStates;
std::vector<Vector> u;
std::vector<Vector> v;
std::vector<Vector> a;
std::vector<ScalarVector> alpha;
std::vector<ScalarVector> weightedNormalStress;
double relativeTime;
double relativeTau;
size_t timeStep;