3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_PYTHON_COMMON_FVASSEMBLER_HH
25#define DUMUX_PYTHON_COMMON_FVASSEMBLER_HH
26
27#include <dune/python/pybind11/pybind11.h>
28#include <dune/python/pybind11/stl.h>
29
30namespace Dumux::Python {
31
32// Python wrapper for the FVAssembler C++ class
33template<class FVAssembler, class... options>
34void registerFVAssembler(pybind11::handle scope, pybind11::class_<FVAssembler, options...> cls)
35{
36 using pybind11::operator""_a;
37
38 using Problem = typename FVAssembler::Problem;
39 using GridGeometry = typename FVAssembler::GridGeometry;
40 using GridVariables = typename FVAssembler::GridVariables;
41 using SolutionVector = typename FVAssembler::ResidualType;
42
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);
48 }));
49
50 // TODO assembler with time loop
51
52 cls.def_property_readonly("numDofs", &FVAssembler::numDofs);
53 cls.def_property_readonly("problem", &FVAssembler::problem);
54 cls.def_property_readonly("gridGeometry", &FVAssembler::gridGeometry);
55 cls.def_property_readonly("gridView", &FVAssembler::gridView);
56 cls.def_property_readonly("jacobian", &FVAssembler::jacobian);
57 cls.def_property_readonly("residual", &FVAssembler::residual);
58 cls.def_property_readonly("prevSol", &FVAssembler::prevSol);
59 cls.def_property_readonly("isStationaryProblem", &FVAssembler::isStationaryProblem);
60 cls.def_property_readonly("gridVariables", [](FVAssembler& self) { return self.gridVariables(); });
61
62 cls.def("assembleResidual", [](FVAssembler& self, const SolutionVector& curSol){
63 self.assembleResidual(curSol);
64 });
65
66 cls.def("assembleJacobianAndResidual", [](FVAssembler& self, const SolutionVector& curSol){
67 self.assembleJacobianAndResidual(curSol);
68 });
69}
70
71} // end namespace Dumux::Python
72
73#endif
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