12#ifndef DUMUX_PYTHON_DISCRETIZATION_GRIDVARIABLES_HH
13#define DUMUX_PYTHON_DISCRETIZATION_GRIDVARIABLES_HH
16#include <dune/istl/bvector.hh>
17#include <dune/python/pybind11/pybind11.h>
18#include <dune/python/common/typeregistry.hh>
25template<
class ElementSolution>
26void registerElementSolution(pybind11::handle scope)
28 using namespace Dune::Python;
30 auto [cls, addedToRegistry] = insertClass<ElementSolution>(
31 scope,
"ElementSolution",
32 GenerateTypeName(Dune::className<ElementSolution>()),
33 IncludeFiles{
"dumux/discretization/elementsolution.hh"}
38 using pybind11::operator
""_a;
40 cls.def(
"__getitem__", [](
const ElementSolution& self, std::size_t i){
42 throw pybind11::index_error();
46 cls.def_property_readonly(
"size", &ElementSolution::size);
54template <
class GV,
class... Options>
57 using pybind11::operator
""_a;
59 using Problem =
typename GV::GridVolumeVariables::Problem;
60 using PrimaryVariables =
typename GV::GridVolumeVariables::VolumeVariables::PrimaryVariables;
61 using GridGeometry =
typename GV::GridGeometry;
62 using Element =
typename GridGeometry::GridView::template Codim<0>::Entity;
63 using SolutionVector = Dune::BlockVector<PrimaryVariables>;
65 cls.def(pybind11::init([](std::shared_ptr<const Problem> problem,
66 std::shared_ptr<const GridGeometry> gridGeometry){
67 return std::make_shared<GV>(problem, gridGeometry);
70 cls.def(
"init", [](GV& self,
const SolutionVector& sol) {
return self.init(sol); });
71 cls.def(
"advanceTimeStep", &GV::advanceTimeStep);
72 cls.def_property_readonly(
"curGridVolVars", [](GV& self) {
return self.curGridVolVars(); });
73 cls.def_property_readonly(
"gridFluxVarsCache", [](GV& self) {
return self.gridFluxVarsCache(); });
74 cls.def_property_readonly(
"prevGridVolVars", [](GV& self) {
return self.prevGridVolVars(); });
75 cls.def_property_readonly(
"gridGeometry", &GV::gridGeometry);
77 cls.def(
"updateAfterGridAdaption", [](GV& self,
const SolutionVector& sol){
78 return self.updateAfterGridAdaption(sol);
81 cls.def(
"resetTimeStep", [](GV& self,
const SolutionVector& sol){
82 return self.resetTimeStep(sol);
85 cls.def(
"update", [](GV& self,
const SolutionVector& sol,
const bool forceFluxCacheUpdate =
false){
86 return self.update(sol, forceFluxCacheUpdate);
89 using ElementSolution = std::decay_t<decltype(elementSolution(std::declval<Element>(),
90 std::declval<SolutionVector>(),
91 std::declval<GridGeometry>()))>;
92 Impl::registerElementSolution<ElementSolution>(scope);
Element solution classes and factory functions.
Definition: python/assembly/fvassembler.hh:18
void registerGridVariables(pybind11::handle scope, pybind11::class_< GV, Options... > cls)
Definition: python/discretization/gridvariables.hh:55