OGS
TwoPhaseFlowWithPPProcess.cpp
Go to the documentation of this file.
1
12
15
16namespace ProcessLib
17{
18namespace TwoPhaseFlowWithPP
19{
21 std::string name,
22 MeshLib::Mesh& mesh,
23 std::unique_ptr<AbstractJacobianAssembler>&& jacobian_assembler,
24 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
25 unsigned const integration_order,
26 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
27 process_variables,
28 TwoPhaseFlowWithPPProcessData&& process_data,
29 SecondaryVariableCollection&& secondary_variables)
30 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
31 integration_order, std::move(process_variables),
32 std::move(secondary_variables)),
33 _process_data(std::move(process_data))
34{
35 DBUG("Create TwoPhaseFlowProcess with PP model.");
36}
37
60
62 const double t, double const dt, std::vector<GlobalVector*> const& x,
63 std::vector<GlobalVector*> const& x_prev, int const process_id,
65{
66 DBUG("Assemble TwoPhaseFlowWithPPProcess.");
67
68 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
70
71 // Call global assembler for each local assembly item.
74 getActiveElementIDs(), dof_table, t, dt, x, x_prev, process_id, &M, &K,
75 &b);
76}
77
79 const double t, double const dt, std::vector<GlobalVector*> const& x,
80 std::vector<GlobalVector*> const& x_prev, int const process_id,
82{
83 DBUG("AssembleWithJacobian TwoPhaseFlowWithPPProcess.");
84
85 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
87
88 // Call global assembler for each local assembly item.
91 _local_assemblers, getActiveElementIDs(), dof_table, t, dt, x, x_prev,
92 process_id, &b, &Jac);
93}
94
95} // namespace TwoPhaseFlowWithPP
96} // namespace ProcessLib
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
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
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:167
SecondaryVariableCollection _secondary_variables
Definition Process.h:370
VectorMatrixAssembler _global_assembler
Definition Process.h:377
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:368
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:208
Handles configuration of several secondary variables from the project file.
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
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 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)