OGS
ProcessLib::HT::HTProcess Class Referencefinal

Detailed Description

HT process

The implementation uses a monolithic approach, i.e., both processes are assembled within one global system of equations.

Process Coupling

The advective term of the heat conduction equation is given by the confined groundwater flow process, i.e., the heat conduction depends on darcy velocity of the groundwater flow process. On the other hand the temperature dependencies of the viscosity and density in the groundwater flow couples the H process to the T process.

Note
At the moment there is not any coupling by source or sink terms, i.e., the coupling is implemented only by density changes due to temperature changes in the buoyance term in the groundwater flow. This coupling schema is referred to as the Boussinesq approximation.

Definition at line 52 of file HTProcess.h.

#include <HTProcess.h>

Inheritance diagram for ProcessLib::HT::HTProcess:
[legend]
Collaboration diagram for ProcessLib::HT::HTProcess:
[legend]

Public Member Functions

 HTProcess (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, HTProcessData &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme, std::unique_ptr< ProcessLib::SurfaceFluxData > &&surfaceflux)
 
Eigen::Vector3d getFlux (std::size_t element_id, MathLib::Point3d const &p, double const t, std::vector< GlobalVector * > const &x) const override
 
void postTimestepConcreteProcess (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double delta_t, int const process_id) override
 
ODESystem interface
bool isLinear () const override
 
- Public Member Functions inherited from ProcessLib::Process
 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)
 
void preTimestep (std::vector< GlobalVector * > const &x, const double t, const double delta_t, const int process_id)
 Preprocessing before starting assembly for new timestep.
 
void postTimestep (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double delta_t, int const process_id)
 Postprocessing after a complete timestep.
 
void postNonLinearSolver (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, double const dt, int const process_id)
 
void preIteration (const unsigned iter, GlobalVector const &x) final
 
void computeSecondaryVariable (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
 compute secondary variables for the coupled equations or for output.
 
NumLib::IterationResult postIteration (GlobalVector const &x) final
 
void initialize (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
 
void setInitialConditions (std::vector< GlobalVector * > &process_solutions, std::vector< GlobalVector * > const &process_solutions_prev, double const t, int const process_id)
 
MathLib::MatrixSpecifications getMatrixSpecifications (const int process_id) const override
 
void updateDeactivatedSubdomains (double const time, const int process_id)
 
virtual bool isMonolithicSchemeUsed () const
 
virtual void extrapolateIntegrationPointValuesToNodes (const double, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > &)
 
void preAssemble (const double t, double const dt, GlobalVector const &x) final
 
void assemble (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) final
 
void assembleWithJacobian (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, GlobalMatrix &Jac) final
 
void preOutput (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)
 
std::vector< NumLib::IndexValueVector< GlobalIndexType > > const * getKnownSolutions (double const t, GlobalVector const &x, int const process_id) const final
 
virtual NumLib::LocalToGlobalIndexMap const & getDOFTable (const int) const
 
MeshLib::MeshgetMesh () const
 
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables (const int process_id) const
 
std::vector< std::size_t > const & getActiveElementIDs () const
 
SecondaryVariableCollection const & getSecondaryVariables () const
 
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const & getIntegrationPointWriters () const
 
virtual void solveReactionEquation (std::vector< GlobalVector * > &, std::vector< GlobalVector * > const &, double const, double const, NumLib::EquationSystem &, int const)
 
- Public Member Functions inherited from ProcessLib::SubmeshAssemblySupport
virtual std::vector< std::string > initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes)
 
virtual ~SubmeshAssemblySupport ()=default
 

Private Member Functions

void initializeConcreteProcess (NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
 Process specific initialization called by initialize().
 
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 assembleWithJacobianConcreteProcess (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, GlobalMatrix &Jac) override
 
std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const override
 

Private Attributes

HTProcessData _process_data
 
std::vector< std::unique_ptr< HTLocalAssemblerInterface > > _local_assemblers
 
std::unique_ptr< ProcessLib::SurfaceFluxData_surfaceflux
 

Additional Inherited Members

- Public Attributes inherited from ProcessLib::Process
std::string const name
 
- Static Public Attributes inherited from ProcessLib::Process
static PROCESSLIB_EXPORT const std::string constant_one_parameter_name = "constant_one"
 
- Protected Member Functions inherited from ProcessLib::Process
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables (int const number_of_processes) const
 
NumLib::ExtrapolatorgetExtrapolator () const
 
NumLib::LocalToGlobalIndexMap const & getSingleComponentDOFTable () const
 
void initializeProcessBoundaryConditionsAndSourceTerms (const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
 
virtual void constructDofTable ()
 
void constructMonolithicProcessDofTable ()
 
void constructDofTableOfSpecifiedProcessStaggeredScheme (const int specified_process_id)
 
std::vector< GlobalIndexTypegetIndicesOfResiduumWithoutInitialCompensation () const override
 
- Protected Attributes inherited from ProcessLib::Process
MeshLib::Mesh_mesh
 
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map
 
SecondaryVariableCollection _secondary_variables
 
std::unique_ptr< ProcessLib::AbstractJacobianAssembler_jacobian_assembler
 
VectorMatrixAssembler _global_assembler
 
const bool _use_monolithic_scheme
 
unsigned const _integration_order
 
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
 
GlobalSparsityPattern _sparsity_pattern
 
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
 
std::vector< BoundaryConditionCollection_boundary_conditions
 

Constructor & Destructor Documentation

◆ HTProcess()

ProcessLib::HT::HTProcess::HTProcess ( 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,
HTProcessData && process_data,
SecondaryVariableCollection && secondary_variables,
bool const use_monolithic_scheme,
std::unique_ptr< ProcessLib::SurfaceFluxData > && surfaceflux )

Definition at line 27 of file HTProcess.cpp.

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 _surfaceflux(std::move(surfaceflux))
44{
45}
std::unique_ptr< ProcessLib::SurfaceFluxData > _surfaceflux
Definition HTProcess.h:114
HTProcessData _process_data
Definition HTProcess.h:110
std::string const name
Definition Process.h:354
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:44

Member Function Documentation

◆ assembleConcreteProcess()

void ProcessLib::HT::HTProcess::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 )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 76 of file HTProcess.cpp.

80{
81 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
83 {
84 DBUG("Assemble HTProcess.");
85 dof_tables.emplace_back(_local_to_global_index_map.get());
86 }
87 else
88 {
90 {
91 DBUG(
92 "Assemble the equations of heat transport process within "
93 "HTProcess.");
94 }
95 else
96 {
97 DBUG(
98 "Assemble the equations of single phase fully saturated "
99 "fluid flow process within HTProcess.");
100 }
101 dof_tables.emplace_back(_local_to_global_index_map.get());
102 dof_tables.emplace_back(_local_to_global_index_map.get());
103 }
104
105 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
106 // Call global assembler for each local assembly item.
109 pv.getActiveElementIDs(), dof_tables, t, dt, x, x_prev, process_id, M,
110 K, b);
111}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
std::vector< std::unique_ptr< HTLocalAssemblerInterface > > _local_assemblers
Definition HTProcess.h:112
std::vector< std::size_t > const & getActiveElementIDs() const
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition Process.h:155
VectorMatrixAssembler _global_assembler
Definition Process.h:367
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:360
const bool _use_monolithic_scheme
Definition Process.h:369
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)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References ProcessLib::Process::_global_assembler, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, _process_data, ProcessLib::Process::_use_monolithic_scheme, ProcessLib::VectorMatrixAssembler::assemble(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), ProcessLib::Process::getProcessVariables(), and ProcessLib::HT::HTProcessData::heat_transport_process_id.

◆ assembleWithJacobianConcreteProcess()

void ProcessLib::HT::HTProcess::assembleWithJacobianConcreteProcess ( 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,
GlobalMatrix & Jac )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 113 of file HTProcess.cpp.

117{
118 DBUG("AssembleWithJacobian HTProcess.");
119
120 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
122 {
123 dof_tables.emplace_back(_local_to_global_index_map.get());
124 }
125 else
126 {
127 dof_tables.emplace_back(_local_to_global_index_map.get());
128 dof_tables.emplace_back(_local_to_global_index_map.get());
129 }
130
131 // Call global assembler for each local assembly item.
132 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
135 _local_assemblers, pv.getActiveElementIDs(), dof_tables, t, dt, x,
136 x_prev, process_id, M, K, b, Jac);
137}
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, 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, GlobalMatrix &Jac)

References ProcessLib::Process::_global_assembler, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, ProcessLib::Process::_use_monolithic_scheme, ProcessLib::VectorMatrixAssembler::assembleWithJacobian(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), and ProcessLib::Process::getProcessVariables().

◆ getDOFTableForExtrapolatorData()

std::tuple< NumLib::LocalToGlobalIndexMap *, bool > ProcessLib::HT::HTProcess::getDOFTableForExtrapolatorData ( ) const
overrideprivatevirtual

Get the address of a LocalToGlobalIndexMap, and the status of its memory. If the LocalToGlobalIndexMap is created as new in this function, the function also returns a true boolean value to let Extrapolator manage the memory by the address returned by this function.

Returns
Address of a LocalToGlobalIndexMap and its memory status.

Reimplemented from ProcessLib::Process.

Definition at line 140 of file HTProcess.cpp.

141{
143 {
144 // For single-variable-single-component processes reuse the existing DOF
145 // table.
146 const bool manage_storage = false;
147 return std::make_tuple(_local_to_global_index_map.get(),
148 manage_storage);
149 }
150
151 // Otherwise construct a new DOF table.
152 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
154
155 const bool manage_storage = true;
156 return std::make_tuple(new NumLib::LocalToGlobalIndexMap(
157 std::move(all_mesh_subsets_single_component),
158 // by location order is needed for output
160 manage_storage);
161}
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:358
@ BY_LOCATION
Ordering data by spatial location.

References ProcessLib::Process::_local_to_global_index_map, ProcessLib::Process::_mesh_subset_all_nodes, ProcessLib::Process::_use_monolithic_scheme, and NumLib::BY_LOCATION.

◆ getFlux()

Eigen::Vector3d ProcessLib::HT::HTProcess::getFlux ( std::size_t element_id,
MathLib::Point3d const & p,
double const t,
std::vector< GlobalVector * > const & x ) const
overridevirtual

Reimplemented from ProcessLib::Process.

Definition at line 163 of file HTProcess.cpp.

167{
168 // fetch local_x from primary variable
169 std::vector<GlobalIndexType> indices_cache;
170 auto const r_c_indices = NumLib::getRowColumnIndices(
171 element_id, *_local_to_global_index_map, indices_cache);
172 std::vector<std::vector<GlobalIndexType>> indices_of_all_coupled_processes{
173 x.size(), r_c_indices.rows};
174 auto const local_x =
175 getCoupledLocalSolutions(x, indices_of_all_coupled_processes);
176
177 return _local_assemblers[element_id]->getFlux(p, t, local_x);
178}
NumLib::LocalToGlobalIndexMap::RowColumnIndices getRowColumnIndices(std::size_t const id, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< GlobalIndexType > &indices)
std::vector< double > getCoupledLocalSolutions(std::vector< GlobalVector * > const &global_solutions, std::vector< std::vector< GlobalIndexType > > const &indices)

References _local_assemblers, ProcessLib::Process::_local_to_global_index_map, ProcessLib::getCoupledLocalSolutions(), and NumLib::getRowColumnIndices().

◆ initializeConcreteProcess()

void ProcessLib::HT::HTProcess::initializeConcreteProcess ( NumLib::LocalToGlobalIndexMap const & dof_table,
MeshLib::Mesh const & mesh,
unsigned const integration_order )
overrideprivatevirtual

Process specific initialization called by initialize().

Implements ProcessLib::Process.

Definition at line 47 of file HTProcess.cpp.

51{
52 int const mesh_space_dimension = _process_data.mesh_space_dimension;
53
55 {
56 ProcessLib::createLocalAssemblers<MonolithicHTFEM>(
57 mesh_space_dimension, mesh.getElements(), dof_table,
59 mesh.isAxiallySymmetric(), _process_data);
60 }
61 else
62 {
63 ProcessLib::createLocalAssemblers<StaggeredHTFEM>(
64 mesh_space_dimension, mesh.getElements(), dof_table,
66 mesh.isAxiallySymmetric(), _process_data);
67 }
68
70 "darcy_velocity",
71 makeExtrapolator(mesh_space_dimension, getExtrapolator(),
74}
virtual std::vector< double > const & getIntPtDarcyVelocity(const double, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const =0
SecondaryVariableCollection _secondary_variables
Definition Process.h:362
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:199
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)

References _local_assemblers, _process_data, ProcessLib::Process::_secondary_variables, ProcessLib::Process::_use_monolithic_scheme, ProcessLib::SecondaryVariableCollection::addSecondaryVariable(), MeshLib::Mesh::getElements(), ProcessLib::Process::getExtrapolator(), ProcessLib::HT::HTLocalAssemblerInterface::getIntPtDarcyVelocity(), MeshLib::Mesh::isAxiallySymmetric(), ProcessLib::makeExtrapolator(), and ProcessLib::HT::HTProcessData::mesh_space_dimension.

◆ isLinear()

bool ProcessLib::HT::HTProcess::isLinear ( ) const
inlineoverride

Definition at line 72 of file HTProcess.h.

72{ return false; }

◆ postTimestepConcreteProcess()

void ProcessLib::HT::HTProcess::postTimestepConcreteProcess ( std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
const double t,
const double delta_t,
int const process_id )
overridevirtual

Reimplemented from ProcessLib::Process.

Definition at line 181 of file HTProcess.cpp.

187{
188 // For the monolithic scheme, process_id is always zero.
189 if (_use_monolithic_scheme && process_id != 0)
190 {
191 OGS_FATAL(
192 "The condition of process_id = 0 must be satisfied for monolithic "
193 "HTProcess, which is a single process.");
194 }
197 {
198 DBUG("This is the thermal part of the staggered HTProcess.");
199 return;
200 }
201 if (!_surfaceflux) // computing the surfaceflux is optional
202 {
203 return;
204 }
205
206 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
207
208 _surfaceflux->integrate(x, t, *this, process_id, _integration_order, _mesh,
210}
#define OGS_FATAL(...)
Definition Error.h:26
MeshLib::Mesh & _mesh
Definition Process.h:357
unsigned const _integration_order
Definition Process.h:374

References ProcessLib::Process::_integration_order, ProcessLib::Process::_mesh, _process_data, _surfaceflux, ProcessLib::Process::_use_monolithic_scheme, DBUG(), ProcessLib::ProcessVariable::getActiveElementIDs(), ProcessLib::Process::getProcessVariables(), ProcessLib::HT::HTProcessData::hydraulic_process_id, and OGS_FATAL.

Member Data Documentation

◆ _local_assemblers

std::vector<std::unique_ptr<HTLocalAssemblerInterface> > ProcessLib::HT::HTProcess::_local_assemblers
private

◆ _process_data

HTProcessData ProcessLib::HT::HTProcess::_process_data
private

◆ _surfaceflux

std::unique_ptr<ProcessLib::SurfaceFluxData> ProcessLib::HT::HTProcess::_surfaceflux
private

Definition at line 114 of file HTProcess.h.

Referenced by postTimestepConcreteProcess().


The documentation for this class was generated from the following files: