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

Detailed Description

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

Definition at line 34 of file SmallDeformationLocalAssemblerMatrix.h.

#include <SmallDeformationLocalAssemblerMatrix.h>

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

Public Types

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

Public Member Functions

 SmallDeformationLocalAssemblerMatrix (SmallDeformationLocalAssemblerMatrix const &)=delete
 
 SmallDeformationLocalAssemblerMatrix (SmallDeformationLocalAssemblerMatrix &&)=delete
 
 SmallDeformationLocalAssemblerMatrix (MeshLib::Element const &e, std::size_t const local_matrix_size, bool is_axially_symmetric, unsigned const integration_order, 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 &, const double, const double, std::vector< double > &, std::vector< double > &, 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. More...
 
std::vector< double > const & getIntPtSigmaXX (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtSigmaYY (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtSigmaZZ (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtSigmaXY (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtSigmaXZ (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtSigmaYZ (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtEpsilonXX (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtEpsilonYY (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtEpsilonZZ (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtEpsilonXY (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtEpsilonXZ (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtEpsilonYZ (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
 
void setInitialConditions (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, bool const use_monolithic_scheme, 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_xdot, 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_xdot, const double dxdot_dx, const double dx_dx, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, 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_dot, 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, double const t, double const dt)
 
void postNonLinearSolver (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, GlobalVector const &xdot, double const t, double const dt, bool const use_monolithic_scheme, int const process_id)
 
virtual std::vector< double > interpolateNodalValuesToIntegrationPoints (std::vector< double > const &)
 
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. More...
 
- Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default
 

Private Member Functions

std::vector< double > const & getIntPtSigma (std::vector< double > &cache, std::size_t const component) const
 
std::vector< double > const & getIntPtEpsilon (std::vector< double > &cache, std::size_t const component) const
 

Private Attributes

SmallDeformationProcessData< DisplacementDim > & _process_data
 
std::vector< IntegrationPointDataMatrix< ShapeMatricesType, BMatricesType, DisplacementDim >, Eigen::aligned_allocator< IntegrationPointDataMatrix< ShapeMatricesType, BMatricesType, DisplacementDim > > > _ip_data
 
IntegrationMethod _integration_method
 
MeshLib::Element const & _element
 
bool const _is_axially_symmetric
 
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
 

Additional Inherited Members

Member Typedef Documentation

◆ BMatricesType

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

Definition at line 43 of file SmallDeformationLocalAssemblerMatrix.h.

◆ NodalDisplacementVectorType

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

Definition at line 47 of file SmallDeformationLocalAssemblerMatrix.h.

◆ NodalForceVectorType

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

Definition at line 46 of file SmallDeformationLocalAssemblerMatrix.h.

◆ NodalMatrixType

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

Definition at line 40 of file SmallDeformationLocalAssemblerMatrix.h.

◆ NodalVectorType

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

Definition at line 41 of file SmallDeformationLocalAssemblerMatrix.h.

◆ ShapeMatrices

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

Definition at line 42 of file SmallDeformationLocalAssemblerMatrix.h.

◆ ShapeMatricesType

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

Definition at line 38 of file SmallDeformationLocalAssemblerMatrix.h.

◆ StiffnessMatrixType

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

Definition at line 45 of file SmallDeformationLocalAssemblerMatrix.h.

Constructor & Destructor Documentation

◆ SmallDeformationLocalAssemblerMatrix() [1/3]

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

◆ SmallDeformationLocalAssemblerMatrix() [2/3]

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

◆ SmallDeformationLocalAssemblerMatrix() [3/3]

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

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

46  : _process_data(process_data),
47  _integration_method(integration_order),
48  _element(e),
49  _is_axially_symmetric(is_axially_symmetric)
50 {
51  unsigned const n_integration_points =
52  _integration_method.getNumberOfPoints();
53 
54  _ip_data.reserve(n_integration_points);
55  _secondary_data.N.resize(n_integration_points);
56 
57  auto const shape_matrices =
59  DisplacementDim>(e, is_axially_symmetric,
61 
63  _process_data.solid_materials, _process_data.material_ids, e.getID());
64 
65  for (unsigned ip = 0; ip < n_integration_points; ip++)
66  {
67  _ip_data.emplace_back(solid_material);
68  auto& ip_data = _ip_data[ip];
69  auto const& sm = shape_matrices[ip];
70  ip_data.N = sm.N;
71  ip_data.dNdx = sm.dNdx;
72  ip_data.integration_weight =
73  _integration_method.getWeightedPoint(ip).getWeight() *
74  sm.integralMeasure * sm.detJ;
75 
76  // Initialize current time step values
77  static const int kelvin_vector_size =
79  ip_data._sigma.setZero(kelvin_vector_size);
80  ip_data._eps.setZero(kelvin_vector_size);
81 
82  // Previous time step values are not initialized and are set later.
83  ip_data._sigma_prev.resize(kelvin_vector_size);
84  ip_data._eps_prev.resize(kelvin_vector_size);
85 
86  ip_data._C.resize(kelvin_vector_size, kelvin_vector_size);
87 
88  _secondary_data.N[ip] = sm.N;
89  }
90 }
std::vector< IntegrationPointDataMatrix< ShapeMatricesType, BMatricesType, DisplacementDim >, Eigen::aligned_allocator< IntegrationPointDataMatrix< ShapeMatricesType, BMatricesType, DisplacementDim > > > _ip_data
MechanicsBase< DisplacementDim > & selectSolidConstitutiveRelation(std::map< int, std::unique_ptr< MechanicsBase< DisplacementDim >>> 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.
Definition: KelvinVector.h:23
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
Definition: SecondaryData.h:28

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

Member Function Documentation

◆ assemble()

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

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

References OGS_FATAL.

◆ assembleWithJacobian()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
void ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::assembleWithJacobian ( double const  t,
double const  dt,
std::vector< double > const &  local_x,
std::vector< double > const &  ,
const double  ,
const double  ,
std::vector< double > &  ,
std::vector< double > &  ,
std::vector< double > &  local_b_data,
std::vector< double > &  local_Jac_data 
)
overridevirtual

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

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

104 {
105  assert(_element.getDimension() == DisplacementDim);
106 
107  auto const local_matrix_size = local_x.size();
108 
109  auto local_Jac = MathLib::createZeroedMatrix<StiffnessMatrixType>(
110  local_Jac_data, local_matrix_size, local_matrix_size);
111 
112  auto local_b = MathLib::createZeroedVector<NodalDisplacementVectorType>(
113  local_b_data, local_matrix_size);
114 
115  unsigned const n_integration_points =
116  _integration_method.getNumberOfPoints();
117 
118  MPL::VariableArray variables;
119  MPL::VariableArray variables_prev;
121  x_position.setElementID(_element.getID());
122 
123  for (unsigned ip = 0; ip < n_integration_points; ip++)
124  {
125  x_position.setIntegrationPoint(ip);
126  auto const& w = _ip_data[ip].integration_weight;
127 
128  auto const& N = _ip_data[ip].N;
129  auto const& dNdx = _ip_data[ip].dNdx;
130  auto const x_coord =
131  NumLib::interpolateXCoordinate<ShapeFunction, ShapeMatricesType>(
132  _element, N);
133  auto const B =
134  LinearBMatrix::computeBMatrix<DisplacementDim,
135  ShapeFunction::NPOINTS,
136  typename BMatricesType::BMatrixType>(
137  dNdx, N, x_coord, _is_axially_symmetric);
138 
139  auto const& eps_prev = _ip_data[ip]._eps_prev;
140  auto const& sigma_prev = _ip_data[ip]._sigma_prev;
141 
142  auto& eps = _ip_data[ip]._eps;
143  auto& sigma = _ip_data[ip]._sigma;
144  auto& state = _ip_data[ip]._material_state_variables;
145 
146  eps.noalias() =
147  B * Eigen::Map<typename BMatricesType::NodalForceVectorType const>(
148  local_x.data(), ShapeFunction::NPOINTS * DisplacementDim);
149 
150  variables[static_cast<int>(
153  eps);
154 
155  variables_prev[static_cast<int>(MaterialPropertyLib::Variable::stress)]
157  sigma_prev);
158  variables_prev[static_cast<int>(
161  eps_prev);
162  variables_prev[static_cast<int>(
164  .emplace<double>(_process_data._reference_temperature);
165 
166  auto&& solution = _ip_data[ip]._solid_material.integrateStress(
167  variables_prev, variables, t, x_position, dt, *state);
168 
169  if (!solution)
170  {
171  OGS_FATAL("Computation of local constitutive relation failed.");
172  }
173 
175  std::tie(sigma, state, C) = std::move(*solution);
176 
177  local_b.noalias() -= B.transpose() * sigma * w;
178  local_Jac.noalias() += B.transpose() * C * B * w;
179  }
180 }
virtual constexpr unsigned getDimension() const =0
Get dimension of the mesh element.
virtual std::size_t getID() const final
Returns the ID of the element.
Definition: Element.h:82
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
Definition: BMatrixPolicy.h:55
std::array< VariableType, static_cast< int >(Variable::number_of_variables)> VariableArray
Definition: VariableType.h:108
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Definition: KelvinVector.h:48
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType
Definition: KelvinVector.h:56
BMatrixType computeBMatrix(DNDX_Type const &dNdx, N_Type const &N, const double radius, const bool is_axially_symmetric)
Fills a B-matrix based on given shape function dN/dx values.
Definition: LinearBMatrix.h:42

References ProcessLib::LinearBMatrix::computeBMatrix(), MaterialPropertyLib::mechanical_strain, OGS_FATAL, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), MaterialPropertyLib::stress, and MaterialPropertyLib::temperature.

◆ computeSecondaryVariableConcreteWithVector()

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

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

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

188 {
189  // Compute average value per element
190  const int n = DisplacementDim == 2 ? 4 : 6;
191  Eigen::VectorXd ele_stress = Eigen::VectorXd::Zero(n);
192  Eigen::VectorXd ele_strain = Eigen::VectorXd::Zero(n);
193 
194  unsigned const n_integration_points =
195  _integration_method.getNumberOfPoints();
196  for (unsigned ip = 0; ip < n_integration_points; ip++)
197  {
198  auto& ip_data = _ip_data[ip];
199 
200  ele_stress += ip_data._sigma;
201  ele_strain += ip_data._eps;
202  }
203  ele_stress /= n_integration_points;
204  ele_strain /= n_integration_points;
205 
206  (*_process_data._mesh_prop_stress_xx)[_element.getID()] = ele_stress[0];
207  (*_process_data._mesh_prop_stress_yy)[_element.getID()] = ele_stress[1];
208  (*_process_data._mesh_prop_stress_zz)[_element.getID()] = ele_stress[2];
209  (*_process_data._mesh_prop_stress_xy)[_element.getID()] = ele_stress[3];
210  if (DisplacementDim == 3)
211  {
212  (*_process_data._mesh_prop_stress_yz)[_element.getID()] = ele_stress[4];
213  (*_process_data._mesh_prop_stress_xz)[_element.getID()] = ele_stress[5];
214  }
215 
216  (*_process_data._mesh_prop_strain_xx)[_element.getID()] = ele_strain[0];
217  (*_process_data._mesh_prop_strain_yy)[_element.getID()] = ele_strain[1];
218  (*_process_data._mesh_prop_strain_zz)[_element.getID()] = ele_strain[2];
219  (*_process_data._mesh_prop_strain_xy)[_element.getID()] = ele_strain[3];
220  if (DisplacementDim == 3)
221  {
222  (*_process_data._mesh_prop_strain_yz)[_element.getID()] = ele_strain[4];
223  (*_process_data._mesh_prop_strain_xz)[_element.getID()] = ele_strain[5];
224  }
225 }

◆ getIntPtEpsilon()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
std::vector<double> const& ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtEpsilon ( std::vector< double > &  cache,
std::size_t const  component 
) const
inlineprivate

Definition at line 242 of file SmallDeformationLocalAssemblerMatrix.h.

244  {
245  cache.clear();
246  cache.reserve(_ip_data.size());
247 
248  for (auto const& ip_data : _ip_data)
249  {
250  if (component < 3)
251  { // xx, yy, zz components
252  cache.push_back(ip_data._eps[component]);
253  }
254  else
255  { // mixed xy, yz, xz components
256  cache.push_back(ip_data._eps[component] / std::sqrt(2));
257  }
258  }
259 
260  return cache;
261  }

References ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::_ip_data.

Referenced by ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtEpsilonXX(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtEpsilonXY(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtEpsilonXZ(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtEpsilonYY(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtEpsilonYZ(), and ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtEpsilonZZ().

◆ getIntPtEpsilonXX()

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

◆ getIntPtEpsilonXY()

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

◆ getIntPtEpsilonXZ()

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

◆ getIntPtEpsilonYY()

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

◆ getIntPtEpsilonYZ()

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

◆ getIntPtEpsilonZZ()

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

◆ getIntPtSigma()

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
std::vector<double> const& ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtSigma ( std::vector< double > &  cache,
std::size_t const  component 
) const
inlineprivate

Definition at line 221 of file SmallDeformationLocalAssemblerMatrix.h.

223  {
224  cache.clear();
225  cache.reserve(_ip_data.size());
226 
227  for (auto const& ip_data : _ip_data)
228  {
229  if (component < 3)
230  { // xx, yy, zz components
231  cache.push_back(ip_data._sigma[component]);
232  }
233  else
234  { // mixed xy, yz, xz components
235  cache.push_back(ip_data._sigma[component] / std::sqrt(2));
236  }
237  }
238 
239  return cache;
240  }

References ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::_ip_data.

Referenced by ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtSigmaXX(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtSigmaXY(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtSigmaXZ(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtSigmaYY(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtSigmaYZ(), and ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::getIntPtSigmaZZ().

◆ getIntPtSigmaXX()

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

◆ getIntPtSigmaXY()

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

◆ getIntPtSigmaXZ()

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

◆ getIntPtSigmaYY()

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

◆ getIntPtSigmaYZ()

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

◆ getIntPtSigmaZZ()

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

◆ getShapeMatrix()

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

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 99 of file SmallDeformationLocalAssemblerMatrix.h.

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

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

◆ preTimestepConcrete()

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

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 83 of file SmallDeformationLocalAssemblerMatrix.h.

86  {
87  unsigned const n_integration_points =
88  _integration_method.getNumberOfPoints();
89 
90  for (unsigned ip = 0; ip < n_integration_points; ip++)
91  {
92  _ip_data[ip].pushBackState();
93  }
94  }

References ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::_integration_method, and ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::_ip_data.

Member Data Documentation

◆ _element

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
MeshLib::Element const& ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::_element
private

Definition at line 272 of file SmallDeformationLocalAssemblerMatrix.h.

◆ _integration_method

◆ _ip_data

◆ _is_axially_symmetric

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

Definition at line 273 of file SmallDeformationLocalAssemblerMatrix.h.

◆ _process_data

template<typename ShapeFunction , typename IntegrationMethod , int DisplacementDim>
SmallDeformationProcessData<DisplacementDim>& ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrix< ShapeFunction, IntegrationMethod, DisplacementDim >::_process_data
private

◆ _secondary_data


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