Loading [MathJax]/jax/output/HTML-CSS/config.js
OGS
EmbeddedAnchor.cpp
Go to the documentation of this file.
1
11#include "EmbeddedAnchor.h"
12
13#include <Eigen/Core>
14#include <cassert>
15#include <numeric>
16#include <range/v3/view/enumerate.hpp>
17#include <unsupported/Eigen/KroneckerProduct>
18
28
29namespace ProcessLib
30{
31// The local indices are, due to the nature of the DOF table, all even
32// for the first node and odd for the second node.
33template <int GlobalDim>
34auto nodeLocalIndices(std::size_t const i)
35{
36 return Eigen::seqN(i, Eigen::fix<GlobalDim>, Eigen::fix<2>);
37}
38
39template <int GlobalDim>
40Eigen::RowVectorXd computeShapeMatrix(
41 MeshLib::Element const& bulk_element,
42 std::span<double const, 3> const nat_coords)
43{
44 using ShapeMatricesType =
45 EigenDynamicShapeMatrixPolicy<void /* dummy */, -1 /* dynamic order */>;
46
47 typename ShapeMatricesType::ShapeMatrices shape_matrix(
48 bulk_element.getDimension(), GlobalDim,
49 bulk_element.getNumberOfNodes());
50
51 bool matched = false;
52
54 [&](auto* tag)
55 {
56 using ElementTrait = std::remove_pointer_t<decltype(tag)>;
57 using ShapeFunction = typename ElementTrait::ShapeFunction;
58 using MeshElement = typename ElementTrait::Element;
59
60 if (!matched &&
61 dynamic_cast<MeshElement const*>(&bulk_element) != nullptr)
62 {
64 ShapeFunction, ShapeMatricesType>(bulk_element);
65
66 fe.template computeShapeFunctions<NumLib::ShapeMatrixType::N>(
67 nat_coords.data(), shape_matrix, GlobalDim,
68 false /* axisymmetric */);
69
70 matched = true;
71 }
72 });
73
74 if (!matched)
75 {
76 OGS_FATAL("Element type '{:s}' is not supported as anchor element.",
78 }
79
80 return shape_matrix.N;
81}
82
83template <int GlobalDim>
84std::tuple<Eigen::VectorXd, Eigen::MatrixXd> assembleLocalBJac(
85 Eigen::Vector<double, GlobalDim> const& f,
86 Eigen::Matrix<double, GlobalDim, GlobalDim> const& Df,
87 std::vector<Eigen::RowVectorXd> const& shape_matrices,
88 std::size_t const num_dof,
89 std::array<std::size_t, 2> const& nodes_per_element)
90{
91 // Signs for the two nodes alternate.
92 constexpr auto even_odd_sign = [](std::size_t const n)
93 { return (n % 2 == 0) ? 1.0 : -1.0; };
94
95 Eigen::VectorXd local_rhs(num_dof);
96 Eigen::MatrixXd local_Jac(num_dof, num_dof);
97 local_rhs.setZero();
98 local_Jac.setZero();
99
100 for (std::size_t node_idx = 0; node_idx < 2; ++node_idx)
101 {
102 // Select the correct block according to the node index.
103 // Replicate shape matrices of the corresponding element to GlobalDim
104 // dimensions. Multiply component-wise with the force vector f.
105 local_rhs.segment(node_idx * GlobalDim * nodes_per_element[0],
106 GlobalDim * nodes_per_element[node_idx]) =
107 even_odd_sign(node_idx) *
108 shape_matrices[node_idx]
109 .transpose()
110 .replicate<GlobalDim, 1>()
111 .cwiseProduct(f.transpose()
112 .replicate(nodes_per_element[node_idx], 1)
113 .reshaped());
114 for (std::size_t node_idx_inner = 0; node_idx_inner < 2;
115 ++node_idx_inner)
116 {
117 Eigen::MatrixXd const& ones = Eigen::MatrixXd::Ones(
118 nodes_per_element[node_idx], nodes_per_element[node_idx_inner]);
119
120 // Select the correct block according to the node index.
121 // Replicate shape matrices of the corresponding element to
122 // GlobalDim dimensions. Multiply component-wise with the derivative
123 // Df of force vector f.
124 local_Jac.block(node_idx * GlobalDim * nodes_per_element[0],
125 node_idx_inner * GlobalDim * nodes_per_element[0],
126 GlobalDim * nodes_per_element[node_idx],
127 GlobalDim * nodes_per_element[node_idx_inner]) =
128 (even_odd_sign(node_idx) *
129 shape_matrices[node_idx]
130 .transpose()
131 .replicate<GlobalDim, 1>() *
132 even_odd_sign(node_idx_inner) *
133 shape_matrices[node_idx_inner].replicate<1, GlobalDim>())
134 .cwiseProduct(kroneckerProduct(Df, ones));
135 }
136 }
137 return {local_rhs, local_Jac};
138}
139
140template <int GlobalDim>
142 MeshLib::Mesh const& bulk_mesh,
143 NumLib::LocalToGlobalIndexMap const& dof_table_bulk,
144 std::size_t const source_term_mesh_id,
145 MeshLib::Mesh const& st_mesh,
146 const int variable_id,
147 ParameterLib::Parameter<double> const& parameter)
148 : dof_table_bulk_(dof_table_bulk),
149 bulk_mesh_(bulk_mesh),
150 source_term_mesh_id_(source_term_mesh_id),
151 st_mesh_(st_mesh),
152 variable_id_(variable_id),
153 component_ids_(
154 []
155 {
156 std::array<int, GlobalDim> arr{};
157 std::iota(arr.begin(), arr.end(), 0);
158 return arr;
159 }()),
160 parameter_(parameter)
161{
162 DBUG("Create EmbeddedAnchor.");
163 std::string_view const bulk_element_ids_string = "bulk_element_ids";
164 bulk_element_ids_ =
165 st_mesh_.getProperties().template getPropertyVector<std::size_t>(
166 bulk_element_ids_string);
167 assert(bulk_element_ids_ != nullptr);
168
169 std::string_view const natural_coordinates_string = "natural_coordinates";
170 natural_coordinates_ =
171 st_mesh_.getProperties().template getPropertyVector<double>(
172 natural_coordinates_string);
173 assert(natural_coordinates_ != nullptr);
174}
175
176template <int GlobalDim>
179 MeshLib::Element const* const anchor_element,
180 std::array<std::size_t, 2>& nodes_per_element,
181 std::vector<Eigen::RowVectorXd>& shape_matrices,
182 std::vector<GlobalIndexType>& global_indices,
183 Eigen::Vector<double, 2 * GlobalDim>& local_x,
184 GlobalVector const& x,
186{
187 for (auto&& [node_idx, anchor_node_id] : ranges::views::enumerate(
188 anchor_element->nodes() | MeshLib::views::ids))
189 {
190 // used pos in integrate function will be of the last node in anchor
191 // element, i,e, only one uniform stiffness is supported
192 auto const bulk_element_id = (*bulk_element_ids_)[anchor_node_id];
193 pos.setElementID(bulk_element_id);
194
195 auto const& bulk_element = *bulk_mesh_.getElement(bulk_element_id);
196 nodes_per_element[node_idx] = bulk_element.nodes().size();
197
198 std::span<double const, 3> const nat_coords(
199 &(*natural_coordinates_)[3 * anchor_node_id], 3);
200 auto N = computeShapeMatrix<GlobalDim>(bulk_element, nat_coords);
201 shape_matrices.push_back(std::move(N));
202
203 for (int component = 0; component < GlobalDim; ++component)
204 {
205 auto const& global_indices_percomponent = getIndices(
206 variable_id_, component, bulk_element_id, dof_table_bulk_);
207
208 Eigen::VectorXd const local_x_element_percomponent =
209 MathLib::toVector(x.get(global_indices_percomponent));
210
211 global_indices.insert(global_indices.cend(),
212 global_indices_percomponent.begin(),
213 global_indices_percomponent.end());
214
215 // this maps the displacement of each component of all nodes
216 // of the element to the desired anchor location in that element
217 // where all even local_x entries are for the first node and all odd
218 // local_x entries are for the second node
220 local_x_element_percomponent,
221 shape_matrices[node_idx],
222 local_x[nodeLocalIndices<GlobalDim>(node_idx)[component]]);
223 }
224 }
225}
226
227template <int GlobalDim>
229 GlobalVector& b,
230 GlobalMatrix* jac) const
231{
232 DBUG("Assemble EmbeddedAnchor.");
233
234 using GlobalDimVector = Eigen::Vector<double, GlobalDim>;
235 using GlobalDimMatrix = Eigen::Matrix<double, GlobalDim, GlobalDim>;
236
237 for (MeshLib::Element const* const anchor_element : st_mesh_.getElements())
238 {
239 std::vector<GlobalIndexType> global_indices;
240 Eigen::Vector<double, 2 * GlobalDim> local_x;
241 std::vector<Eigen::RowVectorXd> shape_matrices;
242 std::array<std::size_t, 2> nodes_per_element;
244 getShapeMatricesAndGlobalIndicesAndDisplacements(
245 anchor_element, nodes_per_element, shape_matrices, global_indices,
246 local_x, x, pos);
247
248 auto node_coords = [anchor_element](int const i)
249 { return anchor_element->getNode(i)->asEigenVector3d(); };
250 GlobalDimVector const l_original =
251 (node_coords(1) - node_coords(0)).template head<GlobalDim>();
252 double const l_original_norm = l_original.norm();
253
254 // Displacement in the two nodes.
255 auto u = [&local_x](int const i)
256 { return local_x(nodeLocalIndices<GlobalDim>(i)); };
257 GlobalDimVector const l = l_original + u(1) - u(0);
258
259 double const K = parameter_(t, pos)[0];
260 GlobalDimVector const f = l_original / l_original_norm * K *
261 (l.norm() - l_original_norm) /
262 l_original_norm;
263
264 GlobalDimMatrix const Df = l_original / l_original_norm * K *
265 l.transpose() / l.norm() / l_original_norm;
266
267 auto const& [local_rhs, local_Jac] = assembleLocalBJac<GlobalDim>(
268 f, Df, shape_matrices, global_indices.size(), nodes_per_element);
269
270 b.add(global_indices, local_rhs);
271 if (jac)
272 {
273 jac->add({global_indices, global_indices}, local_Jac);
274 }
275 }
276}
277
278template class EmbeddedAnchor<2>;
279template class EmbeddedAnchor<3>;
280
281} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
int add(IndexType row, IndexType col, double val)
Definition EigenMatrix.h:87
Global vector based on Eigen vector.
Definition EigenVector.h:26
void add(IndexType rowId, double v)
add entry
Definition EigenVector.h:77
double get(IndexType rowId) const
get entry
Definition EigenVector.h:59
virtual CellType getCellType() const =0
virtual unsigned getNumberOfNodes() const =0
virtual constexpr unsigned getDimension() const =0
Get dimension of the mesh element.
constexpr std::span< Node *const > nodes() const
Span of element's nodes, their pointers actually.
Definition Element.h:74
void setElementID(std::size_t element_id)
void getShapeMatricesAndGlobalIndicesAndDisplacements(MeshLib::Element const *const anchor_element, std::array< std::size_t, 2 > &nodes_per_element, std::vector< Eigen::RowVectorXd > &shape_matrices, std::vector< GlobalIndexType > &global_indices, Eigen::Vector< double, 2 *GlobalDim > &local_x, GlobalVector const &x, ParameterLib::SpatialPosition &pos) const
void integrate(const double t, GlobalVector const &x, GlobalVector &b, GlobalMatrix *jac) const override
EmbeddedAnchor(MeshLib::Mesh const &bulk_mesh, NumLib::LocalToGlobalIndexMap const &dof_table_bulk, std::size_t const source_term_mesh_id, MeshLib::Mesh const &st_mesh, const int variable_id, ParameterLib::Parameter< double > const &parameter)
void foreach(Function &&f)
Definition TMP.h:157
Eigen::Map< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
constexpr ranges::views::view_closure ids
For an element of a range view return its id.
Definition Mesh.h:227
std::string CellType2String(const CellType t)
Given a MeshElemType this returns the appropriate string.
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
NumLib::TemplateIsoparametric< ShapeFunction, ShapeMatricesType > createIsoparametricFiniteElement(MeshLib::Element const &e)
std::tuple< Eigen::VectorXd, Eigen::MatrixXd > assembleLocalBJac(Eigen::Vector< double, GlobalDim > const &f, Eigen::Matrix< double, GlobalDim, GlobalDim > const &Df, std::vector< Eigen::RowVectorXd > const &shape_matrices, std::size_t const num_dof, std::array< std::size_t, 2 > const &nodes_per_element)
Eigen::RowVectorXd computeShapeMatrix(MeshLib::Element const &bulk_element, std::span< double const, 3 > const nat_coords)
auto nodeLocalIndices(std::size_t const i)