OGS
HTFEM.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 <Eigen/Core>
7#include <vector>
8
10#include "HTProcessData.h"
21
22namespace ProcessLib
23{
24namespace HT
25{
26template <typename ShapeFunction, int GlobalDim>
28{
31
34
39
40public:
41 HTFEM(MeshLib::Element const& element,
42 std::size_t const local_matrix_size,
43 NumLib::GenericIntegrationMethod const& integration_method,
44 bool const is_axially_symmetric,
45 HTProcessData const& process_data,
46 const unsigned dof_per_node)
48 _element(element),
49 _process_data(process_data),
50 _integration_method(integration_method)
51 {
52 // This assertion is valid only if all nodal d.o.f. use the same shape
53 // matrices.
54 assert(local_matrix_size == ShapeFunction::NPOINTS * dof_per_node);
55 (void)local_matrix_size;
56 (void)dof_per_node;
57
58 unsigned const n_integration_points =
59 _integration_method.getNumberOfPoints();
60 _ip_data.reserve(n_integration_points);
61
63 pos.setElementID(_element.getID());
64
65 double const aperture_size = _process_data.aperture_size(0.0, pos)[0];
66
67 auto const shape_matrices =
69 GlobalDim>(element, is_axially_symmetric,
71
72 for (unsigned ip = 0; ip < n_integration_points; ip++)
73 {
74 _ip_data.emplace_back(
75 shape_matrices[ip].dNdx,
76 _integration_method.getWeightedPoint(ip).getWeight() *
77 shape_matrices[ip].integralMeasure *
78 shape_matrices[ip].detJ * aperture_size);
79 }
80 }
81
82 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
83 const unsigned integration_point) const override
84 {
85 auto const& N = _process_data.shape_matrix_cache.NsHigherOrder<
86 typename ShapeFunction::MeshElement>()[integration_point];
87
88 // assumes N is stored contiguously in memory
89 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
90 }
91
94 Eigen::Vector3d getFlux(MathLib::Point3d const& pnt_local_coords,
95 double const t,
96 std::vector<double> const& local_x) const override
97 {
98 // Eval shape matrices at given point
99 // Note: Axial symmetry is set to false here, because we only need dNdx
100 // here, which is not affected by axial symmetry.
101 auto const shape_matrices =
103 GlobalDim>(
104 _element, false /*is_axially_symmetric*/,
105 std::array{pnt_local_coords})[0];
106
108 pos.setElementID(this->_element.getID());
109
111
112 // local_x contains the local temperature and pressure values
113 double T_int_pt = 0.0;
114 double p_int_pt = 0.0;
115 NumLib::shapeFunctionInterpolate(local_x, shape_matrices.N, T_int_pt,
116 p_int_pt);
117
118 vars.temperature = T_int_pt;
119 vars.liquid_phase_pressure = p_int_pt;
120
121 auto const& medium =
122 *_process_data.media_map.getMedium(_element.getID());
123 auto const& liquid_phase =
125
126 // TODO (naumov) Temporary value not used by current material models.
127 // Need extension of secondary variables interface.
128 double const dt = std::numeric_limits<double>::quiet_NaN();
129 // fetch permeability, viscosity, density
132 .value(vars, pos, t, dt));
133
134 auto const mu =
136 .template value<double>(vars, pos, t, dt);
137 GlobalDimMatrixType const K_over_mu = K / mu;
138
139 auto const p_nodal_values = Eigen::Map<const NodalVectorType>(
140 &local_x[local_x.size() / 2], ShapeFunction::NPOINTS);
142 -K_over_mu * shape_matrices.dNdx * p_nodal_values;
143
144 if (this->_process_data.has_gravity)
145 {
146 auto const rho_w =
147 liquid_phase
149 .template value<double>(vars, pos, t, dt);
150 auto const b =
151 this->_process_data.projected_specific_body_force_vectors
152 [this->_element.getID()];
153 q += K_over_mu * rho_w * b;
154 }
155
156 Eigen::Vector3d flux;
157 flux.head<GlobalDim>() = q;
158 return flux;
159 }
160
161protected:
164
166 std::vector<IntegrationPointData<GlobalDimNodalMatrixType>> _ip_data;
167
169 MaterialPropertyLib::VariableArray const& vars, const double porosity,
170 const double fluid_density, const double specific_heat_capacity_fluid,
171 ParameterLib::SpatialPosition const& pos, double const t,
172 double const dt)
173 {
174 auto const& medium =
175 *_process_data.media_map.getMedium(this->_element.getID());
176 auto const& solid_phase =
178
179 auto const specific_heat_capacity_solid =
180 solid_phase
181 .property(
183 .template value<double>(vars, pos, t, dt);
184
185 auto const solid_density =
187 .template value<double>(vars, pos, t, dt);
188
189 return solid_density * specific_heat_capacity_solid * (1 - porosity) +
190 fluid_density * specific_heat_capacity_fluid * porosity;
191 }
192
195 const double fluid_density, const double specific_heat_capacity_fluid,
196 const GlobalDimVectorType& velocity,
197 ParameterLib::SpatialPosition const& pos, double const t,
198 double const dt)
199 {
200 auto const& medium =
201 *_process_data.media_map.getMedium(_element.getID());
202
203 auto thermal_conductivity =
205 medium
206 .property(
208 .value(vars, pos, t, dt));
209
210 auto const thermal_dispersivity_transversal =
211 medium
213 thermal_transversal_dispersivity)
214 .template value<double>();
215
216 auto const thermal_dispersivity_longitudinal =
217 medium
219 thermal_longitudinal_dispersivity)
220 .template value<double>();
221
222 // Thermal conductivity is moved outside and zero matrix is passed
223 // instead due to multiplication with fluid's density times specific
224 // heat capacity.
225 return thermal_conductivity +
226 fluid_density * specific_heat_capacity_fluid *
228 _process_data.stabilizer, _element.getID(),
229 GlobalDimMatrixType::Zero(GlobalDim, GlobalDim),
230 velocity, 0 /* phi */, thermal_dispersivity_transversal,
231 thermal_dispersivity_longitudinal);
232 }
233
234 std::vector<double> const& getIntPtDarcyVelocityLocal(
235 const double t, std::vector<double> const& local_x,
236 std::vector<double>& cache) const
237 {
238 std::vector<double> local_p{
239 local_x.data() + pressure_index,
240 local_x.data() + pressure_index + pressure_size};
241 std::vector<double> local_T{
242 local_x.data() + temperature_index,
243 local_x.data() + temperature_index + temperature_size};
244
245 auto const n_integration_points =
246 _integration_method.getNumberOfPoints();
247
248 cache.clear();
249 auto cache_mat = MathLib::createZeroedMatrix<
250 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
251 cache, GlobalDim, n_integration_points);
252
254 pos.setElementID(_element.getID());
255
257
258 auto const p_nodal_values = Eigen::Map<const NodalVectorType>(
259 &local_p[0], ShapeFunction::NPOINTS);
260
261 auto const& medium =
262 *_process_data.media_map.getMedium(_element.getID());
263 auto const& liquid_phase =
265
266 auto const& Ns =
267 _process_data.shape_matrix_cache
268 .NsHigherOrder<typename ShapeFunction::MeshElement>();
269
270 for (unsigned ip = 0; ip < n_integration_points; ++ip)
271 {
272 auto const& ip_data = _ip_data[ip];
273 auto const& dNdx = ip_data.dNdx;
274 auto const& N = Ns[ip];
275
276 double T_int_pt = 0.0;
277 double p_int_pt = 0.0;
278 NumLib::shapeFunctionInterpolate(local_p, N, p_int_pt);
279 NumLib::shapeFunctionInterpolate(local_T, N, T_int_pt);
280
281 vars.temperature = T_int_pt;
282 vars.liquid_phase_pressure = p_int_pt;
283
284 // TODO (naumov) Temporary value not used by current material
285 // models. Need extension of secondary variables interface.
286 double const dt = std::numeric_limits<double>::quiet_NaN();
289 .value(vars, pos, t, dt));
290
291 auto const mu =
292 liquid_phase
294 .template value<double>(vars, pos, t, dt);
295 GlobalDimMatrixType const K_over_mu = K / mu;
296
297 cache_mat.col(ip).noalias() = -K_over_mu * dNdx * p_nodal_values;
298
299 if (_process_data.has_gravity)
300 {
301 auto const rho_w =
302 liquid_phase
304 .template value<double>(vars, pos, t, dt);
305 auto const b =
306 _process_data.projected_specific_body_force_vectors
307 [_element.getID()];
308 // here it is assumed that the vector b is directed 'downwards'
309 cache_mat.col(ip).noalias() += K_over_mu * rho_w * b;
310 }
311 }
312
313 return cache;
314 }
315
316protected:
317 static const int pressure_index = ShapeFunction::NPOINTS;
318 static const int pressure_size = ShapeFunction::NPOINTS;
319 static const int temperature_index = 0;
320 static const int temperature_size = ShapeFunction::NPOINTS;
321};
322
323} // namespace HT
324} // namespace ProcessLib
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:80
void setElementID(std::size_t element_id)
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
Definition HTFEM.h:38
std::vector< double > const & getIntPtDarcyVelocityLocal(const double t, std::vector< double > const &local_x, std::vector< double > &cache) const
Definition HTFEM.h:234
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
Definition HTFEM.h:30
NumLib::GenericIntegrationMethod const & _integration_method
Definition HTFEM.h:165
typename ShapeMatricesType::GlobalDimNodalMatrixType GlobalDimNodalMatrixType
Definition HTFEM.h:36
static const int temperature_index
Definition HTFEM.h:319
typename ShapeMatricesType::NodalVectorType NodalVectorType
Definition HTFEM.h:32
static const int temperature_size
Definition HTFEM.h:320
double getHeatEnergyCoefficient(MaterialPropertyLib::VariableArray const &vars, const double porosity, const double fluid_density, const double specific_heat_capacity_fluid, ParameterLib::SpatialPosition const &pos, double const t, double const dt)
Definition HTFEM.h:168
HTFEM(MeshLib::Element const &element, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, HTProcessData const &process_data, const unsigned dof_per_node)
Definition HTFEM.h:41
static const int pressure_size
Definition HTFEM.h:318
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
Definition HTFEM.h:29
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
Definition HTFEM.h:35
GlobalDimMatrixType getThermalConductivityDispersivity(MaterialPropertyLib::VariableArray const &vars, const double fluid_density, const double specific_heat_capacity_fluid, const GlobalDimVectorType &velocity, ParameterLib::SpatialPosition const &pos, double const t, double const dt)
Definition HTFEM.h:193
Eigen::Vector3d getFlux(MathLib::Point3d const &pnt_local_coords, double const t, std::vector< double > const &local_x) const override
Definition HTFEM.h:94
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
Definition HTFEM.h:82
MeshLib::Element const & _element
Definition HTFEM.h:162
HTProcessData const & _process_data
Definition HTFEM.h:163
std::vector< IntegrationPointData< GlobalDimNodalMatrixType > > _ip_data
Definition HTFEM.h:166
static const int pressure_index
Definition HTFEM.h:317
typename ShapeMatricesType::NodalRowVectorType NodalRowVectorType
Definition HTFEM.h:33
constexpr Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(MaterialPropertyLib::PropertyDataType const &values)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
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)
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > computeShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, PointContainer const &points)
Eigen::MatrixXd computeHydrodynamicDispersion(NumericalStabilization const &stabilizer, std::size_t const element_id, Eigen::MatrixXd const &pore_diffusion_coefficient, Eigen::VectorXd const &velocity, double const porosity, double const solute_dispersivity_transverse, double const solute_dispersivity_longitudinal)
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< GlobalDim, ShapeFunction::NPOINTS > GlobalDimNodalMatrixType
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
VectorType< GlobalDim > GlobalDimVectorType
VectorType< ShapeFunction::NPOINTS > NodalVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType