OGS
ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim > Class Template Reference

Detailed Description

template<typename ShapeFunction, int DisplacementDim>
class ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >

Definition at line 53 of file SmallDeformationNonlocalFEM.h.

#include <SmallDeformationNonlocalFEM.h>

Inheritance diagram for ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >:
[legend]

Public Types

using ShapeMatricesType
 
using NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
 
using NodalVectorType = typename ShapeMatricesType::NodalVectorType
 
using GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
 
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 IpData
 

Public Member Functions

 SmallDeformationNonlocalLocalAssembler (SmallDeformationNonlocalLocalAssembler const &)=delete
 
 SmallDeformationNonlocalLocalAssembler (SmallDeformationNonlocalLocalAssembler &&)=delete
 
 SmallDeformationNonlocalLocalAssembler (MeshLib::Element const &e, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, SmallDeformationNonlocalProcessData< DisplacementDim > &process_data)
 
std::size_t setIPDataInitialConditions (std::string_view const name, double const *values, int const integration_order) override
 
void setIPDataInitialConditionsFromCellData (std::string const &name, std::vector< double > const &value) override
 
double alpha_0 (double const distance2) const
 
void nonlocal (std::size_t const, std::vector< std::unique_ptr< SmallDeformationNonlocalLocalAssemblerInterface< DisplacementDim > > > const &local_assemblers) override
 
Eigen::Vector3d getSingleIntegrationPointCoordinates (int integration_point) const
 
void getIntegrationPointCoordinates (Eigen::Vector3d const &coords, std::vector< double > &distances) const override
 
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 preAssemble (double const t, double const dt, std::vector< double > const &local_x) override
 
void assembleWithJacobian (double const t, double const, std::vector< double > const &local_x, std::vector< double > const &, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
 
void initializeConcrete () override
 
void postTimestepConcrete (Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const) override
 
void computeCrackIntegral (std::size_t mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double &crack_volume) 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 & getNodalValues (std::vector< double > &nodal_values) const override
 
std::vector< double > const & getIntPtFreeEnergyDensity (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtEpsPV (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtEpsPDXX (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtSigma (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtEpsilon (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::size_t setSigma (double const *values)
 
std::vector< double > getSigma () const override
 
void setKappaD (double value)
 
std::vector< double > getKappaD () const override
 
std::vector< double > const & getIntPtDamage (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
unsigned getNumberOfIntegrationPoints () const override
 
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt (int const integration_point) const override
 
- Public Member Functions inherited from ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssemblerInterface< DisplacementDim >
virtual void nonlocal (std::size_t const mesh_item_id, std::vector< std::unique_ptr< SmallDeformationNonlocalLocalAssemblerInterface > > const &local_assemblers)=0
 
- 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 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 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
 
IntegrationPointDataNonlocalInterfacegetIPDataPtr (int const ip) override
 

Private Attributes

SmallDeformationNonlocalProcessData< DisplacementDim > & _process_data
 
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
 
NumLib::GenericIntegrationMethod const & _integration_method
 
MeshLib::Element const & _element
 
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
 
bool const _is_axially_symmetric
 

Static Private Attributes

static const int displacement_size
 

Member Typedef Documentation

◆ BMatricesType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::BMatricesType = BMatrixPolicyType<ShapeFunction, DisplacementDim>

Definition at line 63 of file SmallDeformationNonlocalFEM.h.

◆ BMatrixType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::BMatrixType = typename BMatricesType::BMatrixType

Definition at line 65 of file SmallDeformationNonlocalFEM.h.

◆ GlobalDimVectorType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType

Definition at line 61 of file SmallDeformationNonlocalFEM.h.

◆ IpData

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::IpData
Initial value:
IntegrationPointData<BMatricesType, ShapeMatricesType, DisplacementDim>

Definition at line 70 of file SmallDeformationNonlocalFEM.h.

◆ NodalDisplacementVectorType

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

Definition at line 68 of file SmallDeformationNonlocalFEM.h.

◆ NodalForceVectorType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::NodalForceVectorType = typename BMatricesType::NodalForceVectorType

Definition at line 67 of file SmallDeformationNonlocalFEM.h.

◆ NodalMatrixType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::NodalMatrixType = typename ShapeMatricesType::NodalMatrixType

Definition at line 59 of file SmallDeformationNonlocalFEM.h.

◆ NodalVectorType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::NodalVectorType = typename ShapeMatricesType::NodalVectorType

Definition at line 60 of file SmallDeformationNonlocalFEM.h.

◆ ShapeMatrices

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices

Definition at line 62 of file SmallDeformationNonlocalFEM.h.

◆ ShapeMatricesType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::ShapeMatricesType

◆ StiffnessMatrixType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::StiffnessMatrixType = typename BMatricesType::StiffnessMatrixType

Definition at line 66 of file SmallDeformationNonlocalFEM.h.

Constructor & Destructor Documentation

◆ SmallDeformationNonlocalLocalAssembler() [1/3]

template<typename ShapeFunction , int DisplacementDim>
ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::SmallDeformationNonlocalLocalAssembler ( SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim > const & )
delete

◆ SmallDeformationNonlocalLocalAssembler() [2/3]

template<typename ShapeFunction , int DisplacementDim>
ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::SmallDeformationNonlocalLocalAssembler ( SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim > && )
delete

◆ SmallDeformationNonlocalLocalAssembler() [3/3]

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

Definition at line 78 of file SmallDeformationNonlocalFEM.h.

84 : _process_data(process_data),
85 _integration_method(integration_method),
86 _element(e),
87 _is_axially_symmetric(is_axially_symmetric)
88 {
89 unsigned const n_integration_points =
91
92 _ip_data.reserve(n_integration_points);
93 _secondary_data.N.resize(n_integration_points);
94
95 auto const shape_matrices =
97 DisplacementDim>(e, is_axially_symmetric,
99
100 auto& solid_material =
102 _process_data.solid_materials,
103 _process_data.material_ids,
104 e.getID());
105 auto* ehlers_solid_material = dynamic_cast<
107 &solid_material);
108 if (ehlers_solid_material == nullptr)
109 {
110 OGS_FATAL(
111 "The SmallDeformationNonlocalLocal process supports only "
112 "Ehlers material at the moment. For other materials the "
113 "interface must be extended first.");
114 }
115
116 for (unsigned ip = 0; ip < n_integration_points; ip++)
117 {
118 _ip_data.emplace_back(*ehlers_solid_material);
119 auto& ip_data = _ip_data[ip];
120 auto const& sm = shape_matrices[ip];
121 _ip_data[ip].integration_weight =
123 sm.integralMeasure * sm.detJ;
124
125 ip_data.N = sm.N;
126 ip_data.dNdx = sm.dNdx;
127
128 // Initialize current time step values
129 ip_data.sigma.setZero(
131 DisplacementDim));
133 DisplacementDim));
134
135 // Previous time step values are not initialized and are set later.
136 ip_data.sigma_prev.resize(
138 DisplacementDim));
139 ip_data.eps_prev.resize(
141 DisplacementDim));
142
143 _secondary_data.N[ip] = shape_matrices[ip].N;
144
145 ip_data.coordinates = getSingleIntegrationPointCoordinates(ip);
146 }
147 }
#define OGS_FATAL(...)
Definition Error.h:26
double getWeight() const
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
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::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_integration_method, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_process_data, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_secondary_data, MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getSingleIntegrationPointCoordinates(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), MathLib::KelvinVector::kelvin_vector_dimensions(), ProcessLib::SmallDeformationNonlocal::SecondaryData< ShapeMatrixType >::N, OGS_FATAL, and MaterialLib::Solids::selectSolidConstitutiveRelation().

Member Function Documentation

◆ alpha_0()

template<typename ShapeFunction , int DisplacementDim>
double ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::alpha_0 ( double const distance2) const
inline

Definition at line 193 of file SmallDeformationNonlocalFEM.h.

194 {
195 double const internal_length2 = _process_data.internal_length_squared;
196 return (distance2 > internal_length2)
197 ? 0
198 : (1 - distance2 / (internal_length2)) *
199 (1 - distance2 / (internal_length2));
200 }

References ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_process_data.

Referenced by ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::nonlocal().

◆ assemble()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< 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 312 of file SmallDeformationNonlocalFEM.h.

318 {
319 OGS_FATAL(
320 "SmallDeformationNonlocalLocalAssembler: assembly without jacobian "
321 "is not "
322 "implemented.");
323 }

References OGS_FATAL.

◆ assembleWithJacobian()

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

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 440 of file SmallDeformationNonlocalFEM.h.

445 {
446 auto const local_matrix_size = local_x.size();
447
449 local_Jac_data, local_matrix_size, local_matrix_size);
450
452 local_b_data, local_matrix_size);
453
454 unsigned const n_integration_points =
456
458 x_position.setElementID(_element.getID());
459
460 // Non-local integration.
461 for (unsigned ip = 0; ip < n_integration_points; ip++)
462 {
463 x_position.setIntegrationPoint(ip);
464 auto const& w = _ip_data[ip].integration_weight;
465
466 auto const& N = _ip_data[ip].N;
467 auto const& dNdx = _ip_data[ip].dNdx;
468
469 auto const x_coord =
470 NumLib::interpolateXCoordinate<ShapeFunction,
472 auto const B = LinearBMatrix::computeBMatrix<
473 DisplacementDim, ShapeFunction::NPOINTS,
474 typename BMatricesType::BMatrixType>(dNdx, N, x_coord,
476
477 auto& sigma = _ip_data[ip].sigma;
478 auto& C = _ip_data[ip].C;
479 double& damage = _ip_data[ip].damage;
480
481 {
482 double nonlocal_kappa_d = 0;
483
484 if (_ip_data[ip].active_self || _ip_data[ip].activated)
485 {
486 for (auto const& tuple : _ip_data[ip].non_local_assemblers)
487 {
488 // Get local variable for the integration point l.
489 double const kappa_d_l = tuple.ip_l_pointer->kappa_d;
490 double const a_kl_times_w_l = tuple.alpha_kl_times_w_l;
491 nonlocal_kappa_d += a_kl_times_w_l * kappa_d_l;
492 }
493 }
494
495 auto const& ehlers_material =
497 DisplacementDim> const&>(_ip_data[ip].solid_material);
498
499 //
500 // Overnonlocal formulation
501 //
502 // See (Di Luzio & Bazant 2005, IJSS) for details.
503 // The implementation would go here and would be for a given
504 // gamma_nonlocal:
505 //
506 // Update nonlocal damage with local damage (scaled with 1 -
507 // \gamma_{nonlocal}) for the current integration point and the
508 // nonlocal integral part.
509 // nonlocal_kappa_d = (1. - gamma_nonlocal) * kappa_d +
510 // gamma_nonlocal * nonlocal_kappa_d;
511
512 nonlocal_kappa_d = std::max(0., nonlocal_kappa_d);
513
514 // Update damage based on nonlocal kappa_d
515 {
516 auto const damage_properties =
517 ehlers_material.evaluatedDamageProperties(t,
518 x_position);
519 damage = calculateDamage(nonlocal_kappa_d,
520 damage_properties.alpha_d,
521 damage_properties.beta_d);
522 damage = std::max(0., damage);
523 }
524 sigma = sigma * (1. - damage);
525 }
526
527 local_b.noalias() -= B.transpose() * sigma * w;
528 local_Jac.noalias() += B.transpose() * C * (1. - damage) * B * w;
529 }
530 }
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::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 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.
double calculateDamage(double const kappa_d, double const alpha_d, double const beta_d)
Definition Damage.h:25

References ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_element, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_integration_method, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_is_axially_symmetric, ProcessLib::SmallDeformationNonlocal::calculateDamage(), ProcessLib::LinearBMatrix::computeBMatrix(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), NumLib::interpolateXCoordinate(), ParameterLib::SpatialPosition::setElementID(), and ParameterLib::SpatialPosition::setIntegrationPoint().

◆ computeCrackIntegral()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::computeCrackIntegral ( std::size_t mesh_item_id,
NumLib::LocalToGlobalIndexMap const & dof_table,
GlobalVector const & x,
double & crack_volume )
inlineoverridevirtual

Implements ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssemblerInterface< DisplacementDim >.

Definition at line 557 of file SmallDeformationNonlocalFEM.h.

561 {
562 auto const indices = NumLib::getIndices(mesh_item_id, dof_table);
563 auto local_x = x.get(indices);
564
565 auto u = Eigen::Map<typename BMatricesType::NodalForceVectorType const>(
566 local_x.data(), ShapeFunction::NPOINTS * DisplacementDim);
567
568 int const n_integration_points =
570
571 for (int ip = 0; ip < n_integration_points; ip++)
572 {
573 auto const& dNdx = _ip_data[ip].dNdx;
574 auto const& d = _ip_data[ip].damage;
575 auto const& w = _ip_data[ip].integration_weight;
576
577 double const div_u =
578 Deformation::divergence<DisplacementDim,
579 ShapeFunction::NPOINTS>(u, dNdx);
580 crack_volume += div_u * d * w;
581 }
582 }
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
double divergence(const Eigen::Ref< Eigen::Matrix< double, NPOINTS *DisplacementDim, 1 > const > &u, DNDX_Type const &dNdx)
Divergence of displacement, the volumetric strain.
Definition Divergence.h:19

References ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_integration_method, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::Deformation::divergence(), MathLib::EigenVector::get(), NumLib::getIndices(), and NumLib::GenericIntegrationMethod::getNumberOfPoints().

◆ getIntegrationPointCoordinates()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getIntegrationPointCoordinates ( Eigen::Vector3d const & coords,
std::vector< double > & distances ) const
inlineoverridevirtual

For each of the current element's integration points the squared distance from the current integration point is computed and stored in the given distances cache.

Implements ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssemblerInterface< DisplacementDim >.

Definition at line 296 of file SmallDeformationNonlocalFEM.h.

299 {
300 unsigned const n_integration_points =
302
303 distances.resize(n_integration_points);
304
305 for (unsigned ip = 0; ip < n_integration_points; ip++)
306 {
307 auto const& xyz = _ip_data[ip].coordinates;
308 distances[ip] = (xyz - coords).squaredNorm();
309 }
310 }
constexpr ranges::views::view_closure coords
Definition Mesh.h:232

References ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_integration_method, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, and NumLib::GenericIntegrationMethod::getNumberOfPoints().

◆ getIntPtDamage()

template<typename ShapeFunction , int DisplacementDim>
std::vector< double > const & ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getIntPtDamage ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverridevirtual

◆ getIntPtEpsilon() [1/2]

template<typename ShapeFunction , int DisplacementDim>
std::vector< double > const & ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getIntPtEpsilon ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverridevirtual

◆ getIntPtEpsilon() [2/2]

template<typename ShapeFunction , int DisplacementDim>
std::vector< double > const & ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getIntPtEpsilon ( std::vector< double > & cache,
std::size_t const component ) const
inlineprivate

Definition at line 773 of file SmallDeformationNonlocalFEM.h.

775 {
776 cache.clear();
777 cache.reserve(_ip_data.size());
778
779 for (auto const& ip_data : _ip_data)
780 {
781 if (component < 3) // xx, yy, zz components
782 cache.push_back(ip_data.eps[component]);
783 else // mixed xy, yz, xz components
784 cache.push_back(ip_data.eps[component] / std::sqrt(2));
785 }
786
787 return cache;
788 }

References ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data.

◆ getIntPtEpsPDXX()

template<typename ShapeFunction , int DisplacementDim>
std::vector< double > const & ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getIntPtEpsPDXX ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverridevirtual

Implements ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssemblerInterface< DisplacementDim >.

Definition at line 656 of file SmallDeformationNonlocalFEM.h.

661 {
662 cache.clear();
663 cache.reserve(_ip_data.size());
664
665 transform(cbegin(_ip_data), cend(_ip_data), back_inserter(cache),
666 [](auto const& ip_data) { return *ip_data.eps_p_D_xx; });
667
668 return cache;
669 }

References ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data.

◆ getIntPtEpsPV()

template<typename ShapeFunction , int DisplacementDim>
std::vector< double > const & ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getIntPtEpsPV ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverridevirtual

Implements ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssemblerInterface< DisplacementDim >.

Definition at line 641 of file SmallDeformationNonlocalFEM.h.

646 {
647 cache.clear();
648 cache.reserve(_ip_data.size());
649
650 transform(cbegin(_ip_data), cend(_ip_data), back_inserter(cache),
651 [](auto const& ip_data) { return *ip_data.eps_p_V; });
652
653 return cache;
654 }

References ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data.

◆ getIntPtFreeEnergyDensity()

template<typename ShapeFunction , int DisplacementDim>
std::vector< double > const & ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getIntPtFreeEnergyDensity ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverridevirtual

Implements ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssemblerInterface< DisplacementDim >.

Definition at line 625 of file SmallDeformationNonlocalFEM.h.

630 {
631 cache.clear();
632 cache.reserve(_ip_data.size());
633
634 transform(cbegin(_ip_data), cend(_ip_data), back_inserter(cache),
635 [](auto const& ip_data)
636 { return ip_data.free_energy_density; });
637
638 return cache;
639 }

References ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data.

◆ getIntPtSigma() [1/2]

◆ getIntPtSigma() [2/2]

template<typename ShapeFunction , int DisplacementDim>
std::vector< double > const & ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getIntPtSigma ( std::vector< double > & cache,
std::size_t const component ) const
inlineprivate

Definition at line 752 of file SmallDeformationNonlocalFEM.h.

754 {
755 cache.clear();
756 cache.reserve(_ip_data.size());
757
758 for (auto const& ip_data : _ip_data)
759 {
760 if (component < 3)
761 { // xx, yy, zz components
762 cache.push_back(ip_data.sigma[component]);
763 }
764 else
765 { // mixed xy, yz, xz components
766 cache.push_back(ip_data.sigma[component] / std::sqrt(2));
767 }
768 }
769
770 return cache;
771 }

References ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data.

◆ getIPDataPtr()

◆ getKappaD()

template<typename ShapeFunction , int DisplacementDim>
std::vector< double > ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getKappaD ( ) const
inlineoverridevirtual

Implements ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssemblerInterface< DisplacementDim >.

Definition at line 713 of file SmallDeformationNonlocalFEM.h.

714 {
715 unsigned const n_integration_points =
717
718 std::vector<double> result_values;
719 result_values.resize(n_integration_points);
720
721 for (unsigned ip = 0; ip < n_integration_points; ++ip)
722 {
723 result_values[ip] = _ip_data[ip].kappa_d;
724 }
725
726 return result_values;
727 }

References ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_integration_method, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, and NumLib::GenericIntegrationMethod::getNumberOfPoints().

◆ getMaterialStateVariablesAt()

template<typename ShapeFunction , int DisplacementDim>
MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getMaterialStateVariablesAt ( int const integration_point) const
inlineoverridevirtual

◆ getNodalValues()

template<typename ShapeFunction , int DisplacementDim>
std::vector< double > const & ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getNodalValues ( std::vector< double > & nodal_values) const
inlineoverridevirtual

Implements ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssemblerInterface< DisplacementDim >.

Definition at line 593 of file SmallDeformationNonlocalFEM.h.

595 {
596 nodal_values.clear();
598 nodal_values, ShapeFunction::NPOINTS * DisplacementDim);
599
600 unsigned const n_integration_points =
602
603 for (unsigned ip = 0; ip < n_integration_points; ip++)
604 {
605 auto const& w = _ip_data[ip].integration_weight;
606
607 auto const& N = _ip_data[ip].N;
608 auto const& dNdx = _ip_data[ip].dNdx;
609
610 auto const x_coord =
611 NumLib::interpolateXCoordinate<ShapeFunction,
613 auto const B = LinearBMatrix::computeBMatrix<
614 DisplacementDim, ShapeFunction::NPOINTS,
615 typename BMatricesType::BMatrixType>(dNdx, N, x_coord,
617 auto& sigma = _ip_data[ip].sigma;
618
619 local_b.noalias() += B.transpose() * sigma * w;
620 }
621
622 return nodal_values;
623 }

References ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_element, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_integration_method, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_is_axially_symmetric, ProcessLib::LinearBMatrix::computeBMatrix(), MathLib::createZeroedVector(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), and NumLib::interpolateXCoordinate().

◆ getNumberOfIntegrationPoints()

◆ getShapeMatrix()

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

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 584 of file SmallDeformationNonlocalFEM.h.

586 {
587 auto const& N = _secondary_data.N[integration_point];
588
589 // assumes N is stored contiguously in memory
590 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
591 }

References ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_secondary_data, and ProcessLib::SmallDeformationNonlocal::SecondaryData< ShapeMatrixType >::N.

◆ getSigma()

◆ getSingleIntegrationPointCoordinates()

template<typename ShapeFunction , int DisplacementDim>
Eigen::Vector3d ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getSingleIntegrationPointCoordinates ( int integration_point) const
inline

Definition at line 278 of file SmallDeformationNonlocalFEM.h.

280 {
281 auto const& N = _secondary_data.N[integration_point];
282
283 Eigen::Vector3d xyz = Eigen::Vector3d::Zero(); // Resulting coordinates
284 auto* nodes = _element.getNodes();
285 for (int i = 0; i < N.size(); ++i)
286 {
287 auto const& node_coordinates{nodes[i]->asEigenVector3d()};
288 xyz += node_coordinates * N[i];
289 }
290 return xyz;
291 }
Eigen::Vector3d const & asEigenVector3d() const
Definition Point3d.h:64
virtual Node *const * getNodes() const =0
Get array of element nodes.

References ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_element, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_secondary_data, MathLib::Point3d::asEigenVector3d(), MeshLib::Element::getNodes(), and ProcessLib::SmallDeformationNonlocal::SecondaryData< ShapeMatrixType >::N.

Referenced by ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::SmallDeformationNonlocalLocalAssembler().

◆ initializeConcrete()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::initializeConcrete ( )
inlineoverridevirtual

◆ nonlocal()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::nonlocal ( std::size_t const ,
std::vector< std::unique_ptr< SmallDeformationNonlocalLocalAssemblerInterface< DisplacementDim > > > const & local_assemblers )
inlineoverride

Definition at line 202 of file SmallDeformationNonlocalFEM.h.

207 {
208 auto const search_element_ids = MeshLib::findElementsWithinRadius(
209 _element, _process_data.internal_length_squared);
210
211 unsigned const n_integration_points =
213
214 std::vector<double> distances; // Cache for ip-ip distances.
215 //
216 // For every integration point in this element collect the neighbouring
217 // integration points falling in given radius (internal length) and
218 // compute the alpha_kl weight.
219 //
220 for (unsigned k = 0; k < n_integration_points; k++)
221 {
222 //
223 // Collect the integration points.
224 //
225
226 auto const& xyz = _ip_data[k].coordinates;
227
228 // For all neighbors of element
229 for (auto const search_element_id : search_element_ids)
230 {
231 auto const& la = local_assemblers[search_element_id];
232 la->getIntegrationPointCoordinates(xyz, distances);
233 for (int ip = 0; ip < static_cast<int>(distances.size()); ++ip)
234 {
235 if (distances[ip] >= _process_data.internal_length_squared)
236 {
237 continue;
238 }
239 // save into current ip_k
240 _ip_data[k].non_local_assemblers.push_back(
241 {la->getIPDataPtr(ip),
242 std::numeric_limits<double>::quiet_NaN(),
243 distances[ip]});
244 }
245 }
246 if (_ip_data[k].non_local_assemblers.size() == 0)
247 {
248 OGS_FATAL("no neighbours found!");
249 }
250
251 double a_k_sum_m = 0;
252 for (auto const& tuple : _ip_data[k].non_local_assemblers)
253 {
254 double const distance2_m = tuple.distance2;
255
256 auto const& w_m = tuple.ip_l_pointer->integration_weight;
257
258 a_k_sum_m += w_m * alpha_0(distance2_m);
259 }
260
261 //
262 // Calculate alpha_kl =
263 // alpha_0(|x_k - x_l|) / int_{m \in ip} alpha_0(|x_k - x_m|)
264 //
265 for (auto& tuple : _ip_data[k].non_local_assemblers)
266 {
267 double const distance2_l = tuple.distance2;
268 double const a_kl = alpha_0(distance2_l) / a_k_sum_m;
269
270 // Store the a_kl already multiplied with the integration
271 // weight of that l integration point.
272 auto const w_l = tuple.ip_l_pointer->integration_weight;
273 tuple.alpha_kl_times_w_l = a_kl * w_l;
274 }
275 }
276 }
std::vector< std::size_t > findElementsWithinRadius(Element const &start_element, double const radius_squared)

References ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_element, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_integration_method, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_process_data, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::alpha_0(), MeshLib::findElementsWithinRadius(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), and OGS_FATAL.

◆ postTimestepConcrete()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::postTimestepConcrete ( Eigen::VectorXd const & ,
Eigen::VectorXd const & ,
double const ,
double const ,
int const  )
inlineoverridevirtual

◆ preAssemble()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::preAssemble ( double const t,
double const dt,
std::vector< double > const & local_x )
inlineoverridevirtual

Compute only the local kappa_d.

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 325 of file SmallDeformationNonlocalFEM.h.

327 {
328 auto const n_integration_points =
330
331 MPL::VariableArray variables;
332 MPL::VariableArray variables_prev;
334 x_position.setElementID(_element.getID());
335
336 for (unsigned ip = 0; ip < n_integration_points; ip++)
337 {
338 x_position.setIntegrationPoint(ip);
339
340 auto const& N = _ip_data[ip].N;
341 auto const& dNdx = _ip_data[ip].dNdx;
342
343 auto const x_coord =
344 NumLib::interpolateXCoordinate<ShapeFunction,
346 auto const B = LinearBMatrix::computeBMatrix<
347 DisplacementDim, ShapeFunction::NPOINTS,
348 typename BMatricesType::BMatrixType>(dNdx, N, x_coord,
350 auto const& eps_prev = _ip_data[ip].eps_prev;
351 auto const& sigma_prev = _ip_data[ip].sigma_prev;
352
353 auto& eps = _ip_data[ip].eps;
354 auto& sigma = _ip_data[ip].sigma;
355 auto& C = _ip_data[ip].C;
356 auto& state = _ip_data[ip].material_state_variables;
357 double const& damage_prev = _ip_data[ip].damage_prev;
358
359 eps.noalias() =
360 B *
361 Eigen::Map<typename BMatricesType::NodalForceVectorType const>(
362 local_x.data(), ShapeFunction::NPOINTS * DisplacementDim);
363
364 // sigma is for plastic part only.
365 std::unique_ptr<
367 new_C;
368 std::unique_ptr<typename MaterialLib::Solids::MechanicsBase<
369 DisplacementDim>::MaterialStateVariables>
370 new_state;
371
372 // Compute sigma_eff from damage total stress sigma
374 KelvinVectorType const sigma_eff_prev =
375 sigma_prev /
376 (1. - damage_prev); // damage_prev is in [0,1) range. See
377 // calculateDamage() function.
378
379 variables_prev.stress.emplace<
381 sigma_eff_prev);
382 variables_prev.mechanical_strain.emplace<
384 eps_prev);
385 variables_prev.temperature = _process_data.reference_temperature;
386 variables.mechanical_strain.emplace<
388 variables.temperature = _process_data.reference_temperature;
389
390 auto&& solution = _ip_data[ip].solid_material.integrateStress(
391 variables_prev, variables, t, x_position, dt, *state);
392
393 if (!solution)
394 {
395 OGS_FATAL("Computation of local constitutive relation failed.");
396 }
397
398 std::tie(sigma, state, C) = std::move(*solution);
399
401 {
402 auto const& ehlers_material =
404 DisplacementDim> const&>(_ip_data[ip].solid_material);
405 auto const damage_properties =
406 ehlers_material.evaluatedDamageProperties(t, x_position);
407 auto const material_properties =
408 ehlers_material.evaluatedMaterialProperties(t, x_position);
409
410 // Ehlers material state variables
411 auto& state_vars =
413 DisplacementDim>&>(
414 *_ip_data[ip].material_state_variables);
415
416 double const eps_p_eff_diff =
417 state_vars.eps_p.eff - state_vars.eps_p_prev.eff;
418
420 eps_p_eff_diff, sigma, _ip_data[ip].kappa_d_prev,
421 damage_properties.h_d, material_properties);
422
423 if (!_ip_data[ip].active_self)
424 {
425 _ip_data[ip].active_self |= _ip_data[ip].kappa_d > 0;
426 if (_ip_data[ip].active_self)
427 {
428 for (auto const& tuple :
429 _ip_data[ip].non_local_assemblers)
430 {
431 // Activate the integration point.
432 tuple.ip_l_pointer->activated = true;
433 }
434 }
435 }
436 }
437 }
438 }
DamageProperties evaluatedDamageProperties(double const t, ParameterLib::SpatialPosition const &x) const
Definition Ehlers.h:343
VectorType< _kelvin_vector_size > KelvinVectorType
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
double calculateDamageKappaD(double const eps_p_eff_diff, KelvinVectorType const &sigma, double const kappa_d_prev, double const h_d, MaterialLib::Solids::Ehlers::MaterialProperties const &mp)
Definition Damage.h:39
PlasticStrain< KelvinVector > eps_p
plastic part of the state.
Definition Ehlers.h:244

References ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_element, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_integration_method, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_is_axially_symmetric, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_process_data, ProcessLib::SmallDeformationNonlocal::calculateDamageKappaD(), ProcessLib::LinearBMatrix::computeBMatrix(), MaterialLib::Solids::Ehlers::StateVariables< DisplacementDim >::eps_p, MaterialLib::Solids::Ehlers::SolidEhlers< DisplacementDim >::evaluatedDamageProperties(), MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), NumLib::interpolateXCoordinate(), MaterialPropertyLib::VariableArray::mechanical_strain, OGS_FATAL, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), MaterialPropertyLib::VariableArray::stress, and MaterialPropertyLib::VariableArray::temperature.

◆ setIPDataInitialConditions()

template<typename ShapeFunction , int DisplacementDim>
std::size_t ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::setIPDataInitialConditions ( std::string_view const name,
double const * values,
int const integration_order )
inlineoverridevirtual

Implements ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssemblerInterface< DisplacementDim >.

Definition at line 149 of file SmallDeformationNonlocalFEM.h.

152 {
153 if (integration_order !=
154 static_cast<int>(_integration_method.getIntegrationOrder()))
155 {
156 OGS_FATAL(
157 "Setting integration point initial conditions; The integration "
158 "order of the local assembler for element {:d} is different "
159 "from the integration order in the initial condition.",
160 _element.getID());
161 }
162
163 if (name == "sigma")
164 {
165 return setSigma(values);
166 }
167
168 if (name == "kappa_d")
169 {
172 }
173
174 return 0;
175 }
std::size_t setIntegrationPointScalarData(double const *values, IntegrationPointDataVector &ip_data_vector, MemberType IpData::*const member)

References ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_element, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_integration_method, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_ip_data, MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getIntegrationOrder(), ProcessLib::SmallDeformationNonlocal::IntegrationPointDataNonlocalInterface::kappa_d, OGS_FATAL, ProcessLib::setIntegrationPointScalarData(), and ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::setSigma().

◆ setIPDataInitialConditionsFromCellData()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::setIPDataInitialConditionsFromCellData ( std::string const & name,
std::vector< double > const & value )
inlineoverridevirtual

Implements ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssemblerInterface< DisplacementDim >.

Definition at line 177 of file SmallDeformationNonlocalFEM.h.

179 {
180 if (name == "kappa_d_ip")
181 {
182 if (value.size() != 1)
183 {
184 OGS_FATAL(
185 "CellData for kappa_d initial conditions has wrong number "
186 "of components. 1 expected, got {:d}.",
187 value.size());
188 }
189 setKappaD(value[0]);
190 }
191 }

References OGS_FATAL, and ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::setKappaD().

◆ setKappaD()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::setKappaD ( double value)
inline

◆ setSigma()

Member Data Documentation

◆ _element

◆ _integration_method

template<typename ShapeFunction , int DisplacementDim>
NumLib::GenericIntegrationMethod const& ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::_integration_method
private

Definition at line 800 of file SmallDeformationNonlocalFEM.h.

Referenced by ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::SmallDeformationNonlocalLocalAssembler(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::computeCrackIntegral(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getIntegrationPointCoordinates(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getKappaD(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getNodalValues(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getNumberOfIntegrationPoints(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::initializeConcrete(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::nonlocal(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::postTimestepConcrete(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::preAssemble(), and ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::setIPDataInitialConditions().

◆ _ip_data

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

Definition at line 798 of file SmallDeformationNonlocalFEM.h.

Referenced by ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::SmallDeformationNonlocalLocalAssembler(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::assembleWithJacobian(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::computeCrackIntegral(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getIntegrationPointCoordinates(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getIntPtDamage(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getIntPtEpsilon(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getIntPtEpsilon(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getIntPtEpsPDXX(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getIntPtEpsPV(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getIntPtFreeEnergyDensity(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getIntPtSigma(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getIntPtSigma(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getIPDataPtr(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getKappaD(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getMaterialStateVariablesAt(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getNodalValues(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::getSigma(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::initializeConcrete(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::nonlocal(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::postTimestepConcrete(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::preAssemble(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::setIPDataInitialConditions(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::setKappaD(), and ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, DisplacementDim >::setSigma().

◆ _is_axially_symmetric

◆ _process_data

◆ _secondary_data

◆ displacement_size

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

Definition at line 805 of file SmallDeformationNonlocalFEM.h.


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