OGS
ProjectData.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
4#include "ProjectData.h"
5
6#include <pybind11/eval.h>
7
8#include <algorithm>
9#include <boost/algorithm/string/predicate.hpp>
10#include <cctype>
11#include <range/v3/action/sort.hpp>
12#include <range/v3/action/unique.hpp>
13#include <range/v3/algorithm/contains.hpp>
14#include <range/v3/range/conversion.hpp>
15#include <range/v3/view/adjacent_remove_if.hpp>
16#include <set>
17
18#include "BaseLib/Algorithm.h"
19#include "BaseLib/ConfigTree.h"
20#include "BaseLib/FileTools.h"
21#include "BaseLib/Logging.h"
22#include "BaseLib/StringTools.h"
23#include "GeoLib/GEOObjects.h"
26#include "GeoLib/Raster.h"
27#include "InfoLib/CMakeInfo.h"
31#if defined(USE_LIS)
33#elif defined(USE_PETSC)
35#else
37#endif
41#include "MeshLib/Mesh.h"
47
48// FileIO
53#include "ParameterLib/Utils.h"
55#include "ProcessLib/TimeLoop.h"
56
57#ifdef OGS_BUILD_PROCESS_COMPONENTTRANSPORT
60#endif
61#ifdef OGS_BUILD_PROCESS_STEADYSTATEDIFFUSION
63#endif
64#ifdef OGS_BUILD_PROCESS_HT
66#endif
67#ifdef OGS_BUILD_PROCESS_HEATCONDUCTION
69#endif
70#ifdef OGS_BUILD_PROCESS_HEATTRANSPORTBHE
72#endif
73#ifdef OGS_BUILD_PROCESS_WELLBORESIMULATOR
75#endif
76#ifdef OGS_BUILD_PROCESS_HYDROMECHANICS
78#endif
79#ifdef OGS_BUILD_PROCESS_LARGEDEFORMATION
81#endif
82#ifdef OGS_BUILD_PROCESS_LIE_M
84#endif
85#ifdef OGS_BUILD_PROCESS_LIE_HM
87#endif
88#ifdef OGS_BUILD_PROCESS_LIQUIDFLOW
90#endif
91#ifdef OGS_BUILD_PROCESS_THERMORICHARDSMECHANICS
93#endif
94
95#ifdef OGS_BUILD_PROCESS_PHASEFIELD
97#endif
98#ifdef OGS_BUILD_PROCESS_HMPHASEFIELD
100#endif
101#ifdef OGS_BUILD_PROCESS_RICHARDSCOMPONENTTRANSPORT
103#endif
104#ifdef OGS_BUILD_PROCESS_RICHARDSFLOW
106#endif
107#ifdef OGS_BUILD_PROCESS_RICHARDSMECHANICS
109#endif
110#ifdef OGS_BUILD_PROCESS_SMALLDEFORMATION
112#endif
113#ifdef OGS_BUILD_PROCESS_TH2M
115#endif
116#ifdef OGS_BUILD_PROCESS_THERMALTWOPHASEFLOWWITHPP
118#endif
119#ifdef OGS_BUILD_PROCESS_THERMOHYDROMECHANICS
121#endif
122#ifdef OGS_BUILD_PROCESS_THERMOMECHANICS
124#endif
125#ifdef OGS_BUILD_PROCESS_THERMORICHARDSFLOW
127#endif
128#ifdef OGS_BUILD_PROCESS_TWOPHASEFLOWWITHPP
130#endif
131
132namespace
133{
134void readGeometry(std::string const& fname, GeoLib::GEOObjects& geo_objects,
135 std::string const& dir_first, std::string const& dir_second)
136{
137 DBUG("Reading geometry file '{:s}'.", fname);
138 GeoLib::IO::BoostXmlGmlInterface gml_reader(geo_objects);
139 std::string geometry_file = BaseLib::joinPaths(dir_first, fname);
140 if (!BaseLib::IsFileExisting(geometry_file))
141 {
142 // Fallback to reading gml from prj-file directory
143 geometry_file = BaseLib::joinPaths(dir_second, fname);
144 WARN("File {:s} not found in {:s}! Trying reading from {:s}.", fname,
145 dir_first, dir_second);
146 if (!BaseLib::IsFileExisting(geometry_file))
147 {
148 OGS_FATAL("Could not read geometry file {:s} in {:s}.", fname,
149 dir_second);
150 }
151 }
152 gml_reader.readFile(geometry_file);
153}
154
155std::unique_ptr<MeshLib::Mesh> readSingleMesh(
156 BaseLib::ConfigTree const& mesh_config_parameter,
157 std::string const& directory)
158{
159 std::string const mesh_file = BaseLib::joinPaths(
160 directory, mesh_config_parameter.getValue<std::string>());
161 DBUG("Reading mesh file '{:s}'.", mesh_file);
162
163 auto mesh = std::unique_ptr<MeshLib::Mesh>(MeshLib::IO::readMeshFromFile(
164 mesh_file, true /* compute_element_neighbors */));
165 if (!mesh)
166 {
167 std::filesystem::path abspath{mesh_file};
168 try
169 {
170 abspath = std::filesystem::absolute(mesh_file);
171 }
172 catch (std::filesystem::filesystem_error const& e)
173 {
174 ERR("Determining the absolute path of '{}' failed: {}", mesh_file,
175 e.what());
176 }
177
178 OGS_FATAL("Could not read mesh from '{:s}' file. No mesh added.",
179 abspath.string());
180 }
181
182#ifdef DOXYGEN_DOCU_ONLY
184 mesh_config_parameter.getConfigAttributeOptional<bool>("axially_symmetric");
185#endif // DOXYGEN_DOCU_ONLY
186
187 if (auto const axially_symmetric =
189 mesh_config_parameter.getConfigAttributeOptional<bool>(
190 "axially_symmetric"))
191 {
192 mesh->setAxiallySymmetric(*axially_symmetric);
193 if (mesh->getDimension() == 3 && mesh->isAxiallySymmetric())
194 {
195 OGS_FATAL("3D mesh cannot be axially symmetric.");
196 }
197 }
198
199 return mesh;
200}
201
202std::vector<std::unique_ptr<MeshLib::Mesh>> readMeshes(
203 BaseLib::ConfigTree const& config, std::string const& directory)
204{
205 std::vector<std::unique_ptr<MeshLib::Mesh>> meshes;
206
207 GeoLib::GEOObjects geoObjects;
208
210 auto optional_meshes = config.getConfigSubtreeOptional("meshes");
211 if (optional_meshes)
212 {
213 DBUG("Reading multiple meshes.");
215 auto const configs = optional_meshes->getConfigParameterList("mesh");
216 std::transform(configs.begin(), configs.end(),
217 std::back_inserter(meshes),
218 [&directory](auto const& mesh_config)
219 { return readSingleMesh(mesh_config, directory); });
220 if (auto const geometry_file_name =
222 config.getConfigParameterOptional<std::string>("geometry"))
223 {
224 readGeometry(*geometry_file_name, geoObjects, directory,
225 config.projectDirectory().string());
226 }
227 }
228 else
229 { // Read single mesh with geometry.
230 meshes.push_back(
232 readSingleMesh(config.getConfigParameter("mesh"), directory));
233
234 auto const geometry_file_name =
236 config.getConfigParameter<std::string>("geometry");
237 readGeometry(geometry_file_name, geoObjects, directory,
238 config.projectDirectory().string());
239 }
240
241 if (!geoObjects.getPoints().empty() || !geoObjects.getPolylines().empty() ||
242 !geoObjects.getSurfaces().empty())
243 { // generate meshes from geometries
244 std::unique_ptr<MeshGeoToolsLib::SearchLength> search_length_algorithm =
246 bool const multiple_nodes_allowed = false;
247 auto additional_meshes =
249 geoObjects, *meshes[0], std::move(search_length_algorithm),
250 multiple_nodes_allowed);
251
252 std::move(begin(additional_meshes), end(additional_meshes),
253 std::back_inserter(meshes));
254 }
255
256 auto mesh_names = MeshLib::views::names | ranges::to<std::vector>() |
257 ranges::actions::sort;
258 auto const sorted_names = meshes | mesh_names;
259 auto const unique_names = meshes | mesh_names | ranges::actions::unique;
260 if (unique_names.size() < sorted_names.size())
261 {
262 WARN(
263 "Mesh names aren't unique. From project file read mesh names are:");
264 for (auto const& name : meshes | MeshLib::views::names)
265 {
266 INFO("- {}", name);
267 }
268 }
269
270 auto const zero_mesh_field_data_by_material_ids =
272 config.getConfigParameterOptional<std::vector<int>>(
273 "zero_mesh_field_data_by_material_ids");
274 if (zero_mesh_field_data_by_material_ids)
275 {
276 WARN(
277 "Tag 'zero_mesh_field_data_by_material_ids` is experimental. Its "
278 "name may be changed, or it may be removed due to its "
279 "corresponding feature becomes a single tool. Please use it with "
280 "care!");
282 *meshes[0], *zero_mesh_field_data_by_material_ids);
283 }
284
286
287 return meshes;
288}
289
290std::vector<GeoLib::NamedRaster> readRasters(
291 BaseLib::ConfigTree const& config,
292 GeoLib::MinMaxPoints const& min_max_points)
293{
294 INFO("readRasters ...");
295 std::vector<GeoLib::NamedRaster> named_rasters;
296
298 auto optional_rasters = config.getConfigSubtreeOptional("rasters");
299 if (optional_rasters)
300 {
302 auto const configs = optional_rasters->getConfigSubtreeList("raster");
303 std::transform(
304 configs.begin(), configs.end(), std::back_inserter(named_rasters),
305 [&min_max_points](auto const& raster_config)
306 { return GeoLib::IO::readRaster(raster_config, min_max_points); });
307 }
308 INFO("readRasters done");
309 return named_rasters;
310}
311
312// for debugging raster reading implementation
313// void writeRasters(std::vector<GeoLib::NamedRaster> const& named_rasters,
314// std::string const& output_directory)
315//{
316// for (auto const& named_raster : named_rasters)
317// {
318// #if defined(USE_PETSC)
319// int my_mpi_rank;
320// MPI_Comm_rank(BaseLib::MPI::OGS_COMM_WORLD, &my_mpi_rank);
321// #endif
322// FileIO::AsciiRasterInterface::writeRasterAsASC(
323// named_raster.raster, output_directory + "/" +
324// named_raster.raster_name +
325// #if defined(USE_PETSC)
326// "_" + std::to_string(my_mpi_rank) +
327// #endif
328// ".asc");
329// }
330//}
331
332} // namespace
333
334ProjectData::ProjectData() = default;
335
337 std::string const& output_directory,
338 std::string const& mesh_directory,
339 [[maybe_unused]] std::string const& script_directory)
340 : _mesh_vec(readMeshes(project_config, mesh_directory)),
341 _named_rasters(readRasters(project_config,
342 GeoLib::AABB(_mesh_vec[0]->getNodes().begin(),
343 _mesh_vec[0]->getNodes().end())
344 .getMinMaxPoints()))
345{
346 // for debugging raster reading implementation
347 // writeRasters(_named_rasters, output_directory);
348 if (auto const python_script =
350 project_config.getConfigParameterOptional<std::string>("python_script"))
351 {
352 namespace py = pybind11;
353
354 // Append to python's module search path
355 auto py_path = py::module::import("sys").attr("path");
356 py_path.attr("append")(script_directory); // .prj or -s directory
357
358 auto const script_path =
359 BaseLib::joinPaths(script_directory, *python_script);
360
361 // Evaluate in scope of main module
362 py::object scope = py::module::import("__main__").attr("__dict__");
363 // add (global) variables
364 auto globals = py::dict(scope);
365 globals["ogs_prj_directory"] =
366 project_config.projectDirectory().string();
367 globals["ogs_mesh_directory"] = mesh_directory;
368 globals["ogs_script_directory"] = script_directory;
369 try
370 {
371 py::eval_file(script_path, scope);
372 }
373 catch (py::error_already_set const& e)
374 {
375 OGS_FATAL("Error evaluating python script {}: {}", script_path,
376 e.what());
377 }
378 }
379
381 parseCurves(project_config.getConfigSubtreeOptional("curves"));
382
383 auto parameter_names_for_transformation =
385 parseParameters(project_config.getConfigSubtree("parameters"));
386
387 _local_coordinate_system = ParameterLib::createCoordinateSystem(
389 project_config.getConfigSubtreeOptional("local_coordinate_system"),
390 _parameters);
391
392 for (auto& parameter : _parameters)
393 {
394 if (std::find(begin(parameter_names_for_transformation),
395 end(parameter_names_for_transformation),
396 parameter->name) !=
397 end(parameter_names_for_transformation))
398 {
399 if (!_local_coordinate_system)
400 {
401 OGS_FATAL(
402 "The parameter '{:s}' is using the local coordinate system "
403 "but no local coordinate system was provided.",
404 parameter->name);
405 }
406 parameter->setCoordinateSystem(*_local_coordinate_system);
407 }
408
409 parameter->initialize(_parameters);
410 }
411
413 parseProcessVariables(project_config.getConfigSubtree("process_variables"));
414
416 parseMedia(project_config.getConfigSubtreeOptional("media"));
417
419 parseLinearSolvers(project_config.getConfigSubtree("linear_solvers"));
420
421 auto chemical_solver_interface = parseChemicalSolverInterface(
423 project_config.getConfigSubtreeOptional("chemical_system"),
424 output_directory);
425
427 parseProcesses(project_config.getConfigSubtree("processes"),
428 output_directory,
429 std::move(chemical_solver_interface));
430
432 parseNonlinearSolvers(project_config.getConfigSubtree("nonlinear_solvers"));
433
435 parseTimeLoop(project_config.getConfigSubtree("time_loop"),
436 output_directory);
437}
438
440 BaseLib::ConfigTree const& process_variables_config)
441{
442 DBUG("Parse process variables:");
443
444 std::set<std::string> names;
445
446 for (auto var_config
448 : process_variables_config.getConfigSubtreeList("process_variable"))
449 {
450 // Either the mesh name is given, or the first mesh's name will be
451 // taken. Taking the first mesh's value is deprecated.
452 auto const mesh_name =
454 var_config.getConfigParameter<std::string>("mesh",
455 _mesh_vec[0]->getName());
456
457 auto& mesh = MeshLib::findMeshByName(_mesh_vec, mesh_name);
458
459 auto pv = ProcessLib::ProcessVariable{var_config, mesh, _mesh_vec,
461 if (!names.insert(pv.getName()).second)
462 {
463 OGS_FATAL("A process variable with name `{:s}' already exists.",
464 pv.getName());
465 }
466
467 _process_variables.push_back(std::move(pv));
468 }
469}
470
471std::vector<std::string> ProjectData::parseParameters(
472 BaseLib::ConfigTree const& parameters_config)
473{
474 using namespace ProcessLib;
475
476 std::set<std::string> names;
477 std::vector<std::string> parameter_names_for_transformation;
478
479 DBUG("Reading parameters:");
480 for (auto parameter_config :
482 parameters_config.getConfigSubtreeList("parameter"))
483 {
484 auto p = ParameterLib::createParameter(parameter_config, _mesh_vec,
486 if (!names.insert(p->name).second)
487 {
488 OGS_FATAL("A parameter with name `{:s}' already exists.", p->name);
489 }
490
491 auto const use_local_coordinate_system =
493 parameter_config.getConfigParameterOptional<bool>(
494 "use_local_coordinate_system");
495 if (!!use_local_coordinate_system && *use_local_coordinate_system)
496 {
497 parameter_names_for_transformation.push_back(p->name);
498 }
499
500 _parameters.push_back(std::move(p));
501 }
502
503 _parameters.push_back(
506 _parameters.push_back(
509
510 return parameter_names_for_transformation;
511}
512
514 std::optional<BaseLib::ConfigTree> const& media_config)
515{
516 if (!media_config)
517 {
518 return;
519 }
520
521 DBUG("Reading media:");
522
523 if (_mesh_vec.empty() || _mesh_vec[0] == nullptr)
524 {
525 ERR("A mesh is required to define medium materials.");
526 return;
527 }
528
529 for (auto const& medium_config :
531 media_config->getConfigSubtreeList("medium"))
532 {
533 auto create_medium = [dim = _mesh_vec[0]->getDimension(),
534 &medium_config, this](int const id)
535 {
537 id, _mesh_vec[0]->getDimension(), medium_config, _parameters,
539 _curves);
540 };
541
542 auto const material_id_string =
544 medium_config.getConfigAttribute<std::string>("id", "0");
545
546 std::vector<int> const material_ids_of_this_medium =
547 MaterialLib::parseMaterialIdString(material_id_string,
548 materialIDs(*_mesh_vec[0]));
549
550 for (auto const& id : material_ids_of_this_medium)
551 {
553 id, _media, material_ids_of_this_medium, create_medium);
554 }
555 }
556
557 if (_media.empty())
558 {
559 OGS_FATAL("No entity is found inside <media>.");
560 }
561}
562
563std::unique_ptr<ChemistryLib::ChemicalSolverInterface>
565 std::optional<BaseLib::ConfigTree> const& config,
566 std::string const& output_directory)
567{
568 if (!config)
569 {
570 return nullptr;
571 }
572
573 std::unique_ptr<ChemistryLib::ChemicalSolverInterface>
574 chemical_solver_interface;
575#ifdef OGS_BUILD_PROCESS_COMPONENTTRANSPORT
576 INFO(
577 "Ready for initializing interface to a chemical solver for water "
578 "chemistry calculation.");
579
580 auto const chemical_solver =
582 config->getConfigAttribute<std::string>("chemical_solver");
583
584 if (boost::iequals(chemical_solver, "Phreeqc"))
585 {
586 INFO(
587 "Configuring phreeqc interface for water chemistry calculation "
588 "using file-based approach.");
589
590 chemical_solver_interface = ChemistryLib::createChemicalSolverInterface<
592 *config, output_directory);
593 }
594 else if (boost::iequals(chemical_solver, "PhreeqcKernel"))
595 {
596 OGS_FATAL(
597 "The chemical solver option of PhreeqcKernel is not accessible for "
598 "the time being. Please set 'Phreeqc'' as the chemical solver for "
599 "reactive transport modeling.");
600 }
601 else if (boost::iequals(chemical_solver, "SelfContained"))
602 {
603 INFO(
604 "Use self-contained chemical solver for water chemistry "
605 "calculation.");
606
607 chemical_solver_interface = ChemistryLib::createChemicalSolverInterface<
609 _mesh_vec, _linear_solvers, *config, output_directory);
610 }
611 else
612 {
613 OGS_FATAL(
614 "Unknown chemical solver. Please specify either Phreeqc or "
615 "PhreeqcKernel as the solver for water chemistry calculation "
616 "instead.");
617 }
618#else
619 (void)output_directory;
620
621 OGS_FATAL(
622 "Found the type of the process to be solved is not component transport "
623 "process. Please specify the process type to ComponentTransport. At "
624 "the present, water chemistry calculation is only available for "
625 "component transport process.");
626#endif
627 return chemical_solver_interface;
628}
629
631 BaseLib::ConfigTree const& processes_config,
632 std::string const& output_directory,
633 [[maybe_unused]] std::unique_ptr<ChemistryLib::ChemicalSolverInterface>&&
634 chemical_solver_interface)
635{
636 (void)output_directory; // to avoid compilation warning
637
638 DBUG("Reading processes:");
640 for (auto process_config : processes_config.getConfigSubtreeList("process"))
641 {
642 auto const type =
644 process_config.peekConfigParameter<std::string>("type");
645
646 auto const name =
648 process_config.getConfigParameter<std::string>("name");
649
650 [[maybe_unused]] auto const integration_order =
652 process_config.getConfigParameter<int>("integration_order");
653
654 std::unique_ptr<ProcessLib::Process> process;
655
656 auto jacobian_assembler = ProcessLib::createJacobianAssembler(
658 process_config.getConfigSubtreeOptional("jacobian_assembler"));
659
660#ifdef OGS_BUILD_PROCESS_STEADYSTATEDIFFUSION
661 if (type == "STEADY_STATE_DIFFUSION")
662 {
663 // The existence check of the in the configuration referenced
664 // process variables is checked in the physical process.
665 // TODO at the moment we have only one mesh, later there can be
666 // several meshes. Then we have to assign the referenced mesh
667 // here.
668 process =
670 name, *_mesh_vec[0], std::move(jacobian_assembler),
671 _process_variables, _parameters, integration_order,
672 process_config, _mesh_vec, _media);
673 }
674 else
675#endif
676#ifdef OGS_BUILD_PROCESS_LIQUIDFLOW
677 if (type == "LIQUID_FLOW")
678 {
680 name, *_mesh_vec[0], std::move(jacobian_assembler),
681 _process_variables, _parameters, integration_order,
682 process_config, _mesh_vec, _media);
683 }
684 else
685#endif
686#ifdef OGS_BUILD_PROCESS_TH2M
687 if (type == "TH2M")
688 {
689 switch (_mesh_vec[0]->getDimension())
690 {
691 case 2:
693 name, *_mesh_vec[0], std::move(jacobian_assembler),
695 _local_coordinate_system, integration_order,
696 process_config, _media);
697 break;
698 case 3:
700 name, *_mesh_vec[0], std::move(jacobian_assembler),
702 _local_coordinate_system, integration_order,
703 process_config, _media);
704 break;
705 default:
706 OGS_FATAL("TH2M process does not support given dimension");
707 }
708 }
709 else
710#endif
711#ifdef OGS_BUILD_PROCESS_HEATCONDUCTION
712 if (type == "HEAT_CONDUCTION")
713 {
715 name, *_mesh_vec[0], std::move(jacobian_assembler),
716 _process_variables, _parameters, integration_order,
717 process_config, _media);
718 }
719 else
720#endif
721#ifdef OGS_BUILD_PROCESS_HEATTRANSPORTBHE
722 if (type == "HEAT_TRANSPORT_BHE")
723 {
724 if (_mesh_vec[0]->getDimension() != 3)
725 {
726 OGS_FATAL(
727 "HEAT_TRANSPORT_BHE can only work with a 3-dimensional "
728 "mesh! ");
729 }
730
731 process =
733 name, *_mesh_vec[0], std::move(jacobian_assembler),
734 _process_variables, _parameters, integration_order,
735 process_config, _curves, _media);
736 }
737 else
738#endif
739#ifdef OGS_BUILD_PROCESS_WELLBORESIMULATOR
740 if (type == "WELLBORE_SIMULATOR")
741 {
742 if (_mesh_vec[0]->getDimension() != 1)
743 {
744 OGS_FATAL(
745 "WELLBORE_SIMULATOR can only work with a 1-dimensional "
746 "mesh!");
747 }
748
749 process =
751 name, *_mesh_vec[0], std::move(jacobian_assembler),
752 _process_variables, _parameters, integration_order,
753 process_config, _media);
754 }
755 else
756#endif
757#ifdef OGS_BUILD_PROCESS_HYDROMECHANICS
758 if (type == "HYDRO_MECHANICS")
759 {
760 if (
761 process_config.getConfigParameterOptional<int>("dimension"))
762 {
763 OGS_FATAL(
764 "The 'dimension' tag has been removed in the merge-request "
765 "!4766. The dimension is now taken from the main mesh and "
766 "the tag must be removed. There is a python script in the "
767 "merge-request description for automatic conversion.");
768 }
769 switch (_mesh_vec[0]->getDimension())
770 {
771 case 2:
772 process =
774 2>(name, *_mesh_vec[0],
775 std::move(jacobian_assembler),
777 _local_coordinate_system, integration_order,
778 process_config, _media);
779 break;
780 case 3:
781 process =
783 3>(name, *_mesh_vec[0],
784 std::move(jacobian_assembler),
786 _local_coordinate_system, integration_order,
787 process_config, _media);
788 break;
789 default:
790 OGS_FATAL(
791 "HYDRO_MECHANICS process does not support given "
792 "dimension");
793 }
794 }
795 else
796#endif
797#ifdef OGS_BUILD_PROCESS_LARGEDEFORMATION
798 if (type == "LARGE_DEFORMATION")
799 {
800 switch (_mesh_vec[0]->getDimension())
801 {
802 case 2:
805 name, *_mesh_vec[0], std::move(jacobian_assembler),
807 _local_coordinate_system, integration_order,
808 process_config, _media);
809 break;
810 case 3:
813 name, *_mesh_vec[0], std::move(jacobian_assembler),
815 _local_coordinate_system, integration_order,
816 process_config, _media);
817 break;
818 default:
819 OGS_FATAL(
820 "LARGE_DEFORMATION process does not support given "
821 "dimension");
822 }
823 }
824 else
825#endif
826#ifdef OGS_BUILD_PROCESS_LIE_HM
827 if (type == "HYDRO_MECHANICS_WITH_LIE")
828 {
829 if (
830 process_config.getConfigParameterOptional<int>("dimension"))
831 {
832 OGS_FATAL(
833 "The 'dimension' tag has been removed in the merge-request "
834 "!4766."
835 "The dimension is now taken from the main mesh and the tag "
836 "must be"
837 "removed. There is a python script in the merge-request "
838 "description"
839 "for automatic conversion.");
840 }
841 switch (_mesh_vec[0]->getDimension())
842 {
843 case 2:
846 name, *_mesh_vec[0], std::move(jacobian_assembler),
848 _local_coordinate_system, integration_order,
849 process_config, _media);
850 break;
851 case 3:
854 name, *_mesh_vec[0], std::move(jacobian_assembler),
856 _local_coordinate_system, integration_order,
857 process_config, _media);
858 break;
859 default:
860 OGS_FATAL(
861 "HYDRO_MECHANICS_WITH_LIE process does not support "
862 "given dimension");
863 }
864 }
865 else
866#endif
867#ifdef OGS_BUILD_PROCESS_HT
868 if (type == "HT")
869 {
871 name, *_mesh_vec[0], std::move(jacobian_assembler),
872 _process_variables, _parameters, integration_order,
873 process_config, _mesh_vec, _media);
874 }
875 else
876#endif
877#ifdef OGS_BUILD_PROCESS_COMPONENTTRANSPORT
878 if (type == "ComponentTransport")
879 {
880 process =
882 name, *_mesh_vec[0], std::move(jacobian_assembler),
883 _process_variables, _parameters, integration_order,
884 process_config, _mesh_vec, _media,
885 std::move(chemical_solver_interface));
886 }
887 else
888#endif
889#ifdef OGS_BUILD_PROCESS_PHASEFIELD
890 if (type == "PHASE_FIELD")
891 {
892 switch (_mesh_vec[0]->getDimension())
893 {
894 case 2:
895 process =
897 name, *_mesh_vec[0], std::move(jacobian_assembler),
899 _local_coordinate_system, integration_order,
900 process_config);
901 break;
902 case 3:
903 process =
905 name, *_mesh_vec[0], std::move(jacobian_assembler),
907 _local_coordinate_system, integration_order,
908 process_config);
909 break;
910 }
911 }
912 else
913#endif
914#ifdef OGS_BUILD_PROCESS_HMPHASEFIELD
915 if (type == "HM_PHASE_FIELD")
916 {
917 switch (_mesh_vec[0]->getDimension())
918 {
919 case 2:
920 process =
922 name, *_mesh_vec[0], std::move(jacobian_assembler),
924 _local_coordinate_system, integration_order,
925 process_config, _media);
926 break;
927 case 3:
928 process =
930 name, *_mesh_vec[0], std::move(jacobian_assembler),
932 _local_coordinate_system, integration_order,
933 process_config, _media);
934 break;
935 }
936 }
937 else
938#endif
939#ifdef OGS_BUILD_PROCESS_RICHARDSCOMPONENTTRANSPORT
940 if (type == "RichardsComponentTransport")
941 {
944 name, *_mesh_vec[0], std::move(jacobian_assembler),
945 _process_variables, _parameters, integration_order,
946 process_config, _media);
947 }
948 else
949#endif
950#ifdef OGS_BUILD_PROCESS_SMALLDEFORMATION
951 if (type == "SMALL_DEFORMATION")
952 {
953 switch (_mesh_vec[0]->getDimension())
954 {
955 case 2:
958 name, *_mesh_vec[0], std::move(jacobian_assembler),
960 _local_coordinate_system, integration_order,
961 process_config, _media);
962 break;
963 case 3:
966 name, *_mesh_vec[0], std::move(jacobian_assembler),
968 _local_coordinate_system, integration_order,
969 process_config, _media);
970 break;
971 default:
972 OGS_FATAL(
973 "SMALL_DEFORMATION process does not support given "
974 "dimension");
975 }
976 }
977 else
978#endif
979#ifdef OGS_BUILD_PROCESS_LIE_M
980 if (type == "SMALL_DEFORMATION_WITH_LIE")
981 {
982 if (
983 process_config.getConfigParameterOptional<int>("dimension"))
984 {
985 OGS_FATAL(
986 "The 'dimension' tag has been removed in the merge-request "
987 "!4766."
988 "The dimension is now taken from the main mesh and the tag "
989 "must be"
990 "removed. There is a python script in the merge-request "
991 "description"
992 "for automatic conversion.");
993 }
994 switch (_mesh_vec[0]->getDimension())
995 {
996 case 2:
999 name, *_mesh_vec[0], std::move(jacobian_assembler),
1001 _local_coordinate_system, integration_order,
1002 process_config);
1003 break;
1004 case 3:
1007 name, *_mesh_vec[0], std::move(jacobian_assembler),
1009 _local_coordinate_system, integration_order,
1010 process_config);
1011 break;
1012 default:
1013 OGS_FATAL(
1014 "SMALL_DEFORMATION_WITH_LIE process does not support "
1015 "given dimension");
1016 }
1017 }
1018 else
1019#endif
1020#ifdef OGS_BUILD_PROCESS_THERMOHYDROMECHANICS
1021 if (type == "THERMO_HYDRO_MECHANICS")
1022 {
1023 if (
1024 process_config.getConfigParameterOptional<int>("dimension"))
1025 {
1026 OGS_FATAL(
1027 "The 'dimension' tag has been removed in the merge-request "
1028 "!4766."
1029 "The dimension is now taken from the main mesh and the tag "
1030 "must be"
1031 "removed. There is a python script in the merge-request "
1032 "description"
1033 "for automatic conversion.");
1034 }
1035 switch (_mesh_vec[0]->getDimension())
1036 {
1037 case 2:
1040 name, *_mesh_vec[0], std::move(jacobian_assembler),
1042 _local_coordinate_system, integration_order,
1043 process_config, _media);
1044 break;
1045 case 3:
1048 name, *_mesh_vec[0], std::move(jacobian_assembler),
1050 _local_coordinate_system, integration_order,
1051 process_config, _media);
1052 break;
1053 default:
1054 OGS_FATAL(
1055 "THERMO_HYDRO_MECHANICS process does not support given "
1056 "dimension");
1057 }
1058 }
1059 else
1060#endif
1061#ifdef OGS_BUILD_PROCESS_THERMOMECHANICS
1062 if (type == "THERMO_MECHANICS")
1063 {
1064 switch (_mesh_vec[0]->getDimension())
1065 {
1066 case 2:
1069 name, *_mesh_vec[0], std::move(jacobian_assembler),
1071 _local_coordinate_system, integration_order,
1072 process_config, _media);
1073 break;
1074 case 3:
1077 name, *_mesh_vec[0], std::move(jacobian_assembler),
1079 _local_coordinate_system, integration_order,
1080 process_config, _media);
1081 break;
1082 }
1083 }
1084 else
1085#endif
1086#ifdef OGS_BUILD_PROCESS_RICHARDSFLOW
1087 if (type == "RICHARDS_FLOW")
1088 {
1090 name, *_mesh_vec[0], std::move(jacobian_assembler),
1091 _process_variables, _parameters, integration_order,
1092 process_config, _media);
1093 }
1094 else
1095#endif
1096#ifdef OGS_BUILD_PROCESS_RICHARDSMECHANICS
1097 if (type == "RICHARDS_MECHANICS")
1098 {
1099 if (
1100 process_config.getConfigParameterOptional<int>("dimension"))
1101 {
1102 OGS_FATAL(
1103 "The 'dimension' tag has been removed in the merge-request "
1104 "!4766."
1105 "The dimension is now taken from the main mesh and the tag "
1106 "must be"
1107 "removed. There is a python script in the merge-request "
1108 "description"
1109 "for automatic conversion.");
1110 }
1111 switch (_mesh_vec[0]->getDimension())
1112 {
1113 case 2:
1116 name, *_mesh_vec[0], std::move(jacobian_assembler),
1118 _local_coordinate_system, integration_order,
1119 process_config, _media);
1120 break;
1121 case 3:
1124 name, *_mesh_vec[0], std::move(jacobian_assembler),
1126 _local_coordinate_system, integration_order,
1127 process_config, _media);
1128 break;
1129 }
1130 }
1131 else
1132#endif
1133#ifdef OGS_BUILD_PROCESS_THERMORICHARDSFLOW
1134 if (type == "THERMO_RICHARDS_FLOW")
1135 {
1136 process =
1138 name, *_mesh_vec[0], std::move(jacobian_assembler),
1139 _process_variables, _parameters, integration_order,
1140 process_config, _media);
1141 }
1142 else
1143#endif
1144#ifdef OGS_BUILD_PROCESS_THERMORICHARDSMECHANICS
1145 if (type == "THERMO_RICHARDS_MECHANICS")
1146 {
1147 switch (_mesh_vec[0]->getDimension())
1148 {
1149 case 2:
1152 name, *_mesh_vec[0], std::move(jacobian_assembler),
1154 _local_coordinate_system, integration_order,
1155 process_config, _media);
1156 break;
1157 case 3:
1160 name, *_mesh_vec[0], std::move(jacobian_assembler),
1162 _local_coordinate_system, integration_order,
1163 process_config, _media);
1164 break;
1165 }
1166 }
1167 else
1168#endif
1169
1170#ifdef OGS_BUILD_PROCESS_TWOPHASEFLOWWITHPP
1171 if (type == "TWOPHASE_FLOW_PP")
1172 {
1173 process =
1175 name, *_mesh_vec[0], std::move(jacobian_assembler),
1176 _process_variables, _parameters, integration_order,
1177 process_config, _media);
1178 }
1179 else
1180#endif
1181#ifdef OGS_BUILD_PROCESS_THERMALTWOPHASEFLOWWITHPP
1182 if (type == "THERMAL_TWOPHASE_WITH_PP")
1183 {
1186 name, *_mesh_vec[0], std::move(jacobian_assembler),
1187 _process_variables, _parameters, integration_order,
1188 process_config, _media);
1189 }
1190 else
1191#endif
1192 {
1193 OGS_FATAL("Unknown process type: {:s}", type);
1194 }
1195
1196 if (ranges::contains(_processes, name,
1197 [](std::unique_ptr<ProcessLib::Process> const& p)
1198 { return p->name; }))
1199 {
1200 OGS_FATAL("The process name '{:s}' is not unique.", name);
1201 }
1202 _processes.push_back(std::move(process));
1203 }
1204}
1205
1207 std::string const& output_directory)
1208{
1209 DBUG("Reading time loop configuration.");
1210
1211 bool const compensate_non_equilibrium_initial_residuum = std::any_of(
1212 std::begin(_process_variables),
1213 std::end(_process_variables),
1214 [](auto const& process_variable)
1215 { return process_variable.compensateNonEquilibriumInitialResiduum(); });
1216
1218 config, output_directory, _processes, _nonlinear_solvers, _mesh_vec,
1219 compensate_non_equilibrium_initial_residuum);
1220
1221 if (!_time_loop)
1222 {
1223 OGS_FATAL("Initialization of time loop failed.");
1224 }
1225}
1226
1228{
1229 DBUG("Reading linear solver configuration.");
1230
1232 for (auto conf : config.getConfigSubtreeList("linear_solver"))
1233 {
1235 auto const name = conf.getConfigParameter<std::string>("name");
1236 auto const linear_solver_parser =
1238 auto const solver_options =
1239 linear_solver_parser.parseNameAndOptions("", &conf);
1240
1243 name,
1244 std::make_unique<GlobalLinearSolver>(std::get<0>(solver_options),
1245 std::get<1>(solver_options)),
1246 "The linear solver name is not unique");
1247 }
1248}
1249
1251{
1252 DBUG("Reading non-linear solver configuration.");
1253
1255 for (auto conf : config.getConfigSubtreeList("nonlinear_solver"))
1256 {
1257 auto const ls_name =
1259 conf.getConfigParameter<std::string>("linear_solver");
1260 auto const& linear_solver = BaseLib::getOrError(
1261 _linear_solvers, ls_name,
1262 "A linear solver with the given name does not exist.");
1263
1265 auto const name = conf.getConfigParameter<std::string>("name");
1268 name,
1269 NumLib::createNonlinearSolver(*linear_solver, conf).first,
1270 "The nonlinear solver name is not unique");
1271 }
1272}
1273
1274void ProjectData::parseCurves(std::optional<BaseLib::ConfigTree> const& config)
1275{
1276 if (!config)
1277 {
1278 return;
1279 }
1280
1281 DBUG("Reading curves configuration.");
1282
1284 for (auto conf : config->getConfigSubtreeList("curve"))
1285 {
1287 auto const name = conf.getConfigParameter<std::string>("name");
1289 _curves,
1290 name,
1293 "The curve name is not unique.");
1294 }
1295}
1296
1297MeshLib::Mesh& ProjectData::getMesh(std::string const& mesh_name) const
1298{
1299 return MeshLib::findMeshByName(_mesh_vec, mesh_name);
1300}
1301
1302std::vector<std::string> ProjectData::getMeshNames() const
1303{
1304 return _mesh_vec | MeshLib::views::names | ranges::to<std::vector>;
1305}
Chemical-solver interface used in OGS operator-split reactive transport.
#define OGS_FATAL(...)
Definition Error.h:19
std::vector< std::size_t > getNodes(GeoLib::Point const &pnt, std::vector< MeshLib::Node * > const &nodes, MeshLib::PropertyVector< int > const &mat_ids, std::pair< int, int > const &mat_limits, std::pair< double, double > const &elevation_limits, MeshLib::Mesh const &mesh)
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
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:34
std::optional< ConfigTree > getConfigSubtreeOptional(std::string const &root) const
std::filesystem::path projectDirectory() const
std::optional< T > getConfigParameterOptional(std::string const &param) const
T getConfigParameter(std::string const &param) const
Range< SubtreeIterator > getConfigSubtreeList(std::string const &root) const
std::optional< T > getConfigAttributeOptional(std::string const &attr) const
Container class for geometric objects.
Definition GEOObjects.h:46
std::vector< PolylineVec * > const & getPolylines() const
Read access to polylines w/o using a name.
Definition GEOObjects.h:288
std::vector< PointVec * > const & getPoints() const
Read access to points w/o using a name.
Definition GEOObjects.h:286
std::vector< SurfaceVec * > const & getSurfaces() const
Read access to surfaces w/o using a name.
Definition GEOObjects.h:290
bool readFile(const std::string &fname) override
Reads an xml-file containing OGS geometry.
static PROCESSLIB_EXPORT const std::string constant_one_parameter_name
Definition Process.h:40
std::map< std::string, std::unique_ptr< GlobalLinearSolver > > _linear_solvers
std::optional< ParameterLib::CoordinateSystem > _local_coordinate_system
std::map< std::string, std::unique_ptr< NumLib::NonlinearSolverBase > > _nonlinear_solvers
MeshLib::Mesh & getMesh(std::string const &mesh_name) const
std::unique_ptr< ChemistryLib::ChemicalSolverInterface > parseChemicalSolverInterface(std::optional< BaseLib::ConfigTree > const &config, const std::string &output_directory)
void parseMedia(std::optional< BaseLib::ConfigTree > const &media_config)
Parses media configuration and saves them in an object.
void parseLinearSolvers(BaseLib::ConfigTree const &config)
void parseProcesses(BaseLib::ConfigTree const &processes_config, std::string const &output_directory, std::unique_ptr< ChemistryLib::ChemicalSolverInterface > &&chemical_solver_interface)
std::vector< std::unique_ptr< ParameterLib::ParameterBase > > _parameters
Buffer for each parameter config passed to the process.
std::vector< std::unique_ptr< ProcessLib::Process > > _processes
std::vector< GeoLib::NamedRaster > _named_rasters
void parseNonlinearSolvers(BaseLib::ConfigTree const &config)
std::vector< ProcessLib::ProcessVariable > _process_variables
void parseProcessVariables(BaseLib::ConfigTree const &process_variables_config)
std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation > > _curves
void parseCurves(std::optional< BaseLib::ConfigTree > const &config)
std::vector< std::unique_ptr< MeshLib::Mesh > > _mesh_vec
std::vector< std::string > getMeshNames() const
void parseTimeLoop(BaseLib::ConfigTree const &config, const std::string &output_directory)
Parses the time loop configuration.
std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > _media
std::vector< std::string > parseParameters(BaseLib::ConfigTree const &parameters_config)
std::unique_ptr< ProcessLib::TimeLoop > _time_loop
The time loop used to solve this project's processes.
std::vector< std::unique_ptr< MeshLib::Mesh > > readMeshes(std::vector< std::string > const &filenames)
void insertIfKeyUniqueElseError(Map &map, Key const &key, Value &&value, std::string const &error_message)
Definition Algorithm.h:97
bool IsFileExisting(const std::string &strFilename)
Returns true if given file exists.
Definition FileTools.cpp:23
std::string joinPaths(std::string const &pathA, std::string const &pathB)
OGS_NO_DANGLING Map::mapped_type & getOrError(Map &map, Key const &key, std::string const &error_message)
Definition Algorithm.h:111
std::unique_ptr< ChemicalSolverInterface > createChemicalSolverInterface(std::vector< std::unique_ptr< MeshLib::Mesh > > const &meshes, std::map< std::string, std::unique_ptr< GlobalLinearSolver > > const &linear_solvers, BaseLib::ConfigTree const &config, std::string const &output_directory)
void createMediumForId(int const id, std::map< int, std::shared_ptr< T > > &media, std::vector< int > const &material_ids_of_this_medium, CreateMedium &&create_medium)
std::vector< int > parseMaterialIdString(std::string const &material_id_string, MeshLib::PropertyVector< int > const *const material_ids)
std::unique_ptr< Medium > createMedium(int const material_id, int const geometry_dimension, BaseLib::ConfigTree const &config, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > &parameters, ParameterLib::CoordinateSystem const *const local_coordinate_system, std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation > > const &curves)
std::unique_ptr< CurveType > createPiecewiseLinearCurve(BaseLib::ConfigTree const &config)
std::unique_ptr< MeshGeoToolsLib::SearchLength > createSearchLengthAlgorithm(BaseLib::ConfigTree const &external_config, MeshLib::Mesh const &mesh)
std::vector< std::unique_ptr< MeshLib::Mesh > > constructAdditionalMeshesFromGeoObjects(GeoLib::GEOObjects const &geo_objects, MeshLib::Mesh const &mesh, std::unique_ptr< SearchLength > search_length_algorithm, bool const multiple_nodes_allowed)
MeshLib::Mesh * readMeshFromFile(const std::string &file_name, bool const compute_element_neighbors)
constexpr ranges::views::view_closure names
For an element of a range view return its name.
Definition Mesh.h:222
Mesh & findMeshByName(std::vector< std::unique_ptr< Mesh > > const &meshes, std::string_view const name)
Definition Mesh.cpp:354
void setMeshSpaceDimension(std::vector< std::unique_ptr< Mesh > > const &meshes)
void zeroMeshFieldDataByMaterialIDs(MeshLib::Mesh &mesh, std::vector< int > const &selected_material_ids)
std::pair< std::unique_ptr< NonlinearSolverBase >, NonlinearSolverTag > createNonlinearSolver(GlobalLinearSolver &linear_solver, BaseLib::ConfigTree const &config)
std::optional< ParameterLib::CoordinateSystem > createCoordinateSystem(std::optional< BaseLib::ConfigTree > const &config, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters)
std::unique_ptr< ParameterBase > createParameter(BaseLib::ConfigTree const &config, std::vector< std::unique_ptr< MeshLib::Mesh > > const &meshes, std::vector< GeoLib::NamedRaster > const &named_rasters, std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation > > const &curves)
std::unique_ptr< Process > createComponentTransportProcess(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 > > const &parameters, unsigned const integration_order, BaseLib::ConfigTree const &config, std::vector< std::unique_ptr< MeshLib::Mesh > > const &meshes, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media, std::unique_ptr< ChemistryLib::ChemicalSolverInterface > &&chemical_solver_interface)
template std::unique_ptr< Process > createHMPhaseFieldProcess< 3 >(std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< ProcessVariable > const &variables, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< Process > createHMPhaseFieldProcess< 2 >(std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< ProcessVariable > const &variables, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
std::unique_ptr< Process > createHTProcess(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 > > const &parameters, unsigned const integration_order, BaseLib::ConfigTree const &config, std::vector< std::unique_ptr< MeshLib::Mesh > > const &meshes, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
std::unique_ptr< Process > createHeatConductionProcess(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 > > const &parameters, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
std::unique_ptr< Process > 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)
std::unique_ptr< Process > createHydroMechanicsProcess(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 > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< Process > createHydroMechanicsProcess< 2 >(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 > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< Process > createHydroMechanicsProcess< 3 >(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 > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< Process > createSmallDeformationProcess< 2 >(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 > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config)
template std::unique_ptr< Process > createSmallDeformationProcess< 3 >(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 > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config)
template std::unique_ptr< Process > createLargeDeformationProcess< 2 >(std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< ProcessVariable > const &variables, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< Process > createLargeDeformationProcess< 3 >(std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< ProcessVariable > const &variables, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
std::unique_ptr< Process > createLiquidFlowProcess(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 > > const &parameters, unsigned const integration_order, BaseLib::ConfigTree const &config, std::vector< std::unique_ptr< MeshLib::Mesh > > const &meshes, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< Process > createPhaseFieldProcess< 2 >(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 > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config)
template std::unique_ptr< Process > createPhaseFieldProcess< 3 >(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 > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config)
std::unique_ptr< Process > createRichardsComponentTransportProcess(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 > > const &parameters, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
std::unique_ptr< Process > createRichardsFlowProcess(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 > > const &parameters, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< Process > createRichardsMechanicsProcess< 3 >(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 > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< Process > createRichardsMechanicsProcess< 2 >(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 > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< Process > createSmallDeformationProcess< 3 >(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 > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< Process > createSmallDeformationProcess< 2 >(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 > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
std::unique_ptr< Process > createSteadyStateDiffusion(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 > > const &parameters, unsigned const integration_order, BaseLib::ConfigTree const &config, std::vector< std::unique_ptr< MeshLib::Mesh > > const &meshes, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< Process > createTH2MProcess< 3 >(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 > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< Process > createTH2MProcess< 2 >(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 > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
std::unique_ptr< Process > createThermalTwoPhaseFlowWithPPProcess(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 > > const &parameters, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< Process > createThermoHydroMechanicsProcess< 2 >(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 > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< Process > createThermoHydroMechanicsProcess< 3 >(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 > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< Process > createThermoMechanicsProcess< 3 >(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 > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< Process > createThermoMechanicsProcess< 2 >(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 > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
std::unique_ptr< Process > createThermoRichardsFlowProcess(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 > > const &parameters, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< Process > createThermoRichardsMechanicsProcess< 3 >(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 > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
template std::unique_ptr< Process > createThermoRichardsMechanicsProcess< 2 >(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 > > const &parameters, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
std::unique_ptr< Process > createTwoPhaseFlowWithPPProcess(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 > > const &parameters, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
std::unique_ptr< Process > createWellboreSimulatorProcess(std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< ProcessVariable > const &variables, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, BaseLib::ConfigTree const &config, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
std::unique_ptr< TimeLoop > createTimeLoop(BaseLib::ConfigTree const &config, std::string const &output_directory, const std::vector< std::unique_ptr< Process > > &processes, const std::map< std::string, std::unique_ptr< NumLib::NonlinearSolverBase > > &nonlinear_solvers, std::vector< std::unique_ptr< MeshLib::Mesh > > &meshes, bool const compensate_non_equilibrium_initial_residuum)
Builds a TimeLoop from the given configuration.
std::unique_ptr< AbstractJacobianAssembler > createJacobianAssembler(std::optional< BaseLib::ConfigTree > const &config)
std::vector< GeoLib::NamedRaster > readRasters(BaseLib::ConfigTree const &config, GeoLib::MinMaxPoints const &min_max_points)
std::vector< std::unique_ptr< MeshLib::Mesh > > readMeshes(BaseLib::ConfigTree const &config, std::string const &directory)
void readGeometry(std::string const &fname, GeoLib::GEOObjects &geo_objects, std::string const &dir_first, std::string const &dir_second)
std::unique_ptr< MeshLib::Mesh > readSingleMesh(BaseLib::ConfigTree const &mesh_config_parameter, std::string const &directory)
Single, constant value parameter.
static PROCESSLIB_EXPORT const std::string zero_parameter_name