#ifdef HAVE_CONFIG_H #include "config.h" #endif #include <dune/fufem/boundarypatch.hh> #include <dune/fufem/functions/constantfunction.hh> #include <dune/fufem/assemblers/localassemblers/neumannboundaryassembler.hh> #include <dune/tectonic/globalruinanonlinearity.hh> #include "assemblers.hh" // Assembles Neumann boundary term in f template <class GridView, class LocalVector, class Assembler> void assembleNeumann(GridView const &gridView, Assembler const &assembler, Dune::BitSetVector<1> const &neumannNodes, Dune::BlockVector<LocalVector> &f, Dune::VirtualFunction<double, double> const &neumann, double relativeTime) { // constant sample function on // neumann boundary BoundaryPatch<GridView> const neumannBoundary(gridView, neumannNodes); LocalVector localNeumann(0); neumann.evaluate(relativeTime, localNeumann[0]); ConstantFunction<LocalVector, LocalVector> const fNeumann(localNeumann); NeumannBoundaryAssembler<typename GridView::Grid, LocalVector> neumannBoundaryAssembler(fNeumann); assembler.assembleBoundaryFunctional(neumannBoundaryAssembler, f, neumannBoundary); } // Assembles constant 1-function on frictional boundary in nodalIntegrals template <class GridView, class LocalVector, class Assembler> Dune::shared_ptr<Dune::BlockVector<Dune::FieldVector<double, 1>>> assembleFrictionWeightsal(GridView const &gridView, Assembler const &assembler, Dune::BitSetVector<1> const &frictionalNodes) { using Singleton = Dune::FieldVector<double, 1>; BoundaryPatch<GridView> const frictionalBoundary(gridView, frictionalNodes); ConstantFunction<LocalVector, Singleton> const constantOneFunction(1); NeumannBoundaryAssembler<typename GridView::Grid, Singleton> frictionalBoundaryAssembler(constantOneFunction); auto const nodalIntegrals = Dune::make_shared<Dune::BlockVector<Singleton>>(); assembler.assembleBoundaryFunctional(frictionalBoundaryAssembler, *nodalIntegrals, frictionalBoundary); return nodalIntegrals; } template <class Matrix, class Vector> Dune::shared_ptr<Dune::GlobalNonlinearity<Matrix, Vector>> assembleNonlinearity( Dune::BitSetVector<1> const &frictionalNodes, Dune::BlockVector<Dune::FieldVector<double, 1>> const &nodalIntegrals, FrictionData const &fd) { return Dune::make_shared<Dune::GlobalRuinaNonlinearity<Matrix, Vector>>( frictionalNodes, nodalIntegrals, fd); } #include "assemblers_tmpl.cc"