OGS
TESLocalAssembler-impl.h
Go to the documentation of this file.
1
10#pragma once
11
19#include "TESLocalAssembler.h"
20#include "TESReactionAdaptor.h"
21
22namespace
23{
24enum class MatOutType
25{
26 OGS5,
27 PYTHON
28};
29
31
32// TODO move to some location in the OGS core.
33template <typename Mat>
34void ogs5OutMat(const Mat& mat)
35{
36 for (unsigned r = 0; r < mat.rows(); ++r)
37 {
39 {
41 if (r != 0)
42 {
43 std::printf("\n");
44 }
45 std::printf("|");
46 break;
48 if (r != 0)
49 {
50 std::printf(",\n");
51 }
52 std::printf("[");
53 break;
54 }
55
56 for (unsigned c = 0; c < mat.cols(); ++c)
57 {
59 {
61 std::printf(" %.16e", mat(r, c));
62 break;
64 if (c != 0)
65 {
66 std::printf(",");
67 }
68 std::printf(" %23.16g", mat(r, c));
69 break;
70 }
71 }
72
74 {
76 std::printf(" | ");
77 break;
79 std::printf(" ]");
80 break;
81 }
82 }
83 std::printf("\n");
84}
85
86template <typename Vec>
87void ogs5OutVec(const Vec& vec)
88{
89 for (unsigned r = 0; r < vec.size(); ++r)
90 {
92 {
94 if (r != 0)
95 {
96 std::printf("\n");
97 }
98 std::printf("| %.16e | ", vec[r]);
99 break;
101 if (r != 0)
102 {
103 std::printf(",\n");
104 }
105 std::printf("[ %23.16g ]", vec[r]);
106 break;
107 }
108 }
109 std::printf("\n");
110}
111
112} // anonymous namespace
113
114namespace ProcessLib
115{
116namespace TES
117{
118template <typename ShapeFunction_, int GlobalDim>
120 MeshLib::Element const& e,
121 std::size_t const /*local_matrix_size*/,
122 NumLib::GenericIntegrationMethod const& integration_method,
123 bool const is_axially_symmetric,
124 AssemblyParams const& asm_params)
125 : _element(e),
126 _integration_method(integration_method),
127 _shape_matrices(NumLib::initShapeMatrices<ShapeFunction,
128 ShapeMatricesType, GlobalDim>(
129 e, is_axially_symmetric, _integration_method)),
130 _d(asm_params, _integration_method.getNumberOfPoints(), GlobalDim)
131{
132}
133
134template <typename ShapeFunction_, int GlobalDim>
136 double const /*t*/, double const /*dt*/, std::vector<double> const& local_x,
137 std::vector<double> const& /*local_x_prev*/,
138 std::vector<double>& local_M_data, std::vector<double>& local_K_data,
139 std::vector<double>& local_b_data)
140{
141 auto const local_matrix_size = local_x.size();
142 // This assertion is valid only if all nodal d.o.f. use the same shape
143 // matrices.
144 assert(local_matrix_size == ShapeFunction::NPOINTS * NODAL_DOF);
145
147 local_M_data, local_matrix_size, local_matrix_size);
149 local_K_data, local_matrix_size, local_matrix_size);
151 local_b_data, local_matrix_size);
152
153 unsigned const n_integration_points =
154 _integration_method.getNumberOfPoints();
155
156 _d.preEachAssemble();
157
158 for (unsigned ip = 0; ip < n_integration_points; ip++)
159 {
160 auto const& sm = _shape_matrices[ip];
161 auto const& wp = _integration_method.getWeightedPoint(ip);
162 auto const weight = wp.getWeight();
163
164 _d.assembleIntegrationPoint(ip, local_x, sm, weight, local_M, local_K,
165 local_b);
166 }
167
168 if (_d.getAssemblyParameters().output_element_matrices)
169 {
170 std::puts("### Element: ?");
171
172 std::puts("---Velocity of water");
173 for (auto const& vs : _d.getData().velocity)
174 {
175 std::printf("| ");
176 for (auto v : vs)
177 {
178 std::printf("%23.16e ", v);
179 }
180 std::printf("|\n");
181 }
182
183 std::printf("\n---Mass matrix: \n");
184 ogs5OutMat(local_M);
185 std::printf("\n");
186
187 std::printf("---Laplacian + Advective + Content matrix: \n");
188 ogs5OutMat(local_K);
189 std::printf("\n");
190
191 std::printf("---RHS: \n");
192 ogs5OutVec(local_b);
193 std::printf("\n");
194 }
195}
196
197template <typename ShapeFunction_, int GlobalDim>
198std::vector<double> const&
200 const double /*t*/,
201 std::vector<GlobalVector*> const& /*x*/,
202 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
203 std::vector<double>& /*cache*/) const
204{
205 return _d.getData().solid_density;
206}
207
208template <typename ShapeFunction_, int GlobalDim>
209std::vector<double> const&
211 const double /*t*/,
212 std::vector<GlobalVector*> const& /*x*/,
213 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
214 std::vector<double>& cache) const
215{
216 auto const rho_SR = _d.getData().solid_density;
217 auto const rho_SR_dry = _d.getAssemblyParameters().rho_SR_dry;
218
219 cache.clear();
220 cache.reserve(rho_SR.size());
221
222 transform(cbegin(rho_SR), cend(rho_SR), back_inserter(cache),
223 [&](auto const rho) { return rho / rho_SR_dry - 1.0; });
224
225 return cache;
226}
227
228template <typename ShapeFunction_, int GlobalDim>
229std::vector<double> const&
231 const double /*t*/, std::vector<GlobalVector*> const& /*x*/,
232 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
233 std::vector<double>& cache) const
234{
235 auto const fac = _d.getData().reaction_adaptor->getReactionDampingFactor();
236 auto const num_integration_points = _d.getData().solid_density.size();
237
238 cache.clear();
239 cache.resize(num_integration_points, fac);
240
241 return cache;
242}
243
244template <typename ShapeFunction_, int GlobalDim>
245std::vector<double> const&
247 const double /*t*/,
248 std::vector<GlobalVector*> const& /*x*/,
249 std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
250 std::vector<double>& /*cache*/) const
251{
252 return _d.getData().reaction_rate;
253}
254
255template <typename ShapeFunction_, int GlobalDim>
256std::vector<double> const&
258 const double /*t*/,
259 std::vector<GlobalVector*> const& x,
260 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
261 std::vector<double>& cache) const
262{
263 auto const n_integration_points = _integration_method.getNumberOfPoints();
264
265 constexpr int process_id = 0; // monolithic scheme
266 auto const indices =
267 NumLib::getIndices(_element.getID(), *dof_table[process_id]);
268 assert(!indices.empty());
269 auto const local_x = x[process_id]->get(indices);
270 // local_x is ordered by component, local_x_mat is row major
271 auto const local_x_mat = MathLib::toMatrix<
272 Eigen::Matrix<double, NODAL_DOF, Eigen::Dynamic, Eigen::RowMajor>>(
273 local_x, NODAL_DOF, ShapeFunction_::NPOINTS);
274
275 cache.clear();
276 auto cache_mat = MathLib::createZeroedMatrix<
277 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
278 cache, GlobalDim, n_integration_points);
279
280 for (unsigned i = 0; i < n_integration_points; ++i)
281 {
282 double p, T, x;
283 NumLib::shapeFunctionInterpolate(local_x, _shape_matrices[i].N, p, T,
284 x);
285 const double eta_GR = fluid_viscosity(p, T, x);
286
287 auto const& k = _d.getAssemblyParameters().solid_perm_tensor;
288 cache_mat.col(i).noalias() =
289 k *
290 (_shape_matrices[i].dNdx *
291 local_x_mat.row(COMPONENT_ID_PRESSURE).transpose()) /
292 -eta_GR;
293 }
294
295 return cache;
296}
297
298template <typename ShapeFunction_, int GlobalDim>
300 std::vector<double> const& local_x,
301 std::vector<double> const& local_x_prev_ts)
302{
303 return _d.getReactionAdaptor().checkBounds(local_x, local_x_prev_ts);
304}
305
306} // namespace TES
307} // namespace ProcessLib
std::vector< double > const & getIntPtReactionDampingFactor(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 & getIntPtSolidDensity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
TESLocalAssembler(MeshLib::Element const &e, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, AssemblyParams const &asm_params)
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
std::vector< double > const & getIntPtLoading(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 & getIntPtReactionRate(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
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_b_data) override
bool checkBounds(std::vector< double > const &local_x, std::vector< double > const &local_x_prev_ts) override
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)
Eigen::Map< const Matrix > toMatrix(std::vector< double > const &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
const unsigned NODAL_DOF
double fluid_viscosity(const double p, const double T, const double x)
const unsigned COMPONENT_ID_PRESSURE