OGS
PythonBoundaryCondition.cpp
Go to the documentation of this file.
1
12
13#include <algorithm>
14#include <pybind11/pybind11.h>
15
16#include <iostream>
17
18#include "BaseLib/ConfigTree.h"
19#include "FlushStdoutGuard.h"
24
25namespace
26{
28 std::size_t const nnodes)
29{
30 bc_values.ids.clear();
31 bc_values.values.clear();
32
33 bc_values.ids.reserve(nnodes);
34 bc_values.values.reserve(nnodes);
35}
36
38 NumLib::LocalToGlobalIndexMap const& dof_table,
39 std::vector<std::reference_wrapper<ProcessLib::ProcessVariable>> const& pvs)
40{
41 auto const num_vars_dt = dof_table.getNumberOfVariables();
42 auto const num_vars_pv = pvs.size();
43
44 if (static_cast<std::size_t>(num_vars_dt) != num_vars_pv)
45 {
47 "The number of variables in the d.o.f. table does not match the "
48 "number of process variables: {} != {}.",
49 num_vars_dt, num_vars_pv);
50 }
51
52 for (std::size_t var = 0; var < num_vars_pv; ++var)
53 {
54 auto const num_comp_dt = dof_table.getNumberOfVariableComponents(var);
55 auto const num_comp_pv = pvs[var].get().getNumberOfGlobalComponents();
56
57 if (num_comp_dt != num_comp_pv)
58 {
60 "The number of components of variable #{} in the d.o.f. table "
61 "does not match the number of components of process variable "
62 "#{} ({}): {} != {}.",
63 var, var, pvs[var].get().getName(), num_vars_dt, num_vars_pv);
64 }
65 }
66}
67} // anonymous namespace
68
69namespace ProcessLib
70{
72 PythonBcData&& bc_data, unsigned const integration_order,
73 bool const flush_stdout, unsigned const bulk_mesh_dimension,
74 NumLib::LocalToGlobalIndexMap const& dof_table_bulk)
75 : _bc_data(std::move(bc_data)), _flush_stdout(flush_stdout)
76{
77 checkConsistency(dof_table_bulk,
79
80 std::vector<MeshLib::Node*> const& bc_nodes =
82 MeshLib::MeshSubset bc_mesh_subset(_bc_data.bc_or_st_mesh, bc_nodes);
83
85 dof_table_bulk.deriveBoundaryConstrainedMap(std::move(bc_mesh_subset));
86
89 bulk_mesh_dimension, _bc_data.bc_or_st_mesh.getElements(),
91 NumLib::IntegrationOrder{integration_order},
93}
94
96 const double t, GlobalVector const& x,
98{
99 FlushStdoutGuard guard(_flush_stdout);
100 (void)guard;
101
102 auto const& boundary_mesh = _bc_data.bc_or_st_mesh;
103 auto const& boundary_nodes = boundary_mesh.getNodes();
104 auto const* bc_object = _bc_data.bc_or_st_object;
105
106 initBCValues(bc_values, boundary_nodes.size());
107
108 std::vector<double> primary_variables;
109
110 try
111 {
112 for (auto const* boundary_node : boundary_nodes)
113 {
114 auto const boundary_node_id = boundary_node->getID();
115 auto const dof_idx = getDofIdx(boundary_node_id);
116
117 if (dof_idx == NumLib::MeshComponentMap::nop)
118 {
119 // This d.o.f. has not been found. This can be the case, e.g.,
120 // for Taylor-Hood elements, where the lower order field has no
121 // d.o.f. on higher order nodes.
122 continue;
123 }
124
125 if (dof_idx < 0)
126 {
127 // For the DDC approach (e.g. with PETSc option) a negative
128 // index means that this entry is a ghost entry and should be
129 // dropped.
130 continue;
131 }
132
133 collectPrimaryVariables(primary_variables, *boundary_node, x);
134
135 auto const [apply_bc, bc_value] = bc_object->getDirichletBCValue(
136 t,
137 {(*boundary_node)[0], (*boundary_node)[1], (*boundary_node)[2]},
138 boundary_node_id, primary_variables);
139
140 if (!bc_object->isOverriddenEssential())
141 {
142 DBUG(
143 "Method `getDirichletBCValue' not overridden in Python "
144 "script.");
145 return;
146 }
147
148 if (!apply_bc)
149 {
150 continue;
151 }
152
153 bc_values.ids.emplace_back(dof_idx);
154 bc_values.values.emplace_back(bc_value);
155 }
156 }
157 catch (pybind11::error_already_set const& e)
158 {
159 OGS_FATAL("Error evaluating Dirichlet BC in Python: {}", e.what());
160 }
161}
162
164 std::size_t const boundary_node_id) const
165{
167 MeshLib::MeshItemType::Node, boundary_node_id};
168 return _dof_table_boundary->getGlobalIndex(loc,
170}
171
173 std::size_t const boundary_node_id, int const var, int const comp) const
174{
176 MeshLib::MeshItemType::Node, boundary_node_id};
177 return _dof_table_boundary->getGlobalIndex(loc, var, comp);
178}
179
181 std::vector<double>& primary_variables, MeshLib::Node const& boundary_node,
182 GlobalVector const& x) const
183{
184 primary_variables.clear();
185 auto const num_var = _dof_table_boundary->getNumberOfVariables();
186 auto const boundary_node_id = boundary_node.getID();
187
188 for (int var = 0; var < num_var; ++var)
189 {
190 auto const num_comp =
191 _dof_table_boundary->getNumberOfVariableComponents(var);
192 for (int comp = 0; comp < num_comp; ++comp)
193 {
194 auto const dof_idx = getDofIdx(boundary_node_id, var, comp);
195
196 double const pv_value =
198 ? x[dof_idx]
199 : interpolateToHigherOrderNode(x, var, comp, boundary_node);
200
201 primary_variables.push_back(pv_value);
202 }
203 }
204}
205
207 GlobalVector const& x, int const var, int const comp,
208 MeshLib::Node const& boundary_node) const
209{
210 auto const& boundary_elements =
212 boundary_node.getID());
213
214 if (boundary_elements.size() != 1)
215 {
216 DBUG(
217 "Boundary node {} is associated with {} elements in the boundary "
218 "mesh.",
219 boundary_node.getID(),
220 boundary_elements.size());
221 }
222
223 // Interpolations on all associated elements should return the same. Just
224 // pick any of them.
225 auto const& boundary_element = *boundary_elements.front();
226
227 assert(boundary_element.getNumberOfBaseNodes() <
228 boundary_element.getNumberOfNodes() &&
229 "We expect that the boundary element is a higher order element. "
230 "Otherwise no interpolation should take place.");
231
232 // Search local node id
233 auto const local_node_id_within_boundary_element =
234 getNodeIDinElement(boundary_element, &boundary_node);
235
236 auto const boundary_element_id = boundary_element.getID();
237
238 // Assumption: all boundary elements have a local assembler (in the same
239 // order)
240 auto const& loc_asm = *_local_assemblers[boundary_element_id];
241
242 return loc_asm.interpolate(local_node_id_within_boundary_element,
243 *_dof_table_boundary, x, var, comp);
244}
245
247 const double t, std::vector<GlobalVector*> const& x, int const process_id,
249{
250 FlushStdoutGuard guard(_flush_stdout);
251
252 try
253 {
256 _local_assemblers, *_dof_table_boundary, t, x, process_id, K, b,
257 Jac);
258 }
260 {
261 DBUG("Method `getFlux' not overridden in Python script.");
262 }
263 catch (pybind11::error_already_set const& e)
264 {
265 OGS_FATAL("Error evaluating natural BC in Python: {}", e.what());
266 }
267}
268
269std::unique_ptr<PythonBoundaryCondition> createPythonBoundaryCondition(
270 BaseLib::ConfigTree const& config, MeshLib::Mesh const& boundary_mesh,
271 NumLib::LocalToGlobalIndexMap const& dof_table_bulk,
272 MeshLib::Mesh const& bulk_mesh, int const variable_id,
273 int const component_id, unsigned const integration_order,
274 unsigned const shapefunction_order,
275 std::vector<std::reference_wrapper<ProcessVariable>> const&
276 all_process_variables_for_this_process)
277{
279 config.checkConfigParameter("type", "Python");
280
282 auto const bc_object = config.getConfigParameter<std::string>("bc_object");
284 auto const flush_stdout = config.getConfigParameter("flush_stdout", false);
285
286 // Evaluate Python code in scope of main module
287 pybind11::object scope =
288 pybind11::module::import("__main__").attr("__dict__");
289
290 if (!scope.contains(bc_object))
291 {
292 OGS_FATAL(
293 "Function `{:s}' is not defined in the python script file, or "
294 "there was no python script file specified.",
295 bc_object);
296 }
297
298 auto* bc = scope[bc_object.c_str()]
300
301 if (variable_id >=
302 static_cast<int>(dof_table_bulk.getNumberOfVariables()) ||
303 component_id >=
304 dof_table_bulk.getNumberOfVariableComponents(variable_id))
305 {
306 OGS_FATAL(
307 "Variable id or component id too high. Actual values: ({:d}, "
308 "{:d}), maximum values: ({:d}, {:d}).",
309 variable_id, component_id, dof_table_bulk.getNumberOfVariables(),
310 dof_table_bulk.getNumberOfVariableComponents(variable_id));
311 }
312
313 // In case of partitioned mesh the boundary could be empty, i.e. there is no
314 // boundary condition.
315#ifdef USE_PETSC
316 // This can be extracted to createBoundaryCondition() but then the config
317 // parameters are not read and will cause an error.
318 // TODO (naumov): Add a function to ConfigTree for skipping the tags of the
319 // subtree and move the code up in createBoundaryCondition().
320 if (boundary_mesh.getDimension() == 0 &&
321 boundary_mesh.getNumberOfNodes() == 0 &&
322 boundary_mesh.getNumberOfElements() == 0)
323 {
324 return nullptr;
325 }
326#endif // USE_PETSC
327
328 return std::make_unique<PythonBoundaryCondition>(
330 {bc, dof_table_bulk.getGlobalComponent(variable_id, component_id),
331 boundary_mesh, all_process_variables_for_this_process,
332 shapefunction_order}},
333 integration_order, flush_stdout, bulk_mesh.getDimension(),
334 dof_table_bulk);
335}
336
337} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
GlobalMatrix::IndexType GlobalIndexType
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
std::string getName(std::string const &line)
Returns the name/title from the "Zone"-description.
T getConfigParameter(std::string const &param) const
void checkConfigParameter(std::string const &param, std::string_view const value) const
Global vector based on Eigen vector.
Definition EigenVector.h:25
std::size_t getID() const
A subset of nodes on a single mesh.
Definition MeshSubset.h:26
bool isAxiallySymmetric() const
Definition Mesh.h:137
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition Mesh.h:106
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:109
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition Mesh.h:88
std::size_t getID() const
Get id of the mesh.
Definition Mesh.h:121
std::size_t getNumberOfNodes() const
Get the number of nodes.
Definition Mesh.h:100
std::vector< Element const * > const & getElementsConnectedToNode(std::size_t node_id) const
Definition Mesh.cpp:256
std::size_t getNumberOfElements() const
Get the number of elements.
Definition Mesh.h:97
int getGlobalComponent(int const variable_id, int const component_id) const
int getNumberOfVariableComponents(int variable_id) const
std::unique_ptr< LocalToGlobalIndexMap > deriveBoundaryConstrainedMap(int const variable_id, std::vector< int > const &component_ids, MeshLib::MeshSubset &&new_mesh_subset) const
static constexpr NUMLIB_EXPORT GlobalIndexType const nop
virtual void assemble(std::size_t const id, NumLib::LocalToGlobalIndexMap const &dof_table_boundary, double const t, std::vector< GlobalVector * > const &x, int const process_id, GlobalMatrix *K, GlobalVector &b, GlobalMatrix *Jac)=0
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _dof_table_boundary
Local dof table for the boundary mesh.
double interpolateToHigherOrderNode(GlobalVector const &x, int const var, int const comp, MeshLib::Node const &boundary_node) const
PythonBcData _bc_data
Auxiliary data used by the local assemblers.
void applyNaturalBC(const double t, std::vector< GlobalVector * > const &x, int const process_id, GlobalMatrix *K, GlobalVector &b, GlobalMatrix *Jac) override
std::vector< std::unique_ptr< PythonBoundaryConditionLocalAssemblerInterface > > _local_assemblers
Local assemblers for all elements of the boundary mesh.
void getEssentialBCValues(const double t, const GlobalVector &x, NumLib::IndexValueVector< GlobalIndexType > &bc_values) const override
Writes the values of essential BCs to bc_values.
GlobalIndexType getDofIdx(std::size_t const boundary_node_id) const
void collectPrimaryVariables(std::vector< double > &primary_variables, MeshLib::Node const &boundary_node, GlobalVector const &x) const
PythonBoundaryCondition(PythonBcData &&bc_data, unsigned const integration_order, bool const flush_stdout, unsigned const bulk_mesh_dimension, NumLib::LocalToGlobalIndexMap const &dof_table_bulk)
void createLocalAssemblersPython(const unsigned dimension, std::vector< MeshLib::Element * > const &mesh_elements, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< std::unique_ptr< LocalAssemblerInterface > > &local_assemblers, NumLib::IntegrationOrder const integration_order, ExtraCtorArgs &&... extra_ctor_args)
std::unique_ptr< PythonBoundaryCondition > createPythonBoundaryCondition(BaseLib::ConfigTree const &config, MeshLib::Mesh const &boundary_mesh, NumLib::LocalToGlobalIndexMap const &dof_table_bulk, MeshLib::Mesh const &bulk_mesh, int const variable_id, int const component_id, unsigned const integration_order, unsigned const shapefunction_order, std::vector< std::reference_wrapper< ProcessVariable > > const &all_process_variables_for_this_process)
Creates a new PythonBoundaryCondition object.
void checkConsistency(NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< std::reference_wrapper< ProcessLib::ProcessVariable > > const &pvs)
void initBCValues(NumLib::IndexValueVector< GlobalIndexType > &bc_values, std::size_t const nnodes)
std::vector< IndexType > ids
std::vector< double > values
static void executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)
MeshLib::Mesh const & bc_or_st_mesh
The domain, where this BC or ST will be applied.
Definition BcOrStData.h:46
std::vector< std::reference_wrapper< ProcessVariable > > const & all_process_variables_for_this_process
Definition BcOrStData.h:51
BcOrStPythonSideInterface const *const bc_or_st_object
Python object computing BC or ST values.
Definition BcOrStData.h:39