OGS
HeatConductionFEM.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <vector>
14
29
30namespace ProcessLib
31{
32namespace HeatConduction
33{
34const unsigned NUM_NODAL_DOF = 1;
35
39{
40public:
41 virtual std::vector<double> const& getIntPtHeatFlux(
42 const double t,
43 std::vector<GlobalVector*> const& x,
44 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
45 std::vector<double>& cache) const = 0;
46};
47
48template <typename ShapeFunction, int GlobalDim>
50{
53
55 ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim>;
56
60
61public:
65 MeshLib::Element const& element,
66 std::size_t const local_matrix_size,
67 NumLib::GenericIntegrationMethod const& integration_method,
68 bool is_axially_symmetric,
69 HeatConductionProcessData const& process_data)
70 : _element(element),
71 _process_data(process_data),
72 _integration_method(integration_method),
74 NumLib::initShapeMatrices<ShapeFunction, ShapeMatricesType,
75 GlobalDim>(
76 element, is_axially_symmetric, _integration_method))
77 {
78 // This assertion is valid only if all nodal d.o.f. use the same shape
79 // matrices.
80 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
81 (void)local_matrix_size;
82 }
83
84 void assemble(double const t, double const dt,
85 std::vector<double> const& local_x,
86 std::vector<double> const& /*local_x_prev*/,
87 std::vector<double>& local_M_data,
88 std::vector<double>& local_K_data,
89 std::vector<double>& /*local_b_data*/) override
90 {
91 auto const local_matrix_size = local_x.size();
92 // This assertion is valid only if all nodal d.o.f. use the same shape
93 // matrices.
94 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
95
97 local_M_data, local_matrix_size, local_matrix_size);
99 local_K_data, local_matrix_size, local_matrix_size);
100
101 unsigned const n_integration_points =
103
104 auto const& medium =
107
108 for (unsigned ip = 0; ip < n_integration_points; ip++)
109 {
110 auto const& sm = _shape_matrices[ip];
112 std::nullopt, _element.getID(), ip,
114 NumLib::interpolateCoordinates<ShapeFunction,
116 sm.N))};
117
118 auto const& wp = _integration_method.getWeightedPoint(ip);
119
120 // get the local temperature and put it in the variable array for
121 // access in MPL
122 double T_int_pt = 0.0;
123 NumLib::shapeFunctionInterpolate(local_x, sm.N, T_int_pt);
124 vars.temperature = T_int_pt;
125
127 medium
128 .property(
130 .value(vars, pos, t, dt));
131 auto const specific_heat_capacity =
132 medium
134 specific_heat_capacity)
135 .template value<double>(vars, pos, t, dt);
136 auto const density =
138 .template value<double>(vars, pos, t, dt);
139
140 local_K.noalias() += sm.dNdx.transpose() * k * sm.dNdx * sm.detJ *
141 wp.getWeight() * sm.integralMeasure;
142 local_M.noalias() += sm.N.transpose() * density *
143 specific_heat_capacity * sm.N * sm.detJ *
144 wp.getWeight() * sm.integralMeasure;
145 }
147 {
148 local_M = local_M.colwise().sum().eval().asDiagonal();
149 }
150 }
151
152 void assembleWithJacobian(double const t, double const dt,
153 std::vector<double> const& local_x,
154 std::vector<double> const& local_x_prev,
155 std::vector<double>& local_rhs_data,
156 std::vector<double>& local_Jac_data) override
157 {
158 auto const local_matrix_size = local_x.size();
159 // This assertion is valid only if all nodal d.o.f. use the same shape
160 // matrices.
161 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
162
163 auto x = Eigen::Map<NodalVectorType const>(local_x.data(),
164 local_matrix_size);
165
166 auto x_prev = Eigen::Map<NodalVectorType const>(local_x_prev.data(),
167 local_matrix_size);
168
170 local_Jac_data, local_matrix_size, local_matrix_size);
172 local_rhs_data, local_matrix_size);
173
174 NodalMatrixType laplace =
175 NodalMatrixType::Zero(local_matrix_size, local_matrix_size);
176 NodalMatrixType storage =
177 NodalMatrixType::Zero(local_matrix_size, local_matrix_size);
178
179 unsigned const n_integration_points =
181
184
185 auto const& medium =
188
189 for (unsigned ip = 0; ip < n_integration_points; ip++)
190 {
191 pos.setIntegrationPoint(ip);
192 auto const& sm = _shape_matrices[ip];
193 double const w =
195 sm.integralMeasure;
196
197 // get the local temperature and put it in the variable array for
198 // access in MPL
199 double T_int_pt = 0.0;
200 NumLib::shapeFunctionInterpolate(local_x, sm.N, T_int_pt);
201 vars.temperature = T_int_pt;
202
204 medium
205 .property(
207 .value(vars, pos, t, dt));
208 auto const specific_heat_capacity =
209 medium
211 specific_heat_capacity)
212 .template value<double>(vars, pos, t, dt);
213 auto const density =
215 .template value<double>(vars, pos, t, dt);
216
217 laplace.noalias() += sm.dNdx.transpose() * k * sm.dNdx * w;
218 storage.noalias() +=
219 sm.N.transpose() * density * specific_heat_capacity * sm.N * w;
220 }
222 {
223 storage = storage.colwise().sum().eval().asDiagonal();
224 }
225
226 local_Jac.noalias() += laplace + storage / dt;
227 local_rhs.noalias() -= laplace * x + storage * (x - x_prev) / dt;
228 }
229
230 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
231 const unsigned integration_point) const override
232 {
233 auto const& N = _shape_matrices[integration_point].N;
234
235 // assumes N is stored contiguously in memory
236 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
237 }
238
239 std::vector<double> const& getIntPtHeatFlux(
240 const double t,
241 std::vector<GlobalVector*> const& x,
242 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
243 std::vector<double>& cache) const override
244 {
245 int const process_id = 0; // monolithic case.
246 auto const indices =
247 NumLib::getIndices(_element.getID(), *dof_table[process_id]);
248 assert(!indices.empty());
249 auto const& local_x = x[process_id]->get(indices);
250
251 auto const T_nodal_values = Eigen::Map<const NodalVectorType>(
252 local_x.data(), ShapeFunction::NPOINTS);
253
254 unsigned const n_integration_points =
256
259
260 auto const& medium =
263
264 double const dt = std::numeric_limits<double>::quiet_NaN();
265 cache.clear();
266 auto cache_mat = MathLib::createZeroedMatrix<
267 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
268 cache, GlobalDim, n_integration_points);
269
270 for (unsigned ip = 0; ip < n_integration_points; ip++)
271 {
272 pos.setIntegrationPoint(ip);
273 auto const& sm = _shape_matrices[ip];
274 // get the local temperature and put it in the variable array for
275 // access in MPL
276 vars.temperature = sm.N.dot(T_nodal_values);
277
279 medium
280 .property(
282 .value(vars, pos, t, dt));
283
284 // heat flux only computed for output.
285 cache_mat.col(ip).noalias() = -k * sm.dNdx * T_nodal_values;
286 }
287
288 return cache;
289 }
290
291private:
294
296 std::vector<ShapeMatrices, Eigen::aligned_allocator<ShapeMatrices>>
298};
299
300} // namespace HeatConduction
301} // namespace ProcessLib
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
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
virtual std::vector< double > const & getIntPtHeatFlux(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
HeatConductionProcessData const & _process_data
typename LocalAssemblerTraits::LocalVector NodalVectorType
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
std::vector< double > const & getIntPtHeatFlux(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const 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
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 > &) override
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
LocalAssemblerData(MeshLib::Element const &element, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool is_axially_symmetric, HeatConductionProcessData const &process_data)
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
NumLib::GenericIntegrationMethod const & _integration_method
typename LocalAssemblerTraits::LocalMatrix NodalMatrixType
Eigen::Matrix< double, GlobalDim, GlobalDim > formEigenTensor(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)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
std::array< double, 3 > interpolateCoordinates(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
VectorType< GlobalDim > GlobalDimVectorType
MaterialPropertyLib::MaterialSpatialDistributionMap media_map
bool const mass_lumping
If set mass lumping will be applied to the equation.
Matrix< NNodes *NodalDOF, NNodes *NodalDOF > LocalMatrix