OGS
ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim > Class Template Reference

Detailed Description

template<typename ShapeFunction, int DisplacementDim>
class ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >

Definition at line 81 of file LargeDeformationFEM.h.

#include <LargeDeformationFEM.h>

Inheritance diagram for ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >:
[legend]

Public Types

using ShapeMatricesType
using NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
using GlobalDimNodalMatrixType
using NodalVectorType = typename ShapeMatricesType::NodalVectorType
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
using BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim>
using BMatrixType = typename BMatricesType::BMatrixType
using StiffnessMatrixType = typename BMatricesType::StiffnessMatrixType
using NodalForceVectorType = typename BMatricesType::NodalForceVectorType
using NodalDisplacementVectorType
using GMatricesType = GMatrixPolicyType<ShapeFunction, DisplacementDim>
using GradientVectorType = typename GMatricesType::GradientVectorType
using GradientMatrixType = typename GMatricesType::GradientMatrixType
using VectorTypeForFbar = typename BMatricesType::NodalForceVectorType
using IpData

Public Member Functions

 LargeDeformationLocalAssembler (LargeDeformationLocalAssembler const &)=delete
 LargeDeformationLocalAssembler (LargeDeformationLocalAssembler &&)=delete
 LargeDeformationLocalAssembler (MeshLib::Element const &e, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, LargeDeformationProcessData< DisplacementDim > &process_data)
void initializeConcrete () override
double computeOutputStrainData (double const detF0, BMatrixType const &B, GradientVectorType const &grad_u, Eigen::Ref< Eigen::VectorXd const > const &u, typename ConstitutiveRelations::OutputData< DisplacementDim > &output_data) const
ConstitutiveRelations::ConstitutiveData< DisplacementDim > updateConstitutiveRelations (double const alpha, GradientMatrixType const &G, Eigen::Ref< Eigen::VectorXd const > const &u_prev, ParameterLib::SpatialPosition const &x_position, double const t, double const dt, typename ConstitutiveRelations::ConstitutiveSetting< DisplacementDim > &CS, MaterialPropertyLib::Medium const &medium, typename ConstitutiveRelations::StatefulData< DisplacementDim > &current_state, typename ConstitutiveRelations::StatefulDataPrev< DisplacementDim > const &prev_state, MaterialStateData< DisplacementDim > &material_state, typename ConstitutiveRelations::OutputData< DisplacementDim > &output_data) const
void assemble (double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
std::optional< std::tuple< double, VectorTypeForFbar > > computeFBarVariables (bool const compute_detF0_only, Eigen::VectorXd const &local_x) const
void assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
void postTimestepConcrete (Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, double const t, double const dt, int const) override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point.
Public Member Functions inherited from ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >
 LargeDeformationLocalAssemblerInterface (MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, LargeDeformationProcessData< DisplacementDim > &process_data)
std::size_t setIPDataInitialConditions (std::string_view name, double const *values, int const integration_order)
 Returns number of read integration points.
unsigned getNumberOfIntegrationPoints () const
int getMaterialID () const
std::vector< double > getMaterialStateVariableInternalState (std::function< std::span< double >(typename MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables &)> const &get_values_span, int const &n_components) const
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt (unsigned integration_point) const
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 std::optional< VectorSegmentgetVectorDeformationSegment () const
Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default

Static Public Attributes

static constexpr auto & N_u_op

Static Private Member Functions

static constexpr auto localDOF (std::vector< double > const &x)

Private Attributes

std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data

Static Private Attributes

static const int displacement_size

Additional Inherited Members

Static Public Member Functions inherited from ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >
static auto getReflectionDataForOutput ()
Protected Attributes inherited from ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >
LargeDeformationProcessData< DisplacementDim > & process_data_
std::vector< MaterialStateData< DisplacementDim > > material_states_
std::vector< typename ConstitutiveRelations::StatefulData< DisplacementDim > > current_states_
std::vector< typename ConstitutiveRelations::StatefulDataPrev< DisplacementDim > > prev_states_
std::vector< typename ConstitutiveRelations::OutputData< DisplacementDim > > output_data_
NumLib::GenericIntegrationMethod const & integration_method_
MeshLib::Element const & element_
bool const is_axially_symmetric_
ConstitutiveRelations::SolidConstitutiveRelation< DisplacementDim > const & solid_material_

Member Typedef Documentation

◆ BMatricesType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim>

Definition at line 92 of file LargeDeformationFEM.h.

◆ BMatrixType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::BMatrixType = typename BMatricesType::BMatrixType

Definition at line 94 of file LargeDeformationFEM.h.

◆ GlobalDimNodalMatrixType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::GlobalDimNodalMatrixType
Initial value:
MatrixType< GlobalDim, ShapeFunction::NPOINTS > GlobalDimNodalMatrixType

Definition at line 88 of file LargeDeformationFEM.h.

◆ GMatricesType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::GMatricesType = GMatrixPolicyType<ShapeFunction, DisplacementDim>

Definition at line 100 of file LargeDeformationFEM.h.

◆ GradientMatrixType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::GradientMatrixType = typename GMatricesType::GradientMatrixType

Definition at line 102 of file LargeDeformationFEM.h.

◆ GradientVectorType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::GradientVectorType = typename GMatricesType::GradientVectorType

Definition at line 101 of file LargeDeformationFEM.h.

◆ IpData

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::IpData

◆ NodalDisplacementVectorType

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

Definition at line 97 of file LargeDeformationFEM.h.

◆ NodalForceVectorType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::NodalForceVectorType = typename BMatricesType::NodalForceVectorType

Definition at line 96 of file LargeDeformationFEM.h.

◆ NodalMatrixType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::NodalMatrixType = typename ShapeMatricesType::NodalMatrixType

Definition at line 87 of file LargeDeformationFEM.h.

◆ NodalVectorType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::NodalVectorType = typename ShapeMatricesType::NodalVectorType

Definition at line 90 of file LargeDeformationFEM.h.

◆ ShapeMatrices

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices

Definition at line 91 of file LargeDeformationFEM.h.

◆ ShapeMatricesType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::ShapeMatricesType
Initial value:
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType

Definition at line 85 of file LargeDeformationFEM.h.

◆ StiffnessMatrixType

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::StiffnessMatrixType = typename BMatricesType::StiffnessMatrixType

Definition at line 95 of file LargeDeformationFEM.h.

◆ VectorTypeForFbar

template<typename ShapeFunction, int DisplacementDim>
using ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::VectorTypeForFbar = typename BMatricesType::NodalForceVectorType

Definition at line 104 of file LargeDeformationFEM.h.

Constructor & Destructor Documentation

◆ LargeDeformationLocalAssembler() [1/3]

template<typename ShapeFunction, int DisplacementDim>
ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::LargeDeformationLocalAssembler ( LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim > const & )
delete

◆ LargeDeformationLocalAssembler() [2/3]

template<typename ShapeFunction, int DisplacementDim>
ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::LargeDeformationLocalAssembler ( LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim > && )
delete

◆ LargeDeformationLocalAssembler() [3/3]

template<typename ShapeFunction, int DisplacementDim>
ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::LargeDeformationLocalAssembler ( MeshLib::Element const & e,
std::size_t const ,
NumLib::GenericIntegrationMethod const & integration_method,
bool const is_axially_symmetric,
LargeDeformationProcessData< DisplacementDim > & process_data )
inline

Definition at line 116 of file LargeDeformationFEM.h.

124 {
125 unsigned const n_integration_points =
126 this->integration_method_.getNumberOfPoints();
127
130
131 auto const shape_matrices =
135
136 for (unsigned ip = 0; ip < n_integration_points; ip++)
137 {
138 auto& ip_data = _ip_data[ip];
139 auto const& sm = shape_matrices[ip];
140 _ip_data[ip].integration_weight =
141 this->integration_method_.getWeightedPoint(ip).getWeight() *
142 sm.integralMeasure * sm.detJ;
143
144 ip_data.N_u = sm.N;
145 ip_data.dNdx_u = sm.dNdx;
146
148 }
149 }
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
ShapeMatrixPolicyType< ShapeFunction, DisplacementDim > ShapeMatricesType
LargeDeformationLocalAssemblerInterface(MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, LargeDeformationProcessData< DisplacementDim > &process_data)

References ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::LargeDeformationLocalAssemblerInterface(), _ip_data, _secondary_data, NumLib::initShapeMatrices(), and ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_.

Member Function Documentation

◆ assemble()

template<typename ShapeFunction, int DisplacementDim>
void ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< 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 278 of file LargeDeformationFEM.h.

284 {
285 OGS_FATAL(
286 "LargeDeformationLocalAssembler: assembly without jacobian is not "
287 "implemented.");
288 }
#define OGS_FATAL(...)
Definition Error.h:26

References OGS_FATAL.

◆ assembleWithJacobian()

template<typename ShapeFunction, int DisplacementDim>
void ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian ( double const t,
double const dt,
std::vector< double > const & local_x,
std::vector< double > const & local_x_prev,
std::vector< double > & local_b_data,
std::vector< double > & local_Jac_data )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 317 of file LargeDeformationFEM.h.

322 {
323 auto const local_matrix_size = local_x.size();
324
327
330
331 auto [u] = localDOF(local_x);
332 auto [u_prev] = localDOF(local_x_prev);
333
334 unsigned const n_integration_points =
335 this->integration_method_.getNumberOfPoints();
336
339 auto const& medium =
340 *this->process_data_.media_map.getMedium(this->element_.getID());
341
342 auto const f_bar_variables =
343 computeFBarVariables(/*compute_detF0_only?*/ false, u);
344
345 for (unsigned ip = 0; ip < n_integration_points; ip++)
346 {
347 auto const& w = _ip_data[ip].integration_weight;
348 auto const& N = _ip_data[ip].N_u;
349 auto const& dNdx = _ip_data[ip].dNdx_u;
350
352 std::nullopt, this->element_.getID(),
356 this->element_, N))};
357
358 auto const x_coord = x_position.getCoordinates().value()[0];
359
360 // For the 2D case the 33-component is needed (and the four entries
361 // of the non-symmetric matrix); In 3d there are nine entries.
363 (DisplacementDim == 2 ? 1 : 0),
368
369 GradientVectorType const grad_u = G * u;
370
375
380
381 double const detF0 =
383 double const alpha = computeOutputStrainData(
384 detF0, B_0 + 0.5 * B_N, grad_u, u, this->output_data_[ip]);
385
386 auto const CD = updateConstitutiveRelations(
388 medium, this->current_states_[ip], this->prev_states_[ip],
389 this->material_states_[ip], this->output_data_[ip]);
390
391 BMatrixType B = B_0 + B_N;
392
393 auto const& sigma = this->current_states_[ip].stress_data.sigma;
394 auto const& b = *CD.volumetric_body_force;
395 auto const& C = CD.s_mech_data.stiffness_tensor;
396
397 auto const sigma_geom =
400
401 if (!f_bar_variables)
402 {
403 local_b.noalias() -=
404 (B.transpose() * sigma - N_u_op(N).transpose() * b) * w;
405
406 local_Jac.noalias() += G.transpose() * sigma_geom * G * w;
407
408 local_Jac.noalias() += B.transpose() * C * B * w;
409 }
410 else // with F bar
411 {
412 double const alpha_p2 = alpha * alpha;
413 local_Jac.noalias() +=
414 alpha_p2 * G.transpose() * sigma_geom * G * w;
415
417
418 auto const& F =
419 this->output_data_[ip]
420 .deformation_gradient_data.deformation_gradient;
426
427 auto const q0_q = q0 - q;
428
429 auto const S_q = sigma * (q0_q).transpose();
430
431 local_Jac.noalias() +=
432 2.0 * alpha_p2 * S_q.transpose() * B * w / DisplacementDim;
433
434 auto const twoEplsI =
436 alpha_p2, this->output_data_[ip].eps_data.eps,
438
439 // B bar:
440 BMatrixType const Bbar =
441 alpha_p2 *
442 (B + twoEplsI * (q0_q).transpose() / DisplacementDim);
443
444 local_Jac.noalias() +=
445 2.0 * Bbar.transpose() * S_q * w / DisplacementDim;
446
447 local_Jac.noalias() += Bbar.transpose() * C * Bbar * w;
448
449 local_Jac.noalias() -=
450 alpha_p2 * sigma.dot(twoEplsI) *
451 (q0 * q0.transpose() - q * q.transpose()) * w /
453
454 local_b.noalias() -=
455 (Bbar.transpose() * sigma - N_u_op(N).transpose() * b) * w;
456 }
457 }
458 }
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
std::optional< std::tuple< double, VectorTypeForFbar > > computeFBarVariables(bool const compute_detF0_only, Eigen::VectorXd const &local_x) const
typename GMatricesType::GradientVectorType GradientVectorType
double computeOutputStrainData(double const detF0, BMatrixType const &B, GradientVectorType const &grad_u, Eigen::Ref< Eigen::VectorXd const > const &u, typename ConstitutiveRelations::OutputData< DisplacementDim > &output_data) const
static constexpr auto localDOF(std::vector< double > const &x)
typename BMatricesType::NodalForceVectorType VectorTypeForFbar
ConstitutiveRelations::ConstitutiveData< DisplacementDim > updateConstitutiveRelations(double const alpha, GradientMatrixType const &G, Eigen::Ref< Eigen::VectorXd const > const &u_prev, ParameterLib::SpatialPosition const &x_position, double const t, double const dt, typename ConstitutiveRelations::ConstitutiveSetting< DisplacementDim > &CS, MaterialPropertyLib::Medium const &medium, typename ConstitutiveRelations::StatefulData< DisplacementDim > &current_state, typename ConstitutiveRelations::StatefulDataPrev< DisplacementDim > const &prev_state, MaterialStateData< DisplacementDim > &material_state, typename ConstitutiveRelations::OutputData< DisplacementDim > &output_data) const
typename GMatricesType::GradientMatrixType GradientMatrixType
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)
Eigen::Matrix< double, MPL::tensorSize(DisplacementDim), MPL::tensorSize(DisplacementDim)> computeSigmaGeom(Eigen::Matrix3d const &sigma_tensor)
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > compute2EPlusI(double const alpha_p2, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps_bar, bool const is_axially_symmetric)
std::vector< MaterialStateData< DisplacementDim > > material_states_
std::vector< typename ConstitutiveRelations::StatefulDataPrev< DisplacementDim > > prev_states_
std::vector< typename ConstitutiveRelations::StatefulData< DisplacementDim > > current_states_
std::vector< typename ConstitutiveRelations::OutputData< DisplacementDim > > output_data_

References _ip_data, ProcessLib::NonLinearFbar::compute2EPlusI(), ProcessLib::LinearBMatrix::computeBMatrix(), ProcessLib::NonLinearBMatrix::computeBMatrix(), computeFBarVariables(), ProcessLib::Deformation::computeGMatrix(), computeOutputStrainData(), ProcessLib::NonLinearFbar::computeQVector(), ProcessLib::LargeDeformation::computeSigmaGeom(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::element_, ParameterLib::SpatialPosition::getCoordinates(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateCoordinates(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, MathLib::KelvinVector::kelvinVectorToTensor(), localDOF(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::material_states_, N_u_op, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::output_data_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::process_data_, and updateConstitutiveRelations().

◆ computeFBarVariables()

template<typename ShapeFunction, int DisplacementDim>
std::optional< std::tuple< double, VectorTypeForFbar > > ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::computeFBarVariables ( bool const compute_detF0_only,
Eigen::VectorXd const & local_x ) const
inline

Definition at line 290 of file LargeDeformationFEM.h.

References _ip_data, ProcessLib::NonLinearFbar::computeFBarInitialVariablesAverage(), ProcessLib::NonLinearFbar::computeFBarInitialVariablesElementCenter(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::element_, ProcessLib::NonLinearFbar::ELEMENT_AVERAGE, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, ProcessLib::NonLinearFbar::NONE, and ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::process_data_.

Referenced by assembleWithJacobian(), and postTimestepConcrete().

◆ computeOutputStrainData()

template<typename ShapeFunction, int DisplacementDim>
double ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::computeOutputStrainData ( double const detF0,
BMatrixType const & B,
GradientVectorType const & grad_u,
Eigen::Ref< Eigen::VectorXd const > const & u,
typename ConstitutiveRelations::OutputData< DisplacementDim > & output_data ) const
inline

Definition at line 187 of file LargeDeformationFEM.h.

195 {
196 // Note: Here B=B_{linear}+B_{nonlinear}/2 For Green-Lagrange strain.
197 output_data.eps_data.eps = B * u;
198 output_data.deformation_gradient_data.deformation_gradient =
200 output_data.deformation_gradient_data.volume_ratio =
202 output_data.deformation_gradient_data.deformation_gradient);
203
204 if (this->process_data_.bar_det_f_type ==
206 {
207 return 1.0;
208 }
209
210 double const detF = output_data.deformation_gradient_data.volume_ratio;
211 double const J_ratio = detF0 / detF;
212
213 if (J_ratio < .0)
214 {
215 OGS_FATAL(
216 "det(F0)/det(F) is negative with det(F0) = {:g}, and det(F) = "
217 "{:g} ",
218 detF0, detF);
219 }
220
221 double const alpha =
223
224 output_data.deformation_gradient_data.deformation_gradient
226 output_data.deformation_gradient_data.volume_ratio *=
228
229 double const alpha_p2 = alpha * alpha;
230 output_data.eps_data.eps =
231 alpha_p2 * output_data.eps_data.eps +
232 0.5 * (alpha_p2 - 1) *
235 return alpha;
236 }
double determinant(Eigen::MatrixBase< Derived > const &tensor)
Computes determinant of a vectorized tensor.

References ProcessLib::LargeDeformation::ConstitutiveRelations::OutputData< DisplacementDim >::deformation_gradient_data, MathLib::VectorizedTensor::determinant(), ProcessLib::LargeDeformation::ConstitutiveRelations::OutputData< DisplacementDim >::eps_data, MathLib::VectorizedTensor::identity(), ProcessLib::NonLinearFbar::identityForF(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, ProcessLib::NonLinearFbar::NONE, OGS_FATAL, and ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::process_data_.

Referenced by assembleWithJacobian(), and postTimestepConcrete().

◆ getShapeMatrix()

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

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 533 of file LargeDeformationFEM.h.

535 {
536 auto const& N = _secondary_data.N[integration_point];
537
538 // assumes N is stored contiguously in memory
539 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
540 }

References _secondary_data.

◆ initializeConcrete()

template<typename ShapeFunction, int DisplacementDim>
void ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::initializeConcrete ( )
inlineoverridevirtual

Set initial stress from parameter.

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 151 of file LargeDeformationFEM.h.

152 {
153 unsigned const n_integration_points =
154 this->integration_method_.getNumberOfPoints();
155 for (unsigned ip = 0; ip < n_integration_points; ip++)
156 {
157 auto const& ip_data = _ip_data[ip];
158
160 std::nullopt, this->element_.getID(),
164 this->element_, ip_data.N_u))};
165
167 if (this->process_data_.initial_stress != nullptr)
168 {
169 this->current_states_[ip].stress_data.sigma.noalias() =
171 DisplacementDim>((*this->process_data_.initial_stress)(
173 double>::quiet_NaN() /* time independent */,
174 x_position));
175 }
176
177 double const t = 0; // TODO (naumov) pass t from top
178 auto& material_state = this->material_states_[ip];
179 this->solid_material_.initializeInternalStateVariables(
180 t, x_position, *material_state.material_state_variables);
181
182 material_state.pushBackState();
183 this->prev_states_[ip] = this->current_states_[ip];
184 }
185 }
ConstitutiveRelations::SolidConstitutiveRelation< DisplacementDim > const & solid_material_

References _ip_data, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::element_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateCoordinates(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::material_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::process_data_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::solid_material_, and MathLib::KelvinVector::symmetricTensorToKelvinVector().

◆ localDOF()

template<typename ShapeFunction, int DisplacementDim>
constexpr auto ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::localDOF ( std::vector< double > const & x)
inlinestaticconstexprprivate

Definition at line 543 of file LargeDeformationFEM.h.

References NumLib::localDOF().

Referenced by assembleWithJacobian().

◆ postTimestepConcrete()

template<typename ShapeFunction, int DisplacementDim>
void ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::postTimestepConcrete ( Eigen::VectorXd const & local_x,
Eigen::VectorXd const & local_x_prev,
double const t,
double const dt,
int const  )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 460 of file LargeDeformationFEM.h.

464 {
465 unsigned const n_integration_points =
466 this->integration_method_.getNumberOfPoints();
467
470
471 auto& medium =
472 *this->process_data_.media_map.getMedium(this->element_.getID());
473
474 auto const f_bar_variables =
475 computeFBarVariables(/*compute_detF0_only?*/ true, local_x);
476
477 for (unsigned ip = 0; ip < n_integration_points; ip++)
478 {
479 auto const& N = _ip_data[ip].N_u;
480 auto const& dNdx = _ip_data[ip].dNdx_u;
481
483 std::nullopt, this->element_.getID(),
487 this->element_, N))};
488
489 auto const x_coord = x_position.getCoordinates().value()[0];
490
491 // For the 2D case the 33-component is needed (and the four entries
492 // of the non-symmetric matrix); In 3d there are nine entries.
494 (DisplacementDim == 2 ? 1 : 0),
499
501
502 // B for Green-Lagrange strain.
507
508 B.noalias() += 0.5 * NonLinearBMatrix::computeBMatrix<
511 dNdx, N, grad_u, x_coord,
513
514 double const detF0 =
516 double const alpha = computeOutputStrainData(
517 detF0, B, grad_u, local_x, this->output_data_[ip]);
518
521 medium, this->current_states_[ip], this->prev_states_[ip],
522 this->material_states_[ip], this->output_data_[ip]);
523
524 this->material_states_[ip].pushBackState();
525 }
526
527 for (unsigned ip = 0; ip < n_integration_points; ip++)
528 {
529 this->prev_states_[ip] = this->current_states_[ip];
530 }
531 }

References _ip_data, ProcessLib::LinearBMatrix::computeBMatrix(), ProcessLib::NonLinearBMatrix::computeBMatrix(), computeFBarVariables(), ProcessLib::Deformation::computeGMatrix(), computeOutputStrainData(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::current_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::element_, ParameterLib::SpatialPosition::getCoordinates(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::integration_method_, NumLib::interpolateCoordinates(), ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::is_axially_symmetric_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::material_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::output_data_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::prev_states_, ProcessLib::LargeDeformation::LargeDeformationLocalAssemblerInterface< DisplacementDim >::process_data_, and updateConstitutiveRelations().

◆ updateConstitutiveRelations()

template<typename ShapeFunction, int DisplacementDim>
ConstitutiveRelations::ConstitutiveData< DisplacementDim > ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::updateConstitutiveRelations ( double const alpha,
GradientMatrixType const & G,
Eigen::Ref< Eigen::VectorXd const > const & u_prev,
ParameterLib::SpatialPosition const & x_position,
double const t,
double const dt,
typename ConstitutiveRelations::ConstitutiveSetting< DisplacementDim > & CS,
MaterialPropertyLib::Medium const & medium,
typename ConstitutiveRelations::StatefulData< DisplacementDim > & current_state,
typename ConstitutiveRelations::StatefulDataPrev< DisplacementDim > const & prev_state,
MaterialStateData< DisplacementDim > & material_state,
typename ConstitutiveRelations::OutputData< DisplacementDim > & output_data ) const
inline

Member Data Documentation

◆ _ip_data

template<typename ShapeFunction, int DisplacementDim>
std::vector<IpData, Eigen::aligned_allocator<IpData> > ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data
private

◆ _secondary_data

template<typename ShapeFunction, int DisplacementDim>
SecondaryData<typename ShapeMatrices::ShapeType> ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::_secondary_data
private

Definition at line 552 of file LargeDeformationFEM.h.

Referenced by LargeDeformationLocalAssembler(), and getShapeMatrix().

◆ displacement_size

template<typename ShapeFunction, int DisplacementDim>
const int ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::displacement_size
staticprivate
Initial value:
=
ShapeFunction::NPOINTS * DisplacementDim

Definition at line 554 of file LargeDeformationFEM.h.

◆ N_u_op

template<typename ShapeFunction, int DisplacementDim>
auto& ProcessLib::LargeDeformation::LargeDeformationLocalAssembler< ShapeFunction, DisplacementDim >::N_u_op
staticconstexpr
Initial value:
DisplacementDim, typename ShapeMatricesType::NodalRowVectorType>
constexpr Eigen::CwiseNullaryOp< EigenBlockMatrixViewFunctor< D, M >, typename EigenBlockMatrixViewFunctor< D, M >::Matrix > eigenBlockMatrixView(const Eigen::MatrixBase< M > &matrix)
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType

Definition at line 109 of file LargeDeformationFEM.h.

Referenced by assembleWithJacobian().


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