OGS
ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess Class Referencefinal

Detailed Description

A class to simulate the nonisothermal two-phase flow process in porous media with max. 3 components. The water and air components are usually required on a minimum. The air component is assumed to only exist in the gas phase. The third component (organic contaminant) is optional, and can partition between the liquid and gas phases. Technically, an algorithm for phase appearance/disappearance is not implemented here. However, the liquid phase is allowed to appear/disappear in some exceptional cases, e.g. the heat pipe benchmark (since a pure liquid phase is assumed). For the formulation of global conservation equations and constitutive relationships, see [26].

Definition at line 43 of file ThermalTwoPhaseFlowWithPPProcess.h.

#include <ThermalTwoPhaseFlowWithPPProcess.h>

Inheritance diagram for ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess:
[legend]
Collaboration diagram for ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess:
[legend]

Public Member Functions

 ThermalTwoPhaseFlowWithPPProcess (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, ThermalTwoPhaseFlowWithPPProcessData &&process_data, SecondaryVariableCollection &&secondary_variables)
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 std::vector< 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, GlobalVector &b, GlobalMatrix &Jac) override
void preTimestepConcreteProcess (std::vector< GlobalVector * > const &x, const double t, const double delta_t, const int process_id) override

Private Attributes

ThermalTwoPhaseFlowWithPPProcessData _process_data
std::vector< std::unique_ptr< ThermalTwoPhaseFlowWithPPLocalAssemblerInterface > > _local_assemblers
Assembly::GlobalMatrixOutput _global_output

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

Constructor & Destructor Documentation

◆ ThermalTwoPhaseFlowWithPPProcess()

ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::ThermalTwoPhaseFlowWithPPProcess ( 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,
ThermalTwoPhaseFlowWithPPProcessData && process_data,
SecondaryVariableCollection && secondary_variables )

Definition at line 24 of file ThermalTwoPhaseFlowWithPPProcess.cpp.

34 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
35 integration_order, std::move(process_variables),
36 std::move(secondary_variables)),
37 _process_data(std::move(process_data))
38{
39 DBUG("Create Nonisothermal TwoPhase Flow Process model.");
40}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
std::string const name
Definition Process.h:368
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

References ProcessLib::Process::Process(), _process_data, DBUG(), and ProcessLib::Process::name.

Member Function Documentation

◆ assembleConcreteProcess()

void ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::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 83 of file ThermalTwoPhaseFlowWithPPProcess.cpp.

87{
88 DBUG("Assemble ThermalTwoPhaseFlowWithPPProcess.");
89
90 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
92
93 // Call global assembler for each local assembly item.
96 getActiveElementIDs(), dof_table, t, dt, x, x_prev, process_id, &M, &K,
97 &b);
98
99 _global_output(t, process_id, M, K, b);
100}
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:167
VectorMatrixAssembler _global_assembler
Definition Process.h:383
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:374
std::vector< std::unique_ptr< ThermalTwoPhaseFlowWithPPLocalAssemblerInterface > > _local_assemblers
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, _global_output, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, ProcessLib::VectorMatrixAssembler::assemble(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), and ProcessLib::Process::getActiveElementIDs().

◆ assembleWithJacobianConcreteProcess()

void ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::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 102 of file ThermalTwoPhaseFlowWithPPProcess.cpp.

106{
107 DBUG("AssembleWithJacobian ThermalTwoPhaseFlowWithPPProcess.");
108
109 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
111
112 // Call global assembler for each local assembly item.
115 _local_assemblers, getActiveElementIDs(), dof_table, t, dt, x, x_prev,
116 process_id, &b, &Jac);
117
118 _global_output(t, process_id, b, Jac);
119}
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, GlobalVector *b, GlobalMatrix *Jac)

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

◆ initializeConcreteProcess()

void ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::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 42 of file ThermalTwoPhaseFlowWithPPProcess.cpp.

46{
48 mesh.getDimension(), mesh.getElements(), dof_table, _local_assemblers,
49 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
51
52 _secondary_variables.addSecondaryVariable(
53 "saturation",
55 &ThermalTwoPhaseFlowWithPPLocalAssemblerInterface::
56 getIntPtSaturation));
57
58 _secondary_variables.addSecondaryVariable(
59 "pressure_wetting",
61 &ThermalTwoPhaseFlowWithPPLocalAssemblerInterface::
62 getIntPtWettingPressure));
63
64 _secondary_variables.addSecondaryVariable(
65 "liquid_molar_fraction_contaminant",
67 &ThermalTwoPhaseFlowWithPPLocalAssemblerInterface::
68 getIntPtLiquidMolFracContaminant));
69
70 _secondary_variables.addSecondaryVariable(
71 "gas_molar_fraction_water",
73 &ThermalTwoPhaseFlowWithPPLocalAssemblerInterface::
74 getIntPtGasMolFracWater));
75
76 _secondary_variables.addSecondaryVariable(
77 "gas_molar_fraction_contaminant",
79 &ThermalTwoPhaseFlowWithPPLocalAssemblerInterface::
80 getIntPtGasMolFracContaminant));
81}
SecondaryVariableCollection _secondary_variables
Definition Process.h:376
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:208
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)

References _local_assemblers, _process_data, ProcessLib::Process::_secondary_variables, ProcessLib::createLocalAssemblers(), MeshLib::Mesh::getDimension(), MeshLib::Mesh::getElements(), ProcessLib::Process::getExtrapolator(), MeshLib::Mesh::isAxiallySymmetric(), and ProcessLib::makeExtrapolator().

◆ isLinear()

bool ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::isLinear ( ) const
inlineoverride

Definition at line 58 of file ThermalTwoPhaseFlowWithPPProcess.h.

58{ return false; }

◆ preTimestepConcreteProcess()

void ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::preTimestepConcreteProcess ( std::vector< GlobalVector * > const & x,
const double t,
const double delta_t,
const int process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 121 of file ThermalTwoPhaseFlowWithPPProcess.cpp.

124{
125 DBUG("PreTimestep ThermalTwoPhaseFlowWithPPProcess.");
126
129 getActiveElementIDs(), *_local_to_global_index_map, *x[process_id], t,
130 delta_t);
131}
virtual void preTimestep(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References _local_assemblers, ProcessLib::Process::_local_to_global_index_map, DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), and ProcessLib::LocalAssemblerInterface::preTimestep().

Member Data Documentation

◆ _global_output

Assembly::GlobalMatrixOutput ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::_global_output
private

◆ _local_assemblers

std::vector< std::unique_ptr<ThermalTwoPhaseFlowWithPPLocalAssemblerInterface> > ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::_local_assemblers
private

◆ _process_data

ThermalTwoPhaseFlowWithPPProcessData ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::_process_data
private

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