OGS
TwoPhaseFlowWithPPProcess.cpp
Go to the documentation of this file.
1 
12 
15 
16 namespace ProcessLib
17 {
18 namespace 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  std::map<std::string,
31  std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
32  /*curves*/)
33  : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
34  integration_order, std::move(process_variables),
35  std::move(secondary_variables)),
36  _process_data(std::move(process_data))
37 {
38  DBUG("Create TwoPhaseFlowProcess with PP model.");
39 }
40 
42  NumLib::LocalToGlobalIndexMap const& dof_table,
43  MeshLib::Mesh const& mesh,
44  unsigned const integration_order)
45 {
46  ProcessLib::createLocalAssemblers<TwoPhaseFlowWithPPLocalAssembler>(
47  mesh.getDimension(), mesh.getElements(), dof_table, _local_assemblers,
48  NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
50 
52  "saturation",
56 
58  "pressure_wet",
62 }
63 
65  const double t, double const dt, std::vector<GlobalVector*> const& x,
66  std::vector<GlobalVector*> const& xdot, int const process_id,
68 {
69  DBUG("Assemble TwoPhaseFlowWithPPProcess.");
70 
71  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
72  dof_table = {std::ref(*_local_to_global_index_map)};
73  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
74 
75  // Call global assembler for each local assembly item.
78  pv.getActiveElementIDs(), dof_table, t, dt, x, xdot, process_id, M, K,
79  b);
80 }
81 
83  const double t, double const dt, std::vector<GlobalVector*> const& x,
84  std::vector<GlobalVector*> const& xdot, int const process_id,
86 {
87  DBUG("AssembleWithJacobian TwoPhaseFlowWithPPProcess.");
88 
89  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
90  dof_table = {std::ref(*_local_to_global_index_map)};
91  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
92 
93  // Call global assembler for each local assembly item.
96  _local_assemblers, pv.getActiveElementIDs(), dof_table, t, dt, x, xdot,
97  process_id, M, K, b, Jac);
98 }
99 
100 } // namespace TwoPhaseFlowWithPP
101 } // namespace ProcessLib
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition: Logging.h:29
Global vector based on Eigen vector.
Definition: EigenVector.h:28
bool isAxiallySymmetric() const
Definition: Mesh.h:131
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition: Mesh.h:76
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition: Mesh.h:103
std::vector< std::size_t > const & getActiveElementIDs() const
NumLib::Extrapolator & getExtrapolator() const
Definition: Process.h:184
SecondaryVariableCollection _secondary_variables
Definition: Process.h:335
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition: Process.h:148
VectorMatrixAssembler _global_assembler
Definition: Process.h:337
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition: Process.h:333
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 assembleWithJacobianConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) override
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 &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) 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, std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation >> const &curves)
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, 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< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)
static const double t
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)