OGS
ProcessLib::HeatTransportBHE Namespace Reference

Namespaces

namespace  BHE
namespace  detail

Classes

struct  AlgebraicBCSetting
class  BHEBottomDirichletBoundaryCondition
class  BHEInflowDirichletBoundaryCondition
struct  BHEMeshData
class  HeatTransportBHELocalAssemblerBHE
class  HeatTransportBHELocalAssemblerInterface
class  HeatTransportBHELocalAssemblerSoil
class  HeatTransportBHEProcess
struct  HeatTransportBHEProcessData
struct  IntegrationPointDataBHE
struct  IntegrationPointDataSoil
class  LocalDataInitializer
struct  SecondaryData

Typedefs

using BHECreatorFunc

Functions

BHEMeshData getBHEDataInMesh (MeshLib::Mesh const &mesh)
std::unique_ptr< BHEBottomDirichletBoundaryConditioncreateBHEBottomDirichletBoundaryCondition (std::pair< GlobalIndexType, GlobalIndexType > &&in_out_global_indices)
template<typename BHEUpdateCallback>
std::unique_ptr< BHEInflowDirichletBoundaryCondition< BHEUpdateCallback > > createBHEInflowDirichletBoundaryCondition (std::pair< GlobalIndexType, GlobalIndexType > &&in_out_global_indices, BHEUpdateCallback bhe_update_callback)
void createAndInsertBHE (const std::string &bhe_type, const std::vector< int > &bhe_ids_of_this_bhe, const BaseLib::ConfigTree &bhe_config, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > &parameters, std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation > > const &curves, BHEMeshData const &bhe_mesh_data, std::map< int, BHE::BHETypes > &bhes_map)
std::unique_ptr< ProcesscreateHeatTransportBHEProcess (std::string const &name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< ProcessVariable > const &variables, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > &parameters, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation > > const &curves, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template<template< typename > class LocalAssemblerSoilImplementation, template< typename, typename > class LocalAssemblerBHEImplementation, typename LocalAssemblerInterface, typename... ExtraCtorArgs>
void createLocalAssemblers (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)

Variables

std::map< std::string_view, BHECreatorFuncbheCreators

Typedef Documentation

◆ BHECreatorFunc

Initial value:
std::function<BHE::BHETypes(
std::vector<std::unique_ptr<ParameterLib::ParameterBase>>&,
std::map<std::string,
std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&,
std::vector<MeshLib::Node*> const&)>
std::variant< BHE_1U, BHE_CXA, BHE_CXC, BHE_2U, BHE_1P > BHETypes
Definition BHETypes.h:19

Definition at line 29 of file CreateHeatTransportBHEProcess.cpp.

Function Documentation

◆ createAndInsertBHE()

void ProcessLib::HeatTransportBHE::createAndInsertBHE ( const std::string & bhe_type,
const std::vector< int > & bhe_ids_of_this_bhe,
const BaseLib::ConfigTree & bhe_config,
std::vector< std::unique_ptr< ParameterLib::ParameterBase > > & parameters,
std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation > > const & curves,
BHEMeshData const & bhe_mesh_data,
std::map< int, BHE::BHETypes > & bhes_map )

Definition at line 67 of file CreateHeatTransportBHEProcess.cpp.

77{
78 auto bhe_creator_it = bheCreators.find(bhe_type);
79 if (bhe_creator_it == bheCreators.end())
80 {
81 OGS_FATAL("Unknown BHE type: {:s}", bhe_type);
82 }
83 for (auto const& id : bhe_ids_of_this_bhe)
84 {
85 if (id < 0 || id >= static_cast<int>(bhe_mesh_data.BHE_nodes.size()))
86 {
88 "BHE id {:d} is out of range. The mesh contains {:d} "
89 "BHE(s) (valid ids: 0 to {:d}).",
90 id, bhe_mesh_data.BHE_nodes.size(),
91 static_cast<int>(bhe_mesh_data.BHE_nodes.size()) - 1);
92 }
93 std::pair<std::map<int, BHE::BHETypes>::iterator, bool> result;
94 if (id == bhe_ids_of_this_bhe[0])
95 {
96 auto const& bhe_nodes =
97 bhe_mesh_data.BHE_topology_ordered_nodes[id];
98 result = bhes_map.try_emplace(
99 id, bhe_creator_it->second(bhe_config, parameters, curves,
100 bhe_nodes));
101 }
102 else
103 {
104 // ConfigTree enforces single-read semantics, so we cannot
105 // re-parse the same config block. Grouped BHEs share all
106 // configuration except the borehole geometry, which is rebuilt
107 // per-ID from each BHE's own node set.
108 auto const& first_bhe =
109 bhes_map.find(bhe_ids_of_this_bhe[0])->second;
110 auto const& bhe_nodes =
111 bhe_mesh_data.BHE_topology_ordered_nodes[id];
112 result = bhes_map.try_emplace(
113 id,
114 std::visit(
115 [&](auto const& bhe) -> BHE::BHETypes
116 {
117 return bhe.withGeometry(
118 bhe.borehole_geometry.rebuildForNodes(bhe_nodes));
119 },
120 first_bhe));
121 }
122 if (!result.second)
123 {
124 OGS_FATAL(
125 "BHE with id '{:d}' is already present in the list! Check for "
126 "duplicate definitions of BHE ids.",
127 id);
128 }
129 }
130}
#define OGS_FATAL(...)
Definition Error.h:19
std::map< std::string_view, BHECreatorFunc > bheCreators

References ProcessLib::HeatTransportBHE::BHEMeshData::BHE_nodes, ProcessLib::HeatTransportBHE::BHEMeshData::BHE_topology_ordered_nodes, bheCreators, and OGS_FATAL.

Referenced by createHeatTransportBHEProcess().

◆ createBHEBottomDirichletBoundaryCondition()

std::unique_ptr< BHEBottomDirichletBoundaryCondition > ProcessLib::HeatTransportBHE::createBHEBottomDirichletBoundaryCondition ( std::pair< GlobalIndexType, GlobalIndexType > && in_out_global_indices)

Definition at line 25 of file BHEBottomDirichletBoundaryCondition.cpp.

27{
28 DBUG("Constructing BHEBottomDirichletBoundaryCondition.");
29
30 // In case of partitioned mesh the boundary could be empty, i.e. there is no
31 // boundary condition.
32#ifdef USE_PETSC
33 // For this special boundary condition the boundary condition is not empty
34 // if the global indices are non-negative.
35 if (in_out_global_indices.first < 0 && in_out_global_indices.second < 0)
36 {
37 return nullptr;
38 }
39 // If only one of the global indices (in or out) is negative the
40 // implementation is not valid.
41 if (in_out_global_indices.first < 0 || in_out_global_indices.second < 0)
42 {
44 "The partition cuts the BHE into two independent parts. This "
45 "behaviour is not implemented.");
46 }
47#endif // USE_PETSC
48
49 return std::make_unique<BHEBottomDirichletBoundaryCondition>(
50 std::move(in_out_global_indices));
51}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22

References DBUG(), and OGS_FATAL.

Referenced by ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::createBHEBoundaryConditionTopBottom().

◆ createBHEInflowDirichletBoundaryCondition()

template<typename BHEUpdateCallback>
std::unique_ptr< BHEInflowDirichletBoundaryCondition< BHEUpdateCallback > > ProcessLib::HeatTransportBHE::createBHEInflowDirichletBoundaryCondition ( std::pair< GlobalIndexType, GlobalIndexType > && in_out_global_indices,
BHEUpdateCallback bhe_update_callback )

Definition at line 45 of file BHEInflowDirichletBoundaryCondition.h.

48{
49 DBUG("Constructing BHEInflowDirichletBoundaryCondition.");
50
51 // In case of partitioned mesh the boundary could be empty, i.e. there is no
52 // boundary condition.
53#ifdef USE_PETSC
54 // For this special boundary condition the boundary condition is not empty
55 // if the global indices are non-negative.
56 if (in_out_global_indices.first < 0 && in_out_global_indices.second < 0)
57 {
58 return nullptr;
59 }
60 // If only one of the global indices (in or out) is negative the
61 // implementation is not valid.
62 if (in_out_global_indices.first < 0 || in_out_global_indices.second < 0)
63 {
65 "The partition cuts the BHE into two independent parts. This "
66 "behaviour is not implemented.");
67 }
68#endif // USE_PETSC
69
70 return std::make_unique<
71 BHEInflowDirichletBoundaryCondition<BHEUpdateCallback>>(
72 std::move(in_out_global_indices), bhe_update_callback);
73}

References DBUG(), and OGS_FATAL.

Referenced by ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::createBHEBoundaryConditionTopBottom().

◆ createHeatTransportBHEProcess()

std::unique_ptr< Process > ProcessLib::HeatTransportBHE::createHeatTransportBHEProcess ( std::string const & name,
MeshLib::Mesh & mesh,
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > && jacobian_assembler,
std::vector< ProcessVariable > const & variables,
std::vector< std::unique_ptr< ParameterLib::ParameterBase > > & parameters,
unsigned const integration_order,
BaseLib::ConfigTree const & config,
std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation > > const & curves,
std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & media )
Input File Parameter
prj__processes__process__type

Process Variables

Input File Parameter
prj__processes__process__HEAT_TRANSPORT_BHE__process_variables

Primary process variables as they appear in the global component vector:

Input File Parameter
prj__processes__process__HEAT_TRANSPORT_BHE__process_variables__process_variable

Process Parameters

Input File Parameter
prj__processes__process__HEAT_TRANSPORT_BHE__borehole_heat_exchangers
Input File Parameter
prj__processes__process__HEAT_TRANSPORT_BHE__use_server_communication
Input File Parameter
prj__processes__process__HEAT_TRANSPORT_BHE__mass_lumping
Input File Parameter
prj__processes__process__HEAT_TRANSPORT_BHE__use_algebraic_bc
Input File Parameter
prj__processes__process__HEAT_TRANSPORT_BHE__weighting_factor
Input File Parameter
prj__processes__process__linear
Input File Parameter
prj__processes__process__HEAT_TRANSPORT_BHE__borehole_heat_exchangers__borehole_heat_exchanger
Input File Parameter
prj__processes__process__HEAT_TRANSPORT_BHE__borehole_heat_exchangers__borehole_heat_exchanger__id
Input File Parameter
prj__processes__process__HEAT_TRANSPORT_BHE__borehole_heat_exchangers__borehole_heat_exchanger__type

Python object computing BC values.

Definition at line 132 of file CreateHeatTransportBHEProcess.cpp.

144{
146 config.checkConfigParameter("type", "HEAT_TRANSPORT_BHE");
147
148 DBUG("Create HeatTransportBHE Process.");
149
151
153 auto const pv_config = config.getConfigSubtree("process_variables");
154 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>
155 process_variables;
156
157 // reading primary variables for each
158 // BHE----------------------------------------------------------
160 auto range =
162 pv_config.getConfigParameterList<std::string>("process_variable");
163 std::vector<std::reference_wrapper<ProcessVariable>> per_process_variables;
164
165 for (std::string const& pv_name : range)
166 {
167 if (pv_name != "temperature_soil" &&
168 pv_name.find("temperature_BHE") == std::string::npos)
169 {
170 OGS_FATAL(
171 "Found a process variable name '{}'. It should be "
172 "'temperature_soil' or 'temperature_BHE_X'",
173 pv_name);
174 }
175 auto variable = std::find_if(variables.cbegin(), variables.cend(),
176 [&pv_name](ProcessVariable const& v)
177 { return v.getName() == pv_name; });
178
179 if (variable == variables.end())
180 {
181 OGS_FATAL(
182 "Could not find process variable '{:s}' in the provided "
183 "variables list for config tag <{:s}>.",
184 pv_name, "process_variable");
185 }
186 DBUG("Found process variable '{:s}' for config tag <{:s}>.",
187 variable->getName(), "process_variable");
188
189 per_process_variables.emplace_back(
190 const_cast<ProcessVariable&>(*variable));
191 }
192 process_variables.push_back(std::move(per_process_variables));
193 // end of reading primary variables for each
194 // BHE----------------------------------------------------------
195
197 // reading BHE parameters --------------------------------------------------
198
199 auto bhe_mesh_data = getBHEDataInMesh(mesh);
200
201 auto const& bhe_configs =
203 config.getConfigSubtree("borehole_heat_exchangers");
204
205 auto const using_server_communication =
207 config.getConfigParameter<bool>("use_server_communication", false);
208
209 auto const mass_lumping =
211 config.getConfigParameter<bool>("mass_lumping", false);
212
213 auto const using_algebraic_bc =
215 config.getConfigParameter<bool>("use_algebraic_bc", false);
216
217 auto const weighting_factor =
219 config.getConfigParameter<float>("weighting_factor", 1000.0);
220
221 auto const is_linear =
223 config.getConfigParameter<bool>("linear", false);
224 if (is_linear)
225 {
226 if (!using_algebraic_bc)
227 {
228 WARN(
229 "You specified that the process simulated by OGS is linear. "
230 "With that optimization the process will be assembled only "
231 "once and the non-linear solver will only do iterations per "
232 "time step to fulfill the BHE boundary conditions. No other "
233 "non-linearities will be resolved and OGS will not detect if "
234 "there are any non-linearities. It is your responsibility to "
235 "ensure that the assembled equation systems are linear, "
236 "indeed! There is no safety net!");
237 }
238 else
239 {
240 WARN(
241 "You specified that the process simulated by OGS is linear. "
242 "With that optimization the process will be assembled only "
243 "once and the non-linear solver will do only one iteration per "
244 "time step. No non-linearities will be resolved and OGS will "
245 "not detect if there are any non-linearities. It is your "
246 "responsibility to ensure that the assembled equation systems "
247 "are linear, indeed! There is no safety net!");
248 }
249 }
250
251 std::map<int, BHE::BHETypes> bhes_map;
252
253 int bhe_iterator = 0;
254
255 for (
256 auto const& bhe_config :
258 bhe_configs.getConfigSubtreeList("borehole_heat_exchanger"))
259 {
260 auto const bhe_id_string =
262 bhe_config.getConfigAttribute<std::string>(
263 "id", std::to_string(bhe_iterator));
264
265 std::vector<int> bhe_ids_of_this_bhe;
266
267 if (bhe_id_string == "*")
268 {
269 int const size = static_cast<int>(bhe_mesh_data.BHE_mat_IDs.size());
270 bhe_ids_of_this_bhe.resize(size);
271 std::iota(bhe_ids_of_this_bhe.begin(), bhe_ids_of_this_bhe.end(),
272 0);
273 }
274 else
275 {
276 bhe_ids_of_this_bhe =
278 }
279
280 // read in the parameters
281 const std::string bhe_type =
283 bhe_config.getConfigParameter<std::string>("type");
284
285 createAndInsertBHE(bhe_type, bhe_ids_of_this_bhe, bhe_config,
286 parameters, curves, bhe_mesh_data, bhes_map);
287 bhe_iterator++;
288 }
289
290 if (static_cast<int>(bhes_map.size()) - 1 != bhes_map.rbegin()->first)
291 {
292 OGS_FATAL(
293 "The maximum given BHE id '{:d}' did not match the number of given "
294 "BHE definitions '{:d}'. The BHE ids needs to be defined starting "
295 "from 0, so the maximum BHE id needs to be number of BHE "
296 "definitions minus 1. After all definitions there are no gaps "
297 "allowed between the given ids.",
298 bhes_map.rbegin()->first, bhes_map.size());
299 }
300
301 std::vector<BHE::BHETypes> bhes;
302 bhes.reserve(bhes_map.size());
303 std::ranges::copy(bhes_map | std::views::values, std::back_inserter(bhes));
304 bhe_mesh_data.updateElementSectionIndices(bhes);
305 // end of reading BHE parameters
306 // -------------------------------------------
307
308 auto media_map =
310
311 // find if bhe uses python boundary condition
312 auto const using_tespy =
313 visit([](auto const& bhe) { return bhe.use_python_bcs; }, bhes[0]);
314
316 BHEInflowPythonBoundaryConditionPythonSideInterface* py_object = nullptr;
317 // create a pythonBoundaryCondition object
318 if (using_tespy || using_server_communication)
319 {
320 // Evaluate Python code in scope of main module
321 pybind11::object scope =
322 pybind11::module::import("__main__").attr("__dict__");
323
324 if (!scope.contains("bc_bhe"))
325 OGS_FATAL(
326 "Function 'bc_bhe' is not defined in the python script file, "
327 "or there was no python script file specified.");
328
329 py_object =
330 scope["bc_bhe"]
331 .cast<BHEInflowPythonBoundaryConditionPythonSideInterface*>();
332
333 if (py_object == nullptr)
334 OGS_FATAL(
335 "Not able to access the correct bc pointer from python script "
336 "file specified.");
337
338 // create BHE network dataframe from Python
339 py_object->dataframe_network = py_object->initializeDataContainer();
340 if (!py_object->isOverriddenEssential())
341 {
342 DBUG(
343 "Method `initializeDataContainer' not overridden in Python "
344 "script.");
345 }
346 // clear ogs bc_node_id memory in dataframe
347 std::get<3>(py_object->dataframe_network).clear(); // ogs_bc_node_id
348
349 // here calls the tespyHydroSolver to get the pipe flow velocity in bhe
350 // network
351 /* for 2U type the flowrate initialization process below causes conflict
352 // replace the value in flow velocity Matrix _u
353 auto const tespy_flow_rate = std::get<4>(py_object->dataframe_network);
354 const std::size_t n_bhe = tespy_flow_rate.size();
355 if (bhes.size() != n_bhe)
356 OGS_FATAL(
357 "The number of BHEs defined in OGS and TESPy are not the "
358 "same!");
359
360 for (std::size_t idx_bhe = 0; idx_bhe < n_bhe; idx_bhe++)
361 {
362 // the flow_rate in OGS should be updated from the flow_rate
363 // computed by TESPy.
364 auto update_flow_rate = [&](auto& bhe) {
365 bhe.updateHeatTransferCoefficients(tespy_flow_rate[idx_bhe]);
366 };
367 visit(update_flow_rate, bhes[idx_bhe]);
368 }
369 */
370 }
371
372 HeatTransportBHEProcessData process_data(
373 std::move(media_map), std::move(bhes), py_object, using_tespy,
374 using_server_communication, mass_lumping,
375 {using_algebraic_bc, weighting_factor}, is_linear);
376
377 SecondaryVariableCollection secondary_variables;
378
379 ProcessLib::createSecondaryVariables(config, secondary_variables);
380
381 return std::make_unique<HeatTransportBHEProcess>(
382 std::move(name), mesh, std::move(jacobian_assembler), parameters,
383 integration_order, std::move(process_variables),
384 std::move(process_data), std::move(secondary_variables),
385 std::move(bhe_mesh_data));
386}
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:34
std::vector< int > splitMaterialIdString(std::string const &material_id_string)
MaterialSpatialDistributionMap createMaterialSpatialDistributionMap(std::map< int, std::shared_ptr< Medium > > const &media, MeshLib::Mesh const &mesh)
constexpr int size(int const displacement_dim)
Vectorized tensor size for given displacement dimension.
void createAndInsertBHE(const std::string &bhe_type, const std::vector< int > &bhe_ids_of_this_bhe, const BaseLib::ConfigTree &bhe_config, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > &parameters, std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation > > const &curves, BHEMeshData const &bhe_mesh_data, std::map< int, BHE::BHETypes > &bhes_map)
BHEMeshData getBHEDataInMesh(MeshLib::Mesh const &mesh)
void createSecondaryVariables(BaseLib::ConfigTree const &config, SecondaryVariableCollection &secondary_variables)

References BaseLib::ConfigTree::checkConfigParameter(), createAndInsertBHE(), MaterialPropertyLib::createMaterialSpatialDistributionMap(), ProcessLib::createSecondaryVariables(), DBUG(), getBHEDataInMesh(), BaseLib::ConfigTree::getConfigParameter(), BaseLib::ConfigTree::getConfigParameterList(), BaseLib::ConfigTree::getConfigSubtree(), OGS_FATAL, MaterialLib::splitMaterialIdString(), and WARN().

Referenced by ProjectData::parseProcesses().

◆ createLocalAssemblers()

template<template< typename > class LocalAssemblerSoilImplementation, template< typename, typename > class LocalAssemblerBHEImplementation, typename LocalAssemblerInterface, typename... ExtraCtorArgs>
void ProcessLib::HeatTransportBHE::createLocalAssemblers ( 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 )

Creates local assemblers for each element of the given mesh.

Template Parameters
LocalAssemblerImplementationthe individual local assembler type
LocalAssemblerInterfacethe general local assembler interface
ExtraCtorArgstypes of additional constructor arguments. Those arguments will be passed to the constructor of LocalAssemblerImplementation.

The first two template parameters cannot be deduced from the arguments. Therefore they always have to be provided manually.

Definition at line 68 of file HeatTransportBHE/LocalAssemblers/CreateLocalAssemblers.h.

74{
75 DBUG("Create local assemblers for the HeatTransportBHE process.");
76
77 detail::createLocalAssemblers<LocalAssemblerSoilImplementation,
78 LocalAssemblerBHEImplementation>(
79 dof_table, mesh_elements, local_assemblers, integration_order,
80 std::forward<ExtraCtorArgs>(extra_ctor_args)...);
81}
void createLocalAssemblers(NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< MeshLib::Element * > const &mesh_elements, std::vector< std::unique_ptr< LocalAssemblerInterface > > &local_assemblers, NumLib::IntegrationOrder const integration_order, std::unordered_map< std::size_t, BHE::BHETypes * > const &element_to_bhe_map, ExtraCtorArgs &&... extra_ctor_args)

References createLocalAssemblers(), ProcessLib::HeatTransportBHE::detail::createLocalAssemblers(), and DBUG().

Referenced by createLocalAssemblers(), and ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::initializeConcreteProcess().

◆ getBHEDataInMesh()

BHEMeshData ProcessLib::HeatTransportBHE::getBHEDataInMesh ( MeshLib::Mesh const & mesh)

get data about fracture and matrix elements/nodes from a mesh

Parameters
meshA mesh which includes BHE elements, i.e. 1-dimensional elements. It is assumed that elements forming a BHE have a distinct material ID.

Definition at line 215 of file HeatTransportBHE/BHE/MeshUtils.cpp.

216{
217 std::vector<MeshLib::Element*> const all_bhe_elements =
218 extractOneDimensionalElements(mesh.getElements());
219
220 // finally counting two types of elements
221 // They are (i) soil, and (ii) BHE type of elements
222 DBUG("-> found total {:d} soil elements and {:d} BHE elements",
223 mesh.getNumberOfElements() - all_bhe_elements.size(),
224 all_bhe_elements.size());
225
226 // get BHE material IDs
227 auto const* const opt_material_ids = MeshLib::materialIDs(mesh);
228 if (opt_material_ids == nullptr)
229 {
230 OGS_FATAL("Not able to get material IDs! ");
231 }
232 auto const& material_ids = *opt_material_ids;
233
234 auto const& bhe_material_ids =
235 getUniqueMaterialIds(material_ids, all_bhe_elements);
236 DBUG("-> found {:d} BHE material groups", bhe_material_ids.size());
237
238 // create a vector of BHE elements for each group
239 std::vector<std::vector<MeshLib::Element*>> bhe_elements;
240 bhe_elements.resize(bhe_material_ids.size());
241 for (unsigned bhe_id = 0; bhe_id < bhe_material_ids.size(); bhe_id++)
242 {
243 const auto bhe_mat_id = bhe_material_ids[bhe_id];
244 std::vector<MeshLib::Element*>& vec_elements = bhe_elements[bhe_id];
245 copy_if(begin(all_bhe_elements), end(all_bhe_elements),
246 back_inserter(vec_elements),
247 [&](MeshLib::Element const* const e)
248 { return material_ids[e->getID()] == bhe_mat_id; });
249 DBUG("-> found {:d} elements on the BHE_{:d}", vec_elements.size(),
250 bhe_id);
251 }
252
253 // get a vector of BHE nodes
254 std::vector<std::vector<MeshLib::Node*>> bhe_nodes;
255 bhe_nodes.resize(bhe_material_ids.size());
256 for (unsigned bhe_id = 0; bhe_id < bhe_material_ids.size(); bhe_id++)
257 {
258 std::vector<MeshLib::Node*>& vec_nodes = bhe_nodes[bhe_id];
259 for (MeshLib::Element* e : bhe_elements[bhe_id])
260 {
261 for (unsigned i = 0; i < e->getNumberOfNodes(); i++)
262 {
263 vec_nodes.push_back(const_cast<MeshLib::Node*>(e->getNode(i)));
264 }
265 }
268
269 DBUG("-> found {:d} nodes on the BHE_{:d}", vec_nodes.size(), bhe_id);
270 }
271
272 std::unordered_map<std::size_t, double> bhe_element_distances_from_wellhead;
273 std::vector<std::vector<MeshLib::Node*>> bhe_topology_ordered_nodes(
274 bhe_material_ids.size());
275
276 std::unordered_map<std::size_t, int> bhe_element_section_indices;
277
278 for (unsigned bhe_id = 0; bhe_id < bhe_material_ids.size(); bhe_id++)
279 {
280 auto walk_result =
281 walkChainFromWellhead(bhe_elements[bhe_id], bhe_nodes[bhe_id]);
282
283 bhe_elements[bhe_id] = std::move(walk_result.ordered_elements);
284 bhe_topology_ordered_nodes[bhe_id] =
285 std::move(walk_result.ordered_nodes);
286
287 for (auto const& [element_id, distance] :
288 walk_result.element_distances_from_wellhead)
289 {
290 bhe_element_distances_from_wellhead[element_id] = distance;
291 }
292 }
293
294 return {bhe_material_ids,
295 bhe_elements,
296 bhe_nodes,
297 bhe_topology_ordered_nodes,
298 std::move(bhe_element_distances_from_wellhead),
299 std::move(bhe_element_section_indices)};
300}
virtual unsigned getNumberOfNodes() const =0
virtual const Node * getNode(unsigned idx) const =0
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:80
void makeVectorUnique(std::vector< T > &v)
Definition Algorithm.h:173
bool idsComparator(T const a, T const b)
Definition Mesh.h:199
PropertyVector< int > const * materialIDs(Mesh const &mesh)
Definition Mesh.cpp:258
ChainWalkResult walkChainFromWellhead(std::vector< MeshLib::Element * > const &bhe_elements, std::vector< MeshLib::Node * > const &bhe_nodes)

References DBUG(), MeshLib::Mesh::getElements(), MeshLib::Element::getID(), MeshLib::Mesh::getNumberOfElements(), MeshLib::idsComparator(), BaseLib::makeVectorUnique(), MeshLib::materialIDs(), and OGS_FATAL.

Referenced by createHeatTransportBHEProcess().

Variable Documentation

◆ bheCreators

std::map<std::string_view, BHECreatorFunc> ProcessLib::HeatTransportBHE::bheCreators
Initial value:
= {
{"1U",
[](auto& config, auto& parameters, auto& curves, auto& bhe_nodes)
{
config, parameters, curves, bhe_nodes));
}},
{"2U",
[](auto& config, auto& parameters, auto& curves, auto& bhe_nodes)
{
config, parameters, curves, bhe_nodes));
}},
{"CXA",
[](auto& config, auto& parameters, auto& curves, auto& bhe_nodes)
{
config, parameters, curves, bhe_nodes));
}},
{"CXC",
[](auto& config, auto& parameters, auto& curves, auto& bhe_nodes)
{
config, parameters, curves, bhe_nodes));
}},
{"1P", [](auto& config, auto& parameters, auto& curves, auto& bhe_nodes)
{
config, parameters, curves, bhe_nodes));
}}}
T_BHE createBHECoaxial(BaseLib::ConfigTree const &config, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > &parameters, std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation > > const &curves, std::vector< MeshLib::Node * > const &bhe_nodes)
T_BHE createBHEUType(BaseLib::ConfigTree const &config, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > &parameters, std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation > > const &curves, std::vector< MeshLib::Node * > const &bhe_nodes)
T_BHE createBHE1PType(BaseLib::ConfigTree const &config, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > &parameters, std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation > > const &curves, std::vector< MeshLib::Node * > const &bhe_nodes)

Definition at line 36 of file CreateHeatTransportBHEProcess.cpp.

36 {
37 {"1U",
38 [](auto& config, auto& parameters, auto& curves, auto& bhe_nodes)
39 {
41 config, parameters, curves, bhe_nodes));
42 }},
43 {"2U",
44 [](auto& config, auto& parameters, auto& curves, auto& bhe_nodes)
45 {
46 return BHE::BHE_2U(BHE::createBHEUType<BHE::BHE_2U>(
47 config, parameters, curves, bhe_nodes));
48 }},
49 {"CXA",
50 [](auto& config, auto& parameters, auto& curves, auto& bhe_nodes)
51 {
53 config, parameters, curves, bhe_nodes));
54 }},
55 {"CXC",
56 [](auto& config, auto& parameters, auto& curves, auto& bhe_nodes)
57 {
59 config, parameters, curves, bhe_nodes));
60 }},
61 {"1P", [](auto& config, auto& parameters, auto& curves, auto& bhe_nodes)
62 {
64 config, parameters, curves, bhe_nodes));
65 }}};

Referenced by createAndInsertBHE().