OGS
HeatConductionFEM.h
Go to the documentation of this file.
1 
11 #pragma once
12 
13 #include <vector>
14 
16 #include "MaterialLib/MPL/Medium.h"
26 
27 namespace ProcessLib
28 {
29 namespace HeatConduction
30 {
31 const unsigned NUM_NODAL_DOF = 1;
32 
36 {
37 public:
38  virtual std::vector<double> const& getIntPtHeatFluxX(
39  const double t,
40  std::vector<GlobalVector*> const& x,
41  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
42  std::vector<double>& cache) const = 0;
43 
44  virtual std::vector<double> const& getIntPtHeatFluxY(
45  const double t,
46  std::vector<GlobalVector*> const& x,
47  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
48  std::vector<double>& cache) const = 0;
49 
50  virtual std::vector<double> const& getIntPtHeatFluxZ(
51  const double t,
52  std::vector<GlobalVector*> const& x,
53  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
54  std::vector<double>& cache) const = 0;
55 };
56 
57 template <typename ShapeFunction, typename IntegrationMethod, int GlobalDim>
59 {
62 
64  ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim>;
65 
69 
70 public:
74  std::size_t const local_matrix_size,
75  bool is_axially_symmetric,
76  unsigned const integration_order,
77  HeatConductionProcessData const& process_data)
78  : _element(element),
79  _process_data(process_data),
80  _integration_method(integration_order),
83  GlobalDim>(
84  element, is_axially_symmetric, _integration_method)),
86  GlobalDim,
87  std::vector<double>(_integration_method.getNumberOfPoints()))
88  {
89  // This assertion is valid only if all nodal d.o.f. use the same shape
90  // matrices.
91  assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
92  (void)local_matrix_size;
93  }
94 
95  void assemble(double const t, double const dt,
96  std::vector<double> const& local_x,
97  std::vector<double> const& /*local_xdot*/,
98  std::vector<double>& local_M_data,
99  std::vector<double>& local_K_data,
100  std::vector<double>& /*local_b_data*/) override
101  {
102  auto const local_matrix_size = local_x.size();
103  // This assertion is valid only if all nodal d.o.f. use the same shape
104  // matrices.
105  assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
106 
107  auto local_M = MathLib::createZeroedMatrix<NodalMatrixType>(
108  local_M_data, local_matrix_size, local_matrix_size);
109  auto local_K = MathLib::createZeroedMatrix<NodalMatrixType>(
110  local_K_data, local_matrix_size, local_matrix_size);
111 
112  unsigned const n_integration_points =
113  _integration_method.getNumberOfPoints();
114 
116  pos.setElementID(_element.getID());
117 
118  auto const& medium =
119  *_process_data.media_map->getMedium(_element.getID());
121  vars[static_cast<int>(MaterialPropertyLib::Variable::temperature)] =
122  medium
123  .property(
125  .template value<double>(vars, pos, t, dt);
126 
127  for (unsigned ip = 0; ip < n_integration_points; ip++)
128  {
129  pos.setIntegrationPoint(ip);
130  auto const& sm = _shape_matrices[ip];
131  auto const& wp = _integration_method.getWeightedPoint(ip);
132  auto const k = MaterialPropertyLib::formEigenTensor<GlobalDim>(
133  medium
134  .property(
136  .value(vars, pos, t, dt));
137  auto const heat_capacity =
138  medium
139  .property(
141  .template value<double>(vars, pos, t, dt);
142  auto const density =
143  medium
144  .property(
146  .template value<double>(vars, pos, t, dt);
147 
148  local_K.noalias() += sm.dNdx.transpose() * k * sm.dNdx * sm.detJ *
149  wp.getWeight() * sm.integralMeasure;
150  local_M.noalias() += sm.N.transpose() * density * heat_capacity *
151  sm.N * sm.detJ * wp.getWeight() *
152  sm.integralMeasure;
153  }
155  {
156  local_M = local_M.colwise().sum().eval().asDiagonal();
157  }
158  }
159 
160  void assembleWithJacobian(double const t, double const dt,
161  std::vector<double> const& local_x,
162  std::vector<double> const& local_xdot,
163  const double /*dxdot_dx*/, const double /*dx_dx*/,
164  std::vector<double>& /*local_M_data*/,
165  std::vector<double>& /*local_K_data*/,
166  std::vector<double>& local_rhs_data,
167  std::vector<double>& local_Jac_data) override
168  {
169  auto const local_matrix_size = local_x.size();
170  // This assertion is valid only if all nodal d.o.f. use the same shape
171  // matrices.
172  assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
173 
174  auto x = Eigen::Map<NodalVectorType const>(local_x.data(),
175  local_matrix_size);
176 
177  auto x_dot = Eigen::Map<NodalVectorType const>(local_xdot.data(),
178  local_matrix_size);
179 
180  auto local_Jac = MathLib::createZeroedMatrix<NodalMatrixType>(
181  local_Jac_data, local_matrix_size, local_matrix_size);
182  auto local_rhs = MathLib::createZeroedVector<NodalVectorType>(
183  local_rhs_data, local_matrix_size);
184 
185  NodalMatrixType laplace =
186  NodalMatrixType::Zero(local_matrix_size, local_matrix_size);
188  NodalMatrixType::Zero(local_matrix_size, local_matrix_size);
189 
190  unsigned const n_integration_points =
191  _integration_method.getNumberOfPoints();
192 
194  pos.setElementID(_element.getID());
195 
196  auto const& medium =
197  *_process_data.media_map->getMedium(_element.getID());
199  vars[static_cast<int>(MaterialPropertyLib::Variable::temperature)] =
200  medium
201  .property(
203  .template value<double>(vars, pos, t, dt);
204 
205  for (unsigned ip = 0; ip < n_integration_points; ip++)
206  {
207  pos.setIntegrationPoint(ip);
208  auto const& sm = _shape_matrices[ip];
209  double const w =
210  _integration_method.getWeightedPoint(ip).getWeight() * sm.detJ *
211  sm.integralMeasure;
212 
213  auto const k = MaterialPropertyLib::formEigenTensor<GlobalDim>(
214  medium
215  .property(
217  .value(vars, pos, t, dt));
218  auto const heat_capacity =
219  medium
220  .property(
222  .template value<double>(vars, pos, t, dt);
223  auto const density =
224  medium
225  .property(
227  .template value<double>(vars, pos, t, dt);
228 
229  laplace.noalias() += sm.dNdx.transpose() * k * sm.dNdx * w;
230  storage.noalias() +=
231  sm.N.transpose() * density * heat_capacity * sm.N * w;
232  }
234  {
235  storage = storage.colwise().sum().eval().asDiagonal();
236  }
237 
238  local_Jac.noalias() += laplace + storage / dt;
239  local_rhs.noalias() -= laplace * x + storage * x_dot;
240  }
241 
243  double const t, double const dt, Eigen::VectorXd const& local_x,
244  Eigen::VectorXd const& /*local_x_dot*/) override
245  {
246  unsigned const n_integration_points =
247  _integration_method.getNumberOfPoints();
248 
250  pos.setElementID(_element.getID());
251 
252  auto const& medium =
253  *_process_data.media_map->getMedium(_element.getID());
255  vars[static_cast<int>(MaterialPropertyLib::Variable::temperature)] =
256  medium
257  .property(
259  .template value<double>(vars, pos, t, dt);
260 
261  for (unsigned ip = 0; ip < n_integration_points; ip++)
262  {
263  pos.setIntegrationPoint(ip);
264  auto const& sm = _shape_matrices[ip];
265  auto const k = MaterialPropertyLib::formEigenTensor<GlobalDim>(
266  medium
267  .property(
269  .value(vars, pos, t, dt));
270  // heat flux only computed for output.
271  GlobalDimVectorType const heat_flux = -k * sm.dNdx * local_x;
272 
273  for (int d = 0; d < GlobalDim; ++d)
274  {
275  _heat_fluxes[d][ip] = heat_flux[d];
276  }
277  }
278  }
279 
280  Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
281  const unsigned integration_point) const override
282  {
283  auto const& N = _shape_matrices[integration_point].N;
284 
285  // assumes N is stored contiguously in memory
286  return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
287  }
288 
289  std::vector<double> const& getIntPtHeatFluxX(
290  const double /*t*/,
291  std::vector<GlobalVector*> const& /*x*/,
292  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
293  std::vector<double>& /*cache*/) const override
294  {
295  assert(!_heat_fluxes.empty());
296  return _heat_fluxes[0];
297  }
298 
299  std::vector<double> const& getIntPtHeatFluxY(
300  const double /*t*/,
301  std::vector<GlobalVector*> const& /*x*/,
302  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
303  std::vector<double>& /*cache*/) const override
304  {
305  assert(_heat_fluxes.size() > 1);
306  return _heat_fluxes[1];
307  }
308 
309  std::vector<double> const& getIntPtHeatFluxZ(
310  const double /*t*/,
311  std::vector<GlobalVector*> const& /*x*/,
312  std::vector<NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
313  std::vector<double>& /*cache*/) const override
314  {
315  assert(_heat_fluxes.size() > 2);
316  return _heat_fluxes[2];
317  }
318 
319 private:
322 
323  IntegrationMethod const _integration_method;
324  std::vector<ShapeMatrices, Eigen::aligned_allocator<ShapeMatrices>>
326 
327  std::vector<std::vector<double>> _heat_fluxes;
328 };
329 
330 } // namespace HeatConduction
331 } // namespace ProcessLib
virtual std::size_t getID() const final
Returns the ID of the element.
Definition: Element.h:82
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
virtual std::vector< double > const & getIntPtHeatFluxY(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtHeatFluxX(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtHeatFluxZ(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
std::vector< std::vector< double > > _heat_fluxes
LocalAssemblerData(MeshLib::Element const &element, std::size_t const local_matrix_size, bool is_axially_symmetric, unsigned const integration_order, HeatConductionProcessData const &process_data)
std::vector< double > const & getIntPtHeatFluxY(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override
typename LocalAssemblerTraits::LocalMatrix NodalMatrixType
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
HeatConductionProcessData const & _process_data
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
std::vector< double > const & getIntPtHeatFluxZ(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
typename LocalAssemblerTraits::LocalVector NodalVectorType
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_xdot, const double, const double, std::vector< double > &, std::vector< double > &, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) override
void computeSecondaryVariableConcrete(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &) override
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
std::vector< double > const & getIntPtHeatFluxX(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override
std::array< VariableType, static_cast< int >(Variable::number_of_variables)> VariableArray
Definition: VariableType.h:108
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)
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
VectorType< GlobalDim > GlobalDimVectorType
std::unique_ptr< MaterialPropertyLib::MaterialSpatialDistributionMap > media_map
bool const mass_lumping
If set mass lumping will be applied to the equation.
Matrix< NNodes *NodalDOF, NNodes *NodalDOF > LocalMatrix