77 "extract block data failed");
79#ifdef NEOHOOKEAN_SCALING
94 "Stretch selector is not equal to LOG");
98 "Exponent base is not equal to exp(1)");
106 if (b.blockEnts.find(ent) != b.blockEnts.end()) {
107 return std::make_pair(b.c10, b.K);
112 "Block not found for entity handle. If you mat set "
113 "block, set it to all elements");
119 boost::shared_ptr<PhysicalEquations> physics_ptr)
123 MoFEMErrorCode
doWork(
int side, EntityType type,
124 EntitiesFieldData::EntData &data) {
126 auto neohookean_ptr =
127 boost::dynamic_pointer_cast<HMHNeohookean>(
physicsPtr);
128 *(
scalePtr) = neohookean_ptr->eqScaling;
139 boost::shared_ptr<PhysicalEquations> physics_ptr) {
141 return new OpGetScale(scale_ptr, physics_ptr);
145 using EshelbianPlasticity::OpJacobian::OpJacobian;
152 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
153 boost::shared_ptr<PhysicalEquations> physics_ptr) {
156 data_ptr, physics_ptr);
158 return (
new OpJacobian(tag, eval_rhs, eval_lhs, data_ptr, physics_ptr));
164 PetscOptionsBegin(PETSC_COMM_WORLD,
"neo_hookean_",
"",
"none");
174 CHKERR PetscOptionsScalar(
"-viscosity_alpha_grad_u",
"viscosity",
"",
198 (boost::format(
"%s(.*)") %
"MAT_NEOHOOKEAN").str()
210 for (
auto m : meshset_vec_ptr) {
212 std::vector<double> block_data;
213 CHKERR m->getAttributes(block_data);
214 if (block_data.size() < 2) {
216 "Expected that block has atleast two attributes");
218 auto get_block_ents = [&]() {
224 double c10 = block_data[0];
225 double K = block_data[1];
227 blockData.push_back({c10, K, get_block_ents()});
229 MOFEM_LOG(
"EP", sev) <<
"MatBlock Neo-Hookean c10 = "
231 <<
" K = " <<
blockData.back().K <<
" nb ents. = "
252 ih(
i,
j) = (*t_h_ptr)(
i,
j);
259 enableMinMaxUsingAbs();
292 "gradApproximator not handled");
301 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
302 const double alpha_u);
312 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
313 const double alpha_u) {
326 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
327 boost::shared_ptr<ExternalStrainVec> external_strain_vec_ptr,
328 std::map<std::string, boost::shared_ptr<ScalingMethod>> smv);
339 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
340 boost::shared_ptr<ExternalStrainVec> external_strain_vec_ptr,
341 std::map<std::string, boost::shared_ptr<ScalingMethod>> smv) {
343 external_strain_vec_ptr, smv);
348 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
357 std::string row_field, std::string col_field,
358 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
const double alpha) {
361 row_field, col_field, data_ptr, alpha);
395 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
const double alpha_u)
400 double circumcenter[3],
double *xi,
double *
eta,
407 auto neohookean_ptr =
408 boost::dynamic_pointer_cast<HMHNeohookean>(dataAtPts->physicsPtr);
409 if (!neohookean_ptr) {
411 "Pointer to HMHNeohookean is null");
413 auto [def_c10, def_K] =
414 neohookean_ptr->getMaterialParameters(getFEEntityHandle());
416 double c10 = def_c10 / neohookean_ptr->eqScaling;
417 double alpha_u = alphaU / neohookean_ptr->eqScaling;
418 double K = def_K / neohookean_ptr->eqScaling;
420 auto time_step = getTStimeStep();
421 double alpha_grad_u =
422 neohookean_ptr->alphaGradU / neohookean_ptr->eqScaling;
431 int nb_integration_pts = data.
getN().size1();
432 auto v = getVolume();
433 auto t_w = getFTensor0IntegrationWeight();
434 auto t_approx_P_adjoint_log_du =
435 getFTensor1FromMat<size_symm>(dataAtPts->adjointPdUAtPts);
436 auto t_u = getFTensor2SymmetricFromMat<3>(dataAtPts->stretchTensorAtPts);
437 auto t_grad_h1 = getFTensor2FromMat<3, 3>(dataAtPts->wGradH1AtPts);
439 getFTensor2SymmetricFromMat<3>(dataAtPts->logStretchDotTensorAtPts);
441 getFTensor4DdgFromMat<3, 3, 1>(dataAtPts->diffStretchTensorAtPts);
443 getFTensor2SymmetricFromMat<3>(dataAtPts->logStretchTensorAtPts);
445 getFTensor2FromMat<6, 3>(dataAtPts->gradLogStretchDotTensorAtPts);
447 getFTensor2SymmetricFromMat<3>(dataAtPts->logStretch2H1AtPts);
449 auto t_eigen_vals = getFTensor1FromMat<3>(dataAtPts->eigenVals);
450 auto t_eigen_vecs = getFTensor2FromMat<3, 3>(dataAtPts->eigenVecs);
451 auto &nbUniq = dataAtPts->nbUniq;
464 auto get_ftensor2 = [](
auto &
v) {
466 &
v[0], &
v[1], &
v[2], &
v[3], &
v[4], &
v[5]);
469 int nb_base_functions = data.
getN().size2();
476 for (
int gg = 0; gg != nb_integration_pts; ++gg) {
481 auto t_neohookean_hencky =
487 const double tr = t_log_u(
i,
i);
489 t_P(L) = t_L(
i,
j, L) * (t_neohookean_hencky(
i,
j) +
492 t_viscous_P(L) = alpha_u * (t_L(
i,
j, L) * t_dot_log_u(
i,
j));
495 t_residual(L) = t_approx_P_adjoint_log_du(L) - t_P(L) - t_viscous_P(L);
499 t_grad_residual(L,
i) = alpha_grad_u * t_grad_log_u(L,
i);
500 t_grad_residual(L,
i) *=
a;
502 ++t_approx_P_adjoint_log_du;
507 auto t_nf = getFTensor1FromPtr<size_symm>(&*nF.data().begin());
509 for (; bb != nb_dofs /
size_symm; ++bb) {
510 t_nf(L) -= t_row_base_fun * t_residual(L);
511 t_nf(L) += t_grad_base_fun(
i) * t_grad_residual(L,
i);
516 for (; bb != nb_base_functions; ++bb) {
528 "Not implemented for Neo-Hookean (used ADOL-C)");
542 "gradApproximator not handled");
550 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
551 boost::shared_ptr<ExternalStrainVec> external_strain_vec_ptr,
552 std::map<std::string, boost::shared_ptr<ScalingMethod>> smv)
554 externalStrainVecPtr(external_strain_vec_ptr), scalingMethodsMap{smv} {}
560 auto neohookean_ptr =
561 boost::dynamic_pointer_cast<HMHNeohookean>(dataAtPts->physicsPtr);
562 if (!neohookean_ptr) {
564 "Pointer to HMHNeohookean is null");
575 for (
auto &ext_strain_block : (*externalStrainVecPtr)) {
576 auto block_name =
"(.*)ANALYTICAL_EXTERNALSTRAIN(.*)";
577 std::regex reg_name(block_name);
578 if (std::regex_match(ext_strain_block.blockName, reg_name)) {
580 "Analytical external strain not implemented for Neo-Hookean "
584 if (ext_strain_block.ents.find(fe_ent) != ext_strain_block.ents.end()) {
586 if (scalingMethodsMap.find(ext_strain_block.blockName) !=
587 scalingMethodsMap.end()) {
589 scalingMethodsMap.at(ext_strain_block.blockName)->getScale(time);
592 <<
"No scaling method found for " << ext_strain_block.blockName;
595 auto [def_c10, def_K] =
596 neohookean_ptr->getMaterialParameters(getFEEntityHandle());
597 double c10 = def_c10 / neohookean_ptr->eqScaling;
599 double external_strain_val =
scale * ext_strain_block.val;
600 double K = ext_strain_block.bulkModulusK;
607 int nb_integration_pts = data.
getN().size1();
608 auto v = getVolume();
609 auto t_w = getFTensor0IntegrationWeight();
618 int nb_base_functions = data.
getN().size2();
621 const double tr = 3.0 * external_strain_val;
622 const double sigma_J = K * tr;
624 for (
int gg = 0; gg != nb_integration_pts; ++gg) {
630 t_residual(L) += (t_L(
i,
j, L) *
t_kd(
i,
j)) * sigma_J;
633 auto t_nf = getFTensor1FromPtr<size_symm>(&*nF.data().begin());
635 for (; bb != nb_dofs /
size_symm; ++bb) {
636 t_nf(L) -= t_row_base_fun * t_residual(L);
640 for (; bb != nb_base_functions; ++bb) {
650 std::string row_field, std::string col_field,
651 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
const double alpha)
663 auto neohookean_ptr =
664 boost::dynamic_pointer_cast<HMHNeohookean>(dataAtPts->physicsPtr);
665 if (!neohookean_ptr) {
667 "Pointer to HMHNeohookean is null");
669 auto [def_c10, def_K] =
670 neohookean_ptr->getMaterialParameters(getFEEntityHandle());
672 double c10 = def_c10 / neohookean_ptr->eqScaling;
673 double alpha_u = alphaU / neohookean_ptr->eqScaling;
674 double lambda = def_K / neohookean_ptr->eqScaling;
676 auto time_step = getTStimeStep();
677 double alpha_grad_u =
678 neohookean_ptr->alphaGradU / neohookean_ptr->eqScaling;
689 int nb_integration_pts = row_data.
getN().size1();
690 int row_nb_dofs = row_data.
getIndices().size();
691 int col_nb_dofs = col_data.
getIndices().size();
693 auto get_ftensor2 = [](MatrixDouble &
m,
const int r,
const int c) {
697 &
m(r + 0,
c + 0), &
m(r + 0,
c + 1), &
m(r + 0,
c + 2), &
m(r + 0,
c + 3),
698 &
m(r + 0,
c + 4), &
m(r + 0,
c + 5),
700 &
m(r + 1,
c + 0), &
m(r + 1,
c + 1), &
m(r + 1,
c + 2), &
m(r + 1,
c + 3),
701 &
m(r + 1,
c + 4), &
m(r + 1,
c + 5),
703 &
m(r + 2,
c + 0), &
m(r + 2,
c + 1), &
m(r + 2,
c + 2), &
m(r + 2,
c + 3),
704 &
m(r + 2,
c + 4), &
m(r + 2,
c + 5),
706 &
m(r + 3,
c + 0), &
m(r + 3,
c + 1), &
m(r + 3,
c + 2), &
m(r + 3,
c + 3),
707 &
m(r + 3,
c + 4), &
m(r + 3,
c + 5),
709 &
m(r + 4,
c + 0), &
m(r + 4,
c + 1), &
m(r + 4,
c + 2), &
m(r + 4,
c + 3),
710 &
m(r + 4,
c + 4), &
m(r + 4,
c + 5),
712 &
m(r + 5,
c + 0), &
m(r + 5,
c + 1), &
m(r + 5,
c + 2), &
m(r + 5,
c + 3),
713 &
m(r + 5,
c + 4), &
m(r + 5,
c + 5)
725 auto v = getVolume();
726 auto ts_a = getTSa();
727 auto t_w = getFTensor0IntegrationWeight();
729 int row_nb_base_functions = row_data.
getN().size2();
733 auto t_grad_h1 = getFTensor2FromMat<3, 3>(dataAtPts->wGradH1AtPts);
735 getFTensor4DdgFromMat<3, 3, 1>(dataAtPts->diffStretchTensorAtPts);
737 getFTensor2SymmetricFromMat<3>(dataAtPts->logStretchTensorAtPts);
739 getFTensor2SymmetricFromMat<3>(dataAtPts->logStretch2H1AtPts);
740 auto t_u = getFTensor2SymmetricFromMat<3>(dataAtPts->stretchTensorAtPts);
741 auto t_approx_P_adjoint__dstretch =
742 getFTensor2FromMat<3, 3>(dataAtPts->adjointPdstretchAtPts);
743 auto t_eigen_vals = getFTensor1FromMat<3>(dataAtPts->eigenVals);
744 auto t_eigen_vecs = getFTensor2FromMat<3, 3>(dataAtPts->eigenVecs);
745 auto &nbUniq = dataAtPts->nbUniq;
752 for (
int gg = 0; gg != nb_integration_pts; ++gg) {
757 auto d_neohookean = [c10,
lambda](
auto v) {
762 t_eigen_vals, t_eigen_vecs, neohookean, d_neohookean, t_nb_uniq);
764 const auto tr = t_log_u(
i,
i);
766 t_dP(L,
J) = -t_L(
i,
j, L) * ((t_diff_neohookean(
i,
j,
k,
l) +
768 t_kd_sym(
i,
j) * t_kd_sym(
k,
l)) *
770 t_dP(L,
J) -= (alpha_u * ts_a) *
771 (t_L(
i,
j, L) * (t_diff(
i,
j,
k,
l) * t_L(
k,
l,
J)));
775 t_deltaP(
i,
j) = (t_approx_P_adjoint__dstretch(
i,
j) ||
776 t_approx_P_adjoint__dstretch(
j,
i)) /
781 t_dP(L,
J) += t_L(
i,
j, L) * (t_diff2_uP(
i,
j,
k,
l) * t_L(
k,
l,
J));
783 ++t_approx_P_adjoint__dstretch;
790 for (; rr != row_nb_dofs /
size_symm; ++rr) {
794 auto t_m = get_ftensor2(K, 6 * rr, 0);
795 for (
int cc = 0; cc != col_nb_dofs /
size_symm; ++cc) {
796 double b =
a * t_row_base_fun * t_col_base_fun;
797 double c = (
a * alpha_grad_u * ts_a) *
798 (t_row_grad_fun(
i) * t_col_grad_fun(
i));
799 t_m(L,
J) -= b * t_dP(L,
J);
800 t_m(L,
J) +=
c * t_kd_sym(L,
J);
810 for (; rr != row_nb_base_functions; ++rr) {
822 "Not implemented for Neo-Hookean (used ADOL-C)");
836 "gradApproximator not handled");
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
#define FTENSOR_INDEX(DIM, I)
Kronecker Delta class symmetric.
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
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 getMat(A &&t_val, B &&t_vec, Fun< double > f)
Get the Mat object.
auto getDiffMat(A &&t_val, B &&t_vec, Fun< double > f, Fun< double > d_f, const int nb)
Get the Diff Mat object.
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.
void tetcircumcenter_tp(double a[3], double b[3], double c[3], double d[3], double circumcenter[3], double *xi, double *eta, double *zeta)
ForcesAndSourcesCore::UserDataOperator UserDataOperator
static constexpr auto size_symm
constexpr IntegrationType I
constexpr auto field_name
FTensor::Index< 'm', 3 > m
static enum StretchSelector stretchSelector
static double dynamicTime
static enum RotSelector gradApproximator
static PetscBool dynamicRelaxation
static double exponentBase
static boost::function< double(const double)> f
static boost::function< double(const double)> dd_f
static boost::function< double(const double)> d_f
OpGetScale(boost::shared_ptr< double > scale_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
boost::shared_ptr< double > scalePtr
boost::shared_ptr< PhysicalEquations > physicsPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
MoFEMErrorCode evaluateRhs(EntData &data)
MoFEMErrorCode evaluateLhs(EntData &data)
MoFEMErrorCode integrate(EntData &data)
boost::shared_ptr< ExternalStrainVec > externalStrainVecPtr
std::map< std::string, boost::shared_ptr< ScalingMethod > > scalingMethodsMap
OpSpatialPhysicalExternalStrain(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)
OpSpatialPhysical_du_du(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialPhysical(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha_u)
MoFEMErrorCode integrate(EntData &data)
static double fun_diff_neohookean_bulk(double K, double tr)
Definition of derivative of axiator of Neo-hookean function.
static double fun_d_neohookean(double c10, double v)
Definition of derivative of Neo-hookean function.
static double fun_neohookean_bulk(double K, double tr)
Definition of axiator of Neo-hookean function.
virtual VolUserDataOperator * returnOpSpatialPhysical(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha_u)
std::vector< BlockData > blockData
MoFEMErrorCode extractBlockData(std::vector< const CubitMeshSets * > meshset_vec_ptr, Sev sev)
MoFEMErrorCode extractBlockData(Sev sev)
static constexpr int numberOfDependentVariables
VolUserDataOperator * returnOpSetScale(boost::shared_ptr< double > scale_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
static double fun_neohookean(double c10, double v)
Definition of Neo-hookean function.
MoFEMErrorCode getOptions()
MoFEM::Interface & mField
auto getMaterialParameters(EntityHandle ent)
static constexpr int numberOfActiveVariables
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)
HMHNeohookean(MoFEM::Interface &m_field, const double c10, const double K)
VolUserDataOperator * returnOpSpatialPhysical_du_du(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha)
MoFEMErrorCode recordTape(const int tape, DTensor2Ptr *t_h_ptr)
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)
virtual VolUserDataOperator * returnOpSpatialPhysical(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha_u)
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)
virtual moab::Interface & get_moab()=0
bool sYmm
If true assume that matrix is symmetric structure.
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
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.
EntityHandle getFEEntityHandle() const
Return finite element entity handle.
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
PetscReal ts_t
Current time value.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
double zeta
Viscous hardening.