12 bool delta(
int i,
int j)
14 return (i == j) ? true :
false;
29 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 3, 1>
32 assert(pt.size() ==
size());
33 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 3, 1> res;
37 params_.
lambda_mu(0, 0, 0, pt(0).getValue(), pt(1).getValue(),
size() == 2 ? 0. : pt(2).getValue(), 0, 0, lambda, mu);
52 Eigen::Matrix<double, Eigen::Dynamic, 1> gradient;
61 compute_energy_aux_gradient_fast<3, 2>(data, gradient);
67 compute_energy_aux_gradient_fast<6, 2>(data, gradient);
73 compute_energy_aux_gradient_fast<10, 2>(data, gradient);
79 compute_energy_aux_gradient_fast<Eigen::Dynamic, 2>(data, gradient);
92 compute_energy_aux_gradient_fast<4, 3>(data, gradient);
98 compute_energy_aux_gradient_fast<10, 3>(data, gradient);
104 compute_energy_aux_gradient_fast<20, 3>(data, gradient);
110 compute_energy_aux_gradient_fast<Eigen::Dynamic, 3>(data, gradient);
121 const Eigen::MatrixXd &local_pts,
122 const Eigen::MatrixXd &displacement,
123 Eigen::MatrixXd &tensor)
const
126 assert(displacement.cols() == 1);
128 Eigen::MatrixXd displacement_grad(
size(),
size());
130 for (
long p = 0; p < local_pts.rows(); ++p)
136 const Eigen::MatrixXd def_grad = Eigen::MatrixXd::Identity(
size(),
size()) + displacement_grad;
137 const Eigen::MatrixXd FmT = def_grad.inverse().transpose();
139 const double J = def_grad.determinant();
140 const double tmp1 = mu - lambda * std::log(J);
141 for (
int i = 0, idx = 0; i <
size(); i++)
142 for (
int j = 0; j <
size(); j++)
143 for (
int k = 0; k <
size(); k++)
144 for (
int l = 0; l <
size(); l++)
146 tensor(p, idx) = mu * delta(i, k) * delta(j, l) + tmp1 * FmT(i, l) * FmT(k, j);
150 tensor.row(p) += lambda *
utils::flatten(FmT_vec * FmT_vec.transpose());
179 Eigen::MatrixXd hessian;
187 hessian.resize(6, 6);
189 compute_energy_hessian_aux_fast<3, 2>(data, hessian);
194 hessian.resize(12, 12);
196 compute_energy_hessian_aux_fast<6, 2>(data, hessian);
201 hessian.resize(20, 20);
203 compute_energy_hessian_aux_fast<10, 2>(data, hessian);
210 compute_energy_hessian_aux_fast<Eigen::Dynamic, 2>(data, hessian);
222 hessian.resize(12, 12);
224 compute_energy_hessian_aux_fast<4, 3>(data, hessian);
229 hessian.resize(30, 30);
231 compute_energy_hessian_aux_fast<10, 3>(data, hessian);
236 hessian.resize(60, 60);
238 compute_energy_hessian_aux_fast<20, 3>(data, hessian);
245 compute_energy_hessian_aux_fast<Eigen::Dynamic, 3>(data, hessian);
257 Eigen::MatrixXd &all,
258 const std::function<Eigen::MatrixXd(
const Eigen::MatrixXd &)> &fun)
const
260 const auto &displacement = data.
fun;
262 const auto &bs = data.
bs;
263 const auto &gbs = data.
gbs;
264 const auto el_id = data.
el_id;
265 const auto t = data.
t;
267 Eigen::MatrixXd displacement_grad(
size(),
size());
269 assert(displacement.cols() == 1);
271 all.resize(local_pts.rows(), all_size);
275 const auto I = Eigen::MatrixXd::Identity(
size(),
size());
277 for (
long p = 0; p < local_pts.rows(); ++p)
281 const Eigen::MatrixXd def_grad = I + displacement_grad;
284 all.row(p) = fun(def_grad);
287 const double J = def_grad.determinant();
288 const Eigen::MatrixXd b = def_grad * def_grad.transpose();
293 Eigen::MatrixXd stress_tensor = (lambda * std::log(J) * I + mu * (b - I)) / J;
299 all.row(p) = fun(stress_tensor);
310 return compute_energy_aux<double, 3, 2>(data);
312 return compute_energy_aux<double, 6, 2>(data);
314 return compute_energy_aux<double, 10, 2>(data);
316 return compute_energy_aux<double, Eigen::Dynamic, 2>(data);
325 return compute_energy_aux<double, 4, 3>(data);
327 return compute_energy_aux<double, 10, 3>(data);
329 return compute_energy_aux<double, 20, 3>(data);
331 return compute_energy_aux<double, Eigen::Dynamic, 3>(data);
337 template <
typename T,
int n_basis,
int dim>
340 if constexpr (std::is_same_v<T, double>)
343 local_disp.setZero();
347 for (
size_t ii = 0; ii < bs.global.size(); ++ii)
349 for (
int d = 0; d < dim; ++d)
351 local_disp(i, d) += bs.global[ii].val * data.
x(bs.global[ii].index *
size() + d);
356 Eigen::VectorXd jacs;
360 Eigen::Matrix<T, dim, dim> def_grad(
size(),
size());
364 const int n_pts = data.
da.size();
365 for (
long p = 0; p < n_pts; ++p)
374 const Eigen::Matrix<T, dim, dim> jac_it = data.
vals.
jac_it[p];
376 def_grad = (local_disp.transpose() * grad) * jac_it + Eigen::Matrix<T, dim, dim>::Identity(
size(),
size());
382 const T log_det_j = log(J);
383 const T
val = mu / 2 * (def_grad.squaredNorm() -
size() - 2 * log_det_j) +
384 lambda / 2 * log_det_j * log_det_j;
385 energy +=
val * data.
da(p);
391 typedef Eigen::Matrix<T, Eigen::Dynamic, 1>
AutoDiffVect;
401 const int n_pts = data.
da.size();
402 for (
long p = 0; p < n_pts; ++p)
407 for (
int d = 0; d < dim; ++d)
408 def_grad(d, d) += T(1);
414 const T
val = mu / 2 * ((def_grad.transpose() * def_grad).trace() -
size() - 2 * log_det_j) + lambda / 2 * log_det_j * log_det_j;
416 energy +=
val * data.
da(p);
423 Eigen::Matrix<double, dim, dim>
hat(
const Eigen::Matrix<double, dim, 1> &
x)
426 Eigen::Matrix<double, dim, dim> prod;
423 Eigen::Matrix<double, dim, dim>
hat(
const Eigen::Matrix<double, dim, 1> &
x) {
…}
440 Eigen::Matrix<double, dim, 1>
cross(
const Eigen::Matrix<double, dim, 1> &
x,
const Eigen::Matrix<double, dim, 1> &
y)
443 Eigen::Matrix<double, dim, 1>
z;
446 z(0) =
x(1) *
y(2) -
x(2) *
y(1);
447 z(1) =
x(2) *
y(0) -
x(0) *
y(2);
448 z(2) =
x(0) *
y(1) -
x(1) *
y(0);
440 Eigen::Matrix<double, dim, 1>
cross(
const Eigen::Matrix<double, dim, 1> &
x,
const Eigen::Matrix<double, dim, 1> &
y) {
…}
453 template <
int n_basis,
int dim>
456 assert(data.
x.cols() == 1);
458 const int n_pts = data.
da.size();
461 local_disp.setZero();
465 for (
size_t ii = 0; ii < bs.global.size(); ++ii)
467 for (
int d = 0; d < dim; ++d)
469 local_disp(i, d) += bs.global[ii].val * data.
x(bs.global[ii].index *
size() + d);
474 Eigen::VectorXd jacs;
478 Eigen::Matrix<double, dim, dim> def_grad(
size(),
size());
483 for (
long p = 0; p < n_pts; ++p)
492 Eigen::Matrix<double, dim, dim> jac_it = data.
vals.
jac_it[p];
493 Eigen::Matrix<double, n_basis, dim> delF_delU = grad * jac_it;
496 def_grad = local_disp.transpose() * delF_delU + Eigen::Matrix<double, dim, dim>::Identity(
size(),
size());
498 const double J =
use_robust_jacobian ? jacs(p) * jac_it.determinant() : def_grad.determinant();
499 const double log_det_j = log(J);
501 Eigen::Matrix<double, dim, dim> delJ_delF(
size(),
size());
507 delJ_delF(0, 0) = def_grad(1, 1);
508 delJ_delF(0, 1) = -def_grad(1, 0);
509 delJ_delF(1, 0) = -def_grad(0, 1);
510 delJ_delF(1, 1) = def_grad(0, 0);
516 Eigen::Matrix<double, dim, 1> u(def_grad.rows());
517 Eigen::Matrix<double, dim, 1> v(def_grad.rows());
518 Eigen::Matrix<double, dim, 1> w(def_grad.rows());
524 delJ_delF.col(0) = cross<dim>(v, w);
525 delJ_delF.col(1) = cross<dim>(w, u);
526 delJ_delF.col(2) = cross<dim>(u, v);
532 Eigen::Matrix<double, dim, dim> gradient_temp = mu * def_grad - mu * (1 / J) * delJ_delF + lambda * log_det_j * (1 / J) * delJ_delF;
533 Eigen::Matrix<double, n_basis, dim> gradient = delF_delU * gradient_temp.transpose();
535 double val = mu / 2 * ((def_grad.transpose() * def_grad).trace() -
size() - 2 * log_det_j) + lambda / 2 * log_det_j * log_det_j;
537 G.noalias() += gradient * data.
da(p);
540 Eigen::Matrix<double, dim, n_basis> G_T = G.transpose();
542 constexpr int N = (n_basis == Eigen::Dynamic) ? Eigen::Dynamic : n_basis * dim;
543 Eigen::Matrix<double, N, 1> temp(Eigen::Map<Eigen::Matrix<double, N, 1>>(G_T.data(), G_T.size()));
547 template <
int n_basis,
int dim>
550 assert(data.
x.cols() == 1);
552 constexpr int N = (n_basis == Eigen::Dynamic) ? Eigen::Dynamic : n_basis * dim;
553 const int n_pts = data.
da.size();
556 local_disp.setZero();
560 for (
size_t ii = 0; ii < bs.global.size(); ++ii)
562 local_disp.row(i) += bs.global[ii].val * data.
x.block<dim, 1>(bs.global[ii].index *
size(), 0).transpose();
566 Eigen::VectorXd jacs;
570 Eigen::Matrix<double, dim, dim> def_grad(
size(),
size());
572 for (
long p = 0; p < n_pts; ++p)
581 Eigen::Matrix<double, dim, dim> jac_it = data.
vals.
jac_it[p];
582 const Eigen::Matrix<double, n_basis, dim> delF_delU = grad * jac_it;
585 def_grad = local_disp.transpose() * delF_delU + Eigen::Matrix<double, dim, dim>::Identity(
size(),
size());
587 const double J =
use_robust_jacobian ? jacs(p) * jac_it.determinant() : def_grad.determinant();
588 double log_det_j = log(J);
590 Eigen::Matrix<double, dim, dim> delJ_delF(
size(),
size());
592 Eigen::Matrix<double, dim * dim, dim * dim> del2J_delF2(
size() *
size(),
size() *
size());
593 del2J_delF2.setZero();
597 delJ_delF(0, 0) = def_grad(1, 1);
598 delJ_delF(0, 1) = -def_grad(1, 0);
599 delJ_delF(1, 0) = -def_grad(0, 1);
600 delJ_delF(1, 1) = def_grad(0, 0);
602 del2J_delF2(0, 3) = 1;
603 del2J_delF2(1, 2) = -1;
604 del2J_delF2(2, 1) = -1;
605 del2J_delF2(3, 0) = 1;
607 else if (
size() == 3)
609 Eigen::Matrix<double, dim, 1> u(def_grad.rows());
610 Eigen::Matrix<double, dim, 1> v(def_grad.rows());
611 Eigen::Matrix<double, dim, 1> w(def_grad.rows());
617 delJ_delF.col(0) = cross<dim>(v, w);
618 delJ_delF.col(1) = cross<dim>(w, u);
619 delJ_delF.col(2) = cross<dim>(u, v);
621 del2J_delF2.template block<dim, dim>(0, 6) = hat<dim>(v);
622 del2J_delF2.template block<dim, dim>(6, 0) = -hat<dim>(v);
623 del2J_delF2.template block<dim, dim>(0, 3) = -hat<dim>(w);
624 del2J_delF2.template block<dim, dim>(3, 0) = hat<dim>(w);
625 del2J_delF2.template block<dim, dim>(3, 6) = -hat<dim>(u);
626 del2J_delF2.template block<dim, dim>(6, 3) = hat<dim>(u);
632 Eigen::Matrix<double, dim * dim, dim * dim>
id = Eigen::Matrix<double, dim * dim, dim * dim>::Identity(
size() *
size(),
size() *
size());
634 Eigen::Matrix<double, dim * dim, 1> g_j = Eigen::Map<const Eigen::Matrix<double, dim * dim, 1>>(delJ_delF.data(), delJ_delF.size());
636 Eigen::Matrix<double, dim * dim, dim * dim> hessian_temp = (mu * id) + (((mu + lambda * (1 - log_det_j)) / (J * J)) * (g_j * g_j.transpose())) + (((lambda * log_det_j - mu) / (J)) * del2J_delF2);
638 Eigen::Matrix<double, dim * dim, N> delF_delU_tensor = Eigen::Matrix<double, dim * dim, N>::Zero(jac_it.size(), grad.size());
640 for (
size_t i = 0; i < local_disp.rows(); ++i)
642 for (
size_t j = 0; j < local_disp.cols(); ++j)
644 Eigen::Matrix<double, dim, dim> temp(
size(),
size());
646 temp.row(j) = delF_delU.row(i);
647 Eigen::Matrix<double, dim * dim, 1> temp_flattened(Eigen::Map<Eigen::Matrix<double, dim * dim, 1>>(temp.data(), temp.size()));
648 delF_delU_tensor.col(i *
size() + j) = temp_flattened;
652 Eigen::Matrix<double, N, N> hessian = delF_delU_tensor.transpose() * hessian_temp * delF_delU_tensor;
654 double val = mu / 2 * ((def_grad.transpose() * def_grad).trace() -
size() - 2 * log_det_j) + lambda / 2 * log_det_j * log_det_j;
656 H += hessian * data.
da(p);
662 const Eigen::MatrixXd &mat,
663 Eigen::MatrixXd &stress,
664 Eigen::MatrixXd &result)
const
666 const double t = data.
t;
667 const int el_id = data.
el_id;
668 const Eigen::MatrixXd &local_pts = data.
local_pts;
669 const Eigen::MatrixXd &global_pts = data.
global_pts;
670 const Eigen::MatrixXd &grad_u_i = data.
grad_u_i;
675 Eigen::MatrixXd def_grad = Eigen::MatrixXd::Identity(grad_u_i.rows(), grad_u_i.cols()) + grad_u_i;
676 Eigen::MatrixXd FmT = def_grad.inverse().transpose();
678 stress = mu * (def_grad - FmT) + lambda * std::log(def_grad.determinant()) * FmT;
679 result = mu * mat + FmT * mat.transpose() * FmT * (mu - lambda * std::log(def_grad.determinant())) + lambda * (FmT.array() * mat.array()).sum() * FmT;
684 Eigen::MatrixXd &stress,
685 Eigen::MatrixXd &result)
const
687 const double t = data.
t;
688 const int el_id = data.
el_id;
689 const Eigen::MatrixXd &local_pts = data.
local_pts;
690 const Eigen::MatrixXd &global_pts = data.
global_pts;
691 const Eigen::MatrixXd &grad_u_i = data.
grad_u_i;
696 Eigen::MatrixXd def_grad = Eigen::MatrixXd::Identity(grad_u_i.rows(), grad_u_i.cols()) + grad_u_i;
697 Eigen::MatrixXd FmT = def_grad.inverse().transpose();
699 stress = mu * (def_grad - FmT) + lambda * std::log(def_grad.determinant()) * FmT;
700 result = mu * stress + FmT * stress.transpose() * FmT * (mu - lambda * std::log(def_grad.determinant())) + lambda * (FmT.array() * stress.array()).sum() * FmT;
705 const Eigen::MatrixXd &vect,
706 Eigen::MatrixXd &stress,
707 Eigen::MatrixXd &result)
const
709 const double t = data.
t;
710 const int el_id = data.
el_id;
711 const Eigen::MatrixXd &local_pts = data.
local_pts;
712 const Eigen::MatrixXd &global_pts = data.
global_pts;
713 const Eigen::MatrixXd &grad_u_i = data.
grad_u_i;
718 Eigen::MatrixXd def_grad = Eigen::MatrixXd::Identity(grad_u_i.rows(), grad_u_i.cols()) + grad_u_i;
719 Eigen::MatrixXd FmT = def_grad.inverse().transpose();
721 stress = mu * (def_grad - FmT) + lambda * std::log(def_grad.determinant()) * FmT;
723 if (vect.rows() == 1)
724 for (
int i = 0; i <
size(); ++i)
725 for (
int j = 0; j <
size(); ++j)
726 for (
int l = 0; l <
size(); ++l)
728 result(i *
size() + j, l) += mu * vect(i) * ((j == l) ? 1 : 0);
731 result(i *
size() + j, l) += (mu - lambda * std::log(def_grad.determinant())) * FmT(i, l) * (FmT.col(j).array() * vect.transpose().array()).sum();
732 result(i *
size() + j, l) += lambda * FmT(i, j) * (FmT.col(l).array() * vect.transpose().array()).sum();
735 for (
int i = 0; i <
size(); ++i)
736 for (
int j = 0; j <
size(); ++j)
737 for (
int k = 0; k <
size(); ++k)
739 result(i *
size() + j, k) += mu * vect(j) * ((i == k) ? 1 : 0);
742 result(i *
size() + j, k) += (mu - lambda * std::log(def_grad.determinant())) * FmT(k, i) * (FmT.row(j).array() * vect.transpose().array()).sum();
743 result(i *
size() + j, k) += lambda * FmT(i, j) * (FmT.row(k).array() * vect.transpose().array()).sum();
749 Eigen::MatrixXd &dstress_dmu,
750 Eigen::MatrixXd &dstress_dlambda)
const
752 const double t = data.
t;
753 const int el_id = data.
el_id;
754 const Eigen::MatrixXd &local_pts = data.
local_pts;
755 const Eigen::MatrixXd &global_pts = data.
global_pts;
756 const Eigen::MatrixXd &grad_u_i = data.
grad_u_i;
758 Eigen::MatrixXd def_grad = Eigen::MatrixXd::Identity(grad_u_i.rows(), grad_u_i.cols()) + grad_u_i;
759 Eigen::MatrixXd FmT = def_grad.inverse().transpose();
760 dstress_dmu = def_grad - FmT;
761 dstress_dlambda = std::log(def_grad.determinant()) * FmT;
766 std::map<std::string, ParamFunc> res;
773 params.lambda_mu(uv, p, t, e, lambda, mu);
780 params.lambda_mu(uv, p, t, e, lambda, mu);
786 params.lambda_mu(uv, p, t, e, lambda, mu);
789 return mu * (3.0 * lambda + 2.0 * mu) / (lambda + mu);
791 return 2 * mu * (2.0 * lambda + 2.0 * mu) / (lambda + 2.0 * mu);
797 params.lambda_mu(uv, p, t, e, lambda, mu);
800 return lambda / (2.0 * (lambda + mu));
802 return lambda / (lambda + 2.0 * mu);
ElementAssemblyValues vals
std::string stress() const
stores per element basis values at given quadrature points and geometric mapping
std::vector< AssemblyValues > basis_values
void compute(const int el_index, const bool is_volume, const Eigen::MatrixXd &pts, const basis::ElementBases &basis, const basis::ElementBases &gbasis)
computes the per element values at the local (ref el) points (pts) sets basis_values,...
std::vector< Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > > jac_it
Eigen::VectorXd eval_deformed_jacobian_determinant(const Eigen::VectorXd &disp) const
quadrature::Quadrature quadrature
void lambda_mu(double px, double py, double pz, double x, double y, double z, double t, int el_id, double &lambda, double &mu) const
void add_multimaterial(const int index, const json ¶ms, const bool is_volume, const std::string &stress_unit)
T compute_energy_aux(const NonLinearAssemblerData &data) const
double compute_energy(const NonLinearAssemblerData &data) const override
std::map< std::string, ParamFunc > parameters() const override
void compute_energy_aux_gradient_fast(const NonLinearAssemblerData &data, Eigen::VectorXd &G_flattened) const
VectorNd compute_rhs(const AutodiffHessianPt &pt) const override
void compute_stress_grad_multiply_stress(const OptAssemblerData &data, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const override
void compute_energy_hessian_aux_fast(const NonLinearAssemblerData &data, Eigen::MatrixXd &H) const
void compute_dstress_dmu_dlambda(const OptAssemblerData &data, Eigen::MatrixXd &dstress_dmu, Eigen::MatrixXd &dstress_dlambda) const override
void compute_stress_grad_multiply_vect(const OptAssemblerData &data, const Eigen::MatrixXd &vect, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const override
void compute_stress_grad_multiply_mat(const OptAssemblerData &data, const Eigen::MatrixXd &mat, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const override
void assign_stress_tensor(const OutputData &data, const int all_size, const ElasticityTensorType &type, Eigen::MatrixXd &all, const std::function< Eigen::MatrixXd(const Eigen::MatrixXd &)> &fun) const override
void add_multimaterial(const int index, const json ¶ms, const Units &units) override
Eigen::MatrixXd assemble_hessian(const NonLinearAssemblerData &data) const override
Eigen::VectorXd assemble_gradient(const NonLinearAssemblerData &data) const override
void compute_stiffness_value(const double t, const assembler::ElementAssemblyValues &vals, const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &displacement, Eigen::MatrixXd &tensor) const override
const ElementAssemblyValues & vals
const Eigen::MatrixXd & x
const QuadratureVector & da
const Eigen::MatrixXd & local_pts
const Eigen::MatrixXd & global_pts
const Eigen::MatrixXd & grad_u_i
const basis::ElementBases & bs
const Eigen::MatrixXd & fun
const Eigen::MatrixXd & local_pts
const basis::ElementBases & gbs
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > AutoDiffGradMat
Eigen::Matrix< double, dim, dim > hat(const Eigen::Matrix< double, dim, 1 > &x)
Eigen::Matrix< double, dim, 1 > cross(const Eigen::Matrix< double, dim, 1 > &x, const Eigen::Matrix< double, dim, 1 > &y)
Eigen::Matrix< T, Eigen::Dynamic, 1 > AutoDiffVect
void neo_hookean_3d_function(const AutodiffHessianPt &pt, const double lambda, const double mu, Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > &res)
void neo_hookean_2d_function(const AutodiffHessianPt &pt, const double lambda, const double mu, Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > &res)
Eigen::VectorXd flatten(const Eigen::MatrixXd &X)
Flatten rowwises.
T determinant(const Eigen::Matrix< T, rows, cols, option, maxRow, maxCol > &mat)
Eigen::MatrixXd pk2_from_cauchy(const Eigen::MatrixXd &stress, const Eigen::MatrixXd &F)
void compute_disp_grad_at_quad(const assembler::NonLinearAssemblerData &data, const AutoDiffVect &local_disp, const int p, const int size, AutoDiffGradMat &def_grad)
Eigen::Matrix< AutodiffScalarHessian, Eigen::Dynamic, 1, 0, 3, 1 > AutodiffHessianPt
void compute_diplacement_grad(const int size, const ElementAssemblyValues &vals, const Eigen::MatrixXd &local_pts, const int p, const Eigen::MatrixXd &displacement, Eigen::MatrixXd &displacement_grad)
void get_local_disp(const assembler::NonLinearAssemblerData &data, const int size, AutoDiffVect &local_disp)
Eigen::MatrixXd pk1_from_cauchy(const Eigen::MatrixXd &stress, const Eigen::MatrixXd &F)
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd