Skip to content
Snippets Groups Projects
assemblers.hh 3.21 KiB
Newer Older
#ifndef SRC_ASSEMBLERS_HH
#define SRC_ASSEMBLERS_HH

#include <dune/common/bitsetvector.hh>
#include <dune/common/function.hh>
#include <dune/istl/bcrsmatrix.hh>
#include <dune/istl/bvector.hh>

#include <dune/fufem/assemblers/assembler.hh>
#pragma clang diagnostic push
#pragma clang diagnostic ignored "-Wsign-compare"
#include <dune/fufem/functionspacebases/p0basis.hh>
#pragma clang diagnostic pop
#include <dune/fufem/functionspacebases/p1nodalbasis.hh>
#include <dune/tectonic/globalfriction.hh>
#include <dune/tectonic/globalfrictiondata.hh>
#include "enums.hh"

template <class GridView, int dimension> class MyAssembler {
public:
  using ScalarMatrix = Dune::BCRSMatrix<Dune::FieldMatrix<double, 1, 1>>;
  using Matrix =
      Dune::BCRSMatrix<Dune::FieldMatrix<double, dimension, dimension>>;
  using ScalarVector = Dune::BlockVector<Dune::FieldVector<double, 1>>;
  using Vector = Dune::BlockVector<Dune::FieldVector<double, dimension>>;

  using CellBasis = P0Basis<GridView, double>;
  using VertexBasis = P1NodalBasis<GridView, double>;
Elias Pipping's avatar
Elias Pipping committed
  CellBasis const cellBasis;
  VertexBasis const vertexBasis;

private:
  using Grid = typename GridView::Grid;
  using LocalVector = typename Vector::block_type;
  using LocalScalarVector = typename ScalarVector::block_type;

  using LocalCellBasis = typename CellBasis::LocalFiniteElement;
  using LocalVertexBasis = typename VertexBasis::LocalFiniteElement;

  GridView const &gridView;
  Assembler<CellBasis, CellBasis> cellAssembler;
  Assembler<VertexBasis, VertexBasis> vertexAssembler;

public:
  MyAssembler(GridView const &gridView);

  void assembleFrictionalBoundaryMass(
      BoundaryPatch<GridView> const &frictionalBoundary,
Elias Pipping's avatar
Elias Pipping committed
      ScalarMatrix &frictionalBoundaryMass) const;
  void assembleMass(Dune::VirtualFunction<
                        LocalVector, LocalScalarVector> const &densityFunction,
Elias Pipping's avatar
Elias Pipping committed
                    Matrix &M) const;
Elias Pipping's avatar
Elias Pipping committed
  void assembleElasticity(double E, double nu, Matrix &A) const;
  void assembleViscosity(
      Dune::VirtualFunction<LocalVector, LocalScalarVector> const &
          shearViscosity,
      Dune::VirtualFunction<LocalVector, LocalScalarVector> const &
          bulkViscosity,
Elias Pipping's avatar
Elias Pipping committed
      Matrix &C) const;
  void assembleBodyForce(
      Dune::VirtualFunction<LocalVector, LocalVector> const &gravityField,
Elias Pipping's avatar
Elias Pipping committed
      Vector &f) const;

  void assembleNeumann(BoundaryPatch<GridView> const &neumannBoundary,
                       Vector &f,
                       Dune::VirtualFunction<double, double> const &neumann,
Elias Pipping's avatar
Elias Pipping committed
                       double relativeTime) const;
  void assembleNormalStress(BoundaryPatch<GridView> const &frictionalBoundary,
                            ScalarVector &normalStress, double youngModulus,
Elias Pipping's avatar
Elias Pipping committed
                            double poissonRatio,
                            Vector const &displacement) const;
  std::shared_ptr<GlobalFriction<Matrix, Vector>> assembleFrictionNonlinearity(
      Config::FrictionModel frictionModel,
      BoundaryPatch<GridView> const &frictionalBoundary,
      GlobalFrictionData<dimension> const &frictionInfo,
Elias Pipping's avatar
Elias Pipping committed
      ScalarVector const &normalStress) const;

  void assembleVonMisesStress(double youngModulus, double poissonRatio,
Elias Pipping's avatar
Elias Pipping committed
                              Vector const &u, ScalarVector &stress) const;