40namespace SmallDeformationNonlocal
46template <
typename ShapeMatrixType>
49 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>>
N;
52template <
typename ShapeFunction,
int DisplacementDim>
82 bool const is_axially_symmetric,
89 unsigned const n_integration_points =
92 _ip_data.reserve(n_integration_points);
95 auto const shape_matrices =
97 DisplacementDim>(e, is_axially_symmetric,
100 auto& solid_material =
105 auto* ehlers_solid_material =
dynamic_cast<
108 if (ehlers_solid_material ==
nullptr)
111 "The SmallDeformationNonlocalLocal process supports only "
112 "Ehlers material at the moment. For other materials the "
113 "interface must be extended first.");
116 for (
unsigned ip = 0; ip < n_integration_points; ip++)
118 _ip_data.emplace_back(*ehlers_solid_material);
120 auto const& sm = shape_matrices[ip];
123 sm.integralMeasure * sm.detJ;
126 ip_data.dNdx = sm.dNdx;
129 ip_data.sigma.setZero(
136 ip_data.sigma_prev.resize(
139 ip_data.eps_prev.resize(
150 double const* values,
151 int const integration_order)
override
153 if (integration_order !=
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.",
168 if (name ==
"kappa_d")
178 std::string
const& name, std::vector<double>
const& value)
override
180 if (name ==
"kappa_d_ip")
182 if (value.size() != 1)
185 "CellData for kappa_d initial conditions has wrong number "
186 "of components. 1 expected, got {:d}.",
195 double const internal_length2 =
_process_data.internal_length_squared;
196 return (distance2 > internal_length2)
198 : (1 - distance2 / (internal_length2)) *
199 (1 - distance2 / (internal_length2));
206 DisplacementDim>>>
const& local_assemblers)
override
211 unsigned const n_integration_points =
214 std::vector<double> distances;
220 for (
unsigned k = 0; k < n_integration_points; k++)
226 auto const& xyz =
_ip_data[k].coordinates;
229 for (
auto const search_element_id : search_element_ids)
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)
240 _ip_data[k].non_local_assemblers.push_back(
241 {la->getIPDataPtr(ip),
242 std::numeric_limits<double>::quiet_NaN(),
246 if (
_ip_data[k].non_local_assemblers.size() == 0)
251 double a_k_sum_m = 0;
252 for (
auto const& tuple :
_ip_data[k].non_local_assemblers)
254 double const distance2_m = tuple.distance2;
256 auto const& w_m = tuple.ip_l_pointer->integration_weight;
258 a_k_sum_m += w_m *
alpha_0(distance2_m);
265 for (
auto& tuple :
_ip_data[k].non_local_assemblers)
267 double const distance2_l = tuple.distance2;
268 double const a_kl =
alpha_0(distance2_l) / a_k_sum_m;
272 auto const w_l = tuple.ip_l_pointer->integration_weight;
273 tuple.alpha_kl_times_w_l = a_kl * w_l;
279 int integration_point)
const
283 Eigen::Vector3d xyz = Eigen::Vector3d::Zero();
285 for (
int i = 0; i < N.size(); ++i)
288 xyz += node_coordinates * N[i];
297 Eigen::Vector3d
const& coords,
298 std::vector<double>& distances)
const override
300 unsigned const n_integration_points =
303 distances.resize(n_integration_points);
305 for (
unsigned ip = 0; ip < n_integration_points; ip++)
307 auto const& xyz =
_ip_data[ip].coordinates;
308 distances[ip] = (xyz - coords).squaredNorm();
313 std::vector<double>
const& ,
314 std::vector<double>
const& ,
315 std::vector<double>& ,
316 std::vector<double>& ,
317 std::vector<double>& )
override
320 "SmallDeformationNonlocalLocalAssembler: assembly without jacobian "
326 std::vector<double>
const& local_x)
override
328 auto const n_integration_points =
336 for (
unsigned ip = 0; ip < n_integration_points; ip++)
341 auto const& dNdx =
_ip_data[ip].dNdx;
347 DisplacementDim, ShapeFunction::NPOINTS,
350 auto const& eps_prev =
_ip_data[ip].eps_prev;
351 auto const& sigma_prev =
_ip_data[ip].sigma_prev;
356 auto& state =
_ip_data[ip].material_state_variables;
357 double const& damage_prev =
_ip_data[ip].damage_prev;
361 Eigen::Map<typename BMatricesType::NodalForceVectorType const>(
362 local_x.data(), ShapeFunction::NPOINTS * DisplacementDim);
369 DisplacementDim>::MaterialStateVariables>
374 KelvinVectorType
const sigma_eff_prev =
379 variables_prev.
stress.emplace<
390 auto&& solution =
_ip_data[ip].solid_material.integrateStress(
391 variables_prev, variables, t, x_position, dt, *state);
395 OGS_FATAL(
"Computation of local constitutive relation failed.");
398 std::tie(sigma, state, C) = std::move(*solution);
402 auto const& ehlers_material =
404 DisplacementDim
> const&>(
_ip_data[ip].solid_material);
405 auto const damage_properties =
407 auto const material_properties =
408 ehlers_material.evaluatedMaterialProperties(t, x_position);
414 *
_ip_data[ip].material_state_variables);
416 double const eps_p_eff_diff =
417 state_vars.
eps_p.eff - state_vars.eps_p_prev.eff;
420 eps_p_eff_diff, sigma,
_ip_data[ip].kappa_d_prev,
421 damage_properties.h_d, material_properties);
428 for (
auto const& tuple :
432 tuple.ip_l_pointer->activated =
true;
441 std::vector<double>
const& local_x,
442 std::vector<double>
const& ,
443 std::vector<double>& local_b_data,
444 std::vector<double>& local_Jac_data)
override
446 auto const local_matrix_size = local_x.size();
449 local_Jac_data, local_matrix_size, local_matrix_size);
452 local_b_data, local_matrix_size);
454 unsigned const n_integration_points =
461 for (
unsigned ip = 0; ip < n_integration_points; ip++)
464 auto const& w =
_ip_data[ip].integration_weight;
467 auto const& dNdx =
_ip_data[ip].dNdx;
473 DisplacementDim, ShapeFunction::NPOINTS,
479 double& damage =
_ip_data[ip].damage;
482 double nonlocal_kappa_d = 0;
486 for (
auto const& tuple :
_ip_data[ip].non_local_assemblers)
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;
495 auto const& ehlers_material =
497 DisplacementDim
> const&>(
_ip_data[ip].solid_material);
512 nonlocal_kappa_d = std::max(0., nonlocal_kappa_d);
516 auto const damage_properties =
517 ehlers_material.evaluatedDamageProperties(t,
520 damage_properties.alpha_d,
521 damage_properties.beta_d);
522 damage = std::max(0., damage);
524 sigma = sigma * (1. - damage);
527 local_b.noalias() -= B.transpose() * sigma * w;
528 local_Jac.noalias() += B.transpose() * C * (1. - damage) * B * w;
534 unsigned const n_integration_points =
537 for (
unsigned ip = 0; ip < n_integration_points; ip++)
544 Eigen::VectorXd
const& ,
545 double const ,
double const ,
548 unsigned const n_integration_points =
551 for (
unsigned ip = 0; ip < n_integration_points; ip++)
560 double& crack_volume)
override
563 auto local_x = x.
get(indices);
565 auto u = Eigen::Map<typename BMatricesType::NodalForceVectorType const>(
566 local_x.data(), ShapeFunction::NPOINTS * DisplacementDim);
568 int const n_integration_points =
571 for (
int ip = 0; ip < n_integration_points; ip++)
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;
579 ShapeFunction::NPOINTS>(u, dNdx);
580 crack_volume += div_u * d * w;
585 const unsigned integration_point)
const override
590 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
594 std::vector<double>& nodal_values)
const override
596 nodal_values.clear();
598 nodal_values, ShapeFunction::NPOINTS * DisplacementDim);
600 unsigned const n_integration_points =
603 for (
unsigned ip = 0; ip < n_integration_points; ip++)
605 auto const& w =
_ip_data[ip].integration_weight;
608 auto const& dNdx =
_ip_data[ip].dNdx;
614 DisplacementDim, ShapeFunction::NPOINTS,
619 local_b.noalias() += B.transpose() * sigma * w;
627 std::vector<GlobalVector*>
const& ,
628 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
629 std::vector<double>& cache)
const override
635 [](
auto const& ip_data)
636 {
return ip_data.free_energy_density; });
643 std::vector<GlobalVector*>
const& ,
644 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
645 std::vector<double>& cache)
const override
651 [](
auto const& ip_data) {
return *ip_data.eps_p_V; });
658 std::vector<GlobalVector*>
const& ,
659 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
660 std::vector<double>& cache)
const override
666 [](
auto const& ip_data) {
return *ip_data.eps_p_D_xx; });
673 std::vector<GlobalVector*>
const& ,
674 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
675 std::vector<double>& cache)
const override
683 std::vector<GlobalVector*>
const& ,
684 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
685 std::vector<double>& cache)
const override
710 ip_data.kappa_d = value;
715 unsigned const n_integration_points =
718 std::vector<double> result_values;
719 result_values.resize(n_integration_points);
721 for (
unsigned ip = 0; ip < n_integration_points; ++ip)
723 result_values[ip] =
_ip_data[ip].kappa_d;
726 return result_values;
731 std::vector<GlobalVector*>
const& ,
732 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
733 std::vector<double>& cache)
const override
745 DisplacementDim>::MaterialStateVariables
const&
748 return *
_ip_data[integration_point].material_state_variables;
753 std::size_t
const component)
const
758 for (
auto const& ip_data :
_ip_data)
762 cache.push_back(ip_data.sigma[component]);
766 cache.push_back(ip_data.sigma[component] / std::sqrt(2));
774 std::vector<double>& cache, std::size_t
const component)
const
779 for (
auto const& ip_data :
_ip_data)
782 cache.push_back(ip_data.eps[component]);
784 cache.push_back(ip_data.eps[component] / std::sqrt(2));
798 std::vector<IpData, Eigen::aligned_allocator<IpData>>
_ip_data;
806 ShapeFunction::NPOINTS * DisplacementDim;
DamageProperties evaluatedDamageProperties(double const t, ParameterLib::SpatialPosition const &x) const
KelvinVector mechanical_strain
Global vector based on Eigen vector.
double get(IndexType rowId) const
get entry
Eigen::Vector3d const & asEigenVector3d() const
virtual Node *const * getNodes() const =0
Get array of element nodes.
std::size_t getID() const
Returns the ID of the element.
unsigned getIntegrationOrder() const
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
unsigned getNumberOfPoints() const
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
MatrixType< _number_of_dof, _number_of_dof > StiffnessMatrixType
VectorType< _number_of_dof > NodalForceVectorType
Rhs residual.
VectorType< _kelvin_vector_size > KelvinVectorType
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.
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), kelvin_vector_dimensions(DisplacementDim), Eigen::RowMajor > KelvinMatrixType
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
std::vector< std::size_t > findElementsWithinRadius(Element const &start_element, double const radius_squared)
double interpolateXCoordinate(MeshLib::Element const &e, typename ShapeMatricesType::ShapeMatrices::ShapeType const &N)
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
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)
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.
std::vector< double > const & getIntegrationPointScalarData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)
std::vector< double > const & getIntegrationPointKelvinVectorData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)
std::size_t setIntegrationPointKelvinVectorData(double const *values, IntegrationPointDataVector &ip_data_vector, MemberType IpData::*const member)
std::size_t setIntegrationPointScalarData(double const *values, IntegrationPointDataVector &ip_data_vector, MemberType IpData::*const member)
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
VectorType< GlobalDim > GlobalDimVectorType
VectorType< ShapeFunction::NPOINTS > NodalVectorType
PlasticStrain< KelvinVector > eps_p
plastic part of the state.