OGS
ThermalTwoPhaseFlowWithPPProcess.cpp
Go to the documentation of this file.
1
12
13#include <cassert>
14
19
20namespace ProcessLib
21{
22namespace ThermalTwoPhaseFlowWithPP
23{
25 std::string name,
26 MeshLib::Mesh& mesh,
27 std::unique_ptr<AbstractJacobianAssembler>&& jacobian_assembler,
28 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
29 unsigned const integration_order,
30 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
31 process_variables,
33 SecondaryVariableCollection&& secondary_variables)
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}
41
43 NumLib::LocalToGlobalIndexMap const& dof_table,
44 MeshLib::Mesh const& mesh,
45 unsigned const integration_order)
46{
47 ProcessLib::createLocalAssemblers<ThermalTwoPhaseFlowWithPPLocalAssembler>(
48 mesh.getDimension(), mesh.getElements(), dof_table, _local_assemblers,
49 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
51
53 "saturation",
56 getIntPtSaturation));
57
59 "pressure_wetting",
63
65 "liquid_molar_fraction_contaminant",
69
71 "gas_molar_fraction_water",
75
77 "gas_molar_fraction_contaminant",
81}
82
84 const double t, double const dt, std::vector<GlobalVector*> const& x,
85 std::vector<GlobalVector*> const& x_prev, int const process_id,
87{
88 DBUG("Assemble ThermalTwoPhaseFlowWithPPProcess.");
89
90 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
92 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
93
94 // Call global assembler for each local assembly item.
97 pv.getActiveElementIDs(), dof_table, t, dt, x, x_prev, process_id, M, K,
98 b);
99}
100
102 const double t, double const dt, std::vector<GlobalVector*> const& x,
103 std::vector<GlobalVector*> const& x_prev, int const process_id,
105{
106 DBUG("AssembleWithJacobian ThermalTwoPhaseFlowWithPPProcess.");
107
108 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
110 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
111
112 // Call global assembler for each local assembly item.
115 _local_assemblers, pv.getActiveElementIDs(), dof_table, t, dt, x,
116 x_prev, process_id, M, K, b, Jac);
117}
119 std::vector<GlobalVector*> const& x, double const t, double const delta_t,
120 const int process_id)
121{
122 DBUG("PreTimestep ThermalTwoPhaseFlowWithPPProcess.");
123
124 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
127 pv.getActiveElementIDs(), *_local_to_global_index_map, *x[process_id],
128 t, delta_t);
129}
130
131} // namespace ThermalTwoPhaseFlowWithPP
132} // namespace ProcessLib
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
Definition of the PiecewiseLinearInterpolation class.
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
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)
std::vector< std::size_t > const & getActiveElementIDs() const
SecondaryVariableCollection _secondary_variables
Definition Process.h:362
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
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:199
Handles configuration of several secondary variables from the project file.
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
virtual std::vector< double > const & getIntPtGasMolFracContaminant(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 & getIntPtLiquidMolFracContaminant(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 & getIntPtWettingPressure(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 & getIntPtGasMolFracWater(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
std::vector< std::unique_ptr< ThermalTwoPhaseFlowWithPPLocalAssemblerInterface > > _local_assemblers
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
void preTimestepConcreteProcess(std::vector< GlobalVector * > const &x, const double t, const double delta_t, const int process_id) override
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)
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
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)
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)
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)