OGS
ThermoRichardsFlowProcess.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
8#include "BaseLib/Error.h"
12#include "ProcessLib/Process.h"
16
17namespace ProcessLib
18{
19namespace ThermoRichardsFlow
20{
22 std::string name, MeshLib::Mesh& mesh,
23 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
24 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
25 unsigned const integration_order,
26 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
27 process_variables,
28 ThermoRichardsFlowProcessData&& process_data,
29 SecondaryVariableCollection&& secondary_variables,
30 bool const use_monolithic_scheme, bool const is_linear)
31 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
32 integration_order, std::move(process_variables),
33 std::move(secondary_variables), use_monolithic_scheme),
35 use_monolithic_scheme},
36 _process_data(std::move(process_data))
37{
38 // For numerical Jacobian
39 this->_jacobian_assembler->setNonDeformationComponentIDs({0, 1} /*T, p */);
40
41 // TODO (naumov) remove ip suffix. Probably needs modification of the mesh
42 // properties, s.t. there is no "overlapping" with cell/point data.
43 // See getOrCreateMeshProperty.
44 _integration_point_writer.emplace_back(
45 std::make_unique<MeshLib::IntegrationPointWriter>(
46 "saturation_ip", 1 /*n components*/, integration_order,
47 local_assemblers_, &LocalAssemblerIF::getSaturation));
48
49 _integration_point_writer.emplace_back(
50 std::make_unique<MeshLib::IntegrationPointWriter>(
51 "porosity_ip", 1 /*n components*/, integration_order,
52 local_assemblers_, &LocalAssemblerIF::getPorosity));
53}
54
56 NumLib::LocalToGlobalIndexMap const& dof_table,
57 MeshLib::Mesh const& mesh,
58 unsigned const integration_order)
59{
61 mesh.getDimension(), mesh.getElements(), dof_table, local_assemblers_,
62 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
64
65 auto add_secondary_variable = [&](std::string const& name,
66 int const num_components,
67 auto get_ip_values_function)
68 {
69 _secondary_variables.addSecondaryVariable(
70 name,
71 makeExtrapolator(num_components, getExtrapolator(),
73 std::move(get_ip_values_function)));
74 };
75
76 add_secondary_variable("velocity", mesh.getDimension(),
78
79 add_secondary_variable("saturation", 1,
81
82 add_secondary_variable("porosity", 1, &LocalAssemblerIF::getIntPtPorosity);
83
84 add_secondary_variable("dry_density_solid", 1,
86
88 const_cast<MeshLib::Mesh&>(mesh), "saturation_avg",
90
92 const_cast<MeshLib::Mesh&>(mesh), "porosity_avg",
94
97
98 // Initialize local assemblers after all variables have been set.
102}
103
105 std::vector<GlobalVector*>& x, double const t, int const process_id)
106{
107 if (process_id != 0)
108 {
109 return;
110 }
111 DBUG("SetInitialConditions ThermoRichardsFlowProcess.");
112
115 getDOFTables(x.size()), x, t, process_id);
116}
117
118std::vector<std::vector<std::string>>
120 std::vector<std::reference_wrapper<MeshLib::Mesh>> const& meshes)
121{
122 DBUG("TRF process initializeSubmeshOutput().");
123
124 std::vector<std::vector<std::string>> per_process_residuum_names;
125
126 if (_process_variables.size() == 1) // monolithic
127 {
128 per_process_residuum_names = {{"HeatFlowRate", "MassFlowRate"}};
129 }
130 else // staggered
131 {
132 per_process_residuum_names = {{"HeatFlowRate"}, {"MassFlowRate"}};
133 }
134
136 meshes, per_process_residuum_names);
137
138 return per_process_residuum_names;
139}
140
142 const double t, double const dt, std::vector<GlobalVector*> const& x,
143 std::vector<GlobalVector*> const& x_prev, int const process_id,
145{
146 DBUG("Assemble the equations for ThermoRichardsFlowProcess.");
147
149 process_id, M, K, b);
150}
151
153 const double t, double const dt, std::vector<GlobalVector*> const& x,
154 std::vector<GlobalVector*> const& x_prev, int const process_id,
155 GlobalVector& b, GlobalMatrix& Jac)
156{
157 DBUG(
158 "Assemble the Jacobian of ThermoRichardsFlow for the monolithic "
159 "scheme.");
160
162 t, dt, x, x_prev, process_id, b, Jac);
163}
164
166 std::vector<GlobalVector*> const& /*x*/,
167 const double /*t*/,
168 const double /*dt*/,
169 const int process_id)
170{
171 if (process_id == 0)
172 {
174 }
175}
176
178 std::vector<GlobalVector*> const& x,
179 std::vector<GlobalVector*> const& x_prev, double const t, double const dt,
180 const int process_id)
181{
182 if (process_id != 0)
183 {
184 return;
185 }
186
187 DBUG("PostTimestep ThermoRichardsFlowProcess.");
188
191 getActiveElementIDs(), getDOFTables(x.size()), x, x_prev, t, dt,
192 process_id);
193}
194
196 const double t, const double dt, std::vector<GlobalVector*> const& x,
197 GlobalVector const& x_prev, int const process_id)
198{
199 if (process_id != 0)
200 {
201 return;
202 }
203 DBUG("Compute the secondary variables for ThermoRichardsFlowProcess.");
204
207 getActiveElementIDs(), getDOFTables(x.size()), t, dt, x, x_prev,
208 process_id);
209}
210
211} // namespace ThermoRichardsFlow
212} // namespace ProcessLib
MathLib::EigenMatrix GlobalMatrix
MathLib::EigenVector GlobalVector
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
bool isAxiallySymmetric() const
Definition Mesh.h:128
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:100
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition Mesh.h:79
Properties & getProperties()
Definition Mesh.h:125
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)
virtual void postTimestep(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
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)
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
virtual void setInitialConditions(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, int const process_id)
std::string const name
Definition Process.h:361
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:376
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Definition Process.h:389
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::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
SecondaryVariableCollection _secondary_variables
Definition Process.h:369
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
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:201
Handles configuration of several secondary variables from the project file.
void postTimestepConcreteProcess(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, const int process_id) override
void preTimestepConcreteProcess(std::vector< GlobalVector * > const &, const double, const double, const int process_id) override
void setInitialConditionsConcreteProcess(std::vector< GlobalVector * > &x, double const t, int const) override
void computeSecondaryVariableConcrete(double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id) override
std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes(std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) override
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
ThermoRichardsFlowProcess(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, ThermoRichardsFlowProcessData &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme, bool const is_linear)
std::vector< std::unique_ptr< LocalAssemblerIF > > local_assemblers_
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
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
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)
void createLocalAssemblers(std::vector< MeshLib::Element * > const &mesh_elements, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< std::unique_ptr< LocalAssemblerInterface > > &local_assemblers, ProviderOrOrder const &provider_or_order, ExtraCtorArgs &&... extra_ctor_args)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
void setIPDataInitialConditions(std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const &_integration_point_writer, MeshLib::Properties const &mesh_properties, LocalAssemblersVector &local_assemblers)
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
static void executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)
virtual std::vector< double > const & getIntPtSaturation(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtDryDensitySolid(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtPorosity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0