OGS
ThermalTwoPhaseFlowWithPPProcess.cpp
Go to the documentation of this file.
1 
12 
13 #include <cassert>
14 
18 #include "MeshLib/PropertyVector.h"
21 
22 namespace ProcessLib
23 {
24 namespace ThermalTwoPhaseFlowWithPP
25 {
27  std::string name,
28  MeshLib::Mesh& mesh,
29  std::unique_ptr<AbstractJacobianAssembler>&& jacobian_assembler,
30  std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
31  unsigned const integration_order,
32  std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
33  process_variables,
35  SecondaryVariableCollection&& secondary_variables,
36  BaseLib::ConfigTree const& /*config*/,
37  std::map<std::string,
38  std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
39  /*curves*/)
40  : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
41  integration_order, std::move(process_variables),
42  std::move(secondary_variables)),
43  _process_data(std::move(process_data))
44 {
45  DBUG("Create Nonisothermal TwoPhase Flow Process model.");
46 }
47 
49  NumLib::LocalToGlobalIndexMap const& dof_table,
50  MeshLib::Mesh const& mesh,
51  unsigned const integration_order)
52 {
53  const int monolithic_process_id = 0;
54  ProcessLib::ProcessVariable const& pv =
55  getProcessVariables(monolithic_process_id)[0];
56  ProcessLib::createLocalAssemblers<ThermalTwoPhaseFlowWithPPLocalAssembler>(
57  mesh.getDimension(), mesh.getElements(), dof_table,
59  mesh.isAxiallySymmetric(), integration_order, _process_data);
60 
62  "saturation",
65  getIntPtSaturation));
66 
68  "pressure_wetting",
71  getIntPtWettingPressure));
72 }
73 
75  const double t, double const dt, std::vector<GlobalVector*> const& x,
76  std::vector<GlobalVector*> const& xdot, int const process_id,
78 {
79  DBUG("Assemble ThermalTwoPhaseFlowWithPPProcess.");
80 
81  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
82  dof_table = {std::ref(*_local_to_global_index_map)};
83  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
84 
85  // Call global assembler for each local assembly item.
88  pv.getActiveElementIDs(), dof_table, t, dt, x, xdot, process_id, M, K,
89  b);
90 }
91 
93  const double t, double const dt, std::vector<GlobalVector*> const& x,
94  std::vector<GlobalVector*> const& xdot, const double dxdot_dx,
95  const double dx_dx, int const process_id, GlobalMatrix& M, GlobalMatrix& K,
96  GlobalVector& b, GlobalMatrix& Jac)
97 {
98  DBUG("AssembleWithJacobian ThermalTwoPhaseFlowWithPPProcess.");
99 
100  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
101  dof_table = {std::ref(*_local_to_global_index_map)};
102  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
103 
104  // Call global assembler for each local assembly item.
107  _local_assemblers, pv.getActiveElementIDs(), dof_table, t, dt, x, xdot,
108  dxdot_dx, dx_dx, process_id, M, K, b, Jac);
109 }
111  std::vector<GlobalVector*> const& x, double const t, double const delta_t,
112  const int process_id)
113 {
114  DBUG("PreTimestep ThermalTwoPhaseFlowWithPPProcess.");
115 
116  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
119  pv.getActiveElementIDs(), *_local_to_global_index_map, *x[process_id],
120  t, delta_t);
121 }
122 
123 } // namespace ThermalTwoPhaseFlowWithPP
124 } // namespace ProcessLib
void DBUG(char const *fmt, Args const &... args)
Definition: Logging.h:27
Definition of the PiecewiseLinearInterpolation class.
Global vector based on Eigen vector.
Definition: EigenVector.h:26
bool isAxiallySymmetric() const
Definition: Mesh.h:126
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition: Mesh.h:71
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition: Mesh.h:98
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)
unsigned getShapeFunctionOrder() const
std::vector< std::size_t > const & getActiveElementIDs() const
NumLib::Extrapolator & getExtrapolator() const
Definition: Process.h:185
SecondaryVariableCollection _secondary_variables
Definition: Process.h:331
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition: Process.h:145
VectorMatrixAssembler _global_assembler
Definition: Process.h:333
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition: Process.h:329
Handles configuration of several secondary variables from the project file.
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
std::vector< std::unique_ptr< ThermalTwoPhaseFlowWithPPLocalAssemblerInterface > > _local_assemblers
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, BaseLib::ConfigTree const &config, std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation >> const &curves)
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
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
void assembleWithJacobianConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, const double dxdot_dx, const double dx_dx, 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< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, const double dxdot_dx, const double dx_dx, 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)
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)