OGS
HeatTransportBHEProcess.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
5
6#include <cassert>
7
17
18namespace ProcessLib
19{
20namespace HeatTransportBHE
21{
23 std::string name,
24 MeshLib::Mesh& mesh,
25 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
26 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
27 unsigned const integration_order,
28 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
29 process_variables,
30 HeatTransportBHEProcessData&& process_data,
31 SecondaryVariableCollection&& secondary_variables,
32 BHEMeshData&& bhe_mesh_data)
33 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
34 integration_order, std::move(process_variables),
35 std::move(secondary_variables)),
37 process_data._is_linear,
38 true /*use_monolithic_scheme*/},
39 _process_data(std::move(process_data)),
40 _bheMeshData(std::move(bhe_mesh_data))
41{
42 // For numerical Jacobian assembler
43 this->_jacobian_assembler->setNonDeformationComponentIDs(
44 {0, 1} /* for two temperature variables */);
45
46 if (_bheMeshData.BHE_mat_IDs.size() !=
47 _process_data._vec_BHE_property.size())
48 {
50 "The number of the given BHE properties ({:d}) are not consistent "
51 "with the number of BHE groups in the mesh ({:d}).",
52 _process_data._vec_BHE_property.size(),
53 _bheMeshData.BHE_mat_IDs.size());
54 }
55
56 auto material_ids = MeshLib::materialIDs(mesh);
57 if (material_ids == nullptr)
58 {
59 OGS_FATAL("Not able to get material IDs! ");
60 }
61
62 _process_data._mesh_prop_materialIDs = material_ids;
63
64 // create a map from a material ID to a BHE ID
65 for (int i = 0; i < static_cast<int>(_bheMeshData.BHE_mat_IDs.size()); i++)
66 {
67 // fill in the map structure
68 _process_data._map_materialID_to_BHE_ID[_bheMeshData.BHE_mat_IDs[i]] =
69 i;
70 }
71}
72
74{
75 // Create single component dof in every of the mesh's nodes.
77 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
78
79 //
80 // Soil temperature variable defined on the whole mesh.
81 //
83 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
84 std::vector<MeshLib::MeshSubset> all_mesh_subsets{*_mesh_subset_soil_nodes};
85
86 std::vector<std::vector<MeshLib::Element*> const*> vec_var_elements;
87 vec_var_elements.push_back(&(_mesh.getElements()));
88
89 std::vector<int> vec_n_components{
90 1}; // one component for the soil temperature variable.
91
92 //
93 // BHE nodes with BHE type dependent number of variables.
94 //
95 int const n_BHEs = _process_data._vec_BHE_property.size();
96 assert(n_BHEs == static_cast<int>(_bheMeshData.BHE_mat_IDs.size()));
97 assert(n_BHEs == static_cast<int>(_bheMeshData.BHE_nodes.size()));
98 assert(n_BHEs == static_cast<int>(_bheMeshData.BHE_elements.size()));
99
100 // the BHE nodes need to be cherry-picked from the vector
101 for (int i = 0; i < n_BHEs; i++)
102 {
103 auto const number_of_unknowns =
104 visit([](auto const& bhe) { return bhe.number_of_unknowns; },
105 _process_data._vec_BHE_property[i]);
106 auto const& bhe_nodes = _bheMeshData.BHE_nodes[i];
107 auto const& bhe_elements = _bheMeshData.BHE_elements[i];
108
109 // All the BHE nodes have additional variables.
110 _mesh_subset_BHE_nodes.push_back(
111 std::make_unique<MeshLib::MeshSubset const>(_mesh, bhe_nodes));
112
113 std::generate_n(std::back_inserter(all_mesh_subsets),
114 // Here the number of components equals to the
115 // number of unknowns on the BHE
116 number_of_unknowns,
117 [&ms = _mesh_subset_BHE_nodes.back()]()
118 { return *ms; });
119
120 vec_n_components.push_back(number_of_unknowns);
121 vec_var_elements.push_back(&bhe_elements);
122 }
123
125 std::make_unique<NumLib::LocalToGlobalIndexMap>(
126 std::move(all_mesh_subsets),
127 vec_n_components,
128 vec_var_elements,
130
131 // in case of debugging the dof table, activate the following line
132 // std::cout << *_local_to_global_index_map << "\n";
133}
134
136 NumLib::LocalToGlobalIndexMap const& dof_table,
137 MeshLib::Mesh const& mesh,
138 unsigned const integration_order)
139{
140 // Quick access map to BHE's through element ids.
141 std::unordered_map<std::size_t, BHE::BHETypes*> element_to_bhe_map;
142 int const n_BHEs = _process_data._vec_BHE_property.size();
143 for (int i = 0; i < n_BHEs; i++)
144 {
145 auto const& bhe_elements = _bheMeshData.BHE_elements[i];
146 for (auto const& e : bhe_elements)
147 {
148 element_to_bhe_map[e->getID()] =
149 &_process_data._vec_BHE_property[i];
150 }
151 }
152
153 assert(mesh.getDimension() == 3);
156 mesh.getElements(), dof_table, local_assemblers_,
157 NumLib::IntegrationOrder{integration_order}, element_to_bhe_map,
159
160 // Create BHE boundary conditions for each of the BHEs
162
163 // Store BHE and soil elements to split the assembly and use the matrix
164 // cache in the linear case only for soil elements
165 if (_process_data._is_linear)
166 {
167 _bhes_element_ids = _bheMeshData.BHE_elements | ranges::views::join |
168 MeshLib::views::ids | ranges::to<std::vector>;
169
170 // sort bhe elements if needed
171 if (!std::is_sorted(_bhes_element_ids.begin(), _bhes_element_ids.end()))
172 {
173 std::sort(_bhes_element_ids.begin(), _bhes_element_ids.end());
174 }
175
177 ranges::to<std::vector>();
178
179 // sort soil elements if needed
180 if (!std::is_sorted(_soil_element_ids.begin(), _soil_element_ids.end()))
181 {
182 std::sort(_soil_element_ids.begin(), _soil_element_ids.end());
183 }
184
185 _soil_element_ids = ranges::views::set_difference(_soil_element_ids,
187 ranges::to<std::vector>();
188 }
189
190 if (_process_data._mass_lumping)
191 {
192 std::vector<std::size_t> const bhes_node_ids =
193 _bheMeshData.BHE_nodes | ranges::views::join |
194 ranges::views::transform([](auto const* const node)
195 { return node->getID(); }) |
196 ranges::to<std::vector>;
197
198 // all connected soil elements and also the BHE elements.
199 MeshLib::ElementSearch es{mesh};
200 es.searchByNodeIDs(bhes_node_ids);
201
202 assert(_process_data.mass_lumping_soil_elements.empty());
203 _process_data.mass_lumping_soil_elements.resize(
204 mesh.getNumberOfElements(), false);
205 for (auto const id : es.getSearchedElementIDs())
206 {
207 _process_data.mass_lumping_soil_elements[id] = true;
208 }
209 }
210}
211
212std::vector<std::vector<std::string>>
214 std::vector<std::reference_wrapper<MeshLib::Mesh>> const& meshes)
215{
216 DBUG("T-BHE process initializeSubmeshOutput().");
217
218 namespace R = ranges;
219 namespace RV = ranges::views;
220
221 auto get_residuum_name =
222 [](ProcessLib::ProcessVariable const& pv) -> std::string
223 {
224 std::string const& pv_name = pv.getName();
225 if (pv_name.starts_with("temperature_"))
226 {
227 using namespace std::literals;
228 auto const len_temp = "temperature_"s.size();
229 // TODO is that the correct name?
230 return "HeatFlowRate_" + pv_name.substr(len_temp);
231 }
232
233 OGS_FATAL("Unexpected process variable name '{}'.", pv_name);
234 };
235
236 auto get_residuum_names = [get_residuum_name](auto const& pvs)
237 { return pvs | RV::transform(get_residuum_name); };
238
239 auto residuum_names = _process_variables |
240 RV::transform(get_residuum_names) |
241 R::to<std::vector<std::vector<std::string>>>();
242
244 meshes, residuum_names);
245
246 return residuum_names;
247}
248
250 const double t, double const dt, std::vector<GlobalVector*> const& x,
251 std::vector<GlobalVector*> const& x_prev, int const process_id,
253{
254 DBUG("Assemble HeatTransportBHE process.");
255
256 if (_process_data._is_linear)
257 {
258 if (!getActiveElementIDs().empty())
259 {
260 OGS_FATAL(
261 "Domain Deactivation is currently not implemented with "
262 "linear optimization.");
263 }
264
265 // assemble only soil elements/use cache
266 // note: soil element ids are sorted!
268 t, dt, x, x_prev, process_id, M, K, b, &_soil_element_ids);
269
270 // reset the sparsity pattern for better performance in the BHE assembly
271 auto const& spec = getMatrixSpecifications(process_id);
272 MathLib::setMatrixSparsity(M, *spec.sparsity_pattern);
273 MathLib::setMatrixSparsity(K, *spec.sparsity_pattern);
274
275 // assemble BHE elements
276 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
280 local_assemblers_, _bhes_element_ids, dof_table, t, dt, x, x_prev,
281 process_id, &M, &K, &b);
282 }
283 else
284 {
285 // assemble all elements
287 process_id, M, K, b);
288 }
289
290 // Algebraic BC procedure.
291 if (_process_data._algebraic_BC_Setting._use_algebraic_bc)
292 {
293 algebraicBcConcreteProcess(t, dt, x, x_prev, process_id, M, K, b);
294 }
295
296 //_global_output(t, process_id, M, K, b);
297}
298
300 const double t, double const dt, std::vector<GlobalVector*> const& x,
301 std::vector<GlobalVector*> const& x_prev, int const process_id,
302 GlobalVector& b, GlobalMatrix& Jac)
303{
304 DBUG("AssembleWithJacobian HeatTransportBHE process.");
305
307 t, dt, x, x_prev, process_id, b, Jac);
308}
309
311 double const t, double const dt, std::vector<GlobalVector*> const& x,
312 GlobalVector const& x_prev, int const process_id)
313{
314 DBUG("Compute heat flux for HeatTransportBHE process.");
315
316 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
317 dof_tables.reserve(x.size());
318 std::generate_n(std::back_inserter(dof_tables), x.size(),
319 [&]() { return _local_to_global_index_map.get(); });
320
323 local_assemblers_, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
324 process_id);
325}
326
328 GlobalVector const& x)
329{
330 // if the process use python boundary condition
331 if (_process_data.py_bc_object == nullptr || !(_process_data._use_tespy))
332 {
334 }
335
336 // Here the task is to get current time flowrate and flow temperature from
337 // TESPy and determine whether it converges.
338 auto const Tout_nodes_id =
339 std::get<3>(_process_data.py_bc_object->dataframe_network);
340 const std::size_t n_bc_nodes = Tout_nodes_id.size();
341
342 for (std::size_t i = 0; i < n_bc_nodes; i++)
343 {
344 // read the T_out and store them in dataframe
345 std::get<2>(_process_data.py_bc_object->dataframe_network)[i] =
346 x[Tout_nodes_id[i]];
347 }
348 // Transfer Tin and Tout to TESPy and return the results
349 auto const tespy_result = _process_data.py_bc_object->tespySolver(
350 std::get<0>(_process_data.py_bc_object->dataframe_network), // t
351 std::get<1>(_process_data.py_bc_object->dataframe_network), // T_in
352 std::get<2>(_process_data.py_bc_object->dataframe_network)); // T_out
353 if (!_process_data.py_bc_object->isOverriddenTespy())
354 {
355 DBUG("Method `tespySolver' not overridden in Python script.");
356 }
357
358 // update the Tin and flow rate
359 for (std::size_t i = 0; i < n_bc_nodes; i++)
360 {
361 std::get<1>(_process_data.py_bc_object->dataframe_network)[i] =
362 std::get<2>(tespy_result)[i];
363 std::get<4>(_process_data.py_bc_object->dataframe_network)[i] =
364 std::get<3>(tespy_result)[i];
365 }
366 auto const tespy_has_converged = std::get<1>(tespy_result);
367 if (tespy_has_converged == true)
368 {
370 }
371
373}
374
376 std::vector<GlobalVector*> const& x, const double t, const double dt,
377 int const process_id)
378{
380
381 if (_process_data.py_bc_object == nullptr ||
382 !_process_data._use_server_communication)
383 {
384 return;
385 }
386
387 auto& [time, Tin_value, Tout_value, Tout_nodes_ids, flowrate] =
388 _process_data.py_bc_object->dataframe_network;
389
390 // We found the problem that time != t, but it always equals the last
391 // step. The following line is to correct this, although we do not use
392 // it for server communication.
393 time = t;
394
395 auto const& solution = *x[process_id];
396
397 // Iterate through each BHE
398 const std::size_t n_bc_nodes = Tout_nodes_ids.size();
399 for (std::size_t i = 0; i < n_bc_nodes; i++)
400 {
401 // read the T_out and store them in dataframe
402 Tout_value[i] = solution[Tout_nodes_ids[i]];
403 }
404
405 // Transfer T_out to server_Communication and get back T_in and flowrate
406 auto const server_communication_result =
407 _process_data.py_bc_object->serverCommunicationPreTimestep(
408 t, dt, Tin_value, Tout_value, flowrate);
409 if (!_process_data.py_bc_object
410 ->isOverriddenServerCommunicationPreTimestep())
411 {
412 DBUG("Method `serverCommunication' not overridden in Python script.");
413 }
414
415 auto const& [server_communication_Tin_value,
416 server_communication_flowrate] = server_communication_result;
417
418 std::copy(begin(server_communication_Tin_value),
419 end(server_communication_Tin_value),
420 begin(Tin_value));
421 std::copy(begin(server_communication_flowrate),
422 end(server_communication_flowrate),
423 begin(flowrate));
424}
425
427 std::vector<GlobalVector*> const& x,
428 std::vector<GlobalVector*> const& /*x_prev*/, const double t,
429 const double dt, int const process_id)
430{
431 if (_process_data.py_bc_object == nullptr ||
432 !_process_data._use_server_communication)
433 {
434 return;
435 }
436
437 auto& [time, Tin_value, Tout_value, Tout_nodes_ids, flowrate] =
438 _process_data.py_bc_object->dataframe_network;
439
440 // We found the problem that time != t, but it always equals the last
441 // step. The following line is to correct this, although we do not use
442 // it for server communication.
443 time = t;
444
445 auto const& solution = *x[process_id];
446
447 // Iterate through each BHE
448 const std::size_t n_bc_nodes = Tout_nodes_ids.size();
449 for (std::size_t i = 0; i < n_bc_nodes; i++)
450 {
451 // read the T_out and store them in dataframe
452 Tout_value[i] = solution[Tout_nodes_ids[i]];
453 }
454
455 // Transfer T_out to server_Communication
456 _process_data.py_bc_object->serverCommunicationPostTimestep(
457 t, dt, Tin_value, Tout_value, flowrate);
458 if (!_process_data.py_bc_object
459 ->isOverriddenServerCommunicationPostTimestep())
460 {
461 DBUG("Method `serverCommunication' not overridden in Python script.");
462 }
463}
464
466 [[maybe_unused]] const double t, double const /*dt*/,
467 [[maybe_unused]] std::vector<GlobalVector*> const& x,
468 std::vector<GlobalVector*> const& /*xprev*/, int const /*process_id*/,
469 [[maybe_unused]] GlobalMatrix& M, [[maybe_unused]] GlobalMatrix& K,
470 [[maybe_unused]] GlobalVector& b)
471{
472#ifndef USE_PETSC
473 auto M_normal = M.getRawMatrix();
474 auto K_normal = K.getRawMatrix();
475 auto n_original_rows = K_normal.rows();
476 auto const n_BHE_bottom_pairs = _vec_bottom_BHE_node_indices.size();
477 auto const n_BHE_top_pairs = _vec_top_BHE_node_indices.size();
478
479 // apply weighting factor based on the max value from column wise inner
480 // product and scale it with user defined value
481 const double w_val = _process_data._algebraic_BC_Setting._weighting_factor;
482
483 M_normal.conservativeResize(
484 M_normal.rows() + n_BHE_bottom_pairs + n_BHE_top_pairs,
485 M_normal.cols());
486 K_normal.conservativeResize(
487 K_normal.rows() + n_BHE_bottom_pairs + n_BHE_top_pairs,
488 K_normal.cols());
489
490 for (std::size_t i = 0; i < n_BHE_bottom_pairs; i++)
491 {
492 Eigen::SparseVector<double> M_Plus(M_normal.cols());
493 M_Plus.setZero();
494 M_normal.row(n_original_rows + i) = M_Plus;
495
496 Eigen::SparseVector<double> K_Plus(K_normal.cols());
497 K_Plus.setZero();
498
499 auto const [bhe_idx, first_BHE_bottom_index, second_BHE_bottom_index] =
501
502 K_Plus.insert(first_BHE_bottom_index) = w_val;
503 K_Plus.insert(second_BHE_bottom_index) = -w_val;
504
505 K_normal.row(n_original_rows + i) = K_Plus;
506 }
507
508 auto b_normal = b.getRawVector();
509 Eigen::SparseVector<double> b_Plus(b_normal.rows() + n_BHE_bottom_pairs +
510 n_BHE_top_pairs);
511 b_Plus.setZero();
512
513 // Copy values from the original column vector to the modified one
514 for (int i = 0; i < b_normal.innerSize(); ++i)
515 {
516 b_Plus.insert(i) = b_normal.coeff(i);
517 }
518
519 for (std::size_t i = 0; i < n_BHE_top_pairs; i++)
520 {
521 Eigen::SparseVector<double> M_Plus(M_normal.cols());
522 M_Plus.setZero();
523 M_normal.row(n_original_rows + n_BHE_bottom_pairs + i) = M_Plus;
524
525 Eigen::SparseVector<double> K_Plus(K_normal.cols());
526 K_Plus.setZero();
527
528 auto const [bhe_idx, first_BHE_top_index, second_BHE_top_index] =
530
531 auto first_BHE_top_index_pair = first_BHE_top_index;
532 auto second_BHE_top_index_pair = second_BHE_top_index;
533
534 K_Plus.insert(first_BHE_top_index_pair) =
535 w_val; // for power BC, the inflow node must be positive
536 K_Plus.insert(second_BHE_top_index_pair) =
537 -w_val; // for power BC, the outflow node must be negative
538
539 K_normal.row(n_original_rows + n_BHE_bottom_pairs + i) = K_Plus;
540
541 // get the delta_T value here
542 double const T_out = (*x[0])[second_BHE_top_index_pair];
543
544 auto calculate_delta_T = [&](auto& bhe)
545 {
546 auto const T_in = bhe.updateFlowRateAndTemperature(T_out, t);
547 return T_in - T_out;
548 };
549 auto delta_T = std::visit(calculate_delta_T,
550 _process_data._vec_BHE_property[bhe_idx]);
551
552 b_Plus.insert(n_original_rows + n_BHE_bottom_pairs + i) =
553 delta_T * w_val;
554 }
555
556 M.getRawMatrix() = M_normal;
557 K.getRawMatrix() = K_normal;
558 b.getRawVector() = b_Plus;
559#else
560 OGS_FATAL(
561 "The Algebraic Boundary Condition is not implemented for use with "
562 "PETsc Library! Simulation will be terminated.");
563#endif
564}
565
567 std::vector<std::vector<MeshLib::Node*>> const& all_bhe_nodes)
568{
569 const int process_id = 0;
570 auto& bcs = _boundary_conditions[process_id];
571
572 std::size_t const n_BHEs = _process_data._vec_BHE_property.size();
573
574 // for each BHE
575 for (std::size_t bhe_i = 0; bhe_i < n_BHEs; bhe_i++)
576 {
577 auto const& bhe_nodes = all_bhe_nodes[bhe_i];
578 // find the variable ID
579 // the soil temperature is 0-th variable
580 // the BHE temperature is therefore bhe_i + 1
581 const int variable_id = bhe_i + 1;
582
583 std::vector<MeshLib::Node*> bhe_boundary_nodes;
584
585 // cherry-pick the boundary nodes according to
586 // the number of connected line elements.
587 for (auto const& bhe_node : bhe_nodes)
588 {
589 // Count number of 1d elements connected with every BHE node.
590 auto const& connected_elements =
591 _mesh.getElementsConnectedToNode(*bhe_node);
592 const std::size_t n_line_elements = std::count_if(
593 connected_elements.begin(), connected_elements.end(),
594 [](MeshLib::Element const* elem)
595 { return (elem->getDimension() == 1); });
596
597 if (n_line_elements == 1)
598 {
599 bhe_boundary_nodes.push_back(bhe_node);
600 }
601 }
602
603 if (bhe_boundary_nodes.size() != 2)
604 {
605 OGS_FATAL(
606 "Error!!! The BHE boundary nodes are not correctly found, "
607 "for every single BHE, there should be 2 boundary nodes.");
608 }
609
610 // For 1U, 2U, CXC, CXA type BHE, the node order in the boundary nodes
611 // vector should be rearranged according to its z coordinate in
612 // descending order. In these BHE types, the z coordinate on the top and
613 // bottom node is different. The BHE top node with a higher z coordinate
614 // should be placed at the first, while the BHE bottom node with a lower
615 // z coordinate should be placed at the second. For other horizontal BHE
616 // types e.g. 1P-type BHE, the z coordinate on the top and bottom node
617 // is identical. Thus the node order in the boundary nodes vector can
618 // not be rearranged according to its z coordinate. For these BHE types,
619 // the boundary node order is according to the default node id order in
620 // the model mesh.
621 // for 1P-type BHE
622 if ((*bhe_boundary_nodes[0])[2] == (*bhe_boundary_nodes[1])[2])
623 {
624 INFO(
625 "For 1P-type BHE, the BHE inflow and outflow "
626 "nodes are identified according to their mesh node id in "
627 "ascending order");
628 }
629 // for 1U, 2U, CXC, CXA type BHE
630 else
631 {
632 // swap the boundary nodes if the z coordinate of the
633 // first node is lower than it on the second node
634 if ((*bhe_boundary_nodes[0])[2] < (*bhe_boundary_nodes[1])[2])
635 {
636 std::swap(bhe_boundary_nodes[0], bhe_boundary_nodes[1]);
637 }
638 }
639
640 auto get_global_index =
641 [&](std::size_t const node_id, int const component)
642 {
643 return _local_to_global_index_map->getGlobalIndex(
644 {_mesh.getID(), MeshLib::MeshItemType::Node, node_id},
645 variable_id, component);
646 };
647
648 auto get_global_bhe_bc_indices =
649 [&](std::array<
650 std::pair<std::size_t /*node_id*/, int /*component*/>, 2>
651 nodes_and_components)
652 {
653 return std::make_pair(
654 get_global_index(nodes_and_components[0].first,
655 nodes_and_components[0].second),
656 get_global_index(nodes_and_components[1].first,
657 nodes_and_components[1].second));
658 };
659
660 auto get_global_bhe_bc_indices_with_bhe_idx =
661 [&](std::size_t bhe_idx,
662 std::array<
663 std::pair<std::size_t /*node_id*/, int /*component*/>, 2>
664 nodes_and_components)
665 {
666 return std::make_tuple(
667 bhe_idx,
668 get_global_index(nodes_and_components[0].first,
669 nodes_and_components[0].second),
670 get_global_index(nodes_and_components[1].first,
671 nodes_and_components[1].second));
672 };
673
674 auto createBCs =
675 [&, bc_top_node_id = bhe_boundary_nodes[0]->getID(),
676 bc_bottom_node_id = bhe_boundary_nodes[1]->getID()](auto& bhe)
677 {
678 for (auto const& in_out_component_id :
679 bhe.inflow_outflow_bc_component_ids)
680 {
681 if (bhe.use_python_bcs ||
682 this->_process_data._use_server_communication)
683 // call BHEPythonBoundarycondition
684 {
685 if (this->_process_data
686 .py_bc_object) // the bc object exist
687 {
688 // apply the customized top, inflow BC.
689 bcs.addBoundaryCondition(
691 get_global_bhe_bc_indices(
692 bhe.getBHEInflowDirichletBCNodesAndComponents(
693 bc_top_node_id, bc_bottom_node_id,
694 in_out_component_id.first)),
695 bhe,
696 *(_process_data.py_bc_object)));
697 }
698 else
699 {
700 OGS_FATAL(
701 "The Python Boundary Condition was switched on, "
702 "but the data object does not exist! ");
703 }
704 }
705 else
706 {
707 if (this->_process_data._algebraic_BC_Setting
708 ._use_algebraic_bc &&
709 bhe.isPowerBC())
710 {
711 // for algebraic_bc method, record the pair of indices
712 // in a separate vector
714 get_global_bhe_bc_indices_with_bhe_idx(
715 bhe_i,
716 {{{bc_top_node_id, in_out_component_id.first},
717 {bc_top_node_id,
718 in_out_component_id.second}}}));
719 }
720 else
721 {
722 // Top, inflow, normal case
723 bcs.addBoundaryCondition(
725 get_global_bhe_bc_indices(
726 bhe.getBHEInflowDirichletBCNodesAndComponents(
727 bc_top_node_id, bc_bottom_node_id,
728 in_out_component_id.first)),
729 [&bhe](double const T, double const t)
730 {
731 return bhe.updateFlowRateAndTemperature(T,
732 t);
733 }));
734 }
735 }
736
737 auto const bottom_nodes_and_components =
738 bhe.getBHEBottomDirichletBCNodesAndComponents(
739 bc_bottom_node_id,
740 in_out_component_id.first,
741 in_out_component_id.second);
742
743 if (bottom_nodes_and_components &&
744 !this->_process_data._algebraic_BC_Setting
745 ._use_algebraic_bc)
746 {
747 // Bottom, outflow, all cases | not needed for algebraic_bc
748 // method
749 bcs.addBoundaryCondition(
751 get_global_bhe_bc_indices(
752 {{{bc_bottom_node_id,
753 in_out_component_id.first},
754 {bc_bottom_node_id,
755 in_out_component_id.second}}})));
756 }
757 else if (bottom_nodes_and_components &&
758 this->_process_data._algebraic_BC_Setting
759 ._use_algebraic_bc)
760 {
761 // for algebraic_bc method, record the pair of indices in a
762 // separate vector
764 get_global_bhe_bc_indices_with_bhe_idx(
765 bhe_i,
766 {{{bc_bottom_node_id, in_out_component_id.first},
767 {bc_bottom_node_id,
768 in_out_component_id.second}}}));
769 }
770 }
771 };
772 visit(createBCs, _process_data._vec_BHE_property[bhe_i]);
773 }
774}
775} // namespace HeatTransportBHE
776} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:19
MathLib::EigenMatrix GlobalMatrix
MathLib::EigenVector GlobalVector
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
Element search class.
std::size_t searchByNodeIDs(const std::vector< std::size_t > &nodes)
Marks all elements connecting to any of the given nodes.
bool isAxiallySymmetric() const
Definition Mesh.h:130
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:101
unsigned getDimension() const
Definition Mesh.h:80
std::size_t getNumberOfElements() const
Get the number of elements.
Definition Mesh.h:89
void assemble(double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, std::vector< std::size_t > const *const sorted_element_subset=nullptr, bool const copy_residua_to_mesh=false)
void initializeAssemblyOnSubmeshes(std::vector< std::reference_wrapper< MeshLib::Mesh > > const &submeshes, std::vector< std::vector< std::string > > const &residuum_names)
void assembleWithJacobian(double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac, std::vector< std::size_t > const *const sorted_element_subset=nullptr, bool const copy_residua_to_mesh=false)
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_soil_nodes
void assembleConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) override
std::vector< std::unique_ptr< MeshLib::MeshSubset const > > _mesh_subset_BHE_nodes
std::vector< std::tuple< std::size_t, GlobalIndexType, GlobalIndexType > > _vec_bottom_BHE_node_indices
std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes(std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) override
void computeSecondaryVariableConcrete(double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id) override
HeatTransportBHEProcess(std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, HeatTransportBHEProcessData &&process_data, SecondaryVariableCollection &&secondary_variables, BHEMeshData &&bhe_mesh_data)
NumLib::IterationResult postIterationConcreteProcess(GlobalVector const &x) override
std::vector< std::unique_ptr< HeatTransportBHELocalAssemblerInterface > > local_assemblers_
void preTimestepConcreteProcess(std::vector< GlobalVector * > const &x, const double t, const double dt, int const process_id) override
std::vector< std::tuple< std::size_t, GlobalIndexType, GlobalIndexType > > _vec_top_BHE_node_indices
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
void algebraicBcConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)
void assembleWithJacobianConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac) override
void postTimestepConcreteProcess(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double dt, int const process_id) override
void createBHEBoundaryConditionTopBottom(std::vector< std::vector< MeshLib::Node * > > const &all_bhe_nodes)
virtual void computeSecondaryVariable(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
std::string const name
Definition Process.h:361
std::vector< BoundaryConditionCollection > _boundary_conditions
Definition Process.h:404
Process(std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true)
Definition Process.cpp:37
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:365
MeshLib::Mesh & _mesh
Definition Process.h:364
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
Definition Process.h:399
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:160
VectorMatrixAssembler _global_assembler
Definition Process.h:376
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:367
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > _jacobian_assembler
Definition Process.h:375
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
Definition Process.cpp:203
Handles configuration of several secondary variables from the project file.
void assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix *M, GlobalMatrix *K, GlobalVector *b)
IterationResult
Status flags telling the NonlinearSolver if an iteration succeeded.
void setMatrixSparsity(MATRIX &matrix, SPARSITY_PATTERN const &sparsity_pattern)
constexpr ranges::views::view_closure ids
For an element of a range view return its id.
Definition Mesh.h:218
PropertyVector< int > const * materialIDs(Mesh const &mesh)
Definition Mesh.cpp:258
@ BY_COMPONENT
Ordering data by component type.
std::unique_ptr< BHEInflowDirichletBoundaryCondition< BHEUpdateCallback > > createBHEInflowDirichletBoundaryCondition(std::pair< GlobalIndexType, GlobalIndexType > &&in_out_global_indices, BHEUpdateCallback bhe_update_callback)
std::unique_ptr< BHEBottomDirichletBoundaryCondition > createBHEBottomDirichletBoundaryCondition(std::pair< GlobalIndexType, GlobalIndexType > &&in_out_global_indices)
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)
std::unique_ptr< BHEInflowPythonBoundaryCondition< BHEType > > createBHEInflowPythonBoundaryCondition(std::pair< GlobalIndexType, GlobalIndexType > &&in_out_global_indices, BHEType &bhe, BHEInflowPythonBoundaryConditionPythonSideInterface &py_bc_object)
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)