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

Detailed Description

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

Definition at line 35 of file SmallDeformationLocalAssemblerMatrix.h.

#include <SmallDeformationLocalAssemblerMatrix.h>

Inheritance diagram for ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< 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 StiffnessMatrixType = typename BMatricesType::StiffnessMatrixType
 
using NodalForceVectorType = typename BMatricesType::NodalForceVectorType
 
using NodalDisplacementVectorType
 

Public Member Functions

 SmallDeformationLocalAssemblerMatrix (SmallDeformationLocalAssemblerMatrix const &)=delete
 
 SmallDeformationLocalAssemblerMatrix (SmallDeformationLocalAssemblerMatrix &&)=delete
 
 SmallDeformationLocalAssemblerMatrix (MeshLib::Element const &e, std::size_t const local_matrix_size, 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, std::vector< double > const &local_x, std::vector< double > const &, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) 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)
 
virtual void assembleWithJacobian (double const, double const, Eigen::VectorXd const &, Eigen::VectorXd &, Eigen::MatrixXd &)
 
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.
 
- 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< 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
 

Additional Inherited Members

- Protected Member Functions inherited from ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface

Member Typedef Documentation

◆ BBarMatrixType

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

Definition at line 45 of file SmallDeformationLocalAssemblerMatrix.h.

◆ BMatricesType

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

Definition at line 44 of file SmallDeformationLocalAssemblerMatrix.h.

◆ IntegrationPointDataType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, DisplacementDim >::IntegrationPointDataType
private
Initial value:
IntegrationPointDataMatrix<ShapeMatricesType, BMatricesType,
DisplacementDim>

Definition at line 142 of file SmallDeformationLocalAssemblerMatrix.h.

◆ NodalDisplacementVectorType

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

Definition at line 49 of file SmallDeformationLocalAssemblerMatrix.h.

◆ NodalForceVectorType

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

Definition at line 48 of file SmallDeformationLocalAssemblerMatrix.h.

◆ NodalMatrixType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, DisplacementDim >::NodalMatrixType = typename ShapeMatricesType::NodalMatrixType

Definition at line 41 of file SmallDeformationLocalAssemblerMatrix.h.

◆ NodalVectorType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, DisplacementDim >::NodalVectorType = typename ShapeMatricesType::NodalVectorType

Definition at line 42 of file SmallDeformationLocalAssemblerMatrix.h.

◆ ShapeMatrices

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, DisplacementDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices

Definition at line 43 of file SmallDeformationLocalAssemblerMatrix.h.

◆ ShapeMatricesType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, DisplacementDim >::ShapeMatricesType

◆ StiffnessMatrixType

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

Definition at line 47 of file SmallDeformationLocalAssemblerMatrix.h.

Constructor & Destructor Documentation

◆ SmallDeformationLocalAssemblerMatrix() [1/3]

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

◆ SmallDeformationLocalAssemblerMatrix() [2/3]

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

◆ SmallDeformationLocalAssemblerMatrix() [3/3]

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

Definition at line 38 of file SmallDeformationLocalAssemblerMatrix-impl.h.

45 : _process_data(process_data),
46 _integration_method(integration_method),
47 _element(e),
48 _is_axially_symmetric(is_axially_symmetric)
49{
50 unsigned const n_integration_points =
52
53 _ip_data.reserve(n_integration_points);
54 _secondary_data.N.resize(n_integration_points);
55
56 auto const shape_matrices =
58 DisplacementDim>(e, is_axially_symmetric,
60
62 _process_data.solid_materials, _process_data.material_ids, e.getID());
63
64 for (unsigned ip = 0; ip < n_integration_points; ip++)
65 {
66 _ip_data.emplace_back(solid_material);
67 auto& ip_data = _ip_data[ip];
68 auto const& sm = shape_matrices[ip];
69 ip_data.N_u = sm.N;
70 ip_data.dNdx_u = sm.dNdx;
71 ip_data.integration_weight =
73 sm.integralMeasure * sm.detJ;
74
75 // Initialize current time step values
76 static const int kelvin_vector_size =
78 ip_data._sigma.setZero(kelvin_vector_size);
79 ip_data._eps.setZero(kelvin_vector_size);
80
81 // Previous time step values are not initialized and are set later.
82 ip_data._sigma_prev.resize(kelvin_vector_size);
83 ip_data._eps_prev.resize(kelvin_vector_size);
84
85 ip_data._C.resize(kelvin_vector_size, kelvin_vector_size);
86
87 _secondary_data.N[ip] = sm.N;
88 }
89}
double getWeight() const
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
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.
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N

References ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, DisplacementDim >::_integration_method, ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, DisplacementDim >::_process_data, ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, DisplacementDim >::_secondary_data, MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), MathLib::KelvinVector::kelvin_vector_dimensions(), ProcessLib::LIE::SmallDeformation::SecondaryData< ShapeMatrixType >::N, and MaterialLib::Solids::selectSolidConstitutiveRelation().

Member Function Documentation

◆ assemble()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< 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 64 of file SmallDeformationLocalAssemblerMatrix.h.

70 {
72 "SmallDeformationLocalAssembler: assembly without jacobian is not "
73 "implemented.");
74 }
#define OGS_FATAL(...)
Definition Error.h:26

References OGS_FATAL.

◆ assembleWithJacobian()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, DisplacementDim >::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 )
overridevirtual

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

Definition at line 92 of file SmallDeformationLocalAssemblerMatrix-impl.h.

98{
99 assert(_element.getDimension() == DisplacementDim);
100
101 auto const local_matrix_size = local_x.size();
102
104 local_Jac_data, local_matrix_size, local_matrix_size);
105
107 local_b_data, local_matrix_size);
108
109 unsigned const n_integration_points =
111
112 MPL::VariableArray variables;
113 MPL::VariableArray variables_prev;
115 x_position.setElementID(_element.getID());
116
117 auto const B_dil_bar = getDilatationalBBarMatrix();
118
119 for (unsigned ip = 0; ip < n_integration_points; ip++)
120 {
121 x_position.setIntegrationPoint(ip);
122 auto const& w = _ip_data[ip].integration_weight;
123
124 auto const& N = _ip_data[ip].N_u;
125 auto const& dNdx = _ip_data[ip].dNdx_u;
126 auto const x_coord =
128 _element, N);
129
131 DisplacementDim, ShapeFunction::NPOINTS, BBarMatrixType,
132 typename BMatricesType::BMatrixType>(dNdx, N, B_dil_bar, x_coord,
134
135 auto const& eps_prev = _ip_data[ip]._eps_prev;
136 auto const& sigma_prev = _ip_data[ip]._sigma_prev;
137
138 auto& eps = _ip_data[ip]._eps;
139 auto& sigma = _ip_data[ip]._sigma;
140 auto& state = _ip_data[ip]._material_state_variables;
141
142 eps.noalias() =
143 B * Eigen::Map<typename BMatricesType::NodalForceVectorType const>(
144 local_x.data(), ShapeFunction::NPOINTS * DisplacementDim);
145
146 variables.mechanical_strain
148 eps);
149 variables_prev.stress
151 sigma_prev);
152 variables_prev.mechanical_strain
154 eps_prev);
155 variables_prev.temperature = _process_data.reference_temperature;
156
157 auto&& solution = _ip_data[ip]._solid_material.integrateStress(
158 variables_prev, variables, t, x_position, dt, *state);
159
160 if (!solution)
161 {
162 OGS_FATAL("Computation of local constitutive relation failed.");
163 }
164
166 std::tie(sigma, state, C) = std::move(*solution);
167
168 local_b.noalias() -= B.transpose() * sigma * w;
169 local_Jac.noalias() += B.transpose() * C * B * w;
170 }
171}
virtual constexpr unsigned getDimension() const =0
Get dimension of the mesh element.
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
double interpolateXCoordinate(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
BMatrixType computeBMatrixPossiblyWithBbar(DNDX_Type const &dNdx, N_Type const &N, std::optional< BBarMatrixType > const &B_dil_bar, const double radius, const bool is_axially_symmetric)
Fills a B matrix, or a B bar matrix if required.

References ProcessLib::LinearBMatrix::computeBMatrixPossiblyWithBbar(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), NumLib::interpolateXCoordinate(), MaterialPropertyLib::VariableArray::mechanical_strain, OGS_FATAL, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), MaterialPropertyLib::VariableArray::stress, and MaterialPropertyLib::VariableArray::temperature.

◆ computeSecondaryVariableConcreteWithVector()

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

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

Definition at line 174 of file SmallDeformationLocalAssemblerMatrix-impl.h.

177{
178 // Compute average value per element
180 KV sigma_avg = KV::Zero();
181 auto const e_id = _element.getID();
182
183 unsigned const n_integration_points =
185 for (unsigned ip = 0; ip < n_integration_points; ip++)
186 {
187 sigma_avg += _ip_data[ip]._sigma;
188 }
189 sigma_avg /= n_integration_points;
190
191 Eigen::Map<KV>(
192 &(*_process_data.element_stresses)[e_id * KV::RowsAtCompileTime]) =
194}
Eigen::Matrix< double, 4, 1 > kelvinVectorToSymmetricTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)

References MathLib::KelvinVector::kelvinVectorToSymmetricTensor().

◆ getDilatationalBBarMatrix()

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

Definition at line 154 of file SmallDeformationLocalAssemblerMatrix.h.

155 {
156 if (!(_process_data.use_b_bar))
157 {
158 return std::nullopt;
159 }
160
162 DisplacementDim, ShapeFunction::NPOINTS, ShapeFunction,
166 }
IntegrationPointDataMatrix< ShapeMatricesType, BMatricesType, DisplacementDim > IntegrationPointDataType
BBarMatrixType computeDilatationalBbar(std::vector< IpData, Eigen::aligned_allocator< IpData > > const &ip_data, MeshLib::Element const &element, NumLib::GenericIntegrationMethod const &integration_method, const bool is_axially_symmetric)

References ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, DisplacementDim >::_element, ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, DisplacementDim >::_integration_method, ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, DisplacementDim >::_is_axially_symmetric, ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, DisplacementDim >::_process_data, and ProcessLib::LinearBMatrix::computeDilatationalBbar().

◆ getIntPtEpsilon()

template<typename ShapeFunction , int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< 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

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

Definition at line 210 of file SmallDeformationLocalAssemblerMatrix-impl.h.

216{
219}
std::vector< double > const & getIntegrationPointKelvinVectorData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)

References ProcessLib::getIntegrationPointKelvinVectorData().

◆ getIntPtFractureAperture()

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

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

Definition at line 129 of file SmallDeformationLocalAssemblerMatrix.h.

134 {
135 cache.resize(0);
136 return cache;
137 }

◆ getIntPtFractureStress()

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

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

Definition at line 119 of file SmallDeformationLocalAssemblerMatrix.h.

124 {
125 cache.resize(0);
126 return cache;
127 }

◆ getIntPtSigma()

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

◆ getShapeMatrix()

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

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 98 of file SmallDeformationLocalAssemblerMatrix.h.

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

References ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, DisplacementDim >::_secondary_data, and ProcessLib::LIE::SmallDeformation::SecondaryData< ShapeMatrixType >::N.

◆ preTimestepConcrete()

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

Member Data Documentation

◆ _element

◆ _integration_method

◆ _ip_data

◆ _is_axially_symmetric

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

◆ _process_data

◆ _secondary_data


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