OGS
ThermoRichardsFlowProcess.cpp
Go to the documentation of this file.
1 
12 
13 #include <cassert>
14 
15 #include "BaseLib/Error.h"
16 #include "MeshLib/Elements/Utils.h"
18 #include "ProcessLib/Process.h"
21 #include "ThermoRichardsFlowFEM.h"
22 
23 namespace ProcessLib
24 {
25 namespace ThermoRichardsFlow
26 {
28  std::string name,
29  MeshLib::Mesh& mesh,
30  std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
31  std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
32  unsigned const integration_order,
33  std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
34  process_variables,
35  ThermoRichardsFlowProcessData&& process_data,
36  SecondaryVariableCollection&& secondary_variables,
37  bool const use_monolithic_scheme)
38  : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
39  integration_order, std::move(process_variables),
40  std::move(secondary_variables), use_monolithic_scheme),
41  _process_data(std::move(process_data))
42 {
43  _heat_flux = MeshLib::getOrCreateMeshProperty<double>(
44  mesh, "HeatFlowRate", MeshLib::MeshItemType::Node, 1);
45 
46  _hydraulic_flow = MeshLib::getOrCreateMeshProperty<double>(
47  mesh, "MassFlowRate", MeshLib::MeshItemType::Node, 1);
48 
49  // TODO (naumov) remove ip suffix. Probably needs modification of the mesh
50  // properties, s.t. there is no "overlapping" with cell/point data.
51  // See getOrCreateMeshProperty.
52  _integration_point_writer.emplace_back(
53  std::make_unique<IntegrationPointWriter>(
54  "saturation_ip", 1 /*n components*/, integration_order,
56 
57  _integration_point_writer.emplace_back(
58  std::make_unique<IntegrationPointWriter>(
59  "porosity_ip", 1 /*n components*/, integration_order,
61 }
62 
64  NumLib::LocalToGlobalIndexMap const& dof_table,
65  MeshLib::Mesh const& mesh,
66  unsigned const integration_order)
67 {
68  ProcessLib::createLocalAssemblers<ThermoRichardsFlowLocalAssembler>(
69  mesh.getDimension(), mesh.getElements(), dof_table, _local_assemblers,
70  NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
72 
73  auto add_secondary_variable = [&](std::string const& name,
74  int const num_components,
75  auto get_ip_values_function)
76  {
78  name,
79  makeExtrapolator(num_components, getExtrapolator(),
81  std::move(get_ip_values_function)));
82  };
83 
84  add_secondary_variable("velocity", mesh.getDimension(),
86 
87  add_secondary_variable("saturation", 1,
89 
90  add_secondary_variable("porosity", 1, &LocalAssemblerIF::getIntPtPorosity);
91 
92  add_secondary_variable("dry_density_solid", 1,
94 
95  _process_data.element_saturation = MeshLib::getOrCreateMeshProperty<double>(
96  const_cast<MeshLib::Mesh&>(mesh), "saturation_avg",
98 
99  _process_data.element_porosity = MeshLib::getOrCreateMeshProperty<double>(
100  const_cast<MeshLib::Mesh&>(mesh), "porosity_avg",
102 
105 
106  // Initialize local assemblers after all variables have been set.
110 }
111 
113  std::vector<GlobalVector*>& x, double const t, int const process_id)
114 {
115  if (process_id != 0)
116  {
117  return;
118  }
119  DBUG("SetInitialConditions ThermoRichardsFlowProcess.");
120 
124  process_id);
125 }
126 
128  const double t, double const dt, std::vector<GlobalVector*> const& x,
129  std::vector<GlobalVector*> const& xdot, int const process_id,
131 {
132  DBUG("Assemble the equations for ThermoRichardsFlowProcess.");
133 
134  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
135  dof_table = {std::ref(*_local_to_global_index_map)};
136  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
137 
138  // Call global assembler for each local assembly item.
141  pv.getActiveElementIDs(), dof_table, t, dt, x, xdot, process_id, M, K,
142  b);
143 }
144 
146  const double t, double const dt, std::vector<GlobalVector*> const& x,
147  std::vector<GlobalVector*> const& xdot, int const process_id,
149 {
150  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
151  dof_tables;
152 
153  DBUG(
154  "Assemble the Jacobian of ThermoRichardsFlow for the monolithic "
155  "scheme.");
156  dof_tables.emplace_back(*_local_to_global_index_map);
157 
158  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
159 
162  _local_assemblers, pv.getActiveElementIDs(), dof_tables, t, dt, x, xdot,
163  process_id, M, K, b, Jac);
164 
165  auto copyRhs = [&](int const variable_id, auto& output_vector)
166  {
167  transformVariableFromGlobalVector(b, variable_id, dof_tables[0],
168  output_vector, std::negate<double>());
169  };
170 
171  copyRhs(0, *_heat_flux);
172  copyRhs(1, *_hydraulic_flow);
173 }
174 
176  std::vector<GlobalVector*> const& x, double const t, double const dt,
177  const int process_id)
178 {
179  if (process_id != 0)
180  {
181  return;
182  }
183 
184  DBUG("PostTimestep ThermoRichardsFlowProcess.");
185 
186  auto const dof_tables = getDOFTables(x.size());
187 
188  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
191  pv.getActiveElementIDs(), dof_tables, x, t, dt);
192 }
193 
195  const double t, const double dt, std::vector<GlobalVector*> const& x,
196  GlobalVector const& x_dot, int const process_id)
197 {
198  if (process_id != 0)
199  {
200  return;
201  }
202  DBUG(
203  "Compute the secondary variables for "
204  "ThermoRichardsFlowProcess.");
205  auto const dof_tables = getDOFTables(x.size());
206  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
207 
210  pv.getActiveElementIDs(), dof_tables, t, dt, x, x_dot, process_id);
211 }
212 
213 std::vector<NumLib::LocalToGlobalIndexMap const*>
214 ThermoRichardsFlowProcess::getDOFTables(int const number_of_processes) const
215 {
216  std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
217  dof_tables.reserve(number_of_processes);
218  std::generate_n(std::back_inserter(dof_tables), number_of_processes,
219  [&]() { return _local_to_global_index_map.get(); });
220  return dof_tables;
221 }
222 
223 } // namespace ThermoRichardsFlow
224 } // namespace ProcessLib
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition: Logging.h:29
Global vector based on Eigen vector.
Definition: EigenVector.h:28
bool isAxiallySymmetric() const
Definition: Mesh.h:131
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition: Mesh.h:76
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition: Mesh.h:103
Properties & getProperties()
Definition: Mesh.h:128
virtual void postTimestep(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double const dt)
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
void setInitialConditions(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, bool const use_monolithic_scheme, 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_dot, int const process_id)
std::vector< std::size_t > const & getActiveElementIDs() const
std::string const name
Definition: Process.h:327
NumLib::Extrapolator & getExtrapolator() const
Definition: Process.h:184
std::vector< std::unique_ptr< IntegrationPointWriter > > _integration_point_writer
Definition: Process.h:354
SecondaryVariableCollection _secondary_variables
Definition: Process.h:335
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition: Process.h:148
VectorMatrixAssembler _global_assembler
Definition: Process.h:337
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition: Process.h:333
const bool _use_monolithic_scheme
Definition: Process.h:339
Handles configuration of several secondary variables from the project file.
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
std::vector< std::unique_ptr< LocalAssemblerIF > > _local_assemblers
void postTimestepConcreteProcess(std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
void setInitialConditionsConcreteProcess(std::vector< GlobalVector * > &x, double const t, int const) 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)
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(const int number_of_processes) const
void computeSecondaryVariableConcrete(double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, int const process_id) override
void assembleWithJacobianConcreteProcess(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, GlobalMatrix &Jac) override
void assembleConcreteProcess(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) override
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, 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, GlobalMatrix &Jac)
void assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)
static const double t
void transformVariableFromGlobalVector(GlobalVector const &input_vector, int const variable_id, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, MeshLib::PropertyVector< double > &output_vector, Functor mapFunction)
Definition: DOFTableUtil.h:59
void setIPDataInitialConditions(std::vector< std::unique_ptr< IntegrationPointWriter >> const &_integration_point_writer, MeshLib::Properties const &mesh_properties, LocalAssemblersVector &local_assemblers)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
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 & 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 > getSaturation() 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 & 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 > 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