OGS
ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim > Class Template Reference

Detailed Description

template<typename ShapeFunction, int DisplacementDim>
class ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >

Definition at line 30 of file SmallDeformationLocalAssemblerMatrixNearFracture.h.

#include <SmallDeformationLocalAssemblerMatrixNearFracture.h>

Inheritance diagram for ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >:
[legend]

Public Types

using ShapeMatricesType
using NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
using NodalVectorType = typename ShapeMatricesType::NodalVectorType
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
using BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim>
using BBarMatrixType = typename BMatricesType::BBarMatrixType
using BMatrixType = typename BMatricesType::BMatrixType
using StiffnessMatrixType = typename BMatricesType::StiffnessMatrixType
using NodalForceVectorType = typename BMatricesType::NodalForceVectorType
using NodalDisplacementVectorType

Public Member Functions

 SmallDeformationLocalAssemblerMatrixNearFracture (SmallDeformationLocalAssemblerMatrixNearFracture const &)=delete
 SmallDeformationLocalAssemblerMatrixNearFracture (SmallDeformationLocalAssemblerMatrixNearFracture &&)=delete
 SmallDeformationLocalAssemblerMatrixNearFracture (MeshLib::Element const &e, std::size_t const n_variables, std::size_t const local_matrix_size, std::vector< unsigned > const &dofIndex_to_localIndex, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, SmallDeformationProcessData< DisplacementDim > &process_data)
void assemble (double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
void assembleWithJacobian (double const t, double const dt, Eigen::VectorXd const &local_u, Eigen::VectorXd &local_b, Eigen::MatrixXd &local_J) override
void preTimestepConcrete (std::vector< double > const &, double const, double const) override
void computeSecondaryVariableConcreteWithVector (double const t, Eigen::VectorXd const &local_x) override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point.
std::vector< double > const & getIntPtSigma (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
std::vector< double > const & getIntPtEpsilon (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
std::vector< double > const & getIntPtFractureStress (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
std::vector< double > const & getIntPtFractureAperture (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
Public Member Functions inherited from ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface
 SmallDeformationLocalAssemblerInterface ()
 SmallDeformationLocalAssemblerInterface (std::size_t n_local_size, std::vector< unsigned > dofIndex_to_localIndex)
void assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x_, std::vector< double > const &, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
void computeSecondaryVariableConcrete (double const t, double const, Eigen::VectorXd const &local_x, Eigen::VectorXd const &) override
Public Member Functions inherited from ProcessLib::LocalAssemblerInterface
virtual ~LocalAssemblerInterface ()=default
virtual void setInitialConditions (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, int const process_id)
virtual void initialize (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
virtual void preAssemble (double const, double const, std::vector< double > const &)
virtual void assembleForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
virtual void assembleWithJacobianForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
virtual void computeSecondaryVariable (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
virtual void preTimestep (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)
virtual void postTimestep (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
void postNonLinearSolver (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< double > const &) const
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< std::vector< double > > const &) const
 Fits to staggered scheme.
virtual int getNumberOfVectorElementsForDeformation () const
Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default

Private Types

using IntegrationPointDataType

Private Member Functions

std::optional< BBarMatrixTypegetDilatationalBBarMatrix () const

Private Attributes

SmallDeformationProcessData< DisplacementDim > & _process_data
std::vector< FractureProperty * > _fracture_props
std::vector< JunctionProperty * > _junction_props
std::unordered_map< int, int > _fracID_to_local
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
NumLib::GenericIntegrationMethod const & _integration_method
MeshLib::Element const & _element
bool const _is_axially_symmetric
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data

Member Typedef Documentation

◆ BBarMatrixType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::BBarMatrixType = typename BMatricesType::BBarMatrixType

◆ BMatricesType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim>

◆ BMatrixType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::BMatrixType = typename BMatricesType::BMatrixType

◆ IntegrationPointDataType

◆ NodalDisplacementVectorType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::NodalDisplacementVectorType
Initial value:
VectorType< _number_of_dof > NodalForceVectorType
Rhs residual.

Definition at line 45 of file SmallDeformationLocalAssemblerMatrixNearFracture.h.

◆ NodalForceVectorType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::NodalForceVectorType = typename BMatricesType::NodalForceVectorType

◆ NodalMatrixType

◆ NodalVectorType

◆ ShapeMatrices

◆ ShapeMatricesType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::ShapeMatricesType
Initial value:
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType

Definition at line 34 of file SmallDeformationLocalAssemblerMatrixNearFracture.h.

◆ StiffnessMatrixType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::StiffnessMatrixType = typename BMatricesType::StiffnessMatrixType

Constructor & Destructor Documentation

◆ SmallDeformationLocalAssemblerMatrixNearFracture() [1/3]

◆ SmallDeformationLocalAssemblerMatrixNearFracture() [2/3]

template<typename ShapeFunction, int DisplacementDim>
ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::SmallDeformationLocalAssemblerMatrixNearFracture ( SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim > && )
delete

◆ SmallDeformationLocalAssemblerMatrixNearFracture() [3/3]

template<typename ShapeFunction, int DisplacementDim>
ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::SmallDeformationLocalAssemblerMatrixNearFracture ( MeshLib::Element const & e,
std::size_t const n_variables,
std::size_t const local_matrix_size,
std::vector< unsigned > const & dofIndex_to_localIndex,
NumLib::GenericIntegrationMethod const & integration_method,
bool const is_axially_symmetric,
SmallDeformationProcessData< DisplacementDim > & process_data )

Definition at line 40 of file SmallDeformationLocalAssemblerMatrixNearFracture-impl.h.

54 _element(e),
56{
63
64 unsigned const n_integration_points =
65 _integration_method.getNumberOfPoints();
66
69
71 _process_data.solid_materials, _process_data.material_ids, e.getID());
72
73 for (unsigned ip = 0; ip < n_integration_points; ip++)
74 {
75 _ip_data.emplace_back(solid_material);
76 auto& ip_data = _ip_data[ip];
77 auto const& sm = shape_matrices[ip];
78 ip_data.N_u = sm.N;
79 ip_data.dNdx_u = sm.dNdx;
80 ip_data.integration_weight =
81 _integration_method.getWeightedPoint(ip).getWeight() *
82 sm.integralMeasure * sm.detJ;
83
84 // Initialize current time step values
85 static const int kelvin_vector_size =
87 ip_data._sigma.setZero(kelvin_vector_size);
88 ip_data._eps.setZero(kelvin_vector_size);
89
90 // Previous time step values are not initialized and are set later.
91 ip_data._sigma_prev.resize(kelvin_vector_size);
92 ip_data._eps_prev.resize(kelvin_vector_size);
93
95
96 _secondary_data.N[ip] = sm.N;
97 }
98
99 for (auto fid : process_data.vec_ele_connected_fractureIDs[e.getID()])
100 {
101 _fracID_to_local.insert({fid, _fracture_props.size()});
102 _fracture_props.push_back(&_process_data.fracture_properties[fid]);
103 }
104
105 _junction_props = process_data.vec_ele_connected_junctionIDs[e.getID()] |
107 [&](auto const jid)
108 { return &_process_data.junction_properties[jid]; }) |
110}
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
auto & selectSolidConstitutiveRelation(SolidMaterialsMap const &constitutive_relations, MeshLib::PropertyVector< int > const *const material_ids, std::size_t const element_id)
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices

References ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::SmallDeformationLocalAssemblerInterface(), _element, _fracID_to_local, _fracture_props, _integration_method, _ip_data, _is_axially_symmetric, _junction_props, _process_data, _secondary_data, MeshLib::Element::getID(), NumLib::initShapeMatrices(), MathLib::KelvinVector::kelvin_vector_dimensions(), MaterialLib::Solids::selectSolidConstitutiveRelation(), ProcessLib::LIE::SmallDeformation::SmallDeformationProcessData< DisplacementDim >::vec_ele_connected_fractureIDs, and ProcessLib::LIE::SmallDeformation::SmallDeformationProcessData< DisplacementDim >::vec_ele_connected_junctionIDs.

Member Function Documentation

◆ assemble()

template<typename ShapeFunction, int DisplacementDim>
void ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::assemble ( double const ,
double const ,
std::vector< double > const & ,
std::vector< double > const & ,
std::vector< double > & ,
std::vector< double > & ,
std::vector< double > &  )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 62 of file SmallDeformationLocalAssemblerMatrixNearFracture.h.

68 {
70 "SmallDeformationLocalAssembler: assembly without jacobian is not "
71 "implemented.");
72 }
#define OGS_FATAL(...)
Definition Error.h:19

References OGS_FATAL.

◆ assembleWithJacobian()

template<typename ShapeFunction, int DisplacementDim>
void ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::assembleWithJacobian ( double const t,
double const dt,
Eigen::VectorXd const & local_u,
Eigen::VectorXd & local_b,
Eigen::MatrixXd & local_J )
overridevirtual

Reimplemented from ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface.

Definition at line 115 of file SmallDeformationLocalAssemblerMatrixNearFracture-impl.h.

119{
120 assert(_element.getDimension() == DisplacementDim);
121
123 auto const n_fractures = _fracture_props.size();
124 auto const n_junctions = _junction_props.size();
125 auto const n_enrich_var = n_fractures + n_junctions;
126
127 using BlockVectorType =
129 using BlockMatrixType =
131
132 //--------------------------------------------------------------------------------------
133 // prepare sub vectors, matrices for regular displacement (u) and
134 // displacement jumps (g)
135 //
136 // example with two fractures with one intersection:
137 // |b(u)|
138 // b = |b(g1)|
139 // |b(g2)|
140 // |b(j1)|
141 //
142 // |J(u,u) J(u,g1) J(u,g2) J(u,j1) |
143 // J = |J(g1,u) J(g1,g1) J(g1,g2) J(g1,j1)|
144 // |J(g2,u) J(g2,g1) J(g2,g2) J(g2,j1)|
145 // |J(j1,u) J(j1,g1) J(j1,g2) J(j1,j1)|
146 //--------------------------------------------------------------------------------------
147 auto local_b_u = local_b.segment<N_DOF_PER_VAR>(0);
149 for (unsigned i = 0; i < n_enrich_var; i++)
150 {
151 vec_local_b_g.push_back(
152 local_b.segment<N_DOF_PER_VAR>(N_DOF_PER_VAR * (i + 1)));
153 }
154
155 auto local_J_uu = local_J.block<N_DOF_PER_VAR, N_DOF_PER_VAR>(0, 0);
159 for (unsigned i = 0; i < n_enrich_var; i++)
160 {
162 0, N_DOF_PER_VAR * (i + 1));
163 vec_local_J_ug.push_back(sub_ug);
164
166 N_DOF_PER_VAR * (i + 1), 0);
167 vec_local_J_gu.push_back(sub_gu);
168
169 for (unsigned j = 0; j < n_enrich_var; j++)
170 {
172 N_DOF_PER_VAR * (i + 1), N_DOF_PER_VAR * (j + 1));
173 vec_local_J_gg[i].push_back(sub_gg);
174 }
175 }
176
177 auto const nodal_u = local_u.segment<N_DOF_PER_VAR>(0);
179 for (unsigned i = 0; i < n_enrich_var; i++)
180 {
182 N_DOF_PER_VAR * (i + 1));
183 vec_nodal_g.push_back(sub);
184 }
185
186 //------------------------------------------------
187 // integration
188 //------------------------------------------------
189 unsigned const n_integration_points =
190 _integration_method.getNumberOfPoints();
191
194
196
197 for (unsigned ip = 0; ip < n_integration_points; ip++)
198 {
199 auto& ip_data = _ip_data[ip];
200 auto const& w = _ip_data[ip].integration_weight;
201
202 auto const& N = ip_data.N_u;
203 auto const& dNdx = ip_data.dNdx_u;
204
205 // levelset functions
211
212 // u = u^hat + sum_i(enrich^br_i(x) * [u]_i) + sum_i(enrich^junc_i(x) *
213 // [u]_i)
215 for (unsigned i = 0; i < n_enrich_var; i++)
216 {
218 }
219
221 std::nullopt, this->_element.getID(),
224 this->_element, N))};
225
226 auto const x_coord =
227 x_position.getCoordinates().value()[0]; // r for axisymmetry
232
233 // strain, stress
234 auto const& eps_prev = ip_data._eps_prev;
235 auto const& sigma_prev = ip_data._sigma_prev;
236
237 auto& eps = ip_data._eps;
238 auto& sigma = ip_data._sigma;
239 auto& state = ip_data._material_state_variables;
240
241 eps.noalias() = B * nodal_total_u;
242
243 variables.mechanical_strain
245 eps);
246
247 variables_prev.stress
249 sigma_prev);
250 variables_prev.mechanical_strain
252 eps_prev);
253 variables_prev.temperature = _process_data.reference_temperature;
254
255 auto&& solution = _ip_data[ip]._solid_material.integrateStress(
257
258 if (!solution)
259 {
260 OGS_FATAL("Computation of local constitutive relation failed.");
261 }
262
265
266 // r_u = B^T * Sigma = B^T * C * B * (u+phi*[u])
267 // r_[u] = (phi*B)^T * Sigma = (phi*B)^T * C * B * (u+phi*[u])
268 local_b_u.noalias() -= B.transpose() * sigma * w;
269 for (unsigned i = 0; i < n_enrich_var; i++)
270 {
271 vec_local_b_g[i].noalias() -=
272 levelsets[i] * B.transpose() * sigma * w;
273 }
274
275 // J_uu += B^T * C * B
276 local_J_uu.noalias() += B.transpose() * C * B * w;
277
278 for (unsigned i = 0; i < n_enrich_var; i++)
279 {
280 // J_u[u] += B^T * C * (levelset * B)
281 vec_local_J_ug[i].noalias() +=
282 B.transpose() * C * (levelsets[i] * B) * w;
283
284 // J_[u]u += (levelset * B)^T * C * B
285 vec_local_J_gu[i].noalias() +=
286 (levelsets[i] * B.transpose()) * C * B * w;
287
288 for (unsigned j = 0; j < n_enrich_var; j++)
289 {
290 // J_[u][u] += (levelset * B)^T * C * (levelset * B)
291 vec_local_J_gg[i][j].noalias() +=
292 (levelsets[i] * B.transpose()) * C * (levelsets[j] * B) * w;
293 }
294 }
295 }
296}
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
Eigen::Vector3d computePhysicalCoordinates(MeshLib::Element const &e, Eigen::MatrixBase< Derived > const &shape)
std::vector< double > uGlobalEnrichments(std::vector< FractureProperty * > const &frac_props, std::vector< JunctionProperty * > const &junction_props, std::unordered_map< int, int > const &fracID_to_local, Eigen::Vector3d const &x)

References _element, _fracID_to_local, _fracture_props, _integration_method, _ip_data, _is_axially_symmetric, _junction_props, _process_data, assembleWithJacobian(), ProcessLib::LinearBMatrix::computeBMatrixPossiblyWithBbar(), ProcessLib::LIE::computePhysicalCoordinates(), ParameterLib::SpatialPosition::getCoordinates(), getDilatationalBBarMatrix(), NumLib::interpolateCoordinates(), MaterialPropertyLib::VariableArray::mechanical_strain, OGS_FATAL, MaterialPropertyLib::VariableArray::stress, MaterialPropertyLib::VariableArray::temperature, and ProcessLib::LIE::uGlobalEnrichments().

Referenced by assembleWithJacobian().

◆ computeSecondaryVariableConcreteWithVector()

template<typename ShapeFunction, int DisplacementDim>
void ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::computeSecondaryVariableConcreteWithVector ( double const t,
Eigen::VectorXd const & local_x )
overridevirtual

Implements ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface.

Definition at line 303 of file SmallDeformationLocalAssemblerMatrixNearFracture-impl.h.

306{
307 // Compute average value per element
310 auto const e_id = _element.getID();
311
312 unsigned const n_integration_points =
313 _integration_method.getNumberOfPoints();
314 for (unsigned ip = 0; ip < n_integration_points; ip++)
315 {
316 sigma_avg += _ip_data[ip]._sigma;
317 }
319
321 &(*_process_data.element_stresses)[e_id * KV::RowsAtCompileTime]) =
323}
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)

References _element, _integration_method, _ip_data, _process_data, and MathLib::KelvinVector::kelvinVectorToSymmetricTensor().

◆ getDilatationalBBarMatrix()

template<typename ShapeFunction, int DisplacementDim>
std::optional< BBarMatrixType > ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::getDilatationalBBarMatrix ( ) const
inlineprivate

◆ getIntPtEpsilon()

template<typename ShapeFunction, int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::getIntPtEpsilon ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
overridevirtual

◆ getIntPtFractureAperture()

template<typename ShapeFunction, int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::getIntPtFractureAperture ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverridevirtual

◆ getIntPtFractureStress()

template<typename ShapeFunction, int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::getIntPtFractureStress ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverridevirtual

◆ getIntPtSigma()

◆ getShapeMatrix()

template<typename ShapeFunction, int DisplacementDim>
Eigen::Map< const Eigen::RowVectorXd > ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::getShapeMatrix ( const unsigned integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 95 of file SmallDeformationLocalAssemblerMatrixNearFracture.h.

97 {
98 auto const& N = _secondary_data.N[integration_point];
99
100 // assumes N is stored contiguously in memory
101 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
102 }

References _secondary_data.

◆ preTimestepConcrete()

template<typename ShapeFunction, int DisplacementDim>
void ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::preTimestepConcrete ( std::vector< double > const & ,
double const ,
double const  )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 79 of file SmallDeformationLocalAssemblerMatrixNearFracture.h.

82 {
83 unsigned const n_integration_points =
84 _integration_method.getNumberOfPoints();
85
86 for (unsigned ip = 0; ip < n_integration_points; ip++)
87 {
88 _ip_data[ip].pushBackState();
89 }
90 }

References _integration_method, and _ip_data.

Member Data Documentation

◆ _element

◆ _fracID_to_local

template<typename ShapeFunction, int DisplacementDim>
std::unordered_map<int, int> ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::_fracID_to_local
private

◆ _fracture_props

◆ _integration_method

◆ _ip_data

◆ _is_axially_symmetric

template<typename ShapeFunction, int DisplacementDim>
bool const ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::_is_axially_symmetric
private

◆ _junction_props

◆ _process_data

◆ _secondary_data

template<typename ShapeFunction, int DisplacementDim>
SecondaryData<typename ShapeMatrices::ShapeType> ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::_secondary_data
private

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