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_rhs_data,
74 std::vector<double>& local_Jac_data) override;
75
76 void assemble(double const t, double const dt,
77 std::vector<double> const& local_x,
78 std::vector<double> const& local_x_prev,
79 std::vector<double>& local_M_data,
80 std::vector<double>& local_K_data,
81 std::vector<double>& local_rhs_data) override;
82
83 void initializeConcrete() override
84 {
85 unsigned const n_integration_points =
87
88 for (unsigned ip = 0; ip < n_integration_points; ip++)
89 {
90 auto& ip_data = _ip_data[ip];
91 ip_data.pushBackState();
92 }
93 }
94
95 void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
96 Eigen::VectorXd const& /*local_x_prev*/,
97 double const /*t*/, double const /*dt*/,
98 int const /*process_id*/) override
99 {
100 unsigned const n_integration_points =
102
103 for (unsigned ip = 0; ip < n_integration_points; ip++)
104 {
105 _ip_data[ip].pushBackState();
106 }
107 }
108
110 double const t, double const dt, Eigen::VectorXd const& local_x,
111 Eigen::VectorXd const& local_x_prev) override;
112
113 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
114 const unsigned integration_point) const override
115 {
116 auto const& N = _ip_data[integration_point].N;
117
118 // assumes N is stored contiguously in memory
119 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
120 }
121
122 std::vector<double> const& getIntPtDarcyVelocity(
123 const double t,
124 std::vector<GlobalVector*> const& x,
125 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
126 std::vector<double>& cache) const override;
127
128 std::vector<double> getSaturation() const override;
129 std::vector<double> const& getIntPtSaturation(
130 const double t,
131 std::vector<GlobalVector*> const& x,
132 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
133 std::vector<double>& cache) const override;
134
135 std::vector<double> getPorosity() const override;
136 std::vector<double> const& getIntPtPorosity(
137 const double t,
138 std::vector<GlobalVector*> const& x,
139 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
140 std::vector<double>& cache) const override;
141
142 std::vector<double> const& getIntPtDryDensitySolid(
143 const double t,
144 std::vector<GlobalVector*> const& x,
145 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
146 std::vector<double>& cache) const override;
147
148private:
149 unsigned getNumberOfIntegrationPoints() const override;
150
151private:
153
154 std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
155
159
160 static const int temperature_index = 0;
161 static const int temperature_size = ShapeFunction::NPOINTS;
163 static const int pressure_size = ShapeFunction::NPOINTS;
164};
165
166} // namespace ThermoRichardsFlow
167} // namespace ProcessLib
168
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 > &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