Skip to content
Snippets Groups Projects
program_state.hh 6.58 KiB
Newer Older
#ifndef SRC_PROGRAM_STATE_HH
#define SRC_PROGRAM_STATE_HH

#include <dune/common/parametertree.hh>

podlesny's avatar
.  
podlesny committed
#include <dune/matrix-vector/axpy.hh>

#include <dune/fufem/boundarypatch.hh>
Elias Pipping's avatar
Elias Pipping committed
#include <dune/tnnmg/nonlinearities/zerononlinearity.hh>
#include <dune/tnnmg/problem-classes/blocknonlineartnnmgproblem.hh>
podlesny's avatar
.  
podlesny committed
#include <dune/contact/assemblers/nbodyassembler.hh>

#include <dune/tectonic/body.hh>

#include "assemblers.hh"
#include "matrices.hh"
#include "spatial-solving/solverfactory.hh"
podlesny's avatar
.  
podlesny committed
template <class VectorTEMPLATE, class ScalarVectorTEMPLATE> class BodyState {
  using Vector = VectorTEMPLATE;
  using ScalarVector = ScalarVectorTEMPLATE;

podlesny's avatar
.  
podlesny committed
  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>;

  ProgramState(const std::vector<size_t>& leafVertexCounts)
    : bodyCount_(leafVertexCounts.size()),
      bodies(bodyCount_),
      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),

podlesny's avatar
.  
podlesny committed
      bodies[i] = new BodyState(&u[i], &v[i], &a[i], &alpha[i], &weightedNormalStress[i]);
podlesny's avatar
.  
podlesny committed
    }
  }
podlesny's avatar
.  
podlesny committed
  ~ProgramState() {
    for (size_t i=0; i<bodyCount_; i++) {
      delete bodies[i];
    }
  }

podlesny's avatar
.  
podlesny committed
  size_t size() const {
podlesny's avatar
.  
podlesny committed
      return bodyCount_;
  }


  // Set up initial conditions
  template <class Matrix, class GridView>
  void setupInitialConditions(
podlesny's avatar
.  
podlesny committed
      const Dune::ParameterTree& parset,
      const Dune::Contact::NBodyAssembler<typename GridView::Grid, Vector>& nBodyAssembler,
      std::vector<std::function<void(double, Vector &)>> externalForces,
      const Matrices<Matrix>& matrices,
      const std::vector<std::shared_ptr<MyAssembler<GridView, Vector::block_type::dimension>>>& assemblers,
      const std::vector<Dune::BitSetVector<Vector::block_type::dimension>>& dirichletNodes,
      const std::vector<Dune::BitSetVector<Vector::block_type::dimension>>& noNodes,
      const std::vector<BoundaryPatch<GridView>>& frictionalBoundaries,
      const Body<Vector::block_type::dimension>& body) {

    using LocalVector = typename Vector::block_type;
podlesny's avatar
.  
podlesny committed
    //using LocalMatrix = typename Matrix::block_type;
    auto constexpr dims = LocalVector::dimension;
podlesny's avatar
.  
podlesny committed

    // Solving a linear problem with a multigrid solver
Elias Pipping's avatar
Elias Pipping committed
    auto const solveLinearProblem = [&](
podlesny's avatar
.  
podlesny committed
        const std::vector<Dune::BitSetVector<dims>>& _dirichletNodes, const std::vector<std::shared_ptr<Matrix>>& _matrices,
        const std::vector<Vector>& _rhs, std::vector<Vector>& _x,
Elias Pipping's avatar
Elias Pipping committed
        Dune::ParameterTree const &_localParset) {
podlesny's avatar
.  
podlesny committed
      std::vector<const Matrix*> matrices_ptr(_matrices.size());
      for (size_t i=0; i<matrices_ptr.size(); i++) {
            matrices_ptr[i] = _matrices[i].get();
      }
podlesny's avatar
.  
podlesny committed
      /*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];
podlesny's avatar
.  
podlesny committed
            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);
podlesny's avatar
.  
podlesny committed

podlesny's avatar
.  
podlesny committed
      Vector totalX;
      nBodyAssembler.nodalToTransformed(_x, totalX);
podlesny's avatar
.  
podlesny committed

podlesny's avatar
.  
podlesny committed
      using LinearFactory = SolverFactory<typename GridView::Grid, GlobalFriction<Matrix, Vector>, Matrix, Vector>;
      LinearFactory factory(parset.sub("solver.tnnmg"), nBodyAssembler, _dirichletNodes);

      auto multigridStep = factory.getStep();
      multigridStep->setProblem(bilinearForm, totalX, totalRhs);

      const EnergyNorm<Matrix, Vector> norm(bilinearForm);
podlesny's avatar
.  
podlesny committed

      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();
podlesny's avatar
.  
podlesny committed

podlesny's avatar
.  
podlesny committed
      nBodyAssembler.postprocess(multigridStep->getSol(), _x);
    timeStep = 0;
    relativeTime = 0.0;
    relativeTau = 1e-6;

podlesny's avatar
.  
podlesny committed
    std::vector<Vector> ell0(bodyCount_);
    for (size_t i=0; i<bodyCount_; i++) {
podlesny's avatar
.  
podlesny committed
      // Initial velocity
      v[i] = 0.0;
podlesny's avatar
.  
podlesny committed
      ell0[i].resize(u[i].size());
podlesny's avatar
.  
podlesny committed
      ell0[i] = 0.0;
      // TODO
      //externalForces[i](relativeTime, ell0[i]);
podlesny's avatar
.  
podlesny committed
    }

    // 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
    solveLinearProblem(dirichletNodes, 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
podlesny's avatar
.  
podlesny committed
    std::vector<Vector> accelerationRHS = ell0;
podlesny's avatar
.  
podlesny committed
    for (size_t i=0; i<bodyCount_; i++) {
podlesny's avatar
.  
podlesny committed
      // Initial state
      alpha[i] = parset.get<double>("boundary.friction.initialAlpha");

      // Initial normal stress
      assemblers[i]->assembleWeightedNormalStress(
        frictionalBoundaries[i], weightedNormalStress[i], body.getYoungModulus(),
        body.getPoissonRatio(), u[i]);
podlesny's avatar
.  
podlesny committed
      Dune::MatrixVector::subtractProduct(accelerationRHS[i], *matrices.elasticity[i], u[i]);
    }
podlesny's avatar
.  
podlesny committed
    solveLinearProblem(noNodes, matrices.mass, accelerationRHS, a,
                       parset.sub("a0.solver"));
podlesny's avatar
.  
podlesny committed
private:
  const size_t bodyCount_;

podlesny's avatar
.  
podlesny committed
  std::vector<BodyState* > bodies;
podlesny's avatar
.  
podlesny committed
  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;
podlesny's avatar
.  
podlesny committed

podlesny's avatar
.  
podlesny committed
};