version 3.10-dev
python/assembly/fvassembler.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_PYTHON_COMMON_FVASSEMBLER_HH
13#define DUMUX_PYTHON_COMMON_FVASSEMBLER_HH
14
15#include <dune/python/pybind11/pybind11.h>
16#include <dune/python/pybind11/stl.h>
17
18namespace Dumux::Python {
19
20// Python wrapper for the FVAssembler C++ class
21template<class FVAssembler, class... options>
22void registerFVAssembler(pybind11::handle scope, pybind11::class_<FVAssembler, options...> cls)
23{
24 using pybind11::operator""_a;
25
26 using Problem = typename FVAssembler::Problem;
27 using GridGeometry = typename FVAssembler::GridGeometry;
28 using GridVariables = typename FVAssembler::GridVariables;
29 using SolutionVector = typename FVAssembler::SolutionVector;
30
31 static_assert(std::is_same_v<GridGeometry, typename Problem::GridGeometry>);
32 cls.def(pybind11::init([](std::shared_ptr<const Problem> problem,
33 std::shared_ptr<const GridGeometry> gridGeometry,
34 std::shared_ptr<GridVariables> gridVariables){
35 return std::make_shared<FVAssembler>(problem, gridGeometry, gridVariables);
36 }));
37
38 // TODO assembler with time loop
39
40 cls.def_property_readonly("numDofs", &FVAssembler::numDofs);
41 cls.def_property_readonly("problem", &FVAssembler::problem);
42 cls.def_property_readonly("gridGeometry", &FVAssembler::gridGeometry);
43 cls.def_property_readonly("gridView", &FVAssembler::gridView);
44 cls.def_property_readonly("jacobian", &FVAssembler::jacobian);
45 cls.def_property_readonly("residual", &FVAssembler::residual);
46 cls.def_property_readonly("prevSol", &FVAssembler::prevSol);
47 cls.def_property_readonly("isStationaryProblem", &FVAssembler::isStationaryProblem);
48 cls.def_property_readonly("gridVariables", [](FVAssembler& self) { return self.gridVariables(); });
49
50 cls.def("assembleResidual", [](FVAssembler& self, const SolutionVector& curSol){
51 self.assembleResidual(curSol);
52 });
53
54 cls.def("assembleJacobianAndResidual", [](FVAssembler& self, const SolutionVector& curSol){
55 self.assembleJacobianAndResidual(curSol);
56 });
57}
58
59} // end namespace Dumux::Python
60
61#endif
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: assembly/fvassembler.hh:94
ResidualType & residual()
The residual vector (rhs)
Definition: assembly/fvassembler.hh:294
const GridView & gridView() const
The gridview.
Definition: assembly/fvassembler.hh:278
const Problem & problem() const
The problem.
Definition: assembly/fvassembler.hh:270
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: assembly/fvassembler.hh:318
GetPropType< TypeTag, Properties::SolutionVector > SolutionVector
Definition: assembly/fvassembler.hh:110
std::size_t numDofs() const
Returns the number of degrees of freedom.
Definition: assembly/fvassembler.hh:266
GridVariables & gridVariables()
The global grid variables.
Definition: assembly/fvassembler.hh:282
JacobianMatrix & jacobian()
The jacobian matrix.
Definition: assembly/fvassembler.hh:290
GetPropType< TypeTag, Properties::Problem > Problem
Definition: assembly/fvassembler.hh:116
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: assembly/fvassembler.hh:202
GetPropType< TypeTag, Properties::GridVariables > GridVariables
Definition: assembly/fvassembler.hh:113
const GridGeometry & gridGeometry() const
The global finite volume geometry.
Definition: assembly/fvassembler.hh:274
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:171
GridGeo GridGeometry
Definition: assembly/fvassembler.hh:115
const SolutionVector & prevSol() const
The solution of the previous time step.
Definition: assembly/fvassembler.hh:298
Definition: python/assembly/fvassembler.hh:18
void registerFVAssembler(pybind11::handle scope, pybind11::class_< FVAssembler, options... > cls)
Definition: python/assembly/fvassembler.hh:22