OGS
SteadyStateDiffusionFEM.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include <vector>
7
23
24namespace ProcessLib
25{
27{
28const unsigned NUM_NODAL_DOF = 1;
29
33{
34public:
35 virtual std::vector<double> const& getIntPtDarcyVelocity(
36 const double t,
37 std::vector<GlobalVector*> const& x,
38 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
39 std::vector<double>& cache) const = 0;
40};
41
42template <typename ShapeFunction, int GlobalDim>
44{
47
49 ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim>;
50
54
55public:
59 MeshLib::Element const& element,
60 std::size_t const /*local_matrix_size*/,
61 NumLib::GenericIntegrationMethod const& integration_method,
62 bool const is_axially_symmetric,
63 SteadyStateDiffusionData const& process_data)
64 : _element(element),
65 _process_data(process_data),
66 _integration_method(integration_method),
68 NumLib::initShapeMatrices<ShapeFunction, ShapeMatricesType,
69 GlobalDim>(
70 element, is_axially_symmetric, _integration_method))
71 {
72 }
73
74 void assemble(double const t, double const dt,
75 std::vector<double> const& local_x,
76 std::vector<double> const& /*local_x_prev*/,
77 std::vector<double>& /*local_M_data*/,
78 std::vector<double>& local_K_data,
79 std::vector<double>& /*local_b_data*/) override
80 {
81 auto const local_matrix_size = local_x.size();
82 // This assertion is valid only if all nodal d.o.f. use the same shape
83 // matrices.
84 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
85
87 local_K_data, local_matrix_size, local_matrix_size);
88
89 unsigned const n_integration_points =
90 _integration_method.getNumberOfPoints();
91
93 pos.setElementID(_element.getID());
94
95 auto const& medium =
96 *_process_data.media_map.getMedium(_element.getID());
98 vars.temperature =
99 medium
100 .property(
102 .template value<double>(vars, pos, t, dt);
103
104 for (unsigned ip = 0; ip < n_integration_points; ip++)
105 {
106 auto const& sm = _shape_matrices[ip];
107 auto const& wp = _integration_method.getWeightedPoint(ip);
108
109 pos = {std::nullopt, _element.getID(),
113 _element, sm.N))};
114
115 double p_int_pt = 0.0;
116 NumLib::shapeFunctionInterpolate(local_x, sm.N, p_int_pt);
117 vars.liquid_phase_pressure = p_int_pt;
120 .value(vars, pos, t, dt));
121
122 local_K.noalias() += sm.dNdx.transpose() * k * sm.dNdx * sm.detJ *
123 sm.integralMeasure * wp.getWeight();
124 }
125 }
126
129 Eigen::Vector3d getFlux(MathLib::Point3d const& p_local_coords,
130 double const t,
131 std::vector<double> const& local_x) const override
132 {
133 // TODO (tf) Temporary value not used by current material models. Need
134 // extension of getFlux interface
135 double const dt = std::numeric_limits<double>::quiet_NaN();
136
137 // Eval shape matrices at given point
138 // Note: Axial symmetry is set to false here, because we only need dNdx
139 // here, which is not affected by axial symmetry.
140 auto const shape_matrices =
142 GlobalDim>(
143 _element, false /*is_axially_symmetric*/,
144 std::array{p_local_coords})[0];
145
146 // fetch hydraulic conductivity
148 pos.setElementID(_element.getID());
149 auto const& medium =
150 *_process_data.media_map.getMedium(_element.getID());
151
153 vars.temperature =
154 medium
155 .property(
157 .template value<double>(vars, pos, t, dt);
158 double pressure = 0.0;
159 NumLib::shapeFunctionInterpolate(local_x, shape_matrices.N, pressure);
160 vars.liquid_phase_pressure = pressure;
161
164 .value(vars, pos, t, dt));
165
166 Eigen::Vector3d flux(0.0, 0.0, 0.0);
167 flux.head<GlobalDim>() =
168 -k * shape_matrices.dNdx *
169 Eigen::Map<const NodalVectorType>(local_x.data(), local_x.size());
170
171 return flux;
172 }
173
174 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
175 const unsigned integration_point) const override
176 {
177 auto const& N = _shape_matrices[integration_point].N;
178
179 // assumes N is stored contiguously in memory
180 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
181 }
182
183 std::vector<double> const& getIntPtDarcyVelocity(
184 const double t,
185 std::vector<GlobalVector*> const& x,
186 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
187 std::vector<double>& cache) const override
188 {
189 // TODO (tf) Temporary value not used by current material models. Need
190 // extension of secondary variable interface.
191 double const dt = std::numeric_limits<double>::quiet_NaN();
192
193 auto const n_integration_points =
194 _integration_method.getNumberOfPoints();
195
196 int const process_id = 0; // monolithic scheme
197 auto const indices =
198 NumLib::getIndices(_element.getID(), *dof_table[process_id]);
199 assert(!indices.empty());
200 auto const local_x = x[process_id]->get(indices);
201 auto const local_x_vec =
203 local_x, ShapeFunction::NPOINTS);
204
205 cache.clear();
206 auto cache_mat = MathLib::createZeroedMatrix<
207 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
208 cache, GlobalDim, n_integration_points);
209
211 pos.setElementID(_element.getID());
212
213 auto const& medium =
214 *_process_data.media_map.getMedium(_element.getID());
215
217
218 double pressure = 0.0;
219 for (unsigned i = 0; i < n_integration_points; ++i)
220 {
222 pressure);
223 vars.liquid_phase_pressure = pressure;
224
225 pos = {std::nullopt, _element.getID(),
229 _element, _shape_matrices[i].N))};
230
233 .value(vars, pos, t, dt));
234 // dimensions: (d x 1) = (d x n) * (n x 1)
235 cache_mat.col(i).noalias() =
236 -k * _shape_matrices[i].dNdx * local_x_vec;
237 }
238
239 return cache;
240 }
241
242private:
245
247 std::vector<ShapeMatrices, Eigen::aligned_allocator<ShapeMatrices>>
249};
250
251} // namespace SteadyStateDiffusion
252} // namespace ProcessLib
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
void setElementID(std::size_t element_id)
ProcessLib::LocalAssemblerTraits< ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim > LocalAssemblerTraits
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 > &) override
LocalAssemblerData(MeshLib::Element const &element, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, SteadyStateDiffusionData const &process_data)
typename LocalAssemblerTraits::LocalMatrix NodalMatrixType
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
Eigen::Vector3d getFlux(MathLib::Point3d const &p_local_coords, double const t, std::vector< double > const &local_x) const override
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
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 LocalAssemblerTraits::LocalVector NodalVectorType
NumLib::GenericIntegrationMethod const & _integration_method
virtual 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 =0
constexpr Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
Eigen::Map< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &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)
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > computeShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, PointContainer const &points)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
detail::LocalAssemblerTraitsFixed< ShpPol, NNodes, NodalDOF, Dim > LocalAssemblerTraits
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
VectorType< GlobalDim > GlobalDimVectorType
Matrix< NNodes *NodalDOF, NNodes *NodalDOF > LocalMatrix