OGS
Interpolation.h
Go to the documentation of this file.
1 
11 #pragma once
12 
13 #include <array>
14 #include <cassert>
15 
19 
20 namespace NumLib
21 {
22 namespace detail
23 {
25 template <unsigned DOFOffset, typename NodalValues, typename ShapeMatrix>
26 void shapeFunctionInterpolate(const NodalValues& /*nodal_values*/,
27  const ShapeMatrix& /*shape_matrix_N*/)
28 {
29 }
30 
32 template <unsigned DOFOffset, typename NodalValues, typename ShapeMatrix,
33  typename... ScalarTypes>
34 void shapeFunctionInterpolate(const NodalValues& nodal_values,
35  const ShapeMatrix& shape_matrix_N,
36  double& interpolated_value,
37  ScalarTypes&... interpolated_values)
38 {
39  auto const num_nodes = shape_matrix_N.size();
40  double iv = 0.0;
41 
42  for (auto n = decltype(num_nodes){0}; n < num_nodes; ++n)
43  {
44  iv += nodal_values[DOFOffset * num_nodes + n] * shape_matrix_N[n];
45  }
46 
47  interpolated_value = iv;
48 
49  shapeFunctionInterpolate<DOFOffset + 1>(nodal_values, shape_matrix_N,
50  interpolated_values...);
51 }
52 
53 } // namespace detail
54 
78 template <typename NodalValues, typename ShapeMatrix, typename... ScalarTypes>
79 void shapeFunctionInterpolate(const NodalValues& nodal_values,
80  const ShapeMatrix& shape_matrix_N,
81  double& interpolated_value,
82  ScalarTypes&... interpolated_values)
83 {
84  auto const num_nodal_dof = sizeof...(interpolated_values) + 1;
85  auto const num_nodes = shape_matrix_N.size();
86 
87  assert(num_nodes * num_nodal_dof ==
88  static_cast<std::size_t>(nodal_values.size()));
89  (void)num_nodal_dof;
90  (void)num_nodes; // no warnings when not in debug build
91 
92  detail::shapeFunctionInterpolate<0>(nodal_values, shape_matrix_N,
93  interpolated_value,
94  interpolated_values...);
95 }
96 
104 template <typename LowerOrderShapeFunction, typename HigherOrderMeshElementType,
105  int GlobalDim, typename EigenMatrixType>
107  MeshLib::Element const& element, bool const is_axially_symmetric,
108  Eigen::MatrixBase<EigenMatrixType> const& node_values,
109  MeshLib::PropertyVector<double>& interpolated_values_global_vector)
110 {
111  assert(dynamic_cast<HigherOrderMeshElementType const*>(&element));
112  assert(node_values.cols() == 1); // Scalar quantity only.
113 
114  using SF = LowerOrderShapeFunction;
115  using ShapeMatricesType = ShapeMatrixPolicyType<SF, GlobalDim>;
116 
117  int const number_base_nodes = element.getNumberOfBaseNodes();
118  int const number_all_nodes = element.getNumberOfNodes();
119 
120  // Copy the values for linear nodes.
121  for (int n = 0; n < number_base_nodes; ++n)
122  {
123  std::size_t const global_index = getNodeIndex(element, n);
124  interpolated_values_global_vector[global_index] = node_values[n];
125  }
126 
127  //
128  // Interpolate values for higher order nodes.
129  //
130 
131  int const number_higher_order_nodes = number_all_nodes - number_base_nodes;
132  std::vector<MathLib::Point3d> higher_order_nodes;
133  higher_order_nodes.reserve(number_higher_order_nodes);
134  for (int n = 0; n < number_higher_order_nodes; ++n)
135  {
136  higher_order_nodes.emplace_back(
138  [number_base_nodes + n]);
139  }
140 
141  // Shape matrices evaluated at higher order nodes' coordinates.
142  auto const shape_matrices =
143  computeShapeMatrices<SF, ShapeMatricesType, GlobalDim,
144  ShapeMatrixType::N>(element, is_axially_symmetric,
145  higher_order_nodes);
146 
147  for (int n = 0; n < number_higher_order_nodes; ++n)
148  {
149  std::size_t const global_index =
150  getNodeIndex(element, number_base_nodes + n);
151  interpolated_values_global_vector[global_index] =
152  shape_matrices[n].N * node_values;
153  }
154 }
155 
156 } // namespace NumLib
virtual unsigned getNumberOfNodes() const =0
virtual unsigned getNumberOfBaseNodes() const =0
std::size_t getNodeIndex(Element const &element, unsigned const idx)
Definition: Element.cpp:225
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
Definition: Interpolation.h:26
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > computeShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, PointContainer const &points)
void shapeFunctionInterpolate(const NodalValues &nodal_values, const ShapeMatrix &shape_matrix_N, double &interpolated_value, ScalarTypes &... interpolated_values)
Definition: Interpolation.h:79
void interpolateToHigherOrderNodes(MeshLib::Element const &element, bool const is_axially_symmetric, Eigen::MatrixBase< EigenMatrixType > const &node_values, MeshLib::PropertyVector< double > &interpolated_values_global_vector)