Skip to content
Snippets Groups Projects
program_state.hh 8 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>

podlesny's avatar
.  
podlesny committed
#include <dune/contact/assemblers/nbodyassembler.hh>

podlesny's avatar
.  
podlesny committed
#include <dune/solvers/norms/energynorm.hh>
#include <dune/solvers/solvers/loopsolver.hh>
#include <dune/solvers/iterationsteps/cgstep.hh>

podlesny's avatar
podlesny committed
#include "body/bodydata.hh"
podlesny's avatar
podlesny committed
#include "../assemblers.hh"
podlesny's avatar
podlesny committed
#include "network/contactnetwork.hh"
#include "matrices.hh"
podlesny's avatar
.  
podlesny committed
#include "../spatial-solving/preconditioners/multilevelpatchpreconditioner.hh"
#include "../spatial-solving/makelinearsolver.hh"
#include "../spatial-solving/nonlinearsolver.hh"
podlesny's avatar
podlesny committed
#include "../spatial-solving/solverfactory.hh"
#include "../spatial-solving/functionalfactory.hh"
podlesny's avatar
podlesny committed
#include "../spatial-solving/tnnmg/functional.hh"
#include "../spatial-solving/tnnmg/zerononlinearity.hh"
podlesny's avatar
podlesny committed
#include "../utils/debugutils.hh"
podlesny's avatar
.  
podlesny committed

#include <dune/solvers/iterationsteps/truncatedblockgsstep.hh>
#include <dune/solvers/iterationsteps/multigridstep.hh>

podlesny's avatar
podlesny committed
#include "../spatial-solving/contact/nbodycontacttransfer.hh"
podlesny's avatar
.  
podlesny committed

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>;

private:
    using LocalVector = typename Vector::block_type;
    //using LocalMatrix = typename Matrix::block_type;
    const static int dims = LocalVector::dimension;
podlesny's avatar
.  
podlesny committed
    using BitVector = Dune::BitSetVector<dims>;
podlesny's avatar
.  
podlesny committed
  ProgramState(const std::vector<size_t>& leafVertexCounts)
    : bodyCount_(leafVertexCounts.size()),
      bodyStates(bodyCount_),
podlesny's avatar
.  
podlesny committed
      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]);
podlesny's avatar
.  
podlesny committed
    }
  }
podlesny's avatar
.  
podlesny committed
  ~ProgramState() {
    for (size_t i=0; i<bodyCount_; i++) {
      delete bodyStates[i];
podlesny's avatar
.  
podlesny committed
    }
  }

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

  // Set up initial conditions
podlesny's avatar
.  
podlesny committed
  template <class ContactNetwork>
  void setupInitialConditions(const Dune::ParameterTree& parset, const ContactNetwork& contactNetwork) {
podlesny's avatar
.  
podlesny committed

podlesny's avatar
.  
podlesny committed
    std::cout << "-- setupInitialConditions --" << std::endl;
    std::cout << "----------------------------" << std::endl;

    using Matrix = typename ContactNetwork::Matrix;
    const auto& nBodyAssembler = contactNetwork.nBodyAssembler();
podlesny's avatar
.  
podlesny committed

    auto linearSolver = makeLinearSolver<ContactNetwork, Vector>(parset, contactNetwork);
    auto nonlinearity = ZeroNonlinearity();
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 BitVector& _dirichletNodes, const std::vector<std::shared_ptr<Matrix>>& _matrices,
podlesny's avatar
.  
podlesny committed
        const std::vector<Vector>& _rhs, std::vector<Vector>& _x,
Elias Pipping's avatar
Elias Pipping committed
        Dune::ParameterTree const &_localParset) {
        Vector totalX;
        nBodyAssembler.nodalToTransformed(_x, totalX);

        FunctionalFactory<std::decay_t<decltype(nBodyAssembler)>, decltype(nonlinearity), Matrix, Vector> functionalFactory(nBodyAssembler, nonlinearity, _matrices, _rhs);
        functionalFactory.build();
        auto functional = functionalFactory.functional();

        NonlinearSolver<std::remove_reference_t<decltype(*functional)>, BitVector> solver(parset.sub("solver.tnnmg"), linearSolver, functional, _dirichletNodes);
        solver.solve(_localParset, totalX);

        nBodyAssembler.postprocess(totalX, _x);
podlesny's avatar
.  
podlesny committed
    timeStep = parset.get<size_t>("initialTime.timeStep");
    relativeTime = parset.get<double>("initialTime.relativeTime");
    relativeTau = parset.get<double>("initialTime.relativeTau");
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
podlesny's avatar
podlesny committed
      u[i] = 0.0;
podlesny's avatar
.  
podlesny committed
      v[i] = 0.0;
podlesny's avatar
.  
podlesny committed
      ell0[i].resize(u[i].size());
podlesny's avatar
.  
podlesny committed
      ell0[i] = 0.0;
podlesny's avatar
.  
podlesny committed
      contactNetwork.body(i)->externalForce()(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
podlesny's avatar
.  
podlesny committed
    BitVector dirichletNodes;
podlesny's avatar
.  
podlesny committed
    contactNetwork.totalNodes("dirichlet", dirichletNodes);
podlesny's avatar
podlesny committed

    size_t dof = 0;
    for (size_t i=0; i<bodyCount_; i++) {
        const auto& body = *contactNetwork.body(i);

        if (body.data()->getYoungModulus() == 0.0) {
            for (; dof<body.nVertices(); dof++) {
                dirichletNodes[dof] = true;
            }
        } else {
            dof += body.nVertices();
podlesny's avatar
podlesny committed
        }
podlesny's avatar
podlesny committed

podlesny's avatar
.  
podlesny committed
    std::cout << "solving linear problem for u..." << std::endl;

podlesny's avatar
.  
podlesny committed
    solveLinearProblem(dirichletNodes, contactNetwork.matrices().elasticity, ell0, u,
                       parset.sub("u0.solver"));

podlesny's avatar
podlesny committed
    //print(u, "initial u:");
podlesny's avatar
.  
podlesny committed

    // 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
podlesny's avatar
.  
podlesny committed
      const auto& body = contactNetwork.body(i);
podlesny's avatar
podlesny committed
      /*std::vector<std::shared_ptr<typename ContactNetwork::LeafBody::BoundaryCondition>> frictionBoundaryConditions;
podlesny's avatar
.  
podlesny committed
      body->boundaryConditions("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;
podlesny's avatar
podlesny committed
      }*/

      Dune::MatrixVector::subtractProduct(accelerationRHS[i], *body->matrices().elasticity, u[i]);
podlesny's avatar
.  
podlesny committed
    }
podlesny's avatar
podlesny committed
    for (size_t i=0; i<contactNetwork.nCouplings(); i++) {
      const auto& coupling = contactNetwork.coupling(i);
      const auto& contactCoupling = contactNetwork.nBodyAssembler().getContactCouplings()[i];

      const auto nonmortarIdx = coupling->gridIdx_[0];
      const auto& body = contactNetwork.body(nonmortarIdx);

      ScalarVector frictionBoundaryStress(weightedNormalStress[nonmortarIdx].size());

      body->assembler()->assembleWeightedNormalStress(
        contactCoupling->nonmortarBoundary(), frictionBoundaryStress, body->data()->getYoungModulus(),
        body->data()->getPoissonRatio(), u[nonmortarIdx]);

      weightedNormalStress[nonmortarIdx] += frictionBoundaryStress;
    }

podlesny's avatar
.  
podlesny committed
    std::cout << "solving linear problem for a..." << std::endl;

podlesny's avatar
.  
podlesny committed
    BitVector noNodes(dirichletNodes.size(), false);
podlesny's avatar
.  
podlesny committed
    solveLinearProblem(noNodes, contactNetwork.matrices().mass, accelerationRHS, a,
podlesny's avatar
.  
podlesny committed
                       parset.sub("a0.solver"));
podlesny's avatar
.  
podlesny committed

podlesny's avatar
podlesny committed
    //print(a, "initial a:");
podlesny's avatar
.  
podlesny committed
private:
  const size_t bodyCount_;

  std::vector<BodyState* > bodyStates;
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
};