Loading [MathJax]/extensions/tex2jax.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, ParameterLib::Parameter< double > const &parameter)
 
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_
 
ParameterLib::Parameter< double > const & parameter_
 
MeshLib::PropertyVector< std::size_t > const * bulk_element_ids_ = nullptr
 
MeshLib::PropertyVector< double > const * natural_coordinates_ = 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,
ParameterLib::Parameter< double > const & parameter )
explicit

Definition at line 141 of file EmbeddedAnchor.cpp.

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),
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";
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";
171 st_mesh_.getProperties().template getPropertyVector<double>(
172 natural_coordinates_string);
173 assert(natural_coordinates_ != nullptr);
174}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
Properties & getProperties()
Definition Mesh.h:136
MeshLib::Mesh const & st_mesh_
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_
ParameterLib::Parameter< double > const & parameter_

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

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

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;
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}
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

◆ bulk_element_ids_

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

Definition at line 47 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.

◆ dof_table_bulk_

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

Definition at line 40 of file EmbeddedAnchor.h.

◆ natural_coordinates_

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

Definition at line 48 of file EmbeddedAnchor.h.

◆ parameter_

template<int GlobalDim>
ParameterLib::Parameter<double> const& ProcessLib::EmbeddedAnchor< GlobalDim >::parameter_
private

Definition at line 46 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: