Loading [MathJax]/extensions/MathMenu.js
OGS
ProcessLib::EmbeddedAnchor< GlobalDim > Class Template Referencefinal

Detailed Description

template<int GlobalDim>
class ProcessLib::EmbeddedAnchor< GlobalDim >

Definition at line 18 of file EmbeddedAnchor.h.

#include <EmbeddedAnchor.h>

Inheritance diagram for ProcessLib::EmbeddedAnchor< GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::EmbeddedAnchor< GlobalDim >:
[legend]

Public Member Functions

 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)
 
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
 
- Public Member Functions inherited from ProcessLib::SourceTermBase
virtual ~SourceTermBase ()=default
 

Private Attributes

NumLib::LocalToGlobalIndexMap const & dof_table_bulk_
 
MeshLib::Mesh const & bulk_mesh_
 
std::size_t const source_term_mesh_id_
 
MeshLib::Mesh const & st_mesh_
 
int const variable_id_
 
std::array< int, GlobalDim > const component_ids_
 
MeshLib::PropertyVector< std::size_t > const * bulk_element_ids_ = nullptr
 
MeshLib::PropertyVector< double > const * natural_coordinates_ = nullptr
 
MeshLib::PropertyVector< double > const * maximum_anchor_stress_ = nullptr
 
MeshLib::PropertyVector< double > const * initial_anchor_stress_ = nullptr
 
MeshLib::PropertyVector< double > const * residual_anchor_stress_ = nullptr
 
MeshLib::PropertyVector< double > const * cross_sectional_area_ = nullptr
 
MeshLib::PropertyVector< double > const * anchor_stiffness_ = nullptr
 

Constructor & Destructor Documentation

◆ EmbeddedAnchor()

template<int GlobalDim>
ProcessLib::EmbeddedAnchor< GlobalDim >::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 )
explicit

Definition at line 141 of file EmbeddedAnchor.cpp.

147 : dof_table_bulk_(dof_table_bulk),
148 bulk_mesh_(bulk_mesh),
149 source_term_mesh_id_(source_term_mesh_id),
150 st_mesh_(st_mesh),
151 variable_id_(variable_id),
153 []
154 {
155 std::array<int, GlobalDim> arr{};
156 std::iota(arr.begin(), arr.end(), 0);
157 return arr;
158 }())
159{
160 DBUG("Create EmbeddedAnchor.");
161 std::string_view const bulk_element_ids_string = "bulk_element_ids";
163 st_mesh_.getProperties().template getPropertyVector<std::size_t>(
164 bulk_element_ids_string);
165
166 std::string_view const natural_coordinates_string = "natural_coordinates";
168 st_mesh_.getProperties().template getPropertyVector<double>(
169 natural_coordinates_string);
170
171 std::string_view const maximum_anchor_stress_string =
172 "maximum_anchor_stress";
174 st_mesh_.getProperties().template getPropertyVector<double>(
175 maximum_anchor_stress_string);
176
177 std::string_view const initial_anchor_stress_string =
178 "initial_anchor_stress";
180 st_mesh_.getProperties().template getPropertyVector<double>(
181 initial_anchor_stress_string);
182
183 std::string_view const residual_anchor_stress_string =
184 "residual_anchor_stress";
186 st_mesh_.getProperties().template getPropertyVector<double>(
187 residual_anchor_stress_string);
188
189 std::string_view const anchor_cross_sectional_area_string =
190 "anchor_cross_sectional_area";
192 st_mesh_.getProperties().template getPropertyVector<double>(
193 anchor_cross_sectional_area_string);
194
195 std::string_view const anchor_stiffness_string = "anchor_stiffness";
197 st_mesh_.getProperties().template getPropertyVector<double>(
198 anchor_stiffness_string);
199}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
Properties & getProperties()
Definition Mesh.h:136
MeshLib::PropertyVector< double > const * maximum_anchor_stress_
MeshLib::Mesh const & st_mesh_
MeshLib::PropertyVector< double > const * cross_sectional_area_
MeshLib::PropertyVector< double > const * anchor_stiffness_
MeshLib::PropertyVector< double > const * residual_anchor_stress_
std::array< int, GlobalDim > const component_ids_
MeshLib::PropertyVector< std::size_t > const * bulk_element_ids_
std::size_t const source_term_mesh_id_
NumLib::LocalToGlobalIndexMap const & dof_table_bulk_
MeshLib::PropertyVector< double > const * natural_coordinates_
MeshLib::Mesh const & bulk_mesh_
MeshLib::PropertyVector< double > const * initial_anchor_stress_

Member Function Documentation

◆ getShapeMatricesAndGlobalIndicesAndDisplacements()

template<int GlobalDim>
void ProcessLib::EmbeddedAnchor< GlobalDim >::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

Definition at line 202 of file EmbeddedAnchor.cpp.

211{
212 for (auto&& [node_idx, anchor_node_id] : ranges::views::enumerate(
213 anchor_element->nodes() | MeshLib::views::ids))
214 {
215 // used pos in integrate function will be of the last node in anchor
216 // element, i,e, only one uniform stiffness is supported
217 auto const bulk_element_id = (*bulk_element_ids_)[anchor_node_id];
218 pos.setElementID(bulk_element_id);
219
220 auto const& bulk_element = *bulk_mesh_.getElement(bulk_element_id);
221 nodes_per_element[node_idx] = bulk_element.nodes().size();
222
223 std::span<double const, 3> const nat_coords(
224 &(*natural_coordinates_)[3 * anchor_node_id], 3);
225 auto N = computeShapeMatrix<GlobalDim>(bulk_element, nat_coords);
226 shape_matrices.push_back(std::move(N));
227
228 for (int component = 0; component < GlobalDim; ++component)
229 {
230 auto const& global_indices_percomponent = getIndices(
231 variable_id_, component, bulk_element_id, dof_table_bulk_);
232
233 Eigen::VectorXd const local_x_element_percomponent =
234 MathLib::toVector(x.get(global_indices_percomponent));
235
236 global_indices.insert(global_indices.cend(),
237 global_indices_percomponent.begin(),
238 global_indices_percomponent.end());
239
240 // this maps the displacement of each component of all nodes
241 // of the element to the desired anchor location in that element
242 // where all even local_x entries are for the first node and all odd
243 // local_x entries are for the second node
245 local_x_element_percomponent,
246 shape_matrices[node_idx],
247 local_x[nodeLocalIndices<GlobalDim>(node_idx)[component]]);
248 }
249 }
250}
constexpr std::span< Node *const > nodes() const
Span of element's nodes, their pointers actually.
Definition Element.h:74
const Element * getElement(std::size_t idx) const
Get the element with the given index.
Definition Mesh.h:96
void setElementID(std::size_t element_id)
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
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
Eigen::RowVectorXd computeShapeMatrix(MeshLib::Element const &bulk_element, std::span< double const, 3 > const nat_coords)
auto nodeLocalIndices(std::size_t const i)

References ProcessLib::computeShapeMatrix(), MathLib::EigenVector::get(), MeshLib::views::ids, ProcessLib::nodeLocalIndices(), MeshLib::Element::nodes(), ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), and MathLib::toVector().

◆ integrate()

template<int GlobalDim>
void ProcessLib::EmbeddedAnchor< GlobalDim >::integrate ( const double t,
GlobalVector const & x,
GlobalVector & b,
GlobalMatrix * jac ) const
overridevirtual

Implements ProcessLib::SourceTermBase.

Definition at line 253 of file EmbeddedAnchor.cpp.

257{
258 DBUG("Assemble EmbeddedAnchor.");
259
260 using GlobalDimVector = Eigen::Vector<double, GlobalDim>;
261 using GlobalDimMatrix = Eigen::Matrix<double, GlobalDim, GlobalDim>;
262
263 for (MeshLib::Element const* const anchor_element : st_mesh_.getElements())
264 {
265 auto const anchor_element_id = anchor_element->getID();
266 std::vector<GlobalIndexType> global_indices;
267 Eigen::Vector<double, 2 * GlobalDim> local_x;
268 std::vector<Eigen::RowVectorXd> shape_matrices;
269 std::array<std::size_t, 2> nodes_per_element;
272 anchor_element, nodes_per_element, shape_matrices, global_indices,
273 local_x, x, pos);
274
275 auto node_coords = [anchor_element](int const i)
276 { return anchor_element->getNode(i)->asEigenVector3d(); };
277 GlobalDimVector const l_original =
278 (node_coords(1) - node_coords(0)).template head<GlobalDim>();
279 double const l_original_norm = l_original.norm();
280
281 // Displacement in the two nodes.
282 auto u = [&local_x](int const i)
283 { return local_x(nodeLocalIndices<GlobalDim>(i)); };
284 GlobalDimVector const l = l_original + u(1) - u(0);
285
286 double const K = (*cross_sectional_area_)[anchor_element_id] *
287 (*anchor_stiffness_)[anchor_element_id];
288 double const initial_force =
289 (*cross_sectional_area_)[anchor_element_id] *
290 (*initial_anchor_stress_)[anchor_element_id];
291 double const max_force = (*cross_sectional_area_)[anchor_element_id] *
292 (*maximum_anchor_stress_)[anchor_element_id];
293 double const residual_force =
294 (*cross_sectional_area_)[anchor_element_id] *
295 (*residual_anchor_stress_)[anchor_element_id];
296
297 double const strain = (l.norm() - l_original_norm) / l_original_norm;
298
299 GlobalDimVector const f_friction =
300 residual_force * l_original / l_original_norm;
301 GlobalDimVector const f_elastic =
302 l_original / l_original_norm * (initial_force + K * strain);
303 GlobalDimVector const f =
304 ((f_elastic.norm() < max_force)) ? f_elastic : f_friction;
305
306 GlobalDimMatrix const Df_friction = GlobalDimMatrix::Zero();
307 GlobalDimMatrix const Df_elastic = l_original / l_original_norm * K *
308 l.transpose() / l.norm() /
309 l_original_norm;
310 GlobalDimMatrix const Df =
311 (f_elastic.norm() < max_force) ? Df_elastic : Df_friction;
312
313 auto const& [local_rhs, local_Jac] = assembleLocalBJac<GlobalDim>(
314 f, Df, shape_matrices, global_indices.size(), nodes_per_element);
315
316 b.add(global_indices, local_rhs);
317 if (jac)
318 {
319 jac->add({global_indices, global_indices}, local_Jac);
320 }
321 }
322}
int add(IndexType row, IndexType col, double val)
Definition EigenMatrix.h:87
void add(IndexType rowId, double v)
add entry
Definition EigenVector.h:77
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:111
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
Eigen::Matrix< double, DisplacementDim, DisplacementDim, Eigen::RowMajor > GlobalDimMatrix
Definition Base.h:36
Eigen::Vector< double, DisplacementDim > GlobalDimVector
Definition Base.h:33
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)

References MathLib::EigenMatrix::add(), MathLib::EigenVector::add(), ProcessLib::assembleLocalBJac(), DBUG(), and ProcessLib::nodeLocalIndices().

Member Data Documentation

◆ anchor_stiffness_

template<int GlobalDim>
MeshLib::PropertyVector<double> const* ProcessLib::EmbeddedAnchor< GlobalDim >::anchor_stiffness_ = nullptr
private

Definition at line 52 of file EmbeddedAnchor.h.

◆ bulk_element_ids_

template<int GlobalDim>
MeshLib::PropertyVector<std::size_t> const* ProcessLib::EmbeddedAnchor< GlobalDim >::bulk_element_ids_ = nullptr
private

Definition at line 46 of file EmbeddedAnchor.h.

◆ bulk_mesh_

template<int GlobalDim>
MeshLib::Mesh const& ProcessLib::EmbeddedAnchor< GlobalDim >::bulk_mesh_
private

Definition at line 41 of file EmbeddedAnchor.h.

◆ component_ids_

template<int GlobalDim>
std::array<int, GlobalDim> const ProcessLib::EmbeddedAnchor< GlobalDim >::component_ids_
private

Definition at line 45 of file EmbeddedAnchor.h.

◆ cross_sectional_area_

template<int GlobalDim>
MeshLib::PropertyVector<double> const* ProcessLib::EmbeddedAnchor< GlobalDim >::cross_sectional_area_ = nullptr
private

Definition at line 51 of file EmbeddedAnchor.h.

◆ dof_table_bulk_

template<int GlobalDim>
NumLib::LocalToGlobalIndexMap const& ProcessLib::EmbeddedAnchor< GlobalDim >::dof_table_bulk_
private

Definition at line 40 of file EmbeddedAnchor.h.

◆ initial_anchor_stress_

template<int GlobalDim>
MeshLib::PropertyVector<double> const* ProcessLib::EmbeddedAnchor< GlobalDim >::initial_anchor_stress_ = nullptr
private

Definition at line 49 of file EmbeddedAnchor.h.

◆ maximum_anchor_stress_

template<int GlobalDim>
MeshLib::PropertyVector<double> const* ProcessLib::EmbeddedAnchor< GlobalDim >::maximum_anchor_stress_ = nullptr
private

Definition at line 48 of file EmbeddedAnchor.h.

◆ natural_coordinates_

template<int GlobalDim>
MeshLib::PropertyVector<double> const* ProcessLib::EmbeddedAnchor< GlobalDim >::natural_coordinates_ = nullptr
private

Definition at line 47 of file EmbeddedAnchor.h.

◆ residual_anchor_stress_

template<int GlobalDim>
MeshLib::PropertyVector<double> const* ProcessLib::EmbeddedAnchor< GlobalDim >::residual_anchor_stress_ = nullptr
private

Definition at line 50 of file EmbeddedAnchor.h.

◆ source_term_mesh_id_

template<int GlobalDim>
std::size_t const ProcessLib::EmbeddedAnchor< GlobalDim >::source_term_mesh_id_
private

Definition at line 42 of file EmbeddedAnchor.h.

◆ st_mesh_

template<int GlobalDim>
MeshLib::Mesh const& ProcessLib::EmbeddedAnchor< GlobalDim >::st_mesh_
private

Definition at line 43 of file EmbeddedAnchor.h.

◆ variable_id_

template<int GlobalDim>
int const ProcessLib::EmbeddedAnchor< GlobalDim >::variable_id_
private

Definition at line 44 of file EmbeddedAnchor.h.


The documentation for this class was generated from the following files: