Skip to content
Snippets Groups Projects
program_state.hh 10.3 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/problem-classes/blocknonlineartnnmgproblem.hh>
podlesny's avatar
.  
podlesny committed
#include <dune/contact/assemblers/nbodyassembler.hh>

#include <dune/tectonic/bodydata.hh>
podlesny's avatar
podlesny committed
#include "../assemblers.hh"
podlesny's avatar
.  
podlesny committed
#include "contactnetwork.hh"
#include "matrices.hh"
podlesny's avatar
podlesny committed
#include "../spatial-solving/solverfactory.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
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;

public:
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
  template <class GridType>
podlesny's avatar
.  
podlesny committed
  void setupInitialConditions(const Dune::ParameterTree& parset, const ContactNetwork<GridType, Vector>& contactNetwork) {
    using Matrix = typename ContactNetwork<GridType, Vector>::Matrix;
    const auto& nBodyAssembler = contactNetwork.nBodyAssembler();
podlesny's avatar
.  
podlesny committed

    // Solving a linear problem with a multigrid solver
Elias Pipping's avatar
Elias Pipping committed
    auto const solveLinearProblem = [&](
        const Dune::BitSetVector<dims>& _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) {
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;
podlesny's avatar
.  
podlesny committed
      nBodyAssembler.assembleJacobian(matrices_ptr, bilinearForm);

podlesny's avatar
podlesny committed
      Vector totalRhs, oldTotalRhs;
podlesny's avatar
.  
podlesny committed
      nBodyAssembler.assembleRightHandSide(_rhs, totalRhs);
podlesny's avatar
podlesny committed
      oldTotalRhs = totalRhs;
podlesny's avatar
.  
podlesny committed

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

      // get lower and upper obstacles
      const auto& totalObstacles = nBodyAssembler.totalObstacles_;
      Vector lower(totalObstacles.size());
      Vector upper(totalObstacles.size());

      for (size_t j=0; j<totalObstacles.size(); ++j) {
          const auto& totalObstaclesj = totalObstacles[j];
          auto& lowerj = lower[j];
          auto& upperj = upper[j];
        for (size_t d=0; d<dims; ++d) {
            lowerj[d] = totalObstaclesj[d][0];
            upperj[d] = totalObstaclesj[d][1];
        }
      }

      // print problem
podlesny's avatar
.  
podlesny committed
     /* print(bilinearForm, "bilinearForm");
podlesny's avatar
podlesny committed
      print(totalRhs, "totalRhs");
      print(_dirichletNodes, "ignore");
      print(totalObstacles, "totalObstacles");
      print(lower, "lower");
podlesny's avatar
.  
podlesny committed
      print(upper, "upper");*/
podlesny's avatar
.  
podlesny committed

podlesny's avatar
podlesny committed
      // set up functional and TNMMG solver
podlesny's avatar
.  
podlesny committed
      using Functional = Functional<Matrix&, Vector&, ZeroNonlinearity&, Vector&, Vector&, typename Matrix::field_type>;
podlesny's avatar
.  
podlesny committed
      Functional J(bilinearForm, totalRhs, ZeroNonlinearity(), lower, upper); //TODO
podlesny's avatar
podlesny committed

podlesny's avatar
.  
podlesny committed
    //  using Factory = SolverFactory<Functional, Dune::BitSetVector<dims>>;
    //  Factory factory(parset.sub("solver.tnnmg"), J, _dirichletNodes);
podlesny's avatar
.  
podlesny committed

podlesny's avatar
.  
podlesny committed
   /*   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: ");
podlesny's avatar
podlesny committed
      print(totalX, "totalX: ");
podlesny's avatar
.  
podlesny committed
      print(totalRhs, "totalRhs: ");*/
podlesny's avatar
podlesny committed

podlesny's avatar
.  
podlesny committed
     /* auto tnnmgStep = factory.step();
podlesny's avatar
podlesny committed
      factory.setProblem(totalX);
podlesny's avatar
.  
podlesny committed

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

      LoopSolver<Vector> solver(
podlesny's avatar
podlesny committed
          tnnmgStep.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(tnnmgStep->getSol(), _x);

      Vector res = totalRhs;
      bilinearForm.mmv(tnnmgStep->getSol(), res);
podlesny's avatar
.  
podlesny committed
      std::cout << "TNNMG Res - energy norm: " << norm.operator ()(res) << std::endl; */
podlesny's avatar
podlesny committed

      // TODO: remove after debugging
  /*    using DeformedGridType = typename LevelContactNetwork<GridType, dims>::DeformedGridType;
      using OldLinearFactory = SolverFactoryOld<DeformedGridType, GlobalFriction<Matrix, Vector>, Matrix, Vector>;
      OldLinearFactory oldFactory(parset.sub("solver.tnnmg"), nBodyAssembler, _dirichletNodes);

      auto oldStep = oldFactory.getStep();
      oldStep->setProblem(bilinearForm, oldTotalX, oldTotalRhs);

          LoopSolver<Vector> oldSolver(
              oldStep.get(), _localParset.get<size_t>("maximumIterations"),
              _localParset.get<double>("tolerance"), &norm,
              _localParset.get<Solver::VerbosityMode>("verbosity"),
              false); // absolute error

          oldSolver.preprocess();
          oldSolver.solve();

      Vector oldRes = totalRhs;
      bilinearForm.mmv(oldStep->getSol(), oldRes);
      std::cout << "Old Res - energy norm: " << norm.operator ()(oldRes) << std::endl;*/

podlesny's avatar
.  
podlesny committed
   //   print(tnnmgStep->getSol(), "TNNMG Solution: ");
podlesny's avatar
podlesny committed
   /*   print(oldStep->getSol(), "Old Solution: ");
      auto diff = tnnmgStep->getSol();
      diff -= oldStep->getSol();
      std::cout << "Energy norm: " << norm.operator ()(diff) << std::endl;*/
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
    Dune::BitSetVector<dims> dirichletNodes;
podlesny's avatar
.  
podlesny committed
    contactNetwork.totalNodes("dirichlet", dirichletNodes);
podlesny's avatar
podlesny committed
    /*for (size_t i=0; i<dirichletNodes.size(); i++) {
podlesny's avatar
podlesny committed
        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;
        }
podlesny's avatar
podlesny committed
    }*/
podlesny's avatar
podlesny committed

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:");

    // 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);
      std::vector<std::shared_ptr<typename ContactNetwork<GridType, Vector>::LeafBody::BoundaryCondition>> frictionBoundaryConditions;
      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;
      }

      Dune::MatrixVector::subtractProduct(accelerationRHS[i], *body->matrices().elasticity, u[i]);
podlesny's avatar
.  
podlesny committed
    }
podlesny's avatar
.  
podlesny committed
    Dune::BitSetVector<dims> 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

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