Skip to content
Snippets Groups Projects
sand-wedge.cc 12.2 KiB
Newer Older
#include <Python.h>

#ifdef HAVE_CONFIG_H
#include "config.h"
#endif

#ifndef HAVE_PYTHON
#error Python is required
#endif

#ifdef HAVE_IPOPT
#undef HAVE_IPOPT
#endif

#include <cmath>
#include <exception>
#include <fstream>
#include <iostream>
#include <iomanip>
#include <boost/filesystem/operations.hpp>
#include <boost/filesystem/path.hpp>
#include <boost/format.hpp>

#include <dune/common/bitsetvector.hh>
#include <dune/common/exceptions.hh>
#include <dune/common/fmatrix.hh>
#include <dune/common/function.hh>
#include <dune/common/fvector.hh>
#include <dune/common/parametertree.hh>
#include <dune/common/parametertreeparser.hh>

#include <dune/grid/common/mcmgmapper.hh>
#include <dune/istl/bcrsmatrix.hh>
#include <dune/istl/bvector.hh>

#pragma clang diagnostic push
#pragma clang diagnostic ignored "-Wdeprecated-declarations"
#include <dune/fufem/boundarypatch.hh>
#pragma clang diagnostic pop

#include <dune/fufem/dunepython.hh>
#include <dune/fufem/sharedpointermap.hh>
#include <dune/solvers/norms/energynorm.hh>
#include <dune/solvers/solvers/loopsolver.hh>
#include <dune/solvers/solvers/solver.hh>
#include <dune/tnnmg/problem-classes/convexproblem.hh>

#include <dune/tectonic/myblockproblem.hh>
#include <dune/tectonic/globalfriction.hh>
#include <dune/fufem/hdf5/hdf5file.hh>
#include "adaptivetimestepper.hh"
#include "assemblers.hh"
#include "diameter.hh"
Elias Pipping's avatar
Elias Pipping committed
#include "enumparser.hh"
#include "enums.hh"
#include "gridselector.hh"
#include "hdf5-writer.hh"
#include "hdf5/restart-io.hh"
#include "matrices.hh"
#include "program_state.hh"
#include "sand-wedge-data/mybody.hh"
#include "sand-wedge-data/mygeometry.hh"
#include "sand-wedge-data/myglobalfrictiondata.hh"
#include "sand-wedge-data/mygrid.hh"
#include "sand-wedge-data/weakpatch.hh"
#include "solverfactory.hh"
#include "state.hh"
#include "updaters.hh"
#include "vtk.hh"

size_t const dims = MY_DIM;
namespace fs = boost::filesystem;
void initPython(fs::path const &dataDirectory) {
  Python::start();

  Python::run("import sys");
  Python::run(
      str(boost::format("sys.path.append('%s')") % dataDirectory.string()));
Dune::ParameterTree getParameters(int argc, char *argv[],
                                  fs::path const &dataDirectory) {
  Dune::ParameterTree parset;
  std::string sharedParsetName("parset.cfg");
  fs::path sharedParsetPath(dataDirectory / fs::path(sharedParsetName));
  Dune::ParameterTreeParser::readINITree(sharedParsetPath.string(), parset);

  std::string individualParsetName(str(boost::format("parset-%dD.cfg") % dims));
  fs::path individualParsetPath(dataDirectory / fs::path(individualParsetName));
  Dune::ParameterTreeParser::readINITree(individualParsetPath.string(), parset);

  Dune::ParameterTreeParser::readOptions(argc, argv, parset);

  return parset;
}

int main(int argc, char *argv[]) {
  try {
    auto const dataDirectory =
        fs::system_complete(fs::path(argv[0])).parent_path() /
        fs::path("sand-wedge-data");
    auto const parset = getParameters(argc, argv, dataDirectory);

    MyGeometry::render();
    MyGeometry::write();

Elias Pipping's avatar
Elias Pipping committed
    using GridView = Grid::LeafGridView;
    using MyAssembler = MyAssembler<GridView, dims>;
    using Matrix = MyAssembler::Matrix;
    using Vector = MyAssembler::Vector;
    using LocalVector = Vector::block_type;
    using ScalarMatrix = MyAssembler::ScalarMatrix;
    using ScalarVector = MyAssembler::ScalarVector;

    auto const weakPatch =
        getWeakPatch<LocalVector>(parset.sub("boundary.friction.weakening"));
    // {{{ Set up grid
    GridConstructor<Grid> gridConstructor;
    auto grid = gridConstructor.getGrid();
    refine(*grid, weakPatch,
           parset.get<double>("boundary.friction.smallestDiameter"));

Elias Pipping's avatar
Elias Pipping committed
    double minDiameter = std::numeric_limits<double>::infinity();
    double maxDiameter = 0.0;
    for (auto it = grid->template leafbegin<0>();
         it != grid->template leafend<0>(); ++it) {
      auto const geometry = it->geometry();
      auto const diam = diameter(geometry);
      minDiameter = std::min(minDiameter, diam);
      maxDiameter = std::max(maxDiameter, diam);
    }
    std::cout << "min diameter: " << minDiameter << std::endl;
    std::cout << "max diameter: " << maxDiameter << std::endl;

    auto const leafView = grid->leafGridView();
    auto const leafVertexCount = leafView.size(dims);
Elias Pipping's avatar
Elias Pipping committed
    std::cout << "Number of DOFs: " << leafVertexCount << std::endl;

    auto myFaces = gridConstructor.constructFaces(leafView);

    BoundaryPatch<GridView> const neumannBoundary(leafView);
    BoundaryPatch<GridView> const &frictionalBoundary = myFaces.lower;
    BoundaryPatch<GridView> const &surface = myFaces.upper;

    // Dirichlet Boundary
    Dune::BitSetVector<dims> noNodes(leafVertexCount);
    Dune::BitSetVector<dims> dirichletNodes(leafVertexCount);
    for (size_t i = 0; i < leafVertexCount; ++i) {
      if (myFaces.right.containsVertex(i))
        dirichletNodes[i][0] = true;

      if (myFaces.lower.containsVertex(i))
        dirichletNodes[i][1] = true;
#if MY_DIM == 3
      if (myFaces.front.containsVertex(i) || myFaces.back.containsVertex(i))
        dirichletNodes[i][2] = true;
#endif
    }

    // Set up functions for time-dependent boundary conditions
    using Function = Dune::VirtualFunction<double, double>;
    using FunctionMap = SharedPointerMap<std::string, Function>;
    FunctionMap functions;
    {
      initPython(dataDirectory);
      Python::import("boundaryconditions")
          .get("Functions")
          .toC<typename FunctionMap::Base>(functions);
    }
    auto const &velocityDirichletFunction =
                   functions.get("velocityDirichletCondition"),
               &neumannFunction = functions.get("neumannCondition");

Elias Pipping's avatar
Elias Pipping committed
    MyAssembler const myAssembler(leafView);

    MyBody<dims> const body(parset);

    Matrices<Matrix> matrices;
    myAssembler.assembleElasticity(body.getYoungModulus(),
                                   body.getPoissonRatio(), matrices.elasticity);
    myAssembler.assembleViscosity(body.getShearViscosityField(),
                                  body.getBulkViscosityField(),
                                  matrices.damping);
    myAssembler.assembleMass(body.getDensityField(), matrices.mass);
    ScalarMatrix relativeFrictionalBoundaryMass;
    myAssembler.assembleFrictionalBoundaryMass(frictionalBoundary,
                                               relativeFrictionalBoundaryMass);
    relativeFrictionalBoundaryMass /= frictionalBoundary.area();
    EnergyNorm<ScalarMatrix, ScalarVector> const stateEnergyNorm(
        relativeFrictionalBoundaryMass);

    // Assemble forces
    Vector gravityFunctional;
    myAssembler.assembleBodyForce(body.getGravityField(), gravityFunctional);

    // Problem formulation: right-hand side
Elias Pipping's avatar
Elias Pipping committed
    std::function<void(double, Vector &)> computeExternalForces =
        [&](double _relativeTime, Vector &_ell) {
          myAssembler.assembleNeumann(neumannBoundary, _ell, neumannFunction,
                                      _relativeTime);
          _ell += gravityFunctional;
        };
    ProgramState<Vector, ScalarVector> programState(leafVertexCount);
    auto const firstRestart = parset.get<size_t>("restarts.first");
    auto const restartSpacing = parset.get<size_t>("restarts.spacing");
    auto const restartTemplate = parset.get<std::string>("restarts.template");
    auto const restartDirectory = fs::path(restartTemplate).parent_path();
    HDF5File ioFile("output.h5");
    if (firstRestart != 0)
      RestartIO<ProgramState<Vector, ScalarVector>>(ioFile, leafVertexCount)
          .read(firstRestart, programState);
    else
      programState.setupInitialConditions(parset, computeExternalForces,
                                          matrices, myAssembler, dirichletNodes,
                                          noNodes, frictionalBoundary, body);
    MyGlobalFrictionData<LocalVector> frictionInfo(
        parset.sub("boundary.friction"), weakPatch);
    auto myGlobalFriction = myAssembler.assembleFrictionNonlinearity(
        parset.get<Config::FrictionModel>("boundary.friction.frictionModel"),
        frictionalBoundary, frictionInfo, programState.weightedNormalStress);
    myGlobalFriction->updateAlpha(programState.alpha);
    Vector vertexCoordinates(leafVertexCount);
    {
      Dune::MultipleCodimMultipleGeomTypeMapper<
          GridView, Dune::MCMGVertexLayout> const vertexMapper(leafView);
      for (auto it = leafView.begin<dims>(); it != leafView.end<dims>(); ++it) {
        auto const geometry = it->geometry();
        assert(geometry.corners() == 1);
        vertexCoordinates[vertexMapper.index(*it)] = geometry.corner(0);
    HDF5Writer<ProgramState<Vector, ScalarVector>,
               typename MyAssembler::VertexBasis, GridView>
        hdf5Writer(ioFile, vertexCoordinates, myAssembler.vertexBasis, surface,
                   frictionalBoundary, weakPatch);
    MyVTKWriter<typename MyAssembler::VertexBasis,
                typename MyAssembler::CellBasis> const
        vtkWriter(myAssembler.cellBasis, myAssembler.vertexBasis, "obs");

    IterationRegister iterationCount;
    auto const report = [&](bool initial = false) {
      hdf5Writer.reportSolution(programState, *myGlobalFriction);
      if (!initial)
        hdf5Writer.reportIterations(programState, iterationCount);

      if (!initial and programState.timeStep % restartSpacing == 0)
        hdf5Writer.writeRestart(programState);

      ioFile.flush();
      if (parset.get<bool>("io.writeVTK")) {
        ScalarVector stress;
        myAssembler.assembleVonMisesStress(body.getYoungModulus(),
                                           body.getPoissonRatio(),
                                           programState.u, stress);
        vtkWriter.write(programState.timeStep, programState.u, programState.v,
                        programState.alpha, stress);
      }
    };

    // Set up TNNMG solver
Elias Pipping's avatar
Elias Pipping committed
    using NonlinearFactory = SolverFactory<
        dims,
        MyBlockProblem<ConvexProblem<GlobalFriction<Matrix, Vector>, Matrix>>,
        Grid>;
    NonlinearFactory factory(parset.sub("solver.tnnmg"), *grid, dirichletNodes);
    using MyUpdater = Updaters<RateUpdater<Vector, Matrix, Function, dims>,
                               StateUpdater<ScalarVector, Vector>>;
    MyUpdater current(
        initRateUpdater(parset.get<Config::scheme>("timeSteps.scheme"),
                        velocityDirichletFunction, dirichletNodes, matrices,
                        programState.u, programState.v, programState.a),
        initStateUpdater<ScalarVector, Vector>(
            parset.get<Config::stateModel>("boundary.friction.stateModel"),
            programState.alpha, *frictionalBoundary.getVertices(),
            parset.get<double>("boundary.friction.L"),
            parset.get<double>("boundary.friction.V0")));
    auto const refinementTolerance =
        parset.get<double>("timeSteps.refinementTolerance");
    auto const mustRefine = [&](MyUpdater &coarseUpdater,
                                MyUpdater &fineUpdater) {
      ScalarVector coarseAlpha;
      coarseUpdater.state_->extractAlpha(coarseAlpha);

      ScalarVector fineAlpha;
      fineUpdater.state_->extractAlpha(fineAlpha);

      return stateEnergyNorm.diff(fineAlpha, coarseAlpha) > refinementTolerance;
    };
    AdaptiveTimeStepper<NonlinearFactory, MyUpdater,
                        EnergyNorm<ScalarMatrix, ScalarVector>>
Elias Pipping's avatar
Elias Pipping committed
        adaptiveTimeStepper(factory, parset, myGlobalFriction, current,
                            programState.relativeTime, programState.relativeTau,
Elias Pipping's avatar
Elias Pipping committed
                            computeExternalForces, stateEnergyNorm, mustRefine);
    while (!adaptiveTimeStepper.reachedEnd()) {
      programState.timeStep++;

      iterationCount = adaptiveTimeStepper.advance();
      programState.relativeTime = adaptiveTimeStepper.relativeTime_;
      programState.relativeTau = adaptiveTimeStepper.relativeTau_;
      current.rate_->extractDisplacement(programState.u);
      current.rate_->extractVelocity(programState.v);
      current.rate_->extractAcceleration(programState.a);
      current.state_->extractAlpha(programState.alpha);
      report();
Elias Pipping's avatar
Elias Pipping committed
  } catch (Dune::Exception &e) {
    Dune::derr << "Dune reported error: " << e << std::endl;
Elias Pipping's avatar
Elias Pipping committed
  } catch (std::exception &e) {
    std::cerr << "Standard exception: " << e.what() << std::endl;
  }
}