OGS
ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess Class Referencefinal

Detailed Description

Global assembler for the monolithic scheme of the non-isothermal Richards flow.

Governing equations without vapor diffusion

The energy balance equation is given by

\[ (\rho c_p)^{eff}\dot T - \nabla (\mathbf{k}_T^{eff} \nabla T)+\rho^l c_p^l \nabla T \cdot \mathbf{v}^l = Q_T \]

with \(T\) the temperature, \((\rho c_p)^{eff}\) the effective volumetric heat capacity, \(\mathbf{k}_T^{eff} \) the effective thermal conductivity, \(\rho^l\) the density of liquid, \(c_p^l\) the specific heat capacity of liquid, \(\mathbf{v}^l\) the liquid velocity, and \(Q_T\) the point heat source. The effective volumetric heat can be considered as a composite of the contributions of solid phase and the liquid phase as

\[(\rho c_p)^{eff} = (1-\phi) \rho^s c_p^s + S^l \phi \rho^l c_p^l \]

with \(\phi\) the porosity, \(S^l\) the liquid saturation, \(\rho^s \) the solid density, and \(c_p^s\) the specific heat capacity of solid. Similarly, the effective thermal conductivity is given by

\[\mathbf{k}_T^{eff} = (1-\phi) \mathbf{k}_T^s + S^l \phi k_T^l \mathbf I \]

where \(\mathbf{k}_T^s\) is the thermal conductivity tensor of solid, \( k_T^l\) is the thermal conductivity of liquid, and \(\mathbf I\) is the identity tensor.

The mass balance equation is given by

\begin{eqnarray*}\left(S^l\beta - \phi\frac{\partial S}{\partial p_c}\right) \rho^l\dot p - S \left( \frac{\partial \rho^l}{\partial T} +\rho^l(\alpha_B -S) \alpha_T^s \right)\dot T\\ +\nabla (\rho^l \mathbf{v}^l) + S \alpha_B \rho^l \nabla \cdot \dot {\mathbf u}= Q_H \end{eqnarray*}

where \(p\) is the pore pressure, \(p_c\) is the capillary pressure, which is \(-p\) under the single phase assumption, \(\beta\) is a composite coefficient by the liquid compressibility and solid compressibility, \(\alpha_B\) is the Biot's constant, \(\alpha_T^s\) is the linear thermal expansivity of solid, \(Q_H\) is the point source or sink term, \(H(S-1)\) is the Heaviside function, and \( \mathbf u\) is the displacement. While this process does not contain a fully mechanical coupling, simplified expressions can be given to approximate the latter term under certain stress conditions. The liquid velocity \(\mathbf{v}^l\) is described by the Darcy's law as

\[\mathbf{v}^l=-\frac{{\mathbf k} k_{ref}}{\mu} (\nabla p - \rho^l \mathbf g) \]

with \({\mathbf k}\) the intrinsic permeability, \(k_{ref}\) the relative permeability, \(\mathbf g\) the gravitational force.

Definition at line 75 of file ThermoRichardsFlowProcess.h.

#include <ThermoRichardsFlowProcess.h>

Inheritance diagram for ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess:
[legend]
Collaboration diagram for ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess:
[legend]

Public Member Functions

 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)
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, 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::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables () 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 Eigen::Vector3d getFlux (std::size_t, MathLib::Point3d const &, double const, std::vector< GlobalVector * > const &) const
virtual void solveReactionEquation (std::vector< GlobalVector * > &, std::vector< GlobalVector * > const &, double const, double const, NumLib::EquationSystem &, int const)
bool requiresNormalization () const override
Public Member Functions inherited from ProcessLib::SubmeshAssemblySupport
virtual ~SubmeshAssemblySupport ()=default

Private Types

using LocalAssemblerIF = LocalAssemblerInterface

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 setInitialConditionsConcreteProcess (std::vector< GlobalVector * > &x, double const t, int const) 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
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 preTimestepConcreteProcess (std::vector< GlobalVector * > const &, const double, const double, const int process_id) override
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 computeSecondaryVariableConcrete (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id) override
Private Member Functions inherited from ProcessLib::AssemblyMixin< ThermoRichardsFlowProcess >
void initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &submeshes, std::vector< std::vector< std::string > > const &residuum_names)
void updateActiveElements ()
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 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)
void preOutput (double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)

Private Attributes

std::vector< MeshLib::Node * > _base_nodes
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_base_nodes
ThermoRichardsFlowProcessData _process_data
std::vector< std::unique_ptr< LocalAssemblerIF > > local_assemblers_
GlobalSparsityPattern _sparsity_pattern_with_linear_element

Friends

class AssemblyMixin< ThermoRichardsFlowProcess >

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)
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const
std::vector< GlobalIndexTypegetIndicesOfResiduumWithoutInitialCompensation () const override
void setReleaseNodalForces (GlobalVector const *r_neq, int const process_id) 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
CellAverageData cell_average_data_
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

Member Typedef Documentation

◆ LocalAssemblerIF

Constructor & Destructor Documentation

◆ ThermoRichardsFlowProcess()

ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::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 )

Definition at line 21 of file ThermoRichardsFlowProcess.cpp.

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,
48
49 _integration_point_writer.emplace_back(
50 std::make_unique<MeshLib::IntegrationPointWriter>(
51 "porosity_ip", 1 /*n components*/, integration_order,
53}
std::string const name
Definition Process.h:361
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Definition Process.h:389
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > _jacobian_assembler
Definition Process.h:375
std::vector< std::unique_ptr< LocalAssemblerIF > > local_assemblers_
virtual std::vector< double > getSaturation() const =0
virtual std::vector< double > getPorosity() const =0

References ProcessLib::Process::Process(), ThermoRichardsFlowProcess(), ProcessLib::Process::_jacobian_assembler, and ProcessLib::Process::name.

Referenced by ThermoRichardsFlowProcess(), and AssemblyMixin< ThermoRichardsFlowProcess >.

Member Function Documentation

◆ assembleConcreteProcess()

void ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::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 141 of file ThermoRichardsFlowProcess.cpp.

145{
146 DBUG("Assemble the equations for ThermoRichardsFlowProcess.");
147
149 process_id, M, K, b);
150}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
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)

References ProcessLib::AssemblyMixin< Process >::assemble(), and DBUG().

◆ assembleWithJacobianConcreteProcess()

void ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::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 )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 152 of file ThermoRichardsFlowProcess.cpp.

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}
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)

References ProcessLib::AssemblyMixin< Process >::assembleWithJacobian(), and DBUG().

◆ computeSecondaryVariableConcrete()

void ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::computeSecondaryVariableConcrete ( double const t,
double const dt,
std::vector< GlobalVector * > const & x,
GlobalVector const & x_prev,
int const process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 195 of file ThermoRichardsFlowProcess.cpp.

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}
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)
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:376
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:160
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References ProcessLib::LocalAssemblerInterface::computeSecondaryVariable(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), ProcessLib::Process::getDOFTables(), and local_assemblers_.

◆ initializeAssemblyOnSubmeshes()

std::vector< std::vector< std::string > > ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::initializeAssemblyOnSubmeshes ( std::vector< std::reference_wrapper< MeshLib::Mesh > > const & meshes)
overrideprivatevirtual

Initializes the assembly on submeshes

Parameters
meshesthe submeshes on whom the assembly shall proceed.
Attention
meshes must be a must be a non-overlapping cover of the entire simulation domain (bulk mesh)!
Returns
The names of the residuum vectors that will be assembled for each process: outer vector of size 1 for monolithic schemes and greater for staggered schemes.

Reimplemented from ProcessLib::SubmeshAssemblySupport.

Definition at line 119 of file ThermoRichardsFlowProcess.cpp.

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}
void initializeAssemblyOnSubmeshes(std::vector< std::reference_wrapper< MeshLib::Mesh > > const &submeshes, std::vector< std::vector< std::string > > const &residuum_names)
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
Definition Process.h:399

References ProcessLib::Process::_process_variables, DBUG(), and ProcessLib::AssemblyMixin< Process >::initializeAssemblyOnSubmeshes().

◆ initializeConcreteProcess()

void ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::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 55 of file ThermoRichardsFlowProcess.cpp.

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}
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
SecondaryVariableCollection _secondary_variables
Definition Process.h:369
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:367
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:201
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 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

References ProcessLib::Process::_integration_point_writer, ProcessLib::Process::_local_to_global_index_map, _process_data, ProcessLib::Process::_secondary_variables, MeshLib::Cell, ProcessLib::createLocalAssemblers(), NumLib::SerialExecutor::executeMemberOnDereferenced(), MeshLib::Mesh::getDimension(), MeshLib::Mesh::getElements(), ProcessLib::Process::getExtrapolator(), ProcessLib::ThermoRichardsFlow::LocalAssemblerInterface::getIntPtDarcyVelocity(), ProcessLib::ThermoRichardsFlow::LocalAssemblerInterface::getIntPtDryDensitySolid(), ProcessLib::ThermoRichardsFlow::LocalAssemblerInterface::getIntPtPorosity(), ProcessLib::ThermoRichardsFlow::LocalAssemblerInterface::getIntPtSaturation(), MeshLib::getOrCreateMeshProperty(), MeshLib::Mesh::getProperties(), ProcessLib::LocalAssemblerInterface::initialize(), MeshLib::Mesh::isAxiallySymmetric(), local_assemblers_, ProcessLib::makeExtrapolator(), ProcessLib::Process::name, and ProcessLib::setIPDataInitialConditions().

◆ isLinear()

bool ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::isLinear ( ) const
inlineoverride

◆ postTimestepConcreteProcess()

void ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::postTimestepConcreteProcess ( std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
double const t,
double const dt,
const int process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 177 of file ThermoRichardsFlowProcess.cpp.

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}
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)

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), ProcessLib::Process::getDOFTables(), local_assemblers_, and ProcessLib::LocalAssemblerInterface::postTimestep().

◆ preTimestepConcreteProcess()

void ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::preTimestepConcreteProcess ( std::vector< GlobalVector * > const & ,
const double ,
const double ,
const int process_id )
overrideprivatevirtual

◆ setInitialConditionsConcreteProcess()

void ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::setInitialConditionsConcreteProcess ( std::vector< GlobalVector * > & x,
double const t,
int const process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 104 of file ThermoRichardsFlowProcess.cpp.

106{
107 if (process_id != 0)
108 {
109 return;
110 }
111 DBUG("SetInitialConditions ThermoRichardsFlowProcess.");
112
115 getDOFTables(x.size()), x, t, process_id);
116}
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)

References DBUG(), NumLib::SerialExecutor::executeMemberOnDereferenced(), ProcessLib::Process::getDOFTables(), local_assemblers_, and ProcessLib::LocalAssemblerInterface::setInitialConditions().

◆ AssemblyMixin< ThermoRichardsFlowProcess >

Member Data Documentation

◆ _base_nodes

std::vector<MeshLib::Node*> ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::_base_nodes
private

Definition at line 144 of file ThermoRichardsFlowProcess.h.

◆ _mesh_subset_base_nodes

std::unique_ptr<MeshLib::MeshSubset const> ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::_mesh_subset_base_nodes
private

Definition at line 145 of file ThermoRichardsFlowProcess.h.

◆ _process_data

ThermoRichardsFlowProcessData ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::_process_data
private

Definition at line 146 of file ThermoRichardsFlowProcess.h.

Referenced by initializeConcreteProcess().

◆ _sparsity_pattern_with_linear_element

GlobalSparsityPattern ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::_sparsity_pattern_with_linear_element
private

Sparsity pattern for the flow equation, and it is initialized only if the staggered scheme is used.

Definition at line 152 of file ThermoRichardsFlowProcess.h.

◆ local_assemblers_

std::vector<std::unique_ptr<LocalAssemblerIF> > ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::local_assemblers_
private

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