OGS
HeatTransportBHELocalAssemblerBHE-impl.h
Go to the documentation of this file.
1
11#pragma once
12
16
17namespace ProcessLib
18{
19namespace HeatTransportBHE
20{
21template <typename ShapeFunction, typename BHEType>
24 MeshLib::Element const& e,
25 NumLib::GenericIntegrationMethod const& integration_method,
26 BHEType const& bhe,
27 bool const is_axially_symmetric,
28 HeatTransportBHEProcessData& process_data)
29 : _process_data(process_data),
30 _integration_method(integration_method),
31 _bhe(bhe),
32 _element_id(e.getID())
33{
34 // need to make sure that the BHE elements are one-dimensional
35 assert(e.getDimension() == 1);
36
37 unsigned const n_integration_points =
39
40 _ip_data.reserve(n_integration_points);
41 _secondary_data.N.resize(n_integration_points);
42
43 auto const shape_matrices =
45 3 /* GlobalDim */>(e, is_axially_symmetric,
47
48 // ip data initialization
49 for (unsigned ip = 0; ip < n_integration_points; ip++)
50 {
51 auto const& sm = shape_matrices[ip];
52 // create the class IntegrationPointDataBHE in place
53 _ip_data.push_back(
54 {sm.N, sm.dNdx,
56 sm.integralMeasure * sm.detJ});
57
58 _secondary_data.N[ip] = sm.N;
59 }
60
61 // calculate the element direction vector
62 auto const& p0 = e.getNode(0)->asEigenVector3d();
63 auto const& p1 = e.getNode(1)->asEigenVector3d();
64
65 _element_direction = (p1 - p0).normalized();
66
70 static constexpr int max_num_thermal_exchange_terms = 5;
71 // formulate the local BHE R matrix
72 for (int idx_bhe_unknowns = 0; idx_bhe_unknowns < bhe_unknowns;
73 idx_bhe_unknowns++)
74 {
75 typename ShapeMatricesType::template MatrixType<
77 matBHE_loc_R = ShapeMatricesType::template MatrixType<
81 // Loop over Gauss points
82 for (unsigned ip = 0; ip < n_integration_points; ip++)
83 {
84 auto const& N = _ip_data[ip].N;
85 auto const& w = _ip_data[ip].integration_weight;
86
87 auto const& R = _bhe.thermalResistance(idx_bhe_unknowns);
88 // calculate mass matrix for current unknown
89 matBHE_loc_R += N.transpose() * N * (1 / R) * w;
90 } // end of loop over integration point
91
92 // The following assembly action is according to Diersch (2013) FEFLOW
93 // book please refer to M.127 and M.128 on page 955 and 956
94 // The if check is absolutely necessary because
95 // (i) In the CXA and CXC case, there are 3 exchange terms,
96 // and it is the same as the number of unknowns;
97 // (ii) In the 1U case, there are 4 exchange terms,
98 // and it is again same as the number of unknowns;
99 // (iii) In the 2U case, there are 5 exchange terms,
100 // and it is less than the number of unknowns (8).
101 if (idx_bhe_unknowns < max_num_thermal_exchange_terms)
102 {
103 _bhe.template assembleRMatrices<ShapeFunction::NPOINTS>(
104 idx_bhe_unknowns, matBHE_loc_R, _R_matrix, _R_pi_s_matrix,
106 }
107 } // end of loop over BHE unknowns
108
109 // debugging
110 // std::string sep =
111 // "\n----------------------------------------\n";
112 // Eigen::IOFormat CleanFmt(4, 0, ", ", "\n", "[", "]");
113 // std::cout << "_R_matrix: \n" << sep;
114 // std::cout << _R_matrix.format(CleanFmt) << sep;
115 // std::cout << "_R_s_matrix: \n" << sep;
116 // std::cout << _R_s_matrix.format(CleanFmt) << sep;
117 // std::cout << "_R_pi_s_matrix: \n" << sep;
118 // std::cout << _R_pi_s_matrix.format(CleanFmt) << sep;
119}
120
121template <typename ShapeFunction, typename BHEType>
123 double const /*t*/, double const /*dt*/,
124 std::vector<double> const& /*local_x*/,
125 std::vector<double> const& /*local_x_prev*/,
126 std::vector<double>& local_M_data, std::vector<double>& local_K_data,
127 std::vector<double>& /*local_b_data*/) // local b vector is not touched
128{
129 auto local_M = MathLib::createZeroedMatrix<BheLocalMatrixType>(
130 local_M_data, local_matrix_size, local_matrix_size);
131 auto local_K = MathLib::createZeroedMatrix<BheLocalMatrixType>(
132 local_K_data, local_matrix_size, local_matrix_size);
133
134 unsigned const n_integration_points =
135 _integration_method.getNumberOfPoints();
136
137 auto const& pipe_heat_capacities = _bhe.pipeHeatCapacities();
138 auto const& pipe_heat_conductions = _bhe.pipeHeatConductions();
139 auto const& pipe_advection_vectors =
140 _bhe.pipeAdvectionVectors(_element_direction);
141 auto const& cross_section_areas = _bhe.crossSectionAreas();
142
143 // the mass and conductance matrix terms
144 for (unsigned ip = 0; ip < n_integration_points; ip++)
145 {
146 auto& ip_data = _ip_data[ip];
147
148 auto const& w = ip_data.integration_weight;
149 auto const& N = ip_data.N;
150 auto const& dNdx = ip_data.dNdx;
151
152 // looping over all unknowns.
153 for (int idx_bhe_unknowns = 0; idx_bhe_unknowns < bhe_unknowns;
154 idx_bhe_unknowns++)
155 {
156 // get coefficient of mass from corresponding BHE.
157 auto const& mass_coeff = pipe_heat_capacities[idx_bhe_unknowns];
158 auto const& lambda = pipe_heat_conductions[idx_bhe_unknowns];
159 auto const& advection_vector =
160 pipe_advection_vectors[idx_bhe_unknowns];
161 auto const& A = cross_section_areas[idx_bhe_unknowns];
162
163 int const single_bhe_unknowns_index =
164 bhe_unknowns_index +
165 single_bhe_unknowns_size * idx_bhe_unknowns;
166 // local M
167 local_M
168 .template block<single_bhe_unknowns_size,
169 single_bhe_unknowns_size>(
170 single_bhe_unknowns_index, single_bhe_unknowns_index)
171 .noalias() += N.transpose() * N * mass_coeff * A * w;
172
173 // local K
174 // laplace part
175 local_K
176 .template block<single_bhe_unknowns_size,
177 single_bhe_unknowns_size>(
178 single_bhe_unknowns_index, single_bhe_unknowns_index)
179 .noalias() += dNdx.transpose() * dNdx * lambda * A * w;
180 // advection part
181 local_K
182 .template block<single_bhe_unknowns_size,
183 single_bhe_unknowns_size>(
184 single_bhe_unknowns_index, single_bhe_unknowns_index)
185 .noalias() +=
186 N.transpose() * advection_vector.transpose() * dNdx * A * w;
187 }
188 }
189
190 // add the R matrix to local_K
191 local_K.template block<bhe_unknowns_size, bhe_unknowns_size>(
192 bhe_unknowns_index, bhe_unknowns_index) += _R_matrix;
193
194 // add the R_pi_s matrix to local_K
195 local_K
196 .template block<bhe_unknowns_size, soil_temperature_size>(
197 bhe_unknowns_index, soil_temperature_index)
198 .noalias() += _R_pi_s_matrix;
199 local_K
200 .template block<soil_temperature_size, bhe_unknowns_size>(
201 soil_temperature_index, bhe_unknowns_index)
202 .noalias() += _R_pi_s_matrix.transpose();
203
204 // add the R_s matrix to local_K
205 local_K
206 .template block<soil_temperature_size, soil_temperature_size>(
207 soil_temperature_index, soil_temperature_index)
208 .noalias() += _bhe.number_of_grout_zones * _R_s_matrix;
209
210 // debugging
211 // std::string sep = "\n----------------------------------------\n";
212 // Eigen::IOFormat CleanFmt(4, 0, ", ", "\n", "[", "]");
213 // std::cout << local_K.format(CleanFmt) << sep;
214 // std::cout << local_M.format(CleanFmt) << sep;
215}
216
217template <typename ShapeFunction, typename BHEType>
219 assembleWithJacobian(double const t, double const dt,
220 std::vector<double> const& local_x,
221 std::vector<double> const& local_x_prev,
222 std::vector<double>& local_M_data,
223 std::vector<double>& local_K_data,
224 std::vector<double>& local_rhs_data,
225 std::vector<double>& local_Jac_data)
226{
227 auto const local_matrix_size = local_x.size();
228 // initialize x and x_prev
229 auto x =
230 Eigen::Map<BheLocalVectorType const>(local_x.data(), local_matrix_size);
231 auto x_prev = Eigen::Map<BheLocalVectorType const>(local_x_prev.data(),
232 local_matrix_size);
233 // initialize local_Jac and local_rhs
234 auto local_Jac = MathLib::createZeroedMatrix<BheLocalMatrixType>(
235 local_Jac_data, local_matrix_size, local_matrix_size);
236 auto local_rhs = MathLib::createZeroedVector<BheLocalVectorType>(
237 local_rhs_data, local_matrix_size);
238
239 assemble(t, dt, local_x, local_x_prev, local_M_data, local_K_data,
240 local_rhs_data /*not going to be used*/);
241
242 // convert to matrix
243 auto local_M = MathLib::toMatrix<BheLocalMatrixType>(
244 local_M_data, local_matrix_size, local_matrix_size);
245 auto local_K = MathLib::toMatrix<BheLocalMatrixType>(
246 local_K_data, local_matrix_size, local_matrix_size);
247
248 // Jac matrix and rhs vector operation
249 local_Jac.noalias() += local_K + local_M / dt;
250 local_rhs.noalias() -= local_K * x + local_M * (x - x_prev) / dt;
251
252 local_M.setZero();
253 local_K.setZero();
254}
255} // namespace HeatTransportBHE
256} // namespace ProcessLib
Eigen::Vector3d const & asEigenVector3d() const
Definition Point3d.h:63
double getWeight() const
virtual const Node * getNode(unsigned idx) const =0
virtual constexpr unsigned getDimension() const =0
Get dimension of the mesh element.
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
ShapeMatricesType::template MatrixType< bhe_unknowns_size, bhe_unknowns_size > _R_matrix
ShapeMatricesType::template MatrixType< soil_temperature_size, soil_temperature_size > _R_s_matrix
std::vector< IntegrationPointDataBHE< ShapeMatricesType >, Eigen::aligned_allocator< IntegrationPointDataBHE< ShapeMatricesType > > > _ip_data
ShapeMatricesType::template MatrixType< bhe_unknowns_size, soil_temperature_size > _R_pi_s_matrix
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 > &, std::vector< double > &, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) override
HeatTransportBHELocalAssemblerBHE(HeatTransportBHELocalAssemblerBHE const &)=delete
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