#ifndef SRC_PROGRAM_STATE_HH #define SRC_PROGRAM_STATE_HH #include <dune/common/parametertree.hh> #include <dune/matrix-vector/axpy.hh> #include <dune/fufem/boundarypatch.hh> #include <dune/contact/assemblers/nbodyassembler.hh> #include <dune/solvers/norms/energynorm.hh> #include <dune/solvers/solvers/loopsolver.hh> #include <dune/solvers/iterationsteps/cgstep.hh> #include "body/bodydata.hh" #include "../assemblers.hh" #include "network/contactnetwork.hh" #include "matrices.hh" #include "../spatial-solving/preconditioners/multilevelpatchpreconditioner.hh" #include "../spatial-solving/makelinearsolver.hh" #include "../spatial-solving/nonlinearsolver.hh" #include "../spatial-solving/solverfactory.hh" #include "../spatial-solving/functionalfactory.hh" #include "../spatial-solving/tnnmg/functional.hh" #include "../spatial-solving/tnnmg/zerononlinearity.hh" #include "../utils/debugutils.hh" #include <dune/solvers/iterationsteps/truncatedblockgsstep.hh> #include <dune/solvers/iterationsteps/multigridstep.hh> #include "../spatial-solving/contact/nbodycontacttransfer.hh" template <class VectorTEMPLATE, class ScalarVectorTEMPLATE> class BodyState { public: 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; using BitVector = Dune::BitSetVector<dims>; public: ProgramState(const std::vector<size_t>& leafVertexCounts) : bodyCount_(leafVertexCounts.size()), bodyStates(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), bodyStates[i] = new BodyState(&u[i], &v[i], &a[i], &alpha[i], &weightedNormalStress[i]); } } ~ProgramState() { for (size_t i=0; i<bodyCount_; i++) { delete bodyStates[i]; } } size_t size() const { return bodyCount_; } // Set up initial conditions template <class ContactNetwork> void setupInitialConditions(const Dune::ParameterTree& parset, const ContactNetwork& contactNetwork) { std::cout << "-- setupInitialConditions --" << std::endl; std::cout << "----------------------------" << std::endl; using Matrix = typename ContactNetwork::Matrix; const auto& nBodyAssembler = contactNetwork.nBodyAssembler(); auto linearSolver = makeLinearSolver<ContactNetwork, Vector>(parset, contactNetwork); auto nonlinearity = ZeroNonlinearity(); // Solving a linear problem with a multigrid solver auto const solveLinearProblem = [&]( const BitVector& _dirichletNodes, const std::vector<std::shared_ptr<Matrix>>& _matrices, const std::vector<Vector>& _rhs, std::vector<Vector>& _x, 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); }; timeStep = parset.get<size_t>("initialTime.timeStep"); relativeTime = parset.get<double>("initialTime.relativeTime"); relativeTau = parset.get<double>("initialTime.relativeTau"); std::vector<Vector> ell0(bodyCount_); for (size_t i=0; i<bodyCount_; i++) { // Initial velocity u[i] = 0.0; v[i] = 0.0; ell0[i].resize(u[i].size()); ell0[i] = 0.0; contactNetwork.body(i)->externalForce()(relativeTime, ell0[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 BitVector dirichletNodes; contactNetwork.totalNodes("dirichlet", dirichletNodes); 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(); } } std::cout << "solving linear problem for u..." << std::endl; solveLinearProblem(dirichletNodes, contactNetwork.matrices().elasticity, ell0, u, parset.sub("u0.solver")); //print(u, "initial u:"); // Initial acceleration: Computed in agreement with Ma = ell0 - Au // (without Dirichlet constraints), again assuming dPhi(v = 0) = 0 std::vector<Vector> accelerationRHS = ell0; for (size_t i=0; i<bodyCount_; i++) { // Initial state alpha[i] = parset.get<double>("boundary.friction.initialAlpha"); // Initial normal stress const auto& body = contactNetwork.body(i); /*std::vector<std::shared_ptr<typename ContactNetwork::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]); } 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; } std::cout << "solving linear problem for a..." << std::endl; BitVector noNodes(dirichletNodes.size(), false); solveLinearProblem(noNodes, contactNetwork.matrices().mass, accelerationRHS, a, parset.sub("a0.solver")); //print(a, "initial a:"); } private: const size_t bodyCount_; public: 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; }; #endif