OGS
ThermoRichardsFlowFEM.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <memory>
14#include <vector>
15
28
29namespace ProcessLib
30{
31namespace ThermoRichardsFlow
32{
33namespace MPL = MaterialPropertyLib;
34
35template <typename ShapeFunction, int GlobalDim>
37{
38public:
39 // Note: temperature variable uses the same shape functions as that are used
40 // by pressure variable.
42
45
47
49 delete;
51 delete;
52
54 MeshLib::Element const& e,
55 std::size_t const /*local_matrix_size*/,
56 NumLib::GenericIntegrationMethod const& integration_method,
57 bool const is_axially_symmetric,
58 ThermoRichardsFlowProcessData& process_data);
59
62 std::string_view const name,
63 double const* values,
64 int const integration_order) override;
65
66 void setInitialConditionsConcrete(Eigen::VectorXd const local_x,
67 double const t,
68 int const process_id) override;
69
70 void assembleWithJacobian(double const t, double const dt,
71 std::vector<double> const& local_x,
72 std::vector<double> const& local_x_prev,
73 std::vector<double>& /*local_M_data*/,
74 std::vector<double>& /*local_K_data*/,
75 std::vector<double>& local_rhs_data,
76 std::vector<double>& local_Jac_data) override;
77
78 void assemble(double const t, double const dt,
79 std::vector<double> const& local_x,
80 std::vector<double> const& local_x_prev,
81 std::vector<double>& local_M_data,
82 std::vector<double>& local_K_data,
83 std::vector<double>& local_rhs_data) override;
84
85 void initializeConcrete() override
86 {
87 unsigned const n_integration_points =
89
90 for (unsigned ip = 0; ip < n_integration_points; ip++)
91 {
92 auto& ip_data = _ip_data[ip];
93 ip_data.pushBackState();
94 }
95 }
96
97 void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
98 Eigen::VectorXd const& /*local_x_prev*/,
99 double const /*t*/, double const /*dt*/,
100 int const /*process_id*/) override
101 {
102 unsigned const n_integration_points =
104
105 for (unsigned ip = 0; ip < n_integration_points; ip++)
106 {
107 _ip_data[ip].pushBackState();
108 }
109 }
110
112 double const t, double const dt, Eigen::VectorXd const& local_x,
113 Eigen::VectorXd const& local_x_prev) override;
114
115 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
116 const unsigned integration_point) const override
117 {
118 auto const& N = _ip_data[integration_point].N;
119
120 // assumes N is stored contiguously in memory
121 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
122 }
123
124 std::vector<double> const& getIntPtDarcyVelocity(
125 const double t,
126 std::vector<GlobalVector*> const& x,
127 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
128 std::vector<double>& cache) const override;
129
130 std::vector<double> getSaturation() const override;
131 std::vector<double> const& getIntPtSaturation(
132 const double t,
133 std::vector<GlobalVector*> const& x,
134 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
135 std::vector<double>& cache) const override;
136
137 std::vector<double> getPorosity() const override;
138 std::vector<double> const& getIntPtPorosity(
139 const double t,
140 std::vector<GlobalVector*> const& x,
141 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
142 std::vector<double>& cache) const override;
143
144 std::vector<double> const& getIntPtDryDensitySolid(
145 const double t,
146 std::vector<GlobalVector*> const& x,
147 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
148 std::vector<double>& cache) const override;
149
150private:
151 unsigned getNumberOfIntegrationPoints() const override;
152
153private:
155
156 std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
157
161
162 static const int temperature_index = 0;
163 static const int temperature_size = ShapeFunction::NPOINTS;
165 static const int pressure_size = ShapeFunction::NPOINTS;
166};
167
168} // namespace ThermoRichardsFlow
169} // namespace ProcessLib
170
void postTimestepConcrete(Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const) override
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &, std::vector< double > &, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) override
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 override
std::size_t setIPDataInitialConditions(std::string_view const name, double const *values, int const integration_order) override
ThermoRichardsFlowLocalAssembler(ThermoRichardsFlowLocalAssembler &&)=delete
std::vector< double > const & getIntPtDryDensitySolid(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_rhs_data) override
ThermoRichardsFlowLocalAssembler(ThermoRichardsFlowLocalAssembler const &)=delete
std::vector< double > const & getIntPtPorosity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
void setInitialConditionsConcrete(Eigen::VectorXd const local_x, double const t, int const process_id) override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
void computeSecondaryVariableConcrete(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev) override
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
VectorType< GlobalDim > GlobalDimVectorType