OGS
HeatTransportBHELocalAssemblerSoil-impl.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <vector>
14
24#include "SecondaryData.h"
25
26namespace ProcessLib
27{
28namespace HeatTransportBHE
29{
30template <typename ShapeFunction>
33 MeshLib::Element const& e,
34 NumLib::GenericIntegrationMethod const& integration_method,
35 bool const is_axially_symmetric,
36 HeatTransportBHEProcessData& process_data)
37 : _process_data(process_data),
38 _integration_method(integration_method),
39 _element_id(e.getID())
40{
41 unsigned const n_integration_points =
43
44 _ip_data.reserve(n_integration_points);
45 _secondary_data.N.resize(n_integration_points);
46
49 3 /* GlobalDim */>(e, is_axially_symmetric,
51
53 x_position.setElementID(_element_id);
54
55 // ip data initialization
56 for (unsigned ip = 0; ip < n_integration_points; ip++)
57 {
58 x_position.setIntegrationPoint(ip);
59
60 // create the class IntegrationPointDataBHE in place
61 auto const& sm = _shape_matrices[ip];
62 double const w = _integration_method.getWeightedPoint(ip).getWeight() *
63 sm.integralMeasure * sm.detJ;
64 _ip_data.push_back({sm.N, sm.dNdx, w});
65
66 _secondary_data.N[ip] = sm.N;
67 }
68}
69
70template <typename ShapeFunction>
72 double const t, double const dt, std::vector<double> const& local_x,
73 std::vector<double> const& /*local_x_prev*/,
74 std::vector<double>& local_M_data, std::vector<double>& local_K_data,
75 std::vector<double>& /*local_b_data*/)
76{
77 assert(local_x.size() == ShapeFunction::NPOINTS);
78 (void)local_x; // Avoid unused arg warning.
79
81 local_M_data, ShapeFunction::NPOINTS, ShapeFunction::NPOINTS);
83 local_K_data, ShapeFunction::NPOINTS, ShapeFunction::NPOINTS);
84
86 pos.setElementID(_element_id);
87
88 auto const& medium = *_process_data.media_map.getMedium(_element_id);
89 auto const& solid_phase = medium.phase("Solid");
90 auto const& liquid_phase = medium.phase("AqueousLiquid");
91
93
94 unsigned const n_integration_points =
95 _integration_method.getNumberOfPoints();
96
97 for (unsigned ip = 0; ip < n_integration_points; ip++)
98 {
99 pos.setIntegrationPoint(ip);
100 auto& ip_data = _ip_data[ip];
101 auto const& N = ip_data.N;
102 auto const& dNdx = ip_data.dNdx;
103 auto const& w = ip_data.integration_weight;
104
105 double T_int_pt = 0.0;
106 NumLib::shapeFunctionInterpolate(local_x, N, T_int_pt);
107
108 vars.temperature = T_int_pt;
109
110 // for now only using the solid and liquid phase parameters
111 auto const density_s =
113 .template value<double>(vars, pos, t, dt);
114
115 auto const heat_capacity_s =
116 solid_phase
117 .property(
119 .template value<double>(vars, pos, t, dt);
120
121 auto const density_f =
122 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
123 .template value<double>(vars, pos, t, dt);
124
125 auto const heat_capacity_f =
126 liquid_phase
127 .property(
129 .template value<double>(vars, pos, t, dt);
130
131 auto const porosity =
133 .template value<double>(vars, pos, t, dt);
134
135 auto const velocity =
136 liquid_phase
138 .template value<Eigen::Vector3d>(vars, pos, t, dt);
139
140 // calculate the hydrodynamic thermodispersion tensor
141 auto const thermal_conductivity =
143 medium
144 .property(
146 .value(vars, pos, t, dt));
147
148 auto thermal_conductivity_dispersivity = thermal_conductivity;
149
150 double const velocity_magnitude = velocity.norm();
151
152 if (velocity_magnitude >= std::numeric_limits<double>::epsilon())
153 {
154 auto const thermal_dispersivity_longitudinal =
155 medium
157 thermal_longitudinal_dispersivity)
158 .template value<double>();
159 auto const thermal_dispersivity_transversal =
160 medium
162 thermal_transversal_dispersivity)
163 .template value<double>();
164
165 auto const thermal_dispersivity =
166 density_f * heat_capacity_f *
167 (thermal_dispersivity_transversal * velocity_magnitude *
168 Eigen::Matrix3d::Identity() +
169 (thermal_dispersivity_longitudinal -
170 thermal_dispersivity_transversal) /
171 velocity_magnitude * velocity * velocity.transpose());
172 thermal_conductivity_dispersivity += thermal_dispersivity;
173 }
174
175 // assemble Conductance matrix
176 local_K.noalias() +=
177 (dNdx.transpose() * thermal_conductivity_dispersivity * dNdx +
178 N.transpose() * velocity.transpose() * dNdx * density_f *
179 heat_capacity_f) *
180 w;
181
182 // assemble Mass matrix
183 local_M.noalias() += N.transpose() * N * w *
184 (density_s * heat_capacity_s * (1 - porosity) +
185 density_f * heat_capacity_f * porosity);
186 }
187
188 if (_process_data._mass_lumping)
189 {
190 // only mass lumping at the BHE connected soil elements
191 if (_process_data.mass_lumping_soil_elements[_element_id])
192 {
193 local_M = local_M.colwise().sum().eval().asDiagonal();
194 }
195 }
196
197 // debugging
198 // std::string sep = "\n----------------------------------------\n";
199 // Eigen::IOFormat CleanFmt(4, 0, ", ", "\n", "[", "]");
200 // std::cout << local_K.format(CleanFmt) << sep;
201 // std::cout << local_M.format(CleanFmt) << sep;
202}
203
204template <typename ShapeFunction>
206 double const t, double const dt, std::vector<double> const& local_x,
207 std::vector<double> const& local_x_prev,
208 std::vector<double>& local_rhs_data, std::vector<double>& local_Jac_data)
209{
210 assert(local_x.size() == ShapeFunction::NPOINTS);
211 auto const local_matrix_size = local_x.size();
212 // initialize x and x_prev
213 auto x =
214 Eigen::Map<NodalVectorType const>(local_x.data(), local_matrix_size);
215 auto x_prev = Eigen::Map<NodalVectorType const>(local_x_prev.data(),
216 local_matrix_size);
217 // initialize local_Jac and local_rhs
219 local_Jac_data, local_matrix_size, local_matrix_size);
221 local_rhs_data, local_matrix_size);
222
223 std::vector<double> local_M_data;
224 std::vector<double> local_K_data;
225 assemble(t, dt, local_x, local_x_prev, local_M_data, local_K_data,
226 local_rhs_data /*not going to be used*/);
227
228 // convert to matrix
230 local_M_data, local_matrix_size, local_matrix_size);
232 local_K_data, local_matrix_size, local_matrix_size);
233
234 // Jac matrix and rhs vector operation
235 local_Jac.noalias() += local_K + local_M / dt;
236 local_rhs.noalias() -= local_K * x + local_M * (x - x_prev) / dt;
237
238 local_M.setZero();
239 local_K.setZero();
240}
241
242} // namespace HeatTransportBHE
243} // namespace ProcessLib
double getWeight() const
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
void assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) override
std::vector< IntegrationPointDataSoil< NodalRowVectorType, GlobalDimNodalMatrixType >, Eigen::aligned_allocator< IntegrationPointDataSoil< NodalRowVectorType, GlobalDimNodalMatrixType > > > _ip_data
HeatTransportBHELocalAssemblerSoil(HeatTransportBHELocalAssemblerSoil const &)=delete
template Eigen::Matrix< double, 3, 3 > formEigenTensor< 3 >(MaterialPropertyLib::PropertyDataType const &values)
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< 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< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N