OGS
EmbeddedAnchor.h
Go to the documentation of this file.
1
11#pragma once
12
13#include "SourceTerm.h"
14
15namespace ProcessLib
16{
17template <int GlobalDim>
18class EmbeddedAnchor final : public SourceTermBase
19{
20public:
21 explicit EmbeddedAnchor(MeshLib::Mesh const& bulk_mesh,
22 NumLib::LocalToGlobalIndexMap const& dof_table_bulk,
23 std::size_t const source_term_mesh_id,
24 MeshLib::Mesh const& st_mesh,
25 const int variable_id);
26
28 MeshLib::Element const* const anchor_element,
29 std::array<std::size_t, 2>& nodes_per_element,
30 std::vector<Eigen::RowVectorXd>& shape_matrices,
31 std::vector<GlobalIndexType>& global_indices,
32 Eigen::Vector<double, 2 * GlobalDim>& local_x,
33 GlobalVector const& x,
35
36 void integrate(const double t, GlobalVector const& x, GlobalVector& b,
37 GlobalMatrix* jac) const override;
38
39private:
42 std::size_t const source_term_mesh_id_;
44 int const variable_id_;
45 std::array<int, GlobalDim> const component_ids_;
53};
54
55extern template class EmbeddedAnchor<2>;
56extern template class EmbeddedAnchor<3>;
57} // namespace ProcessLib
Global vector based on Eigen vector.
Definition EigenVector.h:26
MeshLib::PropertyVector< double > const * maximum_anchor_stress_
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
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_
void integrate(const double t, GlobalVector const &x, GlobalVector &b, GlobalMatrix *jac) const override
NumLib::LocalToGlobalIndexMap const & dof_table_bulk_
MeshLib::PropertyVector< double > const * natural_coordinates_
MeshLib::Mesh const & bulk_mesh_
MeshLib::PropertyVector< double > const * initial_anchor_stress_
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)