24#ifndef DUMUX_PYTHON_COMMON_FVASSEMBLER_HH
25#define DUMUX_PYTHON_COMMON_FVASSEMBLER_HH
27#include <dune/python/pybind11/pybind11.h>
28#include <dune/python/pybind11/stl.h>
36 using pybind11::operator
""_a;
43 static_assert(std::is_same_v<GridGeometry, typename Problem::GridGeometry>);
44 cls.def(pybind11::init([](std::shared_ptr<const Problem> problem,
45 std::shared_ptr<const GridGeometry> gridGeometry,
46 std::shared_ptr<GridVariables> gridVariables){
47 return std::make_shared<FVAssembler>(problem, gridGeometry, gridVariables);
62 cls.def(
"assembleResidual", [](
FVAssembler& self,
const SolutionVector& curSol){
66 cls.def(
"assembleJacobianAndResidual", [](
FVAssembler& self,
const SolutionVector& curSol){
Definition: python/assembly/fvassembler.hh:30
void registerFVAssembler(pybind11::handle scope, pybind11::class_< FVAssembler, options... > cls)
Definition: python/assembly/fvassembler.hh:34
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: assembly/fvassembler.hh:119
const GridView & gridView() const
The gridview.
Definition: assembly/fvassembler.hh:342
const Problem & problem() const
The problem.
Definition: assembly/fvassembler.hh:334
SolutionVector ResidualType
Definition: assembly/fvassembler.hh:140
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: assembly/fvassembler.hh:382
std::size_t numDofs() const
Returns the number of degrees of freedom.
Definition: assembly/fvassembler.hh:330
GridVariables & gridVariables()
The global grid variables.
Definition: assembly/fvassembler.hh:346
JacobianMatrix & jacobian()
The jacobian matrix.
Definition: assembly/fvassembler.hh:354
GetPropType< TypeTag, Properties::Problem > Problem
Definition: assembly/fvassembler.hh:137
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: assembly/fvassembler.hh:226
GetPropType< TypeTag, Properties::GridVariables > GridVariables
Definition: assembly/fvassembler.hh:138
const GridGeometry & gridGeometry() const
The global finite volume geometry.
Definition: assembly/fvassembler.hh:338
void assembleJacobianAndResidual(const SolutionVector &curSol, const PartialReassembler *partialReassembler=nullptr)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: assembly/fvassembler.hh:195
GridGeo GridGeometry
Definition: assembly/fvassembler.hh:136
SolutionVector & residual()
The residual vector (rhs)
Definition: assembly/fvassembler.hh:358
const SolutionVector & prevSol() const
The solution of the previous time step.
Definition: assembly/fvassembler.hh:362