41namespace LargeDeformation
47template <
typename ShapeMatrixType>
50 std::vector<ShapeMatrixType, Eigen::aligned_allocator<ShapeMatrixType>>
N;
53template <
int DisplacementDim,
typename ShapeMatricesType>
60 Eigen::Matrix<double, DisplacementDim, DisplacementDim>>;
62 using SigmaGeom = Eigen::Matrix<double,
MPL::tensorSize(DisplacementDim),
64 if constexpr (DisplacementDim == 2)
66 SigmaGeom sigma_geom = SigmaGeom::Zero(5, 5);
67 sigma_geom.template block<4, 4>(0, 0) =
68 sigma_geom_op(sigma_tensor.template block<2, 2>(0, 0).eval());
69 sigma_geom(4, 4) = sigma_tensor(2, 2);
74 if constexpr (DisplacementDim == 3)
76 return sigma_geom_op(sigma_tensor);
80template <
typename ShapeFunction,
int DisplacementDim>
120 bool const is_axially_symmetric,
123 e, integration_method, is_axially_symmetric, process_data)
125 unsigned const n_integration_points =
128 _ip_data.resize(n_integration_points);
131 auto const shape_matrices =
136 for (
unsigned ip = 0; ip < n_integration_points; ip++)
139 auto const& sm = shape_matrices[ip];
142 sm.integralMeasure * sm.detJ;
145 ip_data.dNdx_u = sm.dNdx;
153 unsigned const n_integration_points =
155 for (
unsigned ip = 0; ip < n_integration_points; ip++)
173 double>::quiet_NaN() ,
180 t, x_position, *material_state.material_state_variables);
182 material_state.pushBackState();
192 Eigen::Ref<Eigen::VectorXd const>
const& u,
211 double const J_ratio = detF0 / detF;
216 "det(F0)/det(F) is negative with det(F0) = {:g}, and det(F) = "
222 DisplacementDim == 3 ? std::cbrt(J_ratio) : std::sqrt(J_ratio);
225 .template segment<DisplacementDim * DisplacementDim>(0) *= alpha;
227 std::pow(alpha, DisplacementDim);
229 double const alpha_p2 = alpha * alpha;
231 alpha_p2 * output_data.
eps_data.eps +
232 0.5 * (alpha_p2 - 1) *
241 Eigen::Ref<Eigen::VectorXd const>
const& u_prev,
257 ? (*this->
process_data_.reference_temperature)(t, x_position)[0]
258 : std::numeric_limits<double>::quiet_NaN();
267 models, t, dt, x_position,
271 alpha * (G * u_prev +
273 current_state, prev_state, material_state, tmp, CD);
279 std::vector<double>
const& ,
280 std::vector<double>
const& ,
281 std::vector<double>& ,
282 std::vector<double>& ,
283 std::vector<double>& )
override
286 "LargeDeformationLocalAssembler: assembly without jacobian is not "
291 bool const compute_detF0_only, Eigen::VectorXd
const& local_x)
const
305 _ip_data, compute_detF0_only, local_x,
313 compute_detF0_only, local_x, this->
element_,
318 std::vector<double>
const& local_x,
319 std::vector<double>
const& local_x_prev,
320 std::vector<double>& local_b_data,
321 std::vector<double>& local_Jac_data)
override
323 auto const local_matrix_size = local_x.size();
326 local_Jac_data, local_matrix_size, local_matrix_size);
329 local_b_data, local_matrix_size);
332 auto [u_prev] =
localDOF(local_x_prev);
334 unsigned const n_integration_points =
341 constitutive_setting;
345 auto const f_bar_variables =
348 for (
unsigned ip = 0; ip < n_integration_points; ip++)
350 auto const& w =
_ip_data[ip].integration_weight;
352 auto const& dNdx =
_ip_data[ip].dNdx_u;
362 (DisplacementDim == 2 ? 1 : 0),
363 ShapeFunction::NPOINTS * DisplacementDim);
365 ShapeFunction::NPOINTS>(
371 DisplacementDim, ShapeFunction::NPOINTS,
376 DisplacementDim, ShapeFunction::NPOINTS,
381 f_bar_variables ? std::get<double>(*f_bar_variables) : 0.0;
383 detF0, B_0 + 0.5 * B_N, grad_u, u, this->
output_data_[ip]);
386 alpha, G, u_prev, x_position, t, dt, constitutive_setting,
393 auto const& b = *CD.volumetric_body_force;
394 auto const& C = CD.s_mech_data.stiffness_tensor;
396 auto const sigma_geom =
400 if (!f_bar_variables)
403 (B.transpose() * sigma -
N_u_op(N).transpose() * b) * w;
405 local_Jac.noalias() += G.transpose() * sigma_geom * G * w;
407 local_Jac.noalias() += B.transpose() * C * B * w;
411 double const alpha_p2 = alpha * alpha;
412 local_Jac.noalias() +=
413 alpha_p2 * G.transpose() * sigma_geom * G * w;
415 auto const q0 = std::get<VectorTypeForFbar>(*f_bar_variables);
419 .deformation_gradient_data.deformation_gradient;
426 auto const q0_q = q0 - q;
428 auto const S_q = sigma * (q0_q).transpose();
430 local_Jac.noalias() +=
431 2.0 * alpha_p2 * S_q.transpose() * B * w / DisplacementDim;
433 auto const twoEplsI =
436 this->is_axially_symmetric_);
441 (B + twoEplsI * (q0_q).transpose() / DisplacementDim);
443 local_Jac.noalias() +=
444 2.0 * Bbar.transpose() * S_q * w / DisplacementDim;
446 local_Jac.noalias() += Bbar.transpose() * C * Bbar * w;
448 local_Jac.noalias() -=
449 alpha_p2 * sigma.dot(twoEplsI) *
450 (q0 * q0.transpose() - q * q.transpose()) * w /
454 (Bbar.transpose() * sigma -
N_u_op(N).transpose() * b) * w;
460 Eigen::VectorXd
const& local_x_prev,
461 double const t,
double const dt,
464 unsigned const n_integration_points =
471 constitutive_setting;
476 auto const f_bar_variables =
479 for (
unsigned ip = 0; ip < n_integration_points; ip++)
482 auto const& dNdx =
_ip_data[ip].dNdx_u;
491 (DisplacementDim == 2 ? 1 : 0),
492 ShapeFunction::NPOINTS * DisplacementDim);
494 ShapeFunction::NPOINTS>(
501 DisplacementDim, ShapeFunction::NPOINTS,
506 DisplacementDim, ShapeFunction::NPOINTS,
508 dNdx, N, grad_u, x_coord,
512 f_bar_variables ? std::get<double>(*f_bar_variables) : 0.0;
517 alpha, G, local_x_prev, x_position, t, dt, constitutive_setting,
524 for (
unsigned ip = 0; ip < n_integration_points; ip++)
531 const unsigned integration_point)
const override
536 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
540 static constexpr auto localDOF(std::vector<double>
const& x)
540 static constexpr auto localDOF(std::vector<double>
const& x) {
…}
547 std::vector<IpData, Eigen::aligned_allocator<IpData>>
_ip_data;
552 ShapeFunction::NPOINTS * DisplacementDim;
constexpr double getWeight() const
std::size_t getID() const
Returns the ID of the element.
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
unsigned getNumberOfPoints() const
void setElementID(std::size_t element_id)
MatrixType< _kelvin_vector_size, _number_of_dof > BMatrixType
MatrixType< _number_of_dof, _number_of_dof > StiffnessMatrixType
VectorType< _number_of_dof > NodalForceVectorType
Rhs residual.
VectorTypeFixedSize< MathLib::VectorizedTensor::size(DisplacementDim)> GradientVectorType
MatrixType< MathLib::VectorizedTensor::size(DisplacementDim), _number_of_dof > GradientMatrixType
constexpr int tensorSize(int dim)
See Tensor type for details.
Eigen::Matrix< double, Eigen::MatrixBase< Derived >::RowsAtCompileTime, 1 > symmetricTensorToKelvinVector(Eigen::MatrixBase< Derived > const &v)
Eigen::Matrix< double, 3, 3 > kelvinVectorToTensor(Eigen::Matrix< double, 4, 1, Eigen::ColMajor, 4, 1 > const &v)
double determinant(Eigen::MatrixBase< Derived > const &tensor)
Computes determinant of a vectorized tensor.
constexpr auto identity()
constexpr Eigen::CwiseNullaryOp< EigenBlockMatrixViewFunctor< D, M >, typename EigenBlockMatrixViewFunctor< D, M >::Matrix > eigenBlockMatrixView(const Eigen::MatrixBase< M > &matrix)
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)
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)
auto localDOF(ElementDOFVector const &x)
std::array< double, 3 > interpolateCoordinates(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.
BMatrixType computeBMatrix(DNDX_Type const &dNdx, N_Type const &N, GradientVectorType const &grad_u, const double radius, const bool is_axially_symmetric)
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > identityForF(bool const is_axially_symmetric)
std::tuple< double, VectorTypeForFbar > computeFBarInitialVariablesElementCenter(bool const compute_detF0_only, Eigen::VectorXd const &u, MeshLib::Element const &element, bool const is_axially_symmetric)
VectorTypeForFbar computeQVector(DNDX_Type const &dNdx, GradientVectorType const &F, bool const)
std::tuple< double, VectorTypeForFbar > computeFBarInitialVariablesAverage(std::vector< IpData, Eigen::aligned_allocator< IpData > > const &ip_data, bool const compute_detF0_only, Eigen::VectorXd const &u, NumLib::GenericIntegrationMethod const &integration_method, MeshLib::Element const &element, bool const is_axially_symmetric)
MathLib::KelvinVector::KelvinVectorType< DisplacementDim > compute2EPlusI(double const alpha_p2, MathLib::KelvinVector::KelvinVectorType< DisplacementDim > const &eps_bar, bool const is_axially_symmetric)
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
MatrixType< GlobalDim, ShapeFunction::NPOINTS > GlobalDimNodalMatrixType
VectorType< ShapeFunction::NPOINTS > NodalVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType