OGS
ProcessLib::PythonBoundaryCondition Class Referencefinal

Detailed Description

A boundary condition whose values are computed by a Python script.

Definition at line 57 of file PythonBoundaryCondition.h.

#include <PythonBoundaryCondition.h>

Inheritance diagram for ProcessLib::PythonBoundaryCondition:
[legend]
Collaboration diagram for ProcessLib::PythonBoundaryCondition:
[legend]

Public Member Functions

 PythonBoundaryCondition (PythonBcData &&bc_data, unsigned const integration_order, bool const flush_stdout, unsigned const bulk_mesh_dimension, NumLib::LocalToGlobalIndexMap const &dof_table_bulk)
 
void getEssentialBCValues (const double t, const GlobalVector &x, NumLib::IndexValueVector< GlobalIndexType > &bc_values) const override
 Writes the values of essential BCs to bc_values.
 
void applyNaturalBC (const double t, std::vector< GlobalVector * > const &x, int const process_id, GlobalMatrix &K, GlobalVector &b, GlobalMatrix *Jac) override
 
- Public Member Functions inherited from ProcessLib::BoundaryCondition
virtual void preTimestep (const double, std::vector< GlobalVector * > const &, int const)
 
virtual void postTimestep (const double, std::vector< GlobalVector * > const &, int const)
 
virtual ~BoundaryCondition ()=default
 

Private Member Functions

void collectPrimaryVariables (std::vector< double > &primary_variables, MeshLib::Node const &boundary_node, GlobalVector const &x) const
 
GlobalIndexType getDofIdx (std::size_t const boundary_node_id) const
 
GlobalIndexType getDofIdx (std::size_t const boundary_node_id, int const var, int const comp) const
 
double interpolateToHigherOrderNode (GlobalVector const &x, int const var, int const comp, MeshLib::Node const &boundary_node) const
 

Private Attributes

PythonBcData _bc_data
 Auxiliary data used by the local assemblers.
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_dof_table_boundary
 Local dof table for the boundary mesh.
 
std::vector< std::unique_ptr< PythonBoundaryConditionLocalAssemblerInterface > > _local_assemblers
 Local assemblers for all elements of the boundary mesh.
 
bool const _flush_stdout
 

Constructor & Destructor Documentation

◆ PythonBoundaryCondition()

ProcessLib::PythonBoundaryCondition::PythonBoundaryCondition ( PythonBcData && bc_data,
unsigned const integration_order,
bool const flush_stdout,
unsigned const bulk_mesh_dimension,
NumLib::LocalToGlobalIndexMap const & dof_table_bulk )

Definition at line 70 of file PythonBoundaryCondition.cpp.

74 : _bc_data(std::move(bc_data)), _flush_stdout(flush_stdout)
75{
76 checkConsistency(dof_table_bulk,
78
79 std::vector<MeshLib::Node*> const& bc_nodes =
81 MeshLib::MeshSubset bc_mesh_subset(_bc_data.bc_or_st_mesh, bc_nodes);
82
84 dof_table_bulk.deriveBoundaryConstrainedMap(std::move(bc_mesh_subset));
85
87 PythonBoundaryConditionLocalAssembler>(
88 bulk_mesh_dimension, _bc_data.bc_or_st_mesh.getElements(),
90 NumLib::IntegrationOrder{integration_order},
92}
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
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _dof_table_boundary
Local dof table for the boundary mesh.
PythonBcData _bc_data
Auxiliary data used by the local assemblers.
std::vector< std::unique_ptr< PythonBoundaryConditionLocalAssemblerInterface > > _local_assemblers
Local assemblers for all elements of the boundary mesh.
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)
void checkConsistency(NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< std::reference_wrapper< ProcessLib::ProcessVariable > > const &pvs)
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

References _bc_data, _dof_table_boundary, _local_assemblers, ProcessLib::BoundaryConditionAndSourceTerm::Python::BcOrStData< BcOrStPythonSideInterface >::all_process_variables_for_this_process, ProcessLib::BoundaryConditionAndSourceTerm::Python::BcOrStData< BcOrStPythonSideInterface >::bc_or_st_mesh, ProcessLib::BoundaryConditionAndSourceTerm::createLocalAssemblersPython(), NumLib::LocalToGlobalIndexMap::deriveBoundaryConstrainedMap(), MeshLib::Mesh::getElements(), MeshLib::Mesh::getNodes(), and MeshLib::Mesh::isAxiallySymmetric().

Member Function Documentation

◆ applyNaturalBC()

void ProcessLib::PythonBoundaryCondition::applyNaturalBC ( const double ,
std::vector< GlobalVector * > const & ,
int const ,
GlobalMatrix & ,
GlobalVector & ,
GlobalMatrix *  )
overridevirtual

Applies natural BCs (i.e. non-Dirichlet BCs) to the stiffness matrix K and the vector b.

Reimplemented from ProcessLib::BoundaryCondition.

Definition at line 237 of file PythonBoundaryCondition.cpp.

240{
241 FlushStdoutGuard guard(_flush_stdout);
242
243 try
244 {
247 _local_assemblers, *_dof_table_boundary, t, x, process_id, K, b,
248 Jac);
249 }
250 catch (MethodNotOverriddenInDerivedClassException const& /*e*/)
251 {
252 DBUG("Method `getFlux' not overridden in Python script.");
253 }
254}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
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
static void executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)

References _dof_table_boundary, _flush_stdout, _local_assemblers, ProcessLib::GenericNaturalBoundaryConditionLocalAssemblerInterface::assemble(), DBUG(), and NumLib::SerialExecutor::executeMemberOnDereferenced().

◆ collectPrimaryVariables()

void ProcessLib::PythonBoundaryCondition::collectPrimaryVariables ( std::vector< double > & primary_variables,
MeshLib::Node const & boundary_node,
GlobalVector const & x ) const
private

Collects primary variables at the passed node from the passed GlobalVector to primary_variables.

Primary variables at higher order nodes are interpolated from base nodes if necessary, e.g., for Taylor-Hood elements.

Postcondition
primary_variables will contain the value of each component of each primary variable. Their order is determined by the d.o.f. table.

Definition at line 171 of file PythonBoundaryCondition.cpp.

174{
175 primary_variables.clear();
176 auto const num_var = _dof_table_boundary->getNumberOfVariables();
177 auto const boundary_node_id = boundary_node.getID();
178
179 for (int var = 0; var < num_var; ++var)
180 {
181 auto const num_comp =
182 _dof_table_boundary->getNumberOfVariableComponents(var);
183 for (int comp = 0; comp < num_comp; ++comp)
184 {
185 auto const dof_idx = getDofIdx(boundary_node_id, var, comp);
186
187 double const pv_value =
189 ? x[dof_idx]
190 : interpolateToHigherOrderNode(x, var, comp, boundary_node);
191
192 primary_variables.push_back(pv_value);
193 }
194 }
195}
static constexpr NUMLIB_EXPORT GlobalIndexType const nop
double interpolateToHigherOrderNode(GlobalVector const &x, int const var, int const comp, MeshLib::Node const &boundary_node) const
GlobalIndexType getDofIdx(std::size_t const boundary_node_id) const

References _dof_table_boundary, getDofIdx(), MathLib::Point3dWithID::getID(), interpolateToHigherOrderNode(), and NumLib::MeshComponentMap::nop.

Referenced by getEssentialBCValues().

◆ getDofIdx() [1/2]

GlobalIndexType ProcessLib::PythonBoundaryCondition::getDofIdx ( std::size_t const boundary_node_id) const
private

◆ getDofIdx() [2/2]

GlobalIndexType ProcessLib::PythonBoundaryCondition::getDofIdx ( std::size_t const boundary_node_id,
int const var,
int const comp ) const
private

Get the d.o.f. index at the given boundary_node_id for the given variable and component.

Definition at line 163 of file PythonBoundaryCondition.cpp.

165{
167 MeshLib::MeshItemType::Node, boundary_node_id};
168 return _dof_table_boundary->getGlobalIndex(loc, var, comp);
169}

References _bc_data, _dof_table_boundary, ProcessLib::BoundaryConditionAndSourceTerm::Python::BcOrStData< BcOrStPythonSideInterface >::bc_or_st_mesh, MeshLib::Mesh::getID(), and MeshLib::Node.

◆ getEssentialBCValues()

void ProcessLib::PythonBoundaryCondition::getEssentialBCValues ( const double ,
const GlobalVector & ,
NumLib::IndexValueVector< GlobalIndexType > &  ) const
overridevirtual

Writes the values of essential BCs to bc_values.

Reimplemented from ProcessLib::BoundaryCondition.

Definition at line 94 of file PythonBoundaryCondition.cpp.

97{
98 FlushStdoutGuard guard(_flush_stdout);
99 (void)guard;
100
101 auto const& boundary_mesh = _bc_data.bc_or_st_mesh;
102 auto const& boundary_nodes = boundary_mesh.getNodes();
103 auto const* bc_object = _bc_data.bc_or_st_object;
104
105 initBCValues(bc_values, boundary_nodes.size());
106
107 std::vector<double> primary_variables;
108
109 for (auto const* boundary_node : boundary_nodes)
110 {
111 auto const boundary_node_id = boundary_node->getID();
112 auto const dof_idx = getDofIdx(boundary_node_id);
113
114 if (dof_idx == NumLib::MeshComponentMap::nop)
115 {
116 // This d.o.f. has not been found. This can be the case, e.g., for
117 // Taylor-Hood elements, where the lower order field has no d.o.f.
118 // on higher order nodes.
119 continue;
120 }
121
122 if (dof_idx < 0)
123 {
124 // For the DDC approach (e.g. with PETSc option) a negative
125 // index means that this entry is a ghost entry and should be
126 // dropped.
127 continue;
128 }
129
130 collectPrimaryVariables(primary_variables, *boundary_node, x);
131
132 auto const [apply_bc, bc_value] = bc_object->getDirichletBCValue(
133 t, {(*boundary_node)[0], (*boundary_node)[1], (*boundary_node)[2]},
134 boundary_node_id, primary_variables);
135
136 if (!bc_object->isOverriddenEssential())
137 {
138 DBUG(
139 "Method `getDirichletBCValue' not overridden in Python "
140 "script.");
141 return;
142 }
143
144 if (!apply_bc)
145 {
146 continue;
147 }
148
149 bc_values.ids.emplace_back(dof_idx);
150 bc_values.values.emplace_back(bc_value);
151 }
152}
void collectPrimaryVariables(std::vector< double > &primary_variables, MeshLib::Node const &boundary_node, GlobalVector const &x) const
void initBCValues(NumLib::IndexValueVector< GlobalIndexType > &bc_values, std::size_t const nnodes)
BcOrStPythonSideInterface const *const bc_or_st_object
Python object computing BC or ST values.
Definition BcOrStData.h:39

References _bc_data, _flush_stdout, ProcessLib::BoundaryConditionAndSourceTerm::Python::BcOrStData< BcOrStPythonSideInterface >::bc_or_st_mesh, ProcessLib::BoundaryConditionAndSourceTerm::Python::BcOrStData< BcOrStPythonSideInterface >::bc_or_st_object, collectPrimaryVariables(), DBUG(), getDofIdx(), MeshLib::Mesh::getNodes(), NumLib::IndexValueVector< IndexType >::ids, NumLib::MeshComponentMap::nop, and NumLib::IndexValueVector< IndexType >::values.

◆ interpolateToHigherOrderNode()

double ProcessLib::PythonBoundaryCondition::interpolateToHigherOrderNode ( GlobalVector const & x,
int const var,
int const comp,
MeshLib::Node const & boundary_node ) const
private

Interpolates the given component of the given variable to the given boundary_node.

Definition at line 197 of file PythonBoundaryCondition.cpp.

200{
201 auto const& boundary_elements =
203 boundary_node.getID());
204
205 if (boundary_elements.size() != 1)
206 {
207 DBUG(
208 "Boundary node {} is associated with {} elements in the boundary "
209 "mesh.",
210 boundary_node.getID(),
211 boundary_elements.size());
212 }
213
214 // Interpolations on all associated elements should return the same. Just
215 // pick any of them.
216 auto const& boundary_element = *boundary_elements.front();
217
218 assert(boundary_element.getNumberOfBaseNodes() <
219 boundary_element.getNumberOfNodes() &&
220 "We expect that the boundary element is a higher order element. "
221 "Otherwise no interpolation should take place.");
222
223 // Search local node id
224 auto const local_node_id_within_boundary_element =
225 getNodeIDinElement(boundary_element, &boundary_node);
226
227 auto const boundary_element_id = boundary_element.getID();
228
229 // Assumption: all boundary elements have a local assembler (in the same
230 // order)
231 auto const& loc_asm = *_local_assemblers[boundary_element_id];
232
233 return loc_asm.interpolate(local_node_id_within_boundary_element,
234 *_dof_table_boundary, x, var, comp);
235}
std::vector< Element const * > const & getElementsConnectedToNode(std::size_t node_id) const
Definition Mesh.cpp:256
unsigned getNodeIDinElement(Element const &element, const MeshLib::Node *node)
Returns the position of the given node in the node array of this element.
Definition Element.cpp:206

References _bc_data, _dof_table_boundary, _local_assemblers, ProcessLib::BoundaryConditionAndSourceTerm::Python::BcOrStData< BcOrStPythonSideInterface >::bc_or_st_mesh, DBUG(), MeshLib::Mesh::getElementsConnectedToNode(), and MathLib::Point3dWithID::getID().

Referenced by collectPrimaryVariables().

Member Data Documentation

◆ _bc_data

PythonBcData ProcessLib::PythonBoundaryCondition::_bc_data
private

Auxiliary data used by the local assemblers.

Definition at line 104 of file PythonBoundaryCondition.h.

Referenced by PythonBoundaryCondition(), getDofIdx(), getDofIdx(), getEssentialBCValues(), and interpolateToHigherOrderNode().

◆ _dof_table_boundary

std::unique_ptr<NumLib::LocalToGlobalIndexMap> ProcessLib::PythonBoundaryCondition::_dof_table_boundary
private

Local dof table for the boundary mesh.

Definition at line 107 of file PythonBoundaryCondition.h.

Referenced by PythonBoundaryCondition(), applyNaturalBC(), collectPrimaryVariables(), getDofIdx(), getDofIdx(), and interpolateToHigherOrderNode().

◆ _flush_stdout

bool const ProcessLib::PythonBoundaryCondition::_flush_stdout
private

Whether or not to flush standard output before and after each call to Python code. Ensures right order of output messages and therefore simplifies debugging.

Definition at line 116 of file PythonBoundaryCondition.h.

Referenced by applyNaturalBC(), and getEssentialBCValues().

◆ _local_assemblers

std::vector<std::unique_ptr<PythonBoundaryConditionLocalAssemblerInterface> > ProcessLib::PythonBoundaryCondition::_local_assemblers
private

Local assemblers for all elements of the boundary mesh.

Definition at line 111 of file PythonBoundaryCondition.h.

Referenced by PythonBoundaryCondition(), applyNaturalBC(), and interpolateToHigherOrderNode().


The documentation for this class was generated from the following files: