OGS
ProcessLib::EmbeddedAnchor< GlobalDim > Class Template Referencefinal

Detailed Description

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

Definition at line 11 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 134 of file EmbeddedAnchor.cpp.

146 []
147 {
149 std::iota(arr.begin(), arr.end(), 0);
150 return arr;
151 }())
152{
153 DBUG("Create EmbeddedAnchor.");
154 std::string_view const bulk_element_ids_string = "bulk_element_ids";
156 st_mesh_.getProperties().template getPropertyVector<std::size_t>(
158
159 std::string_view const natural_coordinates_string = "natural_coordinates";
161 st_mesh_.getProperties().template getPropertyVector<double>(
163
165 "maximum_anchor_stress";
167 st_mesh_.getProperties().template getPropertyVector<double>(
169
171 "initial_anchor_stress";
173 st_mesh_.getProperties().template getPropertyVector<double>(
175
177 "residual_anchor_stress";
179 st_mesh_.getProperties().template getPropertyVector<double>(
181
183 "anchor_cross_sectional_area";
185 st_mesh_.getProperties().template getPropertyVector<double>(
187
188 std::string_view const anchor_stiffness_string = "anchor_stiffness";
190 st_mesh_.getProperties().template getPropertyVector<double>(
192}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
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_

References bulk_mesh_, component_ids_, dof_table_bulk_, source_term_mesh_id_, st_mesh_, and variable_id_.

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 195 of file EmbeddedAnchor.cpp.

204{
207 {
208 // used pos in integrate function will be of the last node in anchor
209 // element, i,e, only one uniform stiffness is supported
210 auto const bulk_element_id = (*bulk_element_ids_)[anchor_node_id];
211 pos.setElementID(bulk_element_id);
212
213 auto const& bulk_element = *bulk_mesh_.getElement(bulk_element_id);
214 nodes_per_element[node_idx] = bulk_element.nodes().size();
215
219 shape_matrices.push_back(std::move(N));
220
221 for (int component = 0; component < GlobalDim; ++component)
222 {
225
228
229 global_indices.insert(global_indices.cend(),
232
233 // this maps the displacement of each component of all nodes
234 // of the element to the desired anchor location in that element
235 // where all even local_x entries are for the first node and all odd
236 // local_x entries are for the second node
241 }
242 }
243}
Eigen::Map< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
Eigen::RowVectorXd computeShapeMatrix(MeshLib::Element const &bulk_element, std::span< double const, 3 > const nat_coords)

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

Referenced by integrate().

◆ 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 246 of file EmbeddedAnchor.cpp.

250{
251 DBUG("Assemble EmbeddedAnchor.");
252
255
256 for (MeshLib::Element const* const anchor_element : st_mesh_.getElements())
257 {
258 auto const anchor_element_id = anchor_element->getID();
266 local_x, x, pos);
267
268 auto node_coords = [anchor_element](int const i)
269 { return anchor_element->getNode(i)->asEigenVector3d(); };
271 (node_coords(1) - node_coords(0)).template head<GlobalDim>();
272 double const l_original_norm = l_original.norm();
273
274 // Displacement in the two nodes.
275 auto u = [&local_x](int const i)
277 GlobalDimVector const l = l_original + u(1) - u(0);
278
279 double const K = (*cross_sectional_area_)[anchor_element_id] *
280 (*anchor_stiffness_)[anchor_element_id];
281 double const initial_force =
282 (*cross_sectional_area_)[anchor_element_id] *
283 (*initial_anchor_stress_)[anchor_element_id];
284 double const max_force = (*cross_sectional_area_)[anchor_element_id] *
285 (*maximum_anchor_stress_)[anchor_element_id];
286 double const residual_force =
287 (*cross_sectional_area_)[anchor_element_id] *
288 (*residual_anchor_stress_)[anchor_element_id];
289
290 double const strain = (l.norm() - l_original_norm) / l_original_norm;
291
296 GlobalDimVector const f =
297 ((f_elastic.norm() < max_force)) ? f_elastic : f_friction;
298
301 l.transpose() / l.norm() /
303 GlobalDimMatrix const Df =
305
308
310 if (jac)
311 {
313 }
314 }
315}
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
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(), getShapeMatricesAndGlobalIndicesAndDisplacements(), ProcessLib::nodeLocalIndices(), and st_mesh_.

Member Data Documentation

◆ anchor_stiffness_

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

Definition at line 45 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 39 of file EmbeddedAnchor.h.

◆ bulk_mesh_

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

◆ component_ids_

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

Definition at line 38 of file EmbeddedAnchor.h.

Referenced by EmbeddedAnchor().

◆ cross_sectional_area_

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

Definition at line 44 of file EmbeddedAnchor.h.

◆ dof_table_bulk_

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

◆ initial_anchor_stress_

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

Definition at line 42 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 41 of file EmbeddedAnchor.h.

◆ natural_coordinates_

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

Definition at line 40 of file EmbeddedAnchor.h.

Referenced by getShapeMatricesAndGlobalIndicesAndDisplacements().

◆ residual_anchor_stress_

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

Definition at line 43 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 35 of file EmbeddedAnchor.h.

Referenced by EmbeddedAnchor().

◆ st_mesh_

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

Definition at line 36 of file EmbeddedAnchor.h.

Referenced by EmbeddedAnchor(), and integrate().

◆ variable_id_

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

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