OGS
HTFEM.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <Eigen/Core>
14#include <vector>
15
17#include "HTProcessData.h"
28
29namespace ProcessLib
30{
31namespace HT
32{
33template <typename ShapeFunction, int GlobalDim>
35{
38
41
46
47public:
48 HTFEM(MeshLib::Element const& element,
49 std::size_t const local_matrix_size,
50 NumLib::GenericIntegrationMethod const& integration_method,
51 bool const is_axially_symmetric,
52 HTProcessData const& process_data,
53 const unsigned dof_per_node)
55 _element(element),
56 _process_data(process_data),
57 _integration_method(integration_method)
58 {
59 // This assertion is valid only if all nodal d.o.f. use the same shape
60 // matrices.
61 assert(local_matrix_size == ShapeFunction::NPOINTS * dof_per_node);
62 (void)local_matrix_size;
63 (void)dof_per_node;
64
65 unsigned const n_integration_points =
67 _ip_data.reserve(n_integration_points);
68
71
72 double const aperture_size = _process_data.aperture_size(0.0, pos)[0];
73
74 auto const shape_matrices =
76 GlobalDim>(element, is_axially_symmetric,
78
79 for (unsigned ip = 0; ip < n_integration_points; ip++)
80 {
81 _ip_data.emplace_back(
82 shape_matrices[ip].dNdx,
84 shape_matrices[ip].integralMeasure *
85 shape_matrices[ip].detJ * aperture_size);
86 }
87 }
88
89 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
90 const unsigned integration_point) const override
91 {
93 typename ShapeFunction::MeshElement>()[integration_point];
94
95 // assumes N is stored contiguously in memory
96 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
97 }
98
101 Eigen::Vector3d getFlux(MathLib::Point3d const& pnt_local_coords,
102 double const t,
103 std::vector<double> const& local_x) const override
104 {
105 // Eval shape matrices at given point
106 // Note: Axial symmetry is set to false here, because we only need dNdx
107 // here, which is not affected by axial symmetry.
108 auto const shape_matrices =
110 GlobalDim>(
111 _element, false /*is_axially_symmetric*/,
112 std::array{pnt_local_coords})[0];
113
115 pos.setElementID(this->_element.getID());
116
118
119 // local_x contains the local temperature and pressure values
120 double T_int_pt = 0.0;
121 double p_int_pt = 0.0;
122 NumLib::shapeFunctionInterpolate(local_x, shape_matrices.N, T_int_pt,
123 p_int_pt);
124
125 vars.temperature = T_int_pt;
126 vars.liquid_phase_pressure = p_int_pt;
127
128 auto const& medium =
130 auto const& liquid_phase = medium.phase("AqueousLiquid");
131
132 // TODO (naumov) Temporary value not used by current material models.
133 // Need extension of secondary variables interface.
134 double const dt = std::numeric_limits<double>::quiet_NaN();
135 // fetch permeability, viscosity, density
136 auto const K = MaterialPropertyLib::formEigenTensor<GlobalDim>(
138 .value(vars, pos, t, dt));
139
140 auto const mu =
142 .template value<double>(vars, pos, t, dt);
143 GlobalDimMatrixType const K_over_mu = K / mu;
144
145 auto const p_nodal_values = Eigen::Map<const NodalVectorType>(
146 &local_x[local_x.size() / 2], ShapeFunction::NPOINTS);
148 -K_over_mu * shape_matrices.dNdx * p_nodal_values;
149
150 if (this->_process_data.has_gravity)
151 {
152 auto const rho_w =
153 liquid_phase
155 .template value<double>(vars, pos, t, dt);
156 auto const b =
158 [this->_element.getID()];
159 q += K_over_mu * rho_w * b;
160 }
161
162 Eigen::Vector3d flux;
163 flux.head<GlobalDim>() = q;
164 return flux;
165 }
166
167protected:
170
172 std::vector<IntegrationPointData<GlobalDimNodalMatrixType>> _ip_data;
173
175 MaterialPropertyLib::VariableArray const& vars, const double porosity,
176 const double fluid_density, const double specific_heat_capacity_fluid,
177 ParameterLib::SpatialPosition const& pos, double const t,
178 double const dt)
179 {
180 auto const& medium =
181 *_process_data.media_map.getMedium(this->_element.getID());
182 auto const& solid_phase = medium.phase("Solid");
183
184 auto const specific_heat_capacity_solid =
185 solid_phase
186 .property(
188 .template value<double>(vars, pos, t, dt);
189
190 auto const solid_density =
192 .template value<double>(vars, pos, t, dt);
193
194 return solid_density * specific_heat_capacity_solid * (1 - porosity) +
195 fluid_density * specific_heat_capacity_fluid * porosity;
196 }
197
200 const double fluid_density, const double specific_heat_capacity_fluid,
201 const GlobalDimVectorType& velocity,
202 ParameterLib::SpatialPosition const& pos, double const t,
203 double const dt)
204 {
205 auto const& medium =
207
208 auto thermal_conductivity =
209 MaterialPropertyLib::formEigenTensor<GlobalDim>(
210 medium
211 .property(
213 .value(vars, pos, t, dt));
214
215 auto const thermal_dispersivity_transversal =
216 medium
218 thermal_transversal_dispersivity)
219 .template value<double>();
220
221 auto const thermal_dispersivity_longitudinal =
222 medium
224 thermal_longitudinal_dispersivity)
225 .template value<double>();
226
227 // Thermal conductivity is moved outside and zero matrix is passed
228 // instead due to multiplication with fluid's density times specific
229 // heat capacity.
230 return thermal_conductivity +
231 fluid_density * specific_heat_capacity_fluid *
234 GlobalDimMatrixType::Zero(GlobalDim, GlobalDim),
235 velocity, 0 /* phi */, thermal_dispersivity_transversal,
236 thermal_dispersivity_longitudinal);
237 }
238
239 std::vector<double> const& getIntPtDarcyVelocityLocal(
240 const double t, std::vector<double> const& local_x,
241 std::vector<double>& cache) const
242 {
243 std::vector<double> local_p{
244 local_x.data() + pressure_index,
245 local_x.data() + pressure_index + pressure_size};
246 std::vector<double> local_T{
247 local_x.data() + temperature_index,
248 local_x.data() + temperature_index + temperature_size};
249
250 auto const n_integration_points =
252
253 cache.clear();
254 auto cache_mat = MathLib::createZeroedMatrix<
255 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
256 cache, GlobalDim, n_integration_points);
257
260
262
263 auto const p_nodal_values = Eigen::Map<const NodalVectorType>(
264 &local_p[0], ShapeFunction::NPOINTS);
265
266 auto const& medium =
268 auto const& liquid_phase = medium.phase("AqueousLiquid");
269
270 auto const& Ns =
272 .NsHigherOrder<typename ShapeFunction::MeshElement>();
273
274 for (unsigned ip = 0; ip < n_integration_points; ++ip)
275 {
276 auto const& ip_data = _ip_data[ip];
277 auto const& dNdx = ip_data.dNdx;
278 auto const& N = Ns[ip];
279
280 pos.setIntegrationPoint(ip);
281
282 double T_int_pt = 0.0;
283 double p_int_pt = 0.0;
284 NumLib::shapeFunctionInterpolate(local_p, N, p_int_pt);
285 NumLib::shapeFunctionInterpolate(local_T, N, T_int_pt);
286
287 vars.temperature = T_int_pt;
288 vars.liquid_phase_pressure = p_int_pt;
289
290 // TODO (naumov) Temporary value not used by current material
291 // models. Need extension of secondary variables interface.
292 double const dt = std::numeric_limits<double>::quiet_NaN();
293 auto const K = MaterialPropertyLib::formEigenTensor<GlobalDim>(
295 .value(vars, pos, t, dt));
296
297 auto const mu =
298 liquid_phase
300 .template value<double>(vars, pos, t, dt);
301 GlobalDimMatrixType const K_over_mu = K / mu;
302
303 cache_mat.col(ip).noalias() = -K_over_mu * dNdx * p_nodal_values;
304
306 {
307 auto const rho_w =
308 liquid_phase
310 .template value<double>(vars, pos, t, dt);
311 auto const b =
313 [_element.getID()];
314 // here it is assumed that the vector b is directed 'downwards'
315 cache_mat.col(ip).noalias() += K_over_mu * rho_w * b;
316 }
317 }
318
319 return cache;
320 }
321
322protected:
323 static const int pressure_index = ShapeFunction::NPOINTS;
324 static const int pressure_size = ShapeFunction::NPOINTS;
325 static const int temperature_index = 0;
326 static const int temperature_size = ShapeFunction::NPOINTS;
327};
328
329} // namespace HT
330} // namespace ProcessLib
Phase const & phase(std::size_t index) const
Definition Medium.cpp:33
Property const & property(PropertyType const &p) const
Definition Phase.cpp:53
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
auto const & NsHigherOrder() const
void setElementID(std::size_t element_id)
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
Definition HTFEM.h:45
std::vector< double > const & getIntPtDarcyVelocityLocal(const double t, std::vector< double > const &local_x, std::vector< double > &cache) const
Definition HTFEM.h:239
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
Definition HTFEM.h:37
NumLib::GenericIntegrationMethod const & _integration_method
Definition HTFEM.h:171
typename ShapeMatricesType::GlobalDimNodalMatrixType GlobalDimNodalMatrixType
Definition HTFEM.h:43
static const int temperature_index
Definition HTFEM.h:325
typename ShapeMatricesType::NodalVectorType NodalVectorType
Definition HTFEM.h:39
static const int temperature_size
Definition HTFEM.h:326
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:174
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:48
static const int pressure_size
Definition HTFEM.h:324
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
Definition HTFEM.h:36
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
Definition HTFEM.h:42
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:198
Eigen::Vector3d getFlux(MathLib::Point3d const &pnt_local_coords, double const t, std::vector< double > const &local_x) const override
Definition HTFEM.h:101
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
Definition HTFEM.h:89
MeshLib::Element const & _element
Definition HTFEM.h:168
HTProcessData const & _process_data
Definition HTFEM.h:169
std::vector< IntegrationPointData< GlobalDimNodalMatrixType > > _ip_data
Definition HTFEM.h:172
static const int pressure_index
Definition HTFEM.h:323
typename ShapeMatricesType::NodalRowVectorType NodalRowVectorType
Definition HTFEM.h:40
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
ParameterLib::Parameter< double > const & aperture_size
MaterialPropertyLib::MaterialSpatialDistributionMap media_map
NumLib::NumericalStabilization stabilizer
NumLib::ShapeMatrixCache shape_matrix_cache
caches for each mesh element type the shape matrix
std::vector< Eigen::VectorXd > const projected_specific_body_force_vectors
Projected specific body force vector: R * R^T * b.