OGS
StokesFlowFEM.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <Eigen/LU>
14#include <vector>
15
28
30{
31template <typename ShapeFunctionLiquidVelocity, typename ShapeFunctionPressure,
32 int GlobalDim>
34{
35 static const int liquid_velocity_index = 0;
36 static const int pressure_index =
37 ShapeFunctionLiquidVelocity::NPOINTS * GlobalDim;
38
39 static const int liquid_velocity_size =
40 ShapeFunctionLiquidVelocity::NPOINTS * GlobalDim;
41 static const int pressure_size = ShapeFunctionPressure::NPOINTS;
42
47
49
50public:
52 MeshLib::Element const& element,
53 std::size_t const /*local_matrix_size*/,
54 NumLib::GenericIntegrationMethod const& integration_method,
55 bool const is_axially_symmetric,
56 StokesFlowProcessData const& process_data)
57 : _element(element),
58 _is_axially_symmetric(is_axially_symmetric),
59 _integration_method(integration_method),
60 _process_data(process_data)
61 {
62 unsigned const n_integration_points =
64 _ip_data.reserve(n_integration_points);
65
66 auto const shape_matrices_v =
67 NumLib::initShapeMatrices<ShapeFunctionLiquidVelocity,
69 element, is_axially_symmetric, _integration_method);
70
71 auto const shape_matrices_p =
72 NumLib::initShapeMatrices<ShapeFunctionPressure,
73 ShapeMatrixTypePressure, GlobalDim>(
74 element, is_axially_symmetric, _integration_method);
75
76 for (unsigned ip = 0; ip < n_integration_points; ip++)
77 {
78 _ip_data.emplace_back(
79 shape_matrices_v[ip].N, shape_matrices_v[ip].dNdx,
80 shape_matrices_p[ip].N, shape_matrices_p[ip].dNdx,
82 shape_matrices_v[ip].integralMeasure *
83 shape_matrices_v[ip].detJ);
84
85 auto& ip_data = _ip_data[ip];
86
87 ip_data.N_v_op = ShapeMatrixTypeLiquidVelocity::template MatrixType<
88 GlobalDim, liquid_velocity_size>::Zero(GlobalDim,
90
91 ip_data.dNdx_v_op =
92 ShapeMatrixTypeLiquidVelocity::template MatrixType<
93 GlobalDim * GlobalDim,
94 liquid_velocity_size>::Zero(GlobalDim * GlobalDim,
96
97 for (int i = 0; i < GlobalDim; ++i)
98 {
99 ip_data.N_v_op
100 .template block<1, liquid_velocity_size / GlobalDim>(
101 i, i * liquid_velocity_size / GlobalDim)
102 .noalias() = shape_matrices_v[ip].N;
103
104 ip_data.dNdx_v_op
105 .template block<GlobalDim,
106 liquid_velocity_size / GlobalDim>(
107 i * GlobalDim, i * liquid_velocity_size / GlobalDim)
108 .noalias() = shape_matrices_v[ip].dNdx;
109 }
110 }
111 }
112
113 void assemble(double const t, double const dt,
114 std::vector<double> const& local_x,
115 std::vector<double> const& /*local_x_prev*/,
116 std::vector<double>& /*local_M_data*/,
117 std::vector<double>& local_K_data,
118 std::vector<double>& local_b_data) override
119 {
120 auto const local_matrix_size = liquid_velocity_size + pressure_size;
121 assert(local_x.size() == local_matrix_size);
122
123 auto local_p = Eigen::Map<const NodalVectorType>(
124 &local_x[pressure_index], pressure_size);
125
126 auto local_K = MathLib::createZeroedMatrix<
127 typename ShapeMatrixTypeLiquidVelocity::template MatrixType<
128 local_matrix_size, local_matrix_size>>(
129 local_K_data, local_matrix_size, local_matrix_size);
130 auto local_b = MathLib::createZeroedVector<
131 typename ShapeMatrixTypeLiquidVelocity::template VectorType<
132 local_matrix_size>>(local_b_data, local_matrix_size);
133
134 // Get block matrices
135 // Kvv
136 auto Kvv =
137 local_K.template block<liquid_velocity_size, liquid_velocity_size>(
139 // Kvp
140 auto Kvp = local_K.template block<liquid_velocity_size, pressure_size>(
142 // kpv
143 auto Kpv = local_K.template block<pressure_size, liquid_velocity_size>(
145 // bv
146 auto bv = local_b.template segment<liquid_velocity_size>(
148
149 unsigned const n_integration_points =
151
154
155 auto const& b = _process_data.specific_body_force;
156
158
159 // create unit vector
160 assert(GlobalDim == 2);
161 auto const identity_mat =
162 Eigen::Matrix<double, GlobalDim, GlobalDim>::Identity().eval();
163 auto const unit_vector = Eigen::Map<const Eigen::RowVectorXd>(
164 identity_mat.data(), identity_mat.size())
165 .eval();
166
167 // Get material properties
168 auto const& medium =
170 auto const& phase = medium.phase("AqueousLiquid");
171
172 for (unsigned ip(0); ip < n_integration_points; ++ip)
173 {
174 pos.setIntegrationPoint(ip);
175
176 auto const& ip_data = _ip_data[ip];
177
178 auto const& N_p = ip_data.N_p;
179 auto const& N_v_op = ip_data.N_v_op;
180
181 auto const& dNdx_v = ip_data.dNdx_v;
182 auto const& dNdx_v_op = ip_data.dNdx_v_op;
183 auto const& dNdx_p = ip_data.dNdx_p;
184
185 auto const& w = ip_data.integration_weight;
186
187 double p_int_pt = 0.0;
188
189 NumLib::shapeFunctionInterpolate(local_p, N_p, p_int_pt);
190
191 vars.liquid_phase_pressure = p_int_pt;
192
193 // Use the viscosity model to compute the viscosity
195 .template value<double>(vars, pos, t, dt);
196
197 // Kvv
199 {
200 // permeability
203 .value(vars, pos, t, dt));
204
205 Kvv.noalias() +=
206 w * N_v_op.transpose() * mu * K.inverse() * N_v_op;
207 }
208
209 for (int i = 0; i < GlobalDim; ++i)
210 {
211 Kvv.template block<liquid_velocity_size / GlobalDim,
212 liquid_velocity_size / GlobalDim>(
213 i * liquid_velocity_size / GlobalDim,
214 i * liquid_velocity_size / GlobalDim)
215 .noalias() += w * dNdx_v.transpose() * mu * dNdx_v;
216 }
217
218 // Kvp
219 Kvp.noalias() += w * N_v_op.transpose() * dNdx_p;
220
221 // Kpv
222 Kpv.noalias() += w * N_p.transpose() * unit_vector * dNdx_v_op;
223
224 // bv
225 bv.noalias() += w * N_v_op.transpose() * b;
226 }
227 }
228
229 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
230 const unsigned integration_point) const override
231 {
232 auto const& N_p = _ip_data[integration_point].N_p;
233
234 // assumes N_p is stored contiguously in memory
235 return Eigen::Map<const Eigen::RowVectorXd>(N_p.data(), N_p.size());
236 }
237
239 double const /*t*/,
240 double const /*dt*/,
241 Eigen::VectorXd const& local_x,
242 Eigen::VectorXd const& /*local_x_prev*/) override
243 {
244 auto const local_p =
245 local_x.template segment<pressure_size>(pressure_index);
246
248 ShapeFunctionPressure,
249 typename ShapeFunctionLiquidVelocity::MeshElement, GlobalDim>(
252 }
253
254private:
259
261 ShapeMatrixTypePressure, GlobalDim,
262 ShapeFunctionLiquidVelocity::NPOINTS>,
263 Eigen::aligned_allocator<IntegrationPointData<
265 GlobalDim, ShapeFunctionLiquidVelocity::NPOINTS>>>
267};
268
269} // namespace ProcessLib::StokesFlow
Phase const & phase(std::size_t index) const
Definition Medium.cpp:33
double getWeight() const
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
void computeSecondaryVariableConcrete(double const, double const, Eigen::VectorXd const &local_x, Eigen::VectorXd const &) override
ShapeMatrixPolicyType< ShapeFunctionLiquidVelocity, GlobalDim > ShapeMatrixTypeLiquidVelocity
void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &, std::vector< double > &, std::vector< double > &local_K_data, std::vector< double > &local_b_data) override
StokesFlowProcessData const & _process_data
LocalAssemblerData(MeshLib::Element const &element, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, StokesFlowProcessData const &process_data)
NumLib::GenericIntegrationMethod const & _integration_method
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
typename ShapeMatrixTypePressure::NodalVectorType NodalVectorType
std::vector< IntegrationPointData< ShapeMatrixTypeLiquidVelocity, ShapeMatrixTypePressure, GlobalDim, ShapeFunctionLiquidVelocity::NPOINTS >, Eigen::aligned_allocator< IntegrationPointData< ShapeMatrixTypeLiquidVelocity, ShapeMatrixTypePressure, GlobalDim, ShapeFunctionLiquidVelocity::NPOINTS > > > _ip_data
ShapeMatrixPolicyType< ShapeFunctionPressure, GlobalDim > ShapeMatrixTypePressure
Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
void interpolateToHigherOrderNodes(MeshLib::Element const &element, bool const is_axially_symmetric, Eigen::MatrixBase< EigenMatrixType > const &node_values, MeshLib::PropertyVector< double > &interpolated_values_global_vector)
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
VectorType< ShapeFunction::NPOINTS > NodalVectorType
Eigen::VectorXd const specific_body_force
an external force that applies in the bulk of the fluid, like gravity.
MaterialPropertyLib::MaterialSpatialDistributionMap media_map
MeshLib::PropertyVector< double > * pressure_interpolated