6#include <dune/fufem/boundarypatch.hh> 
    7#include <dune/functions/functionspacebases/lagrangebasis.hh> 
    8#include <dune/functions/functionspacebases/powerbasis.hh> 
    9#include <dune/grid/yaspgrid.hh> 
   10#include <dune/python/common/typeregistry.hh> 
   11#include <dune/python/functions/globalbasis.hh> 
   12#include <dune/python/pybind11/eigen.h> 
   13#include <dune/python/pybind11/functional.h> 
   14#include <dune/python/pybind11/pybind11.h> 
   15#include <dune/python/pybind11/stl.h> 
   25    using pybind11::operator
""_a;
 
   29    using GridView               = 
typename GlobalBasis::GridView;
 
   40            pybind11::keep_alive<1, 2>(), pybind11::keep_alive<1, 3>());
 
   42    using LoadFunction = std::function<Eigen::Vector<double, Traits::worlddim>(
 
   46                [](
const GlobalBasis& 
basis, 
const Element& element, 
const Material& mat,
 
   48            pybind11::keep_alive<1, 2>(), pybind11::keep_alive<1, 3>());
 
   51                              const LoadFunction volumeLoad, 
const BoundaryPatch<GridView>& bp,
 
   52                              const LoadFunction neumannBoundaryLoad) {
 
   55            pybind11::keep_alive<1, 2>(), pybind11::keep_alive<1, 3>(), pybind11::keep_alive<1, 7>());
 
   57    pybind11::module scopedf = pybind11::module::import(
"dune.functions");
 
   59    typedef Dune::Python::LocalViewWrapper<FlatBasis> LocalViewWrapper;
 
   60    auto includes = Dune::Python::IncludeFiles{
"dune/python/functions/globalbasis.hh"};
 
   61    auto lv       = Dune::Python::insertClass<LocalViewWrapper>(
 
   62                  scopedf, 
"LocalViewWrapper",
 
   63                  Dune::Python::GenerateTypeName(
"Dune::Python::LocalViewWrapperWrapper", Dune::MetaType<FlatBasis>()),
 
   66    lv.def(
"bind", &LocalViewWrapper::bind);
 
   67    lv.def(
"unbind", &LocalViewWrapper::unbind);
 
   68    lv.def(
"index", [](
const LocalViewWrapper& localView, 
int index) { 
return localView.index(index); });
 
   69    lv.def(
"__len__", [](LocalViewWrapper& self) -> 
int { 
return self.size(); });
 
   71    Dune::Python::Functions::registerTree<typename LocalViewWrapper::Tree>(lv);
 
   72    lv.def(
"tree", [](
const LocalViewWrapper& view) { 
return view.tree(); });
 
   74    auto basisName = Dune::className<FlatBasis>();
 
   75    auto entry     = Dune::Python::insertClass<FlatBasis>(
 
   76        scopedf, 
"DefaultGlobalBasis", pybind11::buffer_protocol(), Dune::Python::GenerateTypeName(basisName),
 
   77        Dune::Python::IncludeFiles{
"dune/python/functions/globalbasis.hh"});
 
   78    if (entry.second) Dune::Python::registerGlobalBasis(scopedf, entry.first);
 
   83          auto lvWrapped = LocalViewWrapper(self.localView().globalBasis());
 
   86          pybind11::object obj = pybind11::cast(self.localView().element());
 
   90        pybind11::keep_alive<0, 1>());
 
   91    cls.def(
"calculateScalar",
 
   94      return self.calculateVector(req, vec);
 
   99          return self.calculateMatrix(req, mat);
 
  101        pybind11::arg(
"FERequirements"), pybind11::arg(
"elementMatrix").noconvert());
 
  108          self.calculateAt(req, local, resultTypeMap);
 
  114        pybind11::arg(
"resultRequirements"), pybind11::arg(
"local"), pybind11::arg(
"resultType") = 
ResultType::noType);
 
Definition of the LinearElastic class for finite element mechanics computations.
 
ResultType
A strongly typed enum class representing the type of the result request.
Definition: ferequirements.hh:103
 
void init(int argc, char **argv, bool enableFileLogger=true)
Initializes the Ikarus framework.
Definition: init.hh:81
 
Definition: flatassembler.hh:20
 
void registerNonLinearElastic(pybind11::handle scope, pybind11::class_< NonLinearElastic, options... > cls)
Definition: python/finiteelements/nonlinearelastic.hh:24
 
def basis(gv, tree)
Definition: basis.py:10
 
def NonLinearElastic(basis, element, material, volumeLoad=None, bp=None, neumannBoundaryLoad=None)
Definition: finite_elements/__init__.py:149
 
Class representing the requirements for finite element calculations.
Definition: ferequirements.hh:170
 
Class representing a map of result types to result arrays.
Definition: ferequirements.hh:342
 
auto & getSingleResult()
Get the result array for a single result type.
Definition: ferequirements.hh:379
 
ResultArray & getResult(const ResultType &resultType)
Get the result array for a specific result type.
Definition: ferequirements.hh:368
 
Interface classf or materials.
Definition: interface.hh:75
 
NonLinearElastic class represents a non-linear elastic finite element.
Definition: finiteelements/mechanics/nonlinearelastic.hh:48
 
typename Traits::ResultRequirementsType ResultRequirementsType
Definition: finiteelements/mechanics/nonlinearelastic.hh:58
 
typename Traits::Element Element
Definition: finiteelements/mechanics/nonlinearelastic.hh:57
 
typename Traits::FERequirementType FERequirementType
Definition: finiteelements/mechanics/nonlinearelastic.hh:53
 
TraitsFromFE< Basis_, FERequirements_, useEigenRef > Traits
Definition: finiteelements/mechanics/nonlinearelastic.hh:50
 
Material_ Material
Definition: finiteelements/mechanics/nonlinearelastic.hh:60
 
typename Traits::Basis Basis
Definition: finiteelements/mechanics/nonlinearelastic.hh:51
 
typename Traits::FlatBasis FlatBasis
Definition: finiteelements/mechanics/nonlinearelastic.hh:52
 
Definition: resultevaluators.hh:20
 
Definition of the NonLinearElastic class for finite element mechanics computations.
 
Wrapper around Dune-functions global basis.