DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
VariationalReconstruction_bind.hpp
Go to the documentation of this file.
1#pragma once
2
4#include "DNDS/Array_bind.hpp"
5#include <pybind11/functional.h>
6#include <pybind11_json/pybind11_json.hpp>
7#include <pybind11/eigen.h>
8
11
12namespace DNDS::CFV
13{
14 template <int dim = 2>
16
17 template <int dim = 2>
18 // static const int dim = 2;
20 {
21 auto VariationalReconstruction_ = tPy_VariationalReconstruction<dim>(m, fmt::format("VariationalReconstruction_{}", dim).c_str(),
22 tPy_FiniteVolume{m.attr("FiniteVolume")});
24 VariationalReconstruction_
25 .def(py::init([](MPIInfo &mpi, ssp<Geom::UnstructuredMesh> mesh)
26 { return std::make_shared<T>(mpi, mesh); }),
27 py::arg("mpi"), py::arg("mesh"))
28 .def("ConstructMetrics", &T::ConstructMetrics,
30 .def(
31 "ConstructBaseAndWeight",
32 [](T &self, typename T::tFGetBoundaryWeight f)
33 {
34 self.ConstructBaseAndWeight(
35 [f](Geom::t_index id, int order)
36 {
37 py::gil_scoped_acquire scope_gil;
38 return f(id, order);
39 });
40 },
41 py::arg("map_bcId_iOrder_to_bCweight"),
43 .def(
44 "ConstructBaseAndWeight_map",
45 [](T &self, const std::map<std::pair<Geom::t_index, int>, real> &m)
46 {
47 self.ConstructBaseAndWeight(
48 [&](Geom::t_index id, int order)
49 {
50 if (m.count({id, order}))
51 return m.at({id, order});
52 else
53 return 0.0;
54 });
55 },
56 py::arg("map_bcId_iOrder_to_bCweight"),
58 .def("ConstructRecCoeff", &T::ConstructRecCoeff,
60 // TODO: wrap Euler-related calls in EulerSolver inside euler!
61
62 VariationalReconstruction_
63 .def("SetPeriodicTransformations3d", [](T &self, std::array<int, 3> Seq123)
64 { self.SetPeriodicTransformations(Seq123); }, py::arg("Seq123"))
65 .def("SetPeriodicTransformations2d", [](T &self, std::array<int, 2> Seq123)
66 { self.SetPeriodicTransformations(Seq123); }, py::arg("Seq123"))
67 .def("SetPeriodicTransformationsNoOp", [](T &self)
68 { self.SetPeriodicTransformations(); });
69
70 VariationalReconstruction_
71 .def(
72 "ParseSettings", [](T &self, py::object settings)
73 {
74 VRSettings defaultSettings(self.getDim());
75 nlohmann::ordered_json defaultJson;
76 defaultSettings.WriteIntoJson(defaultJson);
77 nlohmann::json settings_json = settings;
78 defaultJson.merge_patch(settings_json);
79 self.parseSettings(defaultJson); },
80 py::arg("settings"));
81 VariationalReconstruction_
82 .def(
83 "GetCellBary", &T::GetCellBary, py::arg("iCell"));
84
85#define DNDS_CFV_VR_PYBIND11_DEFINE_BuildURec(nVarsFixed) \
86 VariationalReconstruction_.def( \
87 ("BuildURec_" + RowSize_To_PySnippet(nVarsFixed)).c_str(), \
88 [](T &self, tURec<nVarsFixed> &u, int nVars, bool buildSon, bool buildTrans) \
89 { self.BuildURec(u, nVars, buildSon, buildTrans); }, \
90 py::arg("u"), py::arg("nVars"), py::arg("buildSon") = true, py::arg("buildTrans") = true, \
91 DNDS_PYBIND11_OSTREAM_GUARD)
92#define DNDS_CFV_VR_PYBIND11_DEFINE_BuildUGrad(nVarsFixed) \
93 VariationalReconstruction_.def( \
94 ("BuildUGrad_" + RowSize_To_PySnippet(nVarsFixed)).c_str(), \
95 [](T &self, tUGrad<nVarsFixed, dim> &u, int nVars, bool buildSon, bool buildTrans) \
96 { self.BuildUGrad(u, nVars, buildSon, buildTrans); }, \
97 py::arg("u"), py::arg("nVars"), py::arg("buildSon") = true, py::arg("buildTrans") = true, \
98 DNDS_PYBIND11_OSTREAM_GUARD)
99
100 // #define DNDS_CFV_VR_PYBIND11_DEFINE_BuildCalls(nVarsFixed) \
101 // { \
102 // DNDS_CFV_VR_PYBIND11_DEFINE_BuildUDof(nVarsFixed); \
103 // DNDS_CFV_VR_PYBIND11_DEFINE_BuildURec(nVarsFixed); \
104 // DNDS_CFV_VR_PYBIND11_DEFINE_BuildUGrad(nVarsFixed); \
105 // }
106
107#define DNDS_CFV_VR_PYBIND11_DEFINE_BuildCalls(nVarsFixed) \
108 { \
109 DNDS_CFV_VR_PYBIND11_DEFINE_BuildURec(nVarsFixed); \
110 DNDS_CFV_VR_PYBIND11_DEFINE_BuildUGrad(nVarsFixed); \
111 }
112 if constexpr (dim == 2)
113 {
120 }
121 else
122 {
128 }
129#undef DNDS_CFV_VR_PYBIND11_DEFINE_BuildUDof
130#undef DNDS_CFV_VR_PYBIND11_DEFINE_BuildURec
131#undef DNDS_CFV_VR_PYBIND11_DEFINE_BuildUGrad
132#undef DNDS_CFV_VR_PYBIND11_DEFINE_BuildCalls
133
134 VariationalReconstruction_
135 .def_property_readonly("matrixAAInvB", &T::get_matrixAAInvB)
136 .def_property_readonly("vectorAInvB", &T::get_vectorAInvB);
137
138 VariationalReconstruction_
139 .def("GetFaceNorm", &T::GetFaceNorm,
140 py::arg("iFace"), py::arg("iG"))
141 .def("GetFaceNormFromCell", &T::GetFaceNormFromCell,
142 py::arg("iFace"), py::arg("iCell"), py::arg("if2c"), py::arg("ig"))
143 .def("GetCellVol", &T::GetCellVol, py::arg("iCell"))
144 .def("GetFaceArea", &T::GetFaceArea, py::arg("iFace"));
145 }
146}
pybind11 bindings for Array / ParArray / ArrayPair.
Shared pybind11 plumbing used by every *_bind.hpp in DNDS (buffer-protocol type check,...
#define DNDS_PYBIND11_OSTREAM_GUARD
#define DNDS_CFV_VR_PYBIND11_DEFINE_BuildCalls(nVarsFixed)
The VR class that provides any information needed in high-order CFV.
py_class_ssp< FiniteVolume > tPy_FiniteVolume
py_class_ssp< VariationalReconstruction< dim > > tPy_VariationalReconstruction
void pybind11_VariationalReconstruction_define(py::module_ &m)
int32_t t_index
Definition Geometric.hpp:6
py::classh< T > py_class_ssp
DNDS_CONSTANT const rowsize DynamicSize
Template parameter flag: "row width is set at runtime but uniform".
Definition Defines.hpp:277
std::shared_ptr< T > ssp
Shortened alias for std::shared_ptr used pervasively in DNDSR.
Definition Defines.hpp:138
double real
Canonical floating-point scalar used throughout DNDSR (double precision).
Definition Defines.hpp:105
A means to translate nlohmann json into c++ primitive data types and back; and stores then during com...
void WriteIntoJson(json &jsonSetting) const
Backward-compatible write (used by Python bindings).
Lightweight bundle of an MPI communicator and the calling rank's coordinates.
Definition MPI.hpp:215
double order
Definition test_ODE.cpp:257