OGS
ThermoRichardsFlowProcess.cpp
Go to the documentation of this file.
1
12
13#include <cassert>
14
15#include "BaseLib/Error.h"
19#include "ProcessLib/Process.h"
23
24namespace ProcessLib
25{
26namespace ThermoRichardsFlow
27{
29 std::string name,
30 MeshLib::Mesh& mesh,
31 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
32 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
33 unsigned const integration_order,
34 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
35 process_variables,
36 ThermoRichardsFlowProcessData&& process_data,
37 SecondaryVariableCollection&& secondary_variables,
38 bool const use_monolithic_scheme)
39 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
40 integration_order, std::move(process_variables),
41 std::move(secondary_variables), use_monolithic_scheme),
42 _process_data(std::move(process_data))
43{
45 mesh, "HeatFlowRate", MeshLib::MeshItemType::Node, 1);
46
48 mesh, "MassFlowRate", MeshLib::MeshItemType::Node, 1);
49
50 // TODO (naumov) remove ip suffix. Probably needs modification of the mesh
51 // properties, s.t. there is no "overlapping" with cell/point data.
52 // See getOrCreateMeshProperty.
53 _integration_point_writer.emplace_back(
54 std::make_unique<MeshLib::IntegrationPointWriter>(
55 "saturation_ip", 1 /*n components*/, integration_order,
57
58 _integration_point_writer.emplace_back(
59 std::make_unique<MeshLib::IntegrationPointWriter>(
60 "porosity_ip", 1 /*n components*/, integration_order,
62}
63
65 NumLib::LocalToGlobalIndexMap const& dof_table,
66 MeshLib::Mesh const& mesh,
67 unsigned const integration_order)
68{
70 mesh.getDimension(), mesh.getElements(), dof_table, _local_assemblers,
71 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
73
74 auto add_secondary_variable = [&](std::string const& name,
75 int const num_components,
76 auto get_ip_values_function)
77 {
79 name,
80 makeExtrapolator(num_components, getExtrapolator(),
82 std::move(get_ip_values_function)));
83 };
84
85 add_secondary_variable("velocity", mesh.getDimension(),
87
88 add_secondary_variable("saturation", 1,
90
91 add_secondary_variable("porosity", 1, &LocalAssemblerIF::getIntPtPorosity);
92
93 add_secondary_variable("dry_density_solid", 1,
95
97 const_cast<MeshLib::Mesh&>(mesh), "saturation_avg",
99
101 const_cast<MeshLib::Mesh&>(mesh), "porosity_avg",
103
106
107 // Initialize local assemblers after all variables have been set.
111}
112
114 std::vector<GlobalVector*>& x, double const t, int const process_id)
115{
116 if (process_id != 0)
117 {
118 return;
119 }
120 DBUG("SetInitialConditions ThermoRichardsFlowProcess.");
121
124 getDOFTables(x.size()), x, t, process_id);
125}
126
128 const double t, double const dt, std::vector<GlobalVector*> const& x,
129 std::vector<GlobalVector*> const& x_prev, int const process_id,
131{
132 DBUG("Assemble the equations for ThermoRichardsFlowProcess.");
133
134 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
136
137 // Call global assembler for each local assembly item.
140 getActiveElementIDs(), dof_table, t, dt, x, x_prev, process_id, &M, &K,
141 &b);
142}
143
145 const double t, double const dt, std::vector<GlobalVector*> const& x,
146 std::vector<GlobalVector*> const& x_prev, int const process_id,
147 GlobalVector& b, GlobalMatrix& Jac)
148{
149 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
150
151 DBUG(
152 "Assemble the Jacobian of ThermoRichardsFlow for the monolithic "
153 "scheme.");
154 dof_tables.emplace_back(_local_to_global_index_map.get());
155
157 dof_tables, t, dt, x, x_prev, process_id, b,
158 Jac);
159
160 auto copyRhs = [&](int const variable_id, auto& output_vector)
161 {
162 transformVariableFromGlobalVector(b, variable_id, *dof_tables[0],
163 output_vector, std::negate<double>());
164 };
165
166 copyRhs(0, *_heat_flux);
167 copyRhs(1, *_hydraulic_flow);
168}
169
171 std::vector<GlobalVector*> const& x,
172 std::vector<GlobalVector*> const& x_prev, double const t, double const dt,
173 const int process_id)
174{
175 if (process_id != 0)
176 {
177 return;
178 }
179
180 DBUG("PostTimestep ThermoRichardsFlowProcess.");
181
184 getActiveElementIDs(), getDOFTables(x.size()), x, x_prev, t, dt,
185 process_id);
186}
187
189 const double t, const double dt, std::vector<GlobalVector*> const& x,
190 GlobalVector const& x_prev, int const process_id)
191{
192 if (process_id != 0)
193 {
194 return;
195 }
196 DBUG("Compute the secondary variables for ThermoRichardsFlowProcess.");
197
200 getActiveElementIDs(), getDOFTables(x.size()), t, dt, x, x_prev,
201 process_id);
202}
203
204} // namespace ThermoRichardsFlow
205} // namespace ProcessLib
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
Global vector based on Eigen vector.
Definition EigenVector.h:25
bool isAxiallySymmetric() const
Definition Mesh.h:137
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:109
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition Mesh.h:88
Properties & getProperties()
Definition Mesh.h:134
void assembleWithJacobian(BaseLib::PolymorphicRandomAccessContainerView< LocalAssemblerInterface > const &local_assemblers, std::vector< std::size_t > const &active_elements, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &xs, std::vector< GlobalVector * > const &x_prevs, int const process_id, GlobalVector &b, GlobalMatrix &Jac)
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:362
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:385
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Definition Process.h:390
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:167
SecondaryVariableCollection _secondary_variables
Definition Process.h:370
VectorMatrixAssembler _global_assembler
Definition Process.h:377
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:368
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:208
Handles configuration of several secondary variables from the project file.
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
void postTimestepConcreteProcess(std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, const int process_id) override
std::vector< std::unique_ptr< LocalAssemblerIF > > _local_assemblers
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
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)
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
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
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)
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 executeSelectedMemberDereferenced(Object &object, 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 > getSaturation() const =0
virtual std::vector< double > getPorosity() 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