OGS
TwoPhaseFlowWithPPProcess.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
5
8
9namespace ProcessLib
10{
11namespace TwoPhaseFlowWithPP
12{
14 std::string name,
15 MeshLib::Mesh& mesh,
16 std::unique_ptr<AbstractJacobianAssembler>&& jacobian_assembler,
17 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
18 unsigned const integration_order,
19 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
20 process_variables,
21 TwoPhaseFlowWithPPProcessData&& process_data,
22 SecondaryVariableCollection&& secondary_variables)
23 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
24 integration_order, std::move(process_variables),
25 std::move(secondary_variables)),
26 _process_data(std::move(process_data))
27{
28 DBUG("Create TwoPhaseFlowProcess with PP model.");
29
30 // For numerical Jacobian
31 this->_jacobian_assembler->setNonDeformationComponentIDs(
32 {0, 1} /* P_g, P_c */);
33}
34
36 std::vector<GlobalVector*>& x, double const t, int const process_id)
37{
38 DBUG("SetInitialConditions ThermoRichardsMechanicsProcess.");
39
42 _local_assemblers, getDOFTables(x.size()), x, t, process_id);
43}
44
67
69 const double t, double const dt, std::vector<GlobalVector*> const& x,
70 std::vector<GlobalVector*> const& x_prev, int const process_id,
72{
73 DBUG("Assemble TwoPhaseFlowWithPPProcess.");
74
75 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
77
78 // Call global assembler for each local assembly item.
81 getActiveElementIDs(), dof_table, t, dt, x, x_prev, process_id, &M, &K,
82 &b);
83}
84
86 const double t, double const dt, std::vector<GlobalVector*> const& x,
87 std::vector<GlobalVector*> const& x_prev, int const process_id,
89{
90 DBUG("AssembleWithJacobian TwoPhaseFlowWithPPProcess.");
91
92 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
94
95 // Call global assembler for each local assembly item.
98 _local_assemblers, getActiveElementIDs(), dof_table, t, dt, x, x_prev,
99 process_id, &b, &Jac);
100}
101
102} // namespace TwoPhaseFlowWithPP
103} // namespace ProcessLib
MathLib::EigenMatrix GlobalMatrix
MathLib::EigenVector GlobalVector
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
bool isAxiallySymmetric() const
Definition Mesh.h:128
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:100
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition Mesh.h:79
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)
std::string const name
Definition Process.h:361
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:376
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:37
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:160
SecondaryVariableCollection _secondary_variables
Definition Process.h:369
VectorMatrixAssembler _global_assembler
Definition Process.h:376
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:367
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > _jacobian_assembler
Definition Process.h:375
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:201
Handles configuration of several secondary variables from the project file.
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 & getIntPtWetPressure(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
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
std::vector< std::unique_ptr< TwoPhaseFlowWithPPLocalAssemblerInterface > > _local_assemblers
TwoPhaseFlowWithPPProcess(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, TwoPhaseFlowWithPPProcessData &&process_data, SecondaryVariableCollection &&secondary_variables)
void setInitialConditionsConcreteProcess(std::vector< GlobalVector * > &x, double const t, int const process_id) 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 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)
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)
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)
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)