OGS
MonolithicHTFEM.h
Go to the documentation of this file.
1
12#pragma once
13
14#include <Eigen/Core>
15#include <typeinfo>
16#include <vector>
17
18#include "HTFEM.h"
19#include "HTProcessData.h"
30
31namespace ProcessLib
32{
33namespace HT
34{
35const unsigned NUM_NODAL_DOF = 2;
36
37template <typename ShapeFunction, int GlobalDim>
38class MonolithicHTFEM : public HTFEM<ShapeFunction, GlobalDim>
39{
42
43 using LocalMatrixType = typename ShapeMatricesType::template MatrixType<
44 NUM_NODAL_DOF * ShapeFunction::NPOINTS,
45 NUM_NODAL_DOF * ShapeFunction::NPOINTS>;
47 typename ShapeMatricesType::template VectorType<NUM_NODAL_DOF *
48 ShapeFunction::NPOINTS>;
49
53
56
57public:
59 std::size_t const local_matrix_size,
60 NumLib::GenericIntegrationMethod const& integration_method,
61 bool is_axially_symmetric,
62 HTProcessData const& process_data)
63 : HTFEM<ShapeFunction, GlobalDim>(
64 element, local_matrix_size, integration_method,
65 is_axially_symmetric, process_data, NUM_NODAL_DOF)
66 {
67 }
68
69 void assemble(double const t, double const dt,
70 std::vector<double> const& local_x,
71 std::vector<double> const& /*local_x_prev*/,
72 std::vector<double>& local_M_data,
73 std::vector<double>& local_K_data,
74 std::vector<double>& local_b_data) override
75 {
76 auto const local_matrix_size = local_x.size();
77 // This assertion is valid only if all nodal d.o.f. use the same shape
78 // matrices.
79 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
80
81 auto local_M = MathLib::createZeroedMatrix<LocalMatrixType>(
82 local_M_data, local_matrix_size, local_matrix_size);
83 auto local_K = MathLib::createZeroedMatrix<LocalMatrixType>(
84 local_K_data, local_matrix_size, local_matrix_size);
85 auto local_b = MathLib::createZeroedVector<LocalVectorType>(
86 local_b_data, local_matrix_size);
87
88 auto KTT = local_K.template block<temperature_size, temperature_size>(
90 auto MTT = local_M.template block<temperature_size, temperature_size>(
92 auto Kpp = local_K.template block<pressure_size, pressure_size>(
94 auto Mpp = local_M.template block<pressure_size, pressure_size>(
96 auto Bp = local_b.template block<pressure_size, 1>(pressure_index, 0);
97
98 auto const& process_data = this->_process_data;
99
101 pos.setElementID(this->_element.getID());
102
103 auto p_nodal_values = Eigen::Map<const NodalVectorType>(
104 &local_x[pressure_index], pressure_size);
105
106 auto const& medium =
107 *process_data.media_map.getMedium(this->_element.getID());
108 auto const& liquid_phase = medium.phase("AqueousLiquid");
109 auto const& solid_phase = medium.phase("Solid");
110
111 auto const& b =
112 process_data
113 .projected_specific_body_force_vectors[this->_element.getID()];
114
116
117 unsigned const n_integration_points =
119
120 std::vector<GlobalDimVectorType> ip_flux_vector;
121 double average_velocity_norm = 0.0;
122 ip_flux_vector.reserve(n_integration_points);
123
124 auto const& Ns =
125 process_data.shape_matrix_cache
126 .template NsHigherOrder<typename ShapeFunction::MeshElement>();
127
128 for (unsigned ip(0); ip < n_integration_points; ip++)
129 {
130 pos.setIntegrationPoint(ip);
131
132 auto const& ip_data = this->_ip_data[ip];
133 auto const& dNdx = ip_data.dNdx;
134 auto const& N = Ns[ip];
135 auto const& w = ip_data.integration_weight;
136
137 double T_int_pt = 0.0;
138 double p_int_pt = 0.0;
139 // Order matters: First T, then P!
140 NumLib::shapeFunctionInterpolate(local_x, N, T_int_pt, p_int_pt);
141
142 vars.temperature = T_int_pt;
143 vars.liquid_phase_pressure = p_int_pt;
144
145 vars.liquid_saturation = 1.0;
146 // \todo the argument to getValue() has to be changed for non
147 // constant storage model
148 auto const specific_storage =
150 .template value<double>(vars, pos, t, dt);
151
152 auto const porosity =
154 .template value<double>(vars, pos, t, dt);
155 vars.porosity = porosity;
156
157 auto const intrinsic_permeability =
158 MaterialPropertyLib::formEigenTensor<GlobalDim>(
159 medium
160 .property(
162 .value(vars, pos, t, dt));
163
164 auto const specific_heat_capacity_fluid =
165 liquid_phase
167 .template value<double>(vars, pos, t, dt);
168
169 // Use the fluid density model to compute the density
170 auto const fluid_density =
171 liquid_phase
173 .template value<double>(vars, pos, t, dt);
174
175 vars.density = fluid_density;
176 // Use the viscosity model to compute the viscosity
177 auto const viscosity =
178 liquid_phase
180 .template value<double>(vars, pos, t, dt);
181 GlobalDimMatrixType K_over_mu = intrinsic_permeability / viscosity;
182
183 GlobalDimVectorType const velocity =
184 process_data.has_gravity
185 ? GlobalDimVectorType(-K_over_mu * (dNdx * p_nodal_values -
186 fluid_density * b))
187 : GlobalDimVectorType(-K_over_mu * dNdx * p_nodal_values);
188
189 // matrix assembly
190 GlobalDimMatrixType const thermal_conductivity_dispersivity =
192 vars, fluid_density, specific_heat_capacity_fluid, velocity,
193 pos, t, dt);
194
195 KTT.noalias() +=
196 dNdx.transpose() * thermal_conductivity_dispersivity * dNdx * w;
197
198 ip_flux_vector.emplace_back(velocity * fluid_density *
199 specific_heat_capacity_fluid);
200 average_velocity_norm += velocity.norm();
201
202 Kpp.noalias() += w * dNdx.transpose() * K_over_mu * dNdx;
203 MTT.noalias() += w *
205 vars, porosity, fluid_density,
206 specific_heat_capacity_fluid, pos, t, dt) *
207 N.transpose() * N;
208 Mpp.noalias() += w * N.transpose() * specific_storage * N;
209 if (process_data.has_gravity)
210 {
211 Bp += w * fluid_density * dNdx.transpose() * K_over_mu * b;
212 }
213 /* with Oberbeck-Boussing assumption density difference only exists
214 * in buoyancy effects */
215 }
216
217 NumLib::assembleAdvectionMatrix<typename ShapeFunction::MeshElement>(
218 process_data.stabilizer, this->_ip_data,
219 process_data.shape_matrix_cache, ip_flux_vector,
220 average_velocity_norm / static_cast<double>(n_integration_points),
221 KTT);
222 }
223
224 std::vector<double> const& getIntPtDarcyVelocity(
225 const double t,
226 std::vector<GlobalVector*> const& x,
227 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
228 std::vector<double>& cache) const override
229 {
230 int const process_id = 0; // monolithic case.
231 auto const indices =
232 NumLib::getIndices(this->_element.getID(), *dof_table[process_id]);
233 assert(!indices.empty());
234 auto const& local_x = x[process_id]->get(indices);
235
236 return this->getIntPtDarcyVelocityLocal(t, local_x, cache);
237 }
238
239private:
240 using HTFEM<ShapeFunction, GlobalDim>::pressure_index;
241 using HTFEM<ShapeFunction, GlobalDim>::pressure_size;
242 using HTFEM<ShapeFunction, GlobalDim>::temperature_index;
243 using HTFEM<ShapeFunction, GlobalDim>::temperature_size;
244};
245
246} // namespace HT
247} // namespace ProcessLib
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
std::vector< double > const & getIntPtDarcyVelocityLocal(const double t, std::vector< double > const &local_x, std::vector< double > &cache) const
Definition HTFEM.h:239
NumLib::GenericIntegrationMethod const & _integration_method
Definition HTFEM.h:171
static const int temperature_index
Definition HTFEM.h:325
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
static const int pressure_size
Definition HTFEM.h:324
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
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::NodalMatrixType NodalMatrixType
typename ShapeMatricesType::template MatrixType< NUM_NODAL_DOF *ShapeFunction::NPOINTS, NUM_NODAL_DOF *ShapeFunction::NPOINTS > LocalMatrixType
typename ShapeMatricesType::template VectorType< NUM_NODAL_DOF * ShapeFunction::NPOINTS > LocalVectorType
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
typename ShapeMatricesType::NodalRowVectorType NodalRowVectorType
MonolithicHTFEM(MeshLib::Element const &element, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool is_axially_symmetric, HTProcessData const &process_data)
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) override
typename ShapeMatricesType::NodalVectorType NodalVectorType
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 ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
const unsigned NUM_NODAL_DOF
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
VectorType< GlobalDim > GlobalDimVectorType
VectorType< ShapeFunction::NPOINTS > NodalVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType