17#include <boost/math/constants/constants.hpp>
19constexpr double third = boost::math::constants::third<double>();
26 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
43 int nb_integration_pts = data.
getN().size1();
46 auto t_approx_P_adjoint_log_du =
47 getFTensor1FromMat<size_symm>(
dataAtPts->adjointPdUAtPts);
48 auto t_P = getFTensor2FromMat<3, 3>(
dataAtPts->PAtPts);
50 getFTensor2SymmetricFromMat<3>(
dataAtPts->logStretchDotTensorAtPts);
52 getFTensor4DdgFromMat<3, 3, 1>(
dataAtPts->diffStretchTensorAtPts);
58 auto get_ftensor2 = [](
auto &
v) {
60 &
v[0], &
v[1], &
v[2], &
v[3], &
v[4], &
v[5]);
63 int nb_base_functions = data.
getN().size2();
65 for (
int gg = 0; gg != nb_integration_pts; ++gg) {
67 auto t_nf = get_ftensor2(
nF);
70 t_Ldiff_u(
i,
j,
L) = t_diff_u(
i,
j,
k,
l) * t_L(
k,
l,
L);
80 a * (t_approx_P_adjoint_log_du(
L) - t_Ldiff_u(
i,
j,
L) * t_P(
i,
j) -
81 t_L(
i,
j,
L) * t_viscous_P(
i,
j));
84 for (; bb != nb_dofs / 6; ++bb) {
85 t_nf(
L) -= t_row_base_fun * t_residual(
L);
89 for (; bb != nb_base_functions; ++bb)
93 ++t_approx_P_adjoint_log_du;
103 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
const double alpha_u) {
109 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
110 boost::shared_ptr<ExternalStrainVec> external_strain_vec_ptr,
111 std::map<std::string, boost::shared_ptr<ScalingMethod>> smv) {
117 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
146 int nb_integration_pts = row_data.
getN().size1();
147 int row_nb_dofs = row_data.
getIndices().size();
148 int col_nb_dofs = col_data.
getIndices().size();
150 auto get_ftensor2 = [](
MatrixDouble &
m,
const int r,
const int c) {
154 &
m(r + 0,
c + 0), &
m(r + 0,
c + 1), &
m(r + 0,
c + 2), &
m(r + 0,
c + 3),
155 &
m(r + 0,
c + 4), &
m(r + 0,
c + 5),
157 &
m(r + 1,
c + 0), &
m(r + 1,
c + 1), &
m(r + 1,
c + 2), &
m(r + 1,
c + 3),
158 &
m(r + 1,
c + 4), &
m(r + 1,
c + 5),
160 &
m(r + 2,
c + 0), &
m(r + 2,
c + 1), &
m(r + 2,
c + 2), &
m(r + 2,
c + 3),
161 &
m(r + 2,
c + 4), &
m(r + 2,
c + 5),
163 &
m(r + 3,
c + 0), &
m(r + 3,
c + 1), &
m(r + 3,
c + 2), &
m(r + 3,
c + 3),
164 &
m(r + 3,
c + 4), &
m(r + 3,
c + 5),
166 &
m(r + 4,
c + 0), &
m(r + 4,
c + 1), &
m(r + 4,
c + 2), &
m(r + 4,
c + 3),
167 &
m(r + 4,
c + 4), &
m(r + 4,
c + 5),
169 &
m(r + 5,
c + 0), &
m(r + 5,
c + 1), &
m(r + 5,
c + 2), &
m(r + 5,
c + 3),
170 &
m(r + 5,
c + 4), &
m(r + 5,
c + 5)
184 auto get_dP = [&]() {
188 auto t_P = getFTensor2FromMat<3, 3>(
dataAtPts->PAtPts);
189 auto r_P_du = getFTensor4FromMat<3, 3, 3, 3>(
dataAtPts->P_du);
190 auto t_approx_P_adjoint__dstretch =
191 getFTensor2FromMat<3, 3>(
dataAtPts->adjointPdstretchAtPts);
193 getFTensor2SymmetricFromMat<3>(
dataAtPts->logStretchDotTensorAtPts);
194 auto t_u = getFTensor2SymmetricFromMat<3>(
dataAtPts->stretchTensorAtPts);
196 getFTensor4DdgFromMat<3, 3, 1>(
dataAtPts->diffStretchTensorAtPts);
197 auto t_eigen_vals = getFTensor1FromMat<3>(
dataAtPts->eigenVals);
198 auto t_eigen_vecs = getFTensor2FromMat<3, 3>(
dataAtPts->eigenVecs);
201 auto t_dP = getFTensor2FromMat<size_symm, size_symm>(
dP);
202 for (
auto gg = 0; gg != nb_integration_pts; ++gg) {
205 t_deltaP(
i,
j) = t_approx_P_adjoint__dstretch(
i,
j) - t_P(
i,
j);
208 t_Ldiff_u(
i,
j,
L) = t_diff_u(
i,
j,
m,
n) * t_L(
m,
n,
L);
210 -t_Ldiff_u(
i,
j,
L) * (r_P_du(
i,
j,
k,
l) * t_Ldiff_u(
k,
l,
J));
213 (t_L(
i,
j,
L) * (t_diff(
i,
j,
k,
l) * t_L(
k,
l,
J)));
218 t_deltaP_sym(
i,
j) = (t_deltaP(
i,
j) || t_deltaP(
j,
i));
219 t_deltaP_sym(
i,
j) /= 2.0;
223 t_dP(
L,
J) += t_L(
i,
j,
L) * (t_diff2_uP2(
i,
j,
k,
l) * t_L(
k,
l,
J));
229 ++t_approx_P_adjoint__dstretch;
237 return getFTensor2FromMat<size_symm, size_symm>(
dP);
240 int row_nb_base_functions = row_data.
getN().size2();
243 auto t_dP = get_dP();
244 for (
int gg = 0; gg != nb_integration_pts; ++gg) {
248 for (; rr != row_nb_dofs / 6; ++rr) {
250 auto t_m = get_ftensor2(
K, 6 * rr, 0);
251 for (
int cc = 0; cc != col_nb_dofs / 6; ++cc) {
252 const double b =
a * t_row_base_fun * t_col_base_fun;
253 t_m(
L,
J) -= b * t_dP(
L,
J);
260 for (; rr != row_nb_base_functions; ++rr) {
271 std::string row_field, std::string col_field,
272 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
const double alpha) {
277 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
278 boost::shared_ptr<double> total_energy_ptr) {
283 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
284 boost::shared_ptr<PhysicalEquations> physics_ptr) {
289 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
290 boost::shared_ptr<PhysicalEquations> physics_ptr) {
311 boost::shared_ptr<double> scale_ptr,
312 boost::shared_ptr<PhysicalEquations> physics_ptr) {
317 OpHMHH(
const int tag,
const bool eval_rhs,
const bool eval_lhs,
318 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
319 boost::shared_ptr<PhysicalEquations> physics_ptr)
320 :
OpJacobian(tag, eval_rhs, eval_lhs, data_ptr, physics_ptr) {}
335 const auto nb_integration_pts = getGaussPts().size2();
336 auto iu = getFTensor2SymmetricFromMat<3>(dataAtPts->stretchTensorAtPts);
337 auto t_grad_h1 = getFTensor2FromMat<3, 3>(dataAtPts->wGradH1AtPts);
339 auto create_data_vec = [nb_integration_pts](
auto &
v) {
340 v.resize(nb_integration_pts,
false);
342 auto create_data_mat = [nb_integration_pts](
auto &
m) {
343 m.resize(9, nb_integration_pts,
false);
346 create_data_mat(dataAtPts->PAtPts);
349 auto r_P = getFTensor2FromMat<3, 3>(dataAtPts->PAtPts);
350 for (
unsigned int gg = 0; gg != nb_integration_pts; ++gg) {
357 physicsPtr->get_h()(
i,
j) = iu(
i,
j);
362 physicsPtr->get_h()(
i,
j) = iu(
i,
k) * t_h1(
k,
j);
365 physicsPtr->get_h()(
i,
j) = iu(
i,
j);
369 int r = ::function(tAg, physicsPtr->dependentVariablesPiola.size(),
370 physicsPtr->activeVariables.size(),
371 &physicsPtr->activeVariables[0],
372 &physicsPtr->dependentVariablesPiola[0]);
375 "ADOL-C function evaluation with error r = %d", r);
380 r_P(
i,
j) = physicsPtr->get_P()(
i,
j);
386 r_P(
k,
l) = physicsPtr->get_P()(
i,
j) * t_h1_du(
i,
j,
k,
l);
389 r_P(
i,
j) = physicsPtr->get_P()(
i,
j);
411 const int number_of_active_variables = physicsPtr->activeVariables.size();
412 const int number_of_dependent_variables =
413 physicsPtr->dependentVariablesPiola.size();
414 std::vector<double *> jac_ptr(number_of_dependent_variables);
415 for (
unsigned int n = 0;
n != number_of_dependent_variables; ++
n) {
418 ->dependentVariablesPiolaDirevatives[
n *
419 number_of_active_variables]);
422 const auto nb_integration_pts = getGaussPts().size2();
424 auto create_data_mat = [nb_integration_pts](
auto &
m) {
425 m.resize(9, nb_integration_pts,
false);
428 dataAtPts->P_du.resize(81, nb_integration_pts,
false);
430 auto iu = getFTensor2SymmetricFromMat<3>(dataAtPts->stretchTensorAtPts);
431 auto t_grad_h1 = getFTensor2FromMat<3, 3>(dataAtPts->wGradH1AtPts);
432 auto r_P_du = getFTensor4FromMat<3, 3, 3, 3>(dataAtPts->P_du);
440 for (
unsigned int gg = 0; gg != nb_integration_pts; ++gg) {
447 physicsPtr->get_h()(
i,
j) = iu(
i,
k) * t_h1(
k,
j);
452 physicsPtr->get_h()(
i,
j) = iu(
i,
k) * t_h1(
k,
j);
455 physicsPtr->get_h()(
i,
j) = iu(
i,
j);
460 int r = ::jacobian(tAg, number_of_dependent_variables,
461 number_of_active_variables,
462 &physicsPtr->activeVariables[0], &jac_ptr[0]);
465 "ADOL-C function evaluation with error");
469 t_P_dh_tmp(
i,
j, N0,
k) = physicsPtr->get_P_dh0()(
i,
j,
k);
470 t_P_dh_tmp(
i,
j, N1,
k) = physicsPtr->get_P_dh1()(
i,
j,
k);
471 t_P_dh_tmp(
i,
j, N2,
k) = physicsPtr->get_P_dh2()(
i,
j,
k);
478 (t_P_dh_tmp(
i,
j, o, p) * t_h1_du(o, p,
m,
n)) * t_h1_du(
i,
j,
k,
l);
484 (t_P_dh_tmp(
i,
j, o, p) * t_h1_du(o, p,
m,
n)) * t_h1_du(
i,
j,
k,
l);
488 r_P_du(
i,
j,
m,
n) = t_P_dh_tmp(
i,
j,
m,
n);
500 const int tag,
const bool eval_rhs,
const bool eval_lhs,
501 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
502 boost::shared_ptr<PhysicalEquations> physics_ptr) {
503 return (
new OpHMHH(tag, eval_rhs, eval_lhs, data_ptr, physics_ptr));
516 const int tape,
const double lambda,
const double mu,
517 const double sigma_y) {
519 physicalEquations = boost::make_shared<HMHStVenantKirchhoff>(
lambda,
mu);
520 CHKERR physicalEquations->recordTape(tape,
nullptr);
525 const int tape,
const double alpha,
const double beta,
const double lambda,
526 const double sigma_y) {
529 boost::make_shared<HMHPMooneyRivlinWriggersEq63>(alpha, beta,
lambda);
530 CHKERR physicalEquations->recordTape(tape,
nullptr);
537 physicalEquations = boost::make_shared<HMHNeohookean>(mField, c10,
K);
538 CHKERR physicalEquations->recordTape(tape,
nullptr);
544 physicalEquations = boost::make_shared<HMHHencky>(mField,
E, nu);
Auxilary functions for Eshelbian plasticity.
Eshelbian plasticity interface.
Implementation StVenantKirchhoff.
Implementation of NeoHookean material.
Implement of MooneyRivlinWriggersEq63.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_OPERATION_UNSUCCESSFUL
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
FTensor::Index< 'i', SPACE_DIM > i
const double c
speed of light (cm/ns)
const double v
phase velocity of light in medium (cm/ns)
const double n
refractive index of diffusive medium
FTensor::Index< 'J', DIM1 > J
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
auto getDiffDiffMat(A &&t_val, B &&t_vec, Fun< double > f, Fun< double > d_f, Fun< double > dd_f, C &&t_S, const int nb)
Get the Diff Diff Mat object.
ForcesAndSourcesCore::UserDataOperator UserDataOperator
static constexpr auto size_symm
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
constexpr auto field_name
FTensor::Index< 'm', 3 > m
static enum StretchSelector stretchSelector
MoFEMErrorCode addMaterial_Hencky(double E, double nu)
MoFEMErrorCode addMaterial_HMHHStVenantKirchhoff(const int tape, const double lambda, const double mu, const double sigma_y)
MoFEMErrorCode addMaterial_HMHMooneyRivlin(const int tape, const double alpha, const double beta, const double lambda, const double sigma_y)
static enum RotSelector gradApproximator
static boost::function< double(const double)> f
static boost::function< double(const double)> dd_f
static boost::function< double(const double)> d_f
MoFEMErrorCode addMaterial_HMHNeohookean(const int tape, const double c10, const double K)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< double > scalePtr
OpGetScale(boost::shared_ptr< double > scale_ptr)
OpHMHH(const int tag, const bool eval_rhs, const bool eval_lhs, boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
MoFEMErrorCode evaluateLhs(EntData &data)
MoFEMErrorCode evaluateRhs(EntData &data)
MoFEMErrorCode integratePiola(EntData &row_data, EntData &col_data)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialPhysical_du_du(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha)
OpSpatialPhysical(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha_u)
MoFEMErrorCode integrate(EntData &data)
virtual VolUserDataOperator * returnOpCalculateVarStretchFromStress(boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
virtual VolUserDataOperator * returnOpCalculateEnergy(boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< double > total_energy_ptr)
virtual VolUserDataOperator * returnOpSpatialPhysical(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha_u)
virtual VolUserDataOperator * returnOpSpatialPhysicalExternalStrain(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< ExternalStrainVec > external_strain_vec_ptr, std::map< std::string, boost::shared_ptr< ScalingMethod > > smv)
virtual VolUserDataOperator * returnOpCalculateStretchFromStress(boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
virtual VolUserDataOperator * returnOpSetScale(boost::shared_ptr< double > scale_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
virtual VolUserDataOperator * returnOpSpatialPhysical_du_du(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha)
virtual UserDataOperator * returnOpJacobian(const int tag, const bool eval_rhs, const bool eval_lhs, boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
bool sYmm
If true assume that matrix is symmetric structure.
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorInt & getIndices() const
Get global indices of degrees of freedom on entity.
auto getFTensor0IntegrationWeight()
Get integration weights.
@ OPROW
operator doWork function is executed on FE rows
@ OPROWCOL
operator doWork is executed on FE rows &columns
@ OPSPACE
operator do Work is execute on space data
double getVolume() const
element volume (linear geometry)
MatrixDouble K
local tangent matrix
VectorDouble nF
local right hand side vector
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
data at integration pts