v0.13.0
VEC-3: Nonlinear dynamic elastic
Note
Prerequisites of this tutorial include MSH-1: Create a 2D mesh from Gmsh


Note
Intended learning outcome:
  • general structure of a program developed using MoFEM
  • idea of Simple Interface in MoFEM and how to use it
  • second order equation in time
  • idea of domain element in MoFEM and how to use it
  • Use default forms integrals
  • how to push the developed UDOs to the Pipeline
  • utilisation of tools to convert outputs (MOAB) and visualise them (Paraview)
/**
* \file dynamic_elastic.cpp
* \example dynamic_elastic.cpp
*
* Plane stress elastic dynamic problem
*
*/
/* This file is part of MoFEM.
* MoFEM is free software: you can redistribute it and/or modify it under
* the terms of the GNU Lesser General Public License as published by the
* Free Software Foundation, either version 3 of the License, or (at your
* option) any later version.
*
* MoFEM is distributed in the hope that it will be useful, but WITHOUT
* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
* FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
* License for more details.
*
* You should have received a copy of the GNU Lesser General Public
* License along with MoFEM. If not, see <http://www.gnu.org/licenses/>. */
#include <MoFEM.hpp>
using namespace MoFEM;
template <int DIM> struct ElementsAndOps {};
template <> struct ElementsAndOps<2> {
using DomainEle = FaceElementForcesAndSourcesCoreBase;
};
template <> struct ElementsAndOps<3> {
using DomainEle = VolumeElementForcesAndSourcesCoreBase;
};
constexpr int SPACE_DIM =
EXECUTABLE_DIMENSION; //< Space dimension of problem, mesh
using OpK = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::BiLinearForm<
GAUSS>::OpGradTensorGrad<1, SPACE_DIM, SPACE_DIM, 1>;
using OpMass = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::BiLinearForm<
using OpInternalForce = FormsIntegrators<DomainEleOp>::Assembly<
using OpBodyForce = FormsIntegrators<DomainEleOp>::Assembly<PETSC>::LinearForm<
GAUSS>::OpSource<1, SPACE_DIM>;
using OpInertiaForce = FormsIntegrators<DomainEleOp>::Assembly<
constexpr bool is_quasi_static = false;
constexpr double rho = 1;
constexpr double omega = 2.4;
constexpr double young_modulus = 1;
constexpr double poisson_ratio = 0.25;
constexpr double bulk_modulus_K = young_modulus / (3 * (1 - 2 * poisson_ratio));
constexpr double shear_modulus_G = young_modulus / (2 * (1 + poisson_ratio));
#include <HenckyOps.hpp>
using namespace Tutorial;
using namespace HenckyOps;
static double *ts_time_ptr;
static double *ts_aa_ptr;
struct Example {
Example(MoFEM::Interface &m_field) : mField(m_field) {}
MoFEMErrorCode runProblem();
private:
MoFEMErrorCode readMesh();
MoFEMErrorCode setupProblem();
MoFEMErrorCode createCommonData();
MoFEMErrorCode boundaryCondition();
MoFEMErrorCode assembleSystem();
MoFEMErrorCode solveSystem();
MoFEMErrorCode outputResults();
MoFEMErrorCode checkResults();
boost::shared_ptr<MatrixDouble> matGradPtr;
boost::shared_ptr<MatrixDouble> matStrainPtr;
boost::shared_ptr<MatrixDouble> matStressPtr;
boost::shared_ptr<MatrixDouble> matAccelerationPtr;
boost::shared_ptr<MatrixDouble> matInertiaPtr;
boost::shared_ptr<MatrixDouble> matDPtr;
boost::shared_ptr<MatrixDouble> matTangentPtr;
};
//! [Create common data]
auto set_matrial_stiffness = [&]() {
constexpr double A =
(SPACE_DIM == 2) ? 2 * shear_modulus_G /
: 1;
auto t_D = getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, 0>(*matDPtr);
t_D(i, j, k, l) = 2 * shear_modulus_G * ((t_kd(i, k) ^ t_kd(j, l)) / 4.) +
A * (bulk_modulus_K - (2. / 3.) * shear_modulus_G) *
t_kd(i, j) * t_kd(k, l);
};
matGradPtr = boost::make_shared<MatrixDouble>();
matStrainPtr = boost::make_shared<MatrixDouble>();
matStressPtr = boost::make_shared<MatrixDouble>();
matAccelerationPtr = boost::make_shared<MatrixDouble>();
matInertiaPtr = boost::make_shared<MatrixDouble>();
matDPtr = boost::make_shared<MatrixDouble>();
matTangentPtr = boost::make_shared<MatrixDouble>();
constexpr auto size_symm = (SPACE_DIM * (SPACE_DIM + 1)) / 2;
matDPtr->resize(size_symm * size_symm, 1);
CHKERR set_matrial_stiffness();
}
//! [Create common data]
//! [Run problem]
CHKERR readMesh();
CHKERR setupProblem();
CHKERR createCommonData();
CHKERR boundaryCondition();
CHKERR assembleSystem();
CHKERR solveSystem();
CHKERR outputResults();
CHKERR checkResults();
}
//! [Run problem]
//! [Read mesh]
auto simple = mField.getInterface<Simple>();
CHKERR simple->getOptions();
CHKERR simple->loadFile();
}
//! [Read mesh]
//! [Set up problem]
Simple *simple = mField.getInterface<Simple>();
// Add field
int order = 2;
CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
CHKERR simple->setFieldOrder("U", order);
CHKERR simple->setUp();
}
//! [Set up problem]
//! [Boundary condition]
auto simple = mField.getInterface<Simple>();
auto bc_mng = mField.getInterface<BcManager>();
CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "FIX_X",
"U", 0, 0);
CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "FIX_Y",
"U", 1, 1);
CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "FIX_Z",
"U", 2, 2);
CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "FIX_ALL",
"U", 0, 3);
}
//! [Boundary condition]
//! [Push operators to pipeline]
auto *simple = mField.getInterface<Simple>();
auto *pipeline_mng = mField.getInterface<PipelineManager>();
if (SPACE_DIM == 2) {
auto det_ptr = boost::make_shared<VectorDouble>();
auto jac_ptr = boost::make_shared<MatrixDouble>();
auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
pipeline_mng->getOpDomainLhsPipeline().push_back(
new OpCalculateHOJacForFace(jac_ptr));
pipeline_mng->getOpDomainLhsPipeline().push_back(
new OpInvertMatrix<2>(jac_ptr, det_ptr, inv_jac_ptr));
pipeline_mng->getOpDomainLhsPipeline().push_back(
new OpSetInvJacH1ForFace(inv_jac_ptr));
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculateHOJacForFace(jac_ptr));
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpInvertMatrix<2>(jac_ptr, det_ptr, inv_jac_ptr));
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpSetInvJacH1ForFace(inv_jac_ptr));
}
// Get pointer to U_tt shift in domain element
auto get_rho = [this](const double, const double, const double) {
auto *pipeline_mng = mField.getInterface<PipelineManager>();
auto &fe_domain_lhs = pipeline_mng->getDomainLhsFE();
return rho * fe_domain_lhs->ts_aa;
};
auto get_body_force = [this](const double, const double, const double) {
auto *pipeline_mng = mField.getInterface<PipelineManager>();
auto fe_domain_rhs = pipeline_mng->getDomainRhsFE();
t_source(i) = 0.;
t_source(0) = 0.1;
t_source(1) = 1.;
const auto time = fe_domain_rhs->ts_t;
t_source(i) *= sin(time * omega * M_PI);
return t_source;
};
auto henky_common_data_ptr = boost::make_shared<HenckyOps::CommonData>();
henky_common_data_ptr->matGradPtr = matGradPtr;
henky_common_data_ptr->matDPtr = matDPtr;
pipeline_mng->getOpDomainLhsPipeline().push_back(
matGradPtr));
pipeline_mng->getOpDomainLhsPipeline().push_back(
new OpCalculateEigenVals<SPACE_DIM>("U", henky_common_data_ptr));
pipeline_mng->getOpDomainLhsPipeline().push_back(
new OpCalculateLogC<SPACE_DIM>("U", henky_common_data_ptr));
pipeline_mng->getOpDomainLhsPipeline().push_back(
new OpCalculateLogC_dC<SPACE_DIM>("U", henky_common_data_ptr));
pipeline_mng->getOpDomainLhsPipeline().push_back(
new OpCalculateHenckyStress<SPACE_DIM>("U", henky_common_data_ptr));
pipeline_mng->getOpDomainLhsPipeline().push_back(
new OpCalculatePiolaStress<SPACE_DIM>("U", henky_common_data_ptr));
pipeline_mng->getOpDomainLhsPipeline().push_back(
new OpHenckyTangent<SPACE_DIM>("U", henky_common_data_ptr));
pipeline_mng->getOpDomainLhsPipeline().push_back(
new OpK("U", "U", henky_common_data_ptr->getMatTangent()));
pipeline_mng->getOpDomainLhsPipeline().push_back(
new OpMass("U", "U", get_rho));
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpBodyForce("U", get_body_force));
pipeline_mng->getOpDomainRhsPipeline().push_back(
matGradPtr));
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculateEigenVals<SPACE_DIM>("U", henky_common_data_ptr));
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculateLogC<SPACE_DIM>("U", henky_common_data_ptr));
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculateLogC_dC<SPACE_DIM>("U", henky_common_data_ptr));
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculateHenckyStress<SPACE_DIM>("U", henky_common_data_ptr));
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculatePiolaStress<SPACE_DIM>("U", henky_common_data_ptr));
pipeline_mng->getOpDomainRhsPipeline().push_back(new OpInternalForce(
"U", henky_common_data_ptr->getMatFirstPiolaStress()));
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpCalculateVectorFieldValuesDotDot<SPACE_DIM>("U",
matAccelerationPtr));
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpScaleMatrix("U", rho, matAccelerationPtr, matInertiaPtr));
pipeline_mng->getOpDomainRhsPipeline().push_back(new OpInertiaForce(
"U", matInertiaPtr, [](double, double, double) { return 1.; }));
}
auto integration_rule = [](int, int, int approx_order) {
return 2 * approx_order;
};
CHKERR pipeline_mng->setDomainRhsIntegrationRule(integration_rule);
CHKERR pipeline_mng->setDomainLhsIntegrationRule(integration_rule);
}
//! [Push operators to pipeline]
/**
* @brief Monitor solution
*
* This functions is called by TS solver at the end of each step. It is used
* to output results to the hard drive.
*/
struct Monitor : public FEMethod {
Monitor(SmartPetscObj<DM> dm, boost::shared_ptr<PostProcEle> post_proc)
: dM(dm), postProc(post_proc){};
MoFEMErrorCode postProcess() {
constexpr int save_every_nth_step = 1;
if (ts_step % save_every_nth_step == 0) {
CHKERR DMoFEMLoopFiniteElements(dM, "dFE", postProc);
CHKERR postProc->writeFile(
"out_step_" + boost::lexical_cast<std::string>(ts_step) + ".h5m");
}
}
private:
SmartPetscObj<DM> dM;
boost::shared_ptr<PostProcEle> postProc;
};
//! [Solve]
auto *simple = mField.getInterface<Simple>();
auto *pipeline_mng = mField.getInterface<PipelineManager>();
auto dm = simple->getDM();
ts = pipeline_mng->createTSIM();
else
ts = pipeline_mng->createTSIM2();
// Setup postprocessing
auto post_proc_fe = boost::make_shared<PostProcEle>(mField);
post_proc_fe->generateReferenceElementMesh();
if (SPACE_DIM == 2) {
auto det_ptr = boost::make_shared<VectorDouble>();
auto jac_ptr = boost::make_shared<MatrixDouble>();
auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
post_proc_fe->getOpPtrVector().push_back(
new OpCalculateHOJacForFace(jac_ptr));
post_proc_fe->getOpPtrVector().push_back(
new OpInvertMatrix<2>(jac_ptr, det_ptr, inv_jac_ptr));
post_proc_fe->getOpPtrVector().push_back(
new OpSetInvJacH1ForFace(inv_jac_ptr));
}
post_proc_fe->getOpPtrVector().push_back(
matGradPtr));
post_proc_fe->getOpPtrVector().push_back(
new OpSymmetrizeTensor<SPACE_DIM>("U", matGradPtr, matStrainPtr));
post_proc_fe->getOpPtrVector().push_back(
new OpTensorTimesSymmetricTensor<SPACE_DIM, SPACE_DIM>(
"U", matStrainPtr, matStressPtr, matDPtr));
post_proc_fe->getOpPtrVector().push_back(new OpPostProcElastic<SPACE_DIM>(
"U", post_proc_fe->postProcMesh, post_proc_fe->mapGaussPts, matStrainPtr,
matStressPtr));
post_proc_fe->addFieldValuesPostProc("U");
// Add monitor to time solver
boost::shared_ptr<FEMethod> null_fe;
auto monitor_ptr = boost::make_shared<Monitor>(dm, post_proc_fe);
CHKERR DMMoFEMTSSetMonitor(dm, ts, simple->getDomainFEName(), null_fe,
null_fe, monitor_ptr);
double ftime = 1;
// CHKERR TSSetMaxTime(ts, ftime);
CHKERR TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);
auto T = smartCreateDMVector(simple->getDM());
CHKERR DMoFEMMeshToLocalVector(simple->getDM(), T, INSERT_VALUES,
SCATTER_FORWARD);
CHKERR TSSetSolution(ts, T);
CHKERR TSSetFromOptions(ts);
} else {
auto TT = smartVectorDuplicate(T);
CHKERR TS2SetSolution(ts, T, TT);
CHKERR TSSetFromOptions(ts);
}
CHKERR TSSolve(ts, NULL);
CHKERR TSGetTime(ts, &ftime);
PetscInt steps, snesfails, rejects, nonlinits, linits;
#if PETSC_VERSION_GE(3, 8, 0)
CHKERR TSGetStepNumber(ts, &steps);
#else
CHKERR TSGetTimeStepNumber(ts, &steps);
#endif
CHKERR TSGetSNESFailures(ts, &snesfails);
CHKERR TSGetStepRejections(ts, &rejects);
CHKERR TSGetSNESIterations(ts, &nonlinits);
CHKERR TSGetKSPIterations(ts, &linits);
MOFEM_LOG_C("EXAMPLE", Sev::inform,
"steps %D (%D rejected, %D SNES fails), ftime %g, nonlinits "
"%D, linits %D\n",
steps, rejects, snesfails, ftime, nonlinits, linits);
}
//! [Solve]
//! [Postprocess results]
PetscBool test_flg = PETSC_FALSE;
CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-test", &test_flg, PETSC_NULL);
if (test_flg) {
auto *simple = mField.getInterface<Simple>();
auto T = smartCreateDMVector(simple->getDM());
CHKERR DMoFEMMeshToLocalVector(simple->getDM(), T, INSERT_VALUES,
SCATTER_FORWARD);
double nrm2;
CHKERR VecNorm(T, NORM_2, &nrm2);
MOFEM_LOG("EXAMPLE", Sev::inform) << "Regression norm " << nrm2;
constexpr double regression_value = 1.09572;
if (fabs(nrm2 - regression_value) > 1e-2)
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"Regression test faileed; wrong norm value.");
}
}
//! [Postprocess results]
//! [Check]
}
//! [Check]
static char help[] = "...\n\n";
int main(int argc, char *argv[]) {
// Initialisation of MoFEM/PETSc and MOAB data structures
const char param_file[] = "param_file.petsc";
// Add logging channel for example
auto core_log = logging::core::get();
core_log->add_sink(
LogManager::setLog("EXAMPLE");
MOFEM_LOG_TAG("EXAMPLE", "example");
try {
//! [Register MoFEM discrete manager in PETSc]
DMType dm_name = "DMMOFEM";
//! [Register MoFEM discrete manager in PETSc
//! [Create MoAB]
moab::Core mb_instance; ///< mesh database
moab::Interface &moab = mb_instance; ///< mesh database interface
//! [Create MoAB]
//! [Create MoFEM]
MoFEM::Core core(moab); ///< finite element database
MoFEM::Interface &m_field = core; ///< finite element database insterface
//! [Create MoFEM]
//! [Example]
Example ex(m_field);
CHKERR ex.runProblem();
//! [Example]
}
}
std::string param_file
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:314
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69
int main(int argc, char *argv[])
static char help[]
EntitiesFieldData::EntData EntData
constexpr int SPACE_DIM
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
Kronecker Delta class symmetric.
constexpr double rho
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesTensor< 1, SPACE_DIM, SPACE_DIM > OpInternalForce
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpMass
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradTensorGrad< 1, SPACE_DIM, SPACE_DIM, 1 > OpK
static double * ts_aa_ptr
constexpr double poisson_ratio
constexpr double shear_modulus_G
constexpr double bulk_modulus_K
constexpr double omega
static double * ts_time_ptr
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpSource< 1, SPACE_DIM > OpBodyForce
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesVector< 1, SPACE_DIM, 1 > OpInertiaForce
constexpr bool is_quasi_static
constexpr double young_modulus
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition: definitions.h:77
@ H1
continuous field
Definition: definitions.h:98
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:53
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
#define CHKERR
Inline error check.
Definition: definitions.h:548
auto integration_rule
constexpr auto t_kd
auto smartCreateDMVector
Get smart vector from DM.
Definition: DMMoFEM.hpp:956
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition: DMMMoFEM.cpp:481
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMMoFEM.cpp:59
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition: DMMMoFEM.cpp:544
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:364
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:311
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:342
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
const double T
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:67
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:21
PetscErrorCode DMMoFEMTSSetMonitor(DM dm, TS ts, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
Set Monitor To TS solver.
Definition: DMMMoFEM.cpp:994
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
OpCalculateHOJacForFaceImpl< 2 > OpCalculateHOJacForFace
CoreTmp< 0 > Core
Definition: Core.hpp:1096
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1965
SmartPetscObj< Vec > smartVectorDuplicate(SmartPetscObj< Vec > &vec)
Create duplicate vector of smart vector.
double A
int save_every_nth_step
#define EXECUTABLE_DIMENSION
Definition: plastic.cpp:24
static constexpr int approx_order
[Example]
MoFEMErrorCode boundaryCondition()
[Set up problem]
MoFEMErrorCode assembleSystem()
[Boundary condition]
MoFEMErrorCode readMesh()
[Run problem]
MoFEMErrorCode checkResults()
[Postprocess results]
MoFEMErrorCode solveSystem()
[Solve]
MoFEMErrorCode createCommonData()
[Create common data]
MoFEMErrorCode runProblem()
[Create common data]
MoFEMErrorCode setupProblem()
[Read mesh]
MoFEMErrorCode outputResults()
[Solve]
Core (interface) class.
Definition: Core.hpp:92
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:85
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:125
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
Definition: LogManager.cpp:279
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
Definition: LogManager.cpp:323
VolumeElementForcesAndSourcesCoreBase::UserDataOperator UserDataOperator
[Push operators to pipeline]
Postprocess on face.
Post processing.
/* This file is part of MoFEM.
* MoFEM is free software: you can redistribute it and/or modify it under
* the terms of the GNU Lesser General Public License as published by the
* Free Software Foundation, either version 3 of the License, or (at your
* option) any later version.
*
* MoFEM is distributed in the hope that it will be useful, but WITHOUT
* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
* FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
* License for more details.
*
* You should have received a copy of the GNU Lesser General Public
* License along with MoFEM. If not, see <http://www.gnu.org/licenses/>. */
namespace HenckyOps {
constexpr double eps = std::numeric_limits<float>::epsilon();
auto f = [](double v) { return 0.5 * std::log(static_cast<long double>(v)); };
auto d_f = [](double v) { return 0.5 / static_cast<long double>(v); };
auto dd_f = [](double v) {
return -0.5 / (static_cast<long double>(v) * static_cast<long double>(v));
};
inline auto is_eq(const double &a, const double &b) {
return std::abs(a - b) < eps;
};
template <int DIM> inline auto get_uniq_nb(double *ptr) {
std::array<double, DIM> tmp;
std::copy(ptr, &ptr[DIM], tmp.begin());
std::sort(tmp.begin(), tmp.end());
return std::distance(tmp.begin(), std::unique(tmp.begin(), tmp.end(), is_eq));
};
template <int DIM>
std::array<std::pair<double, size_t>, DIM> tmp;
auto is_eq_pair = [](const auto &a, const auto &b) {
return a.first < b.first;
};
for (size_t i = 0; i != DIM; ++i)
tmp[i] = std::make_pair(eig(i), i);
std::sort(tmp.begin(), tmp.end(), is_eq_pair);
int i, j, k;
if (is_eq_pair(tmp[0], tmp[1])) {
i = tmp[0].second;
j = tmp[2].second;
k = tmp[1].second;
} else {
i = tmp[0].second;
j = tmp[1].second;
k = tmp[2].second;
}
eigen_vec(i, 0), eigen_vec(i, 1), eigen_vec(i, 2),
eigen_vec(j, 0), eigen_vec(j, 1), eigen_vec(j, 2),
eigen_vec(k, 0), eigen_vec(k, 1), eigen_vec(k, 2)};
FTensor::Tensor1<double, 3> eig_c{eig(i), eig(j), eig(k)};
{
eig(i) = eig_c(i);
eigen_vec(i, j) = eigen_vec_c(i, j);
}
};
struct CommonData : public boost::enable_shared_from_this<CommonData> {
boost::shared_ptr<MatrixDouble> matGradPtr;
boost::shared_ptr<MatrixDouble> matDPtr;
boost::shared_ptr<MatrixDouble> matLogCPlastic;
inline auto getMatFirstPiolaStress() {
return boost::shared_ptr<MatrixDouble>(shared_from_this(),
}
inline auto getMatHenckyStress() {
return boost::shared_ptr<MatrixDouble>(shared_from_this(),
}
inline auto getMatLogC() {
return boost::shared_ptr<MatrixDouble>(shared_from_this(), &matLogC);
}
inline auto getMatTangent() {
return boost::shared_ptr<MatrixDouble>(shared_from_this(), &matTangent);
}
};
template <int DIM> struct OpCalculateEigenVals : public DomainEleOp {
OpCalculateEigenVals(const std::string field_name,
boost::shared_ptr<CommonData> common_data)
: DomainEleOp(field_name, DomainEleOp::OPROW),
commonDataPtr(common_data) {
std::fill(&doEntities[MBEDGE], &doEntities[MBMAXTYPE], false);
}
MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
// const size_t nb_gauss_pts = matGradPtr->size2();
const size_t nb_gauss_pts = getGaussPts().size2();
auto t_grad = getFTensor2FromMat<DIM, DIM>(*(commonDataPtr->matGradPtr));
commonDataPtr->matEigVal.resize(DIM, nb_gauss_pts, false);
commonDataPtr->matEigVec.resize(DIM * DIM, nb_gauss_pts, false);
auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
F(i, j) = t_grad(i, j) + t_kd(i, j);
C(i, j) = F(k, i) ^ F(k, j);
for (int ii = 0; ii != DIM; ii++)
for (int jj = 0; jj != DIM; jj++)
eigen_vec(ii, jj) = C(ii, jj);
CHKERR computeEigenValuesSymmetric<DIM>(eigen_vec, eig);
// rare case when two eigen values are equal
auto nb_uniq = get_uniq_nb<DIM>(&eig(0));
if (DIM == 3 && nb_uniq == 2)
sort_eigen_vals<DIM>(eig, eigen_vec);
t_eig_val(i) = eig(i);
t_eig_vec(i, j) = eigen_vec(i, j);
++t_grad;
++t_eig_val;
++t_eig_vec;
}
}
private:
boost::shared_ptr<CommonData> commonDataPtr;
};
template <int DIM> struct OpCalculateLogC : public DomainEleOp {
OpCalculateLogC(const std::string field_name,
boost::shared_ptr<CommonData> common_data)
: DomainEleOp(field_name, DomainEleOp::OPROW),
commonDataPtr(common_data) {
std::fill(&doEntities[MBEDGE], &doEntities[MBMAXTYPE], false);
}
MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
// const size_t nb_gauss_pts = matGradPtr->size2();
const size_t nb_gauss_pts = getGaussPts().size2();
constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
commonDataPtr->matLogC.resize(size_symm, nb_gauss_pts, false);
auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
// rare case when two eigen values are equal
auto nb_uniq = get_uniq_nb<DIM>(&(t_eig_val(0)));
eig(i) = t_eig_val(i);
eigen_vec(i, j) = t_eig_vec(i, j);
auto logC = EigenMatrix::getMat(eig, eigen_vec, f);
t_logC(i, j) = logC(i, j);
++t_eig_val;
++t_eig_vec;
++t_logC;
}
}
private:
boost::shared_ptr<CommonData> commonDataPtr;
};
template <int DIM> struct OpCalculateLogC_dC : public DomainEleOp {
OpCalculateLogC_dC(const std::string field_name,
boost::shared_ptr<CommonData> common_data)
: DomainEleOp(field_name, DomainEleOp::OPROW),
commonDataPtr(common_data) {
std::fill(&doEntities[MBEDGE], &doEntities[MBMAXTYPE], false);
}
MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
// const size_t nb_gauss_pts = matGradPtr->size2();
const size_t nb_gauss_pts = getGaussPts().size2();
constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
commonDataPtr->matLogCdC.resize(size_symm * size_symm, nb_gauss_pts, false);
auto t_logC_dC = getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->matLogCdC);
auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
eig(i) = t_eig_val(i);
eigen_vec(i, j) = t_eig_vec(i, j);
// rare case when two eigen values are equal
auto nb_uniq = get_uniq_nb<DIM>(&eig(0));
auto dlogC_dC = EigenMatrix::getDiffMat(eig, eigen_vec, f, d_f, nb_uniq);
dlogC_dC(i, j, k, l) *= 2;
t_logC_dC(i, j, k, l) = dlogC_dC(i, j, k, l);
++t_logC_dC;
++t_eig_val;
++t_eig_vec;
}
}
private:
boost::shared_ptr<CommonData> commonDataPtr;
};
template <int DIM> struct OpCalculateHenckyStress : public DomainEleOp {
OpCalculateHenckyStress(const std::string field_name,
boost::shared_ptr<CommonData> common_data)
: DomainEleOp(field_name, DomainEleOp::OPROW),
commonDataPtr(common_data) {
std::fill(&doEntities[MBEDGE], &doEntities[MBMAXTYPE], false);
}
MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
// const size_t nb_gauss_pts = matGradPtr->size2();
const size_t nb_gauss_pts = getGaussPts().size2();
auto t_D = getFTensor4DdgFromMat<DIM, DIM, 0>(*commonDataPtr->matDPtr);
auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
commonDataPtr->matHenckyStress.resize(size_symm, nb_gauss_pts, false);
auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
t_T(i, j) = t_D(i, j, k, l) * t_logC(k, l);
++t_logC;
++t_T;
++t_D;
}
}
private:
boost::shared_ptr<CommonData> commonDataPtr;
};
template <int DIM> struct OpCalculateHenckyPlasticStress : public DomainEleOp {
OpCalculateHenckyPlasticStress(const std::string field_name,
boost::shared_ptr<CommonData> common_data,
boost::shared_ptr<MatrixDouble> mat_D_ptr,
const double scale = 1)
: DomainEleOp(field_name, DomainEleOp::OPROW), commonDataPtr(common_data),
scaleStress(scale), matDPtr(mat_D_ptr) {
std::fill(&doEntities[MBEDGE], &doEntities[MBMAXTYPE], false);
matLogCPlastic = commonDataPtr->matLogCPlastic;
}
MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
// const size_t nb_gauss_pts = matGradPtr->size2();
const size_t nb_gauss_pts = getGaussPts().size2();
auto t_D = getFTensor4DdgFromMat<DIM, DIM, 0>(*matDPtr);
auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
auto t_logCPlastic = getFTensor2SymmetricFromMat<DIM>(*matLogCPlastic);
constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
commonDataPtr->matHenckyStress.resize(size_symm, nb_gauss_pts, false);
auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
t_T(i, j) = t_D(i, j, k, l) * (t_logC(k, l) - t_logCPlastic(k, l));
t_T(i, j) /= scaleStress;
++t_logC;
++t_T;
++t_D;
++t_logCPlastic;
}
}
private:
boost::shared_ptr<CommonData> commonDataPtr;
boost::shared_ptr<MatrixDouble> matDPtr;
boost::shared_ptr<MatrixDouble> matLogCPlastic;
const double scaleStress;
};
template <int DIM> struct OpCalculatePiolaStress : public DomainEleOp {
OpCalculatePiolaStress(const std::string field_name,
boost::shared_ptr<CommonData> common_data)
: DomainEleOp(field_name, DomainEleOp::OPROW),
commonDataPtr(common_data) {
std::fill(&doEntities[MBEDGE], &doEntities[MBMAXTYPE], false);
}
MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
// const size_t nb_gauss_pts = matGradPtr->size2();
const size_t nb_gauss_pts = getGaussPts().size2();
auto t_D = getFTensor4DdgFromMat<DIM, DIM, 0>(*commonDataPtr->matDPtr);
auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
auto t_logC_dC = getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->matLogCdC);
constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
commonDataPtr->matFirstPiolaStress.resize(DIM * DIM, nb_gauss_pts, false);
commonDataPtr->matSecondPiolaStress.resize(size_symm, nb_gauss_pts, false);
auto t_P = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matFirstPiolaStress);
auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
auto t_S =
getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matSecondPiolaStress);
auto t_grad = getFTensor2FromMat<DIM, DIM>(*(commonDataPtr->matGradPtr));
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
t_F(i, j) = t_grad(i, j) + t_kd(i, j);
t_S(k, l) = t_T(i, j) * t_logC_dC(i, j, k, l);
t_P(i, l) = t_F(i, k) * t_S(k, l);
++t_grad;
++t_logC;
++t_logC_dC;
++t_P;
++t_T;
++t_S;
++t_D;
}
}
private:
boost::shared_ptr<CommonData> commonDataPtr;
};
template <int DIM> struct OpHenckyTangent : public DomainEleOp {
OpHenckyTangent(const std::string field_name,
boost::shared_ptr<CommonData> common_data,
boost::shared_ptr<MatrixDouble> mat_D_ptr = nullptr)
: DomainEleOp(field_name, DomainEleOp::OPROW),
commonDataPtr(common_data) {
std::fill(&doEntities[MBEDGE], &doEntities[MBMAXTYPE], false);
if (mat_D_ptr)
matDPtr = mat_D_ptr;
else
matDPtr = commonDataPtr->matDPtr;
}
MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
// const size_t nb_gauss_pts = matGradPtr->size2();
const size_t nb_gauss_pts = getGaussPts().size2();
constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
commonDataPtr->matTangent.resize(DIM * DIM * DIM * DIM, nb_gauss_pts);
auto dP_dF =
getFTensor4FromMat<DIM, DIM, DIM, DIM, 1>(commonDataPtr->matTangent);
auto t_D = getFTensor4DdgFromMat<DIM, DIM, 0>(*matDPtr);
auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
auto t_S =
getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matSecondPiolaStress);
auto t_grad = getFTensor2FromMat<DIM, DIM>(*(commonDataPtr->matGradPtr));
auto t_logC_dC = getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->matLogCdC);
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
t_F(i, j) = t_grad(i, j) + t_kd(i, j);
eig(i) = t_eig_val(i);
eigen_vec(i, j) = t_eig_vec(i, j);
T(i, j) = t_T(i, j);
// rare case when two eigen values are equal
auto nb_uniq = get_uniq_nb<DIM>(&eig(0));
dC_dF(i, j, k, l) = (t_kd(i, l) * t_F(k, j)) + (t_kd(j, l) * t_F(k, i));
auto TL =
EigenMatrix::getDiffDiffMat(eig, eigen_vec, f, d_f, dd_f, T, nb_uniq);
TL(i, j, k, l) *= 4;
P_D_P_plus_TL(i, j, k, l) =
TL(i, j, k, l) +
(t_logC_dC(i, j, o, p) * t_D(o, p, m, n)) * t_logC_dC(m, n, k, l);
P_D_P_plus_TL(i, j, k, l) *= 0.5;
dP_dF(i, j, m, n) = t_kd(i, m) * (t_kd(k, n) * t_S(k, j));
dP_dF(i, j, m, n) +=
t_F(i, k) * (P_D_P_plus_TL(k, j, o, p) * dC_dF(o, p, m, n));
++dP_dF;
++t_grad;
++t_eig_val;
++t_eig_vec;
++t_logC_dC;
++t_S;
++t_T;
++t_D;
}
}
private:
boost::shared_ptr<CommonData> commonDataPtr;
boost::shared_ptr<MatrixDouble> matDPtr;
};
template <int DIM> struct OpPostProcHencky : public DomainEleOp {
OpPostProcHencky(const std::string field_name,
moab::Interface &post_proc_mesh,
std::vector<EntityHandle> &map_gauss_pts,
boost::shared_ptr<CommonData> common_data_ptr);
MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
private:
std::vector<EntityHandle> &mapGaussPts;
boost::shared_ptr<CommonData> commonDataPtr;
};
template <int DIM>
const std::string field_name, moab::Interface &post_proc_mesh,
std::vector<EntityHandle> &map_gauss_pts,
boost::shared_ptr<CommonData> common_data_ptr)
: DomainEleOp(field_name, DomainEleOp::OPROW), postProcMesh(post_proc_mesh),
mapGaussPts(map_gauss_pts), commonDataPtr(common_data_ptr) {
// Opetor is only executed for vertices
std::fill(&doEntities[MBEDGE], &doEntities[MBMAXTYPE], false);
}
//! [Postprocessing]
template <int DIM>
MoFEMErrorCode OpPostProcHencky<DIM>::doWork(int side, EntityType type,
EntData &data) {
std::array<double, 9> def;
std::fill(def.begin(), def.end(), 0);
auto get_tag = [&](const std::string name, size_t size) {
Tag th;
CHKERR postProcMesh.tag_get_handle(name.c_str(), size, MB_TYPE_DOUBLE, th,
MB_TAG_CREAT | MB_TAG_SPARSE,
def.data());
return th;
};
MatrixDouble3by3 mat(3, 3);
auto set_matrix = [&](auto &t) -> MatrixDouble3by3 & {
mat.clear();
for (size_t r = 0; r != SPACE_DIM; ++r)
for (size_t c = 0; c != SPACE_DIM; ++c)
mat(r, c) = t(r, c);
return mat;
};
auto set_tag = [&](auto th, auto gg, MatrixDouble3by3 &mat) {
return postProcMesh.tag_set_data(th, &mapGaussPts[gg], 1,
&*mat.data().begin());
};
auto th_grad = get_tag("GRAD", 9);
auto th_stress = get_tag("STRESS", 9);
auto t_piola =
getFTensor2FromMat<DIM, DIM>(commonDataPtr->matFirstPiolaStress);
auto t_grad = getFTensor2FromMat<DIM, DIM>(*commonDataPtr->matGradPtr);
for (int gg = 0; gg != commonDataPtr->matGradPtr->size2(); ++gg) {
F(i, j) = t_grad(i, j) + kronecker_delta(i, j);
double inv_t_detF;
if (DIM == 2)
determinantTensor2by2(F, inv_t_detF);
else
determinantTensor3by3(F, inv_t_detF);
inv_t_detF = 1. / inv_t_detF;
cauchy_stress(i, j) = inv_t_detF * (t_piola(i, k) ^ F(j, k));
CHKERR set_tag(th_grad, gg, set_matrix(t_grad));
CHKERR set_tag(th_stress, gg, set_matrix(cauchy_stress));
++t_piola;
++t_grad;
}
}
} // namespace HenckyOps
static Index< 'n', 3 > n
static Index< 'o', 3 > o
static Index< 'p', 3 > p
constexpr double a
Kronecker Delta class.
double v
phase velocity of light in medium (cm/ns)
FTensor::Ddg< double, 3, 3 > getDiffMat(Val< double, 3 > &t_val, Vec< double, 3 > &t_vec, Fun< double > f, Fun< double > d_f, const int nb)
Get the Diff Mat object.
FTensor::Tensor2_symmetric< double, 3 > getMat(Val< double, 3 > &t_val, Vec< double, 3 > &t_vec, Fun< double > f)
Get the Mat object.
FTensor::Ddg< double, 3, 3 > getDiffDiffMat(Val< double, 3 > &t_val, Vec< double, 3 > &t_vec, Fun< double > f, Fun< double > d_f, Fun< double > dd_f, FTensor::Tensor2< double, 3, 3 > &t_S, const int nb)
Tensor2_Expr< Kronecker_Delta< T >, T, Dim0, Dim1, i, j > kronecker_delta(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
Rank 2.
auto get_uniq_nb(double *ptr)
Definition: HenckyOps.hpp:29
auto sort_eigen_vals(FTensor::Tensor1< double, DIM > &eig, FTensor::Tensor2< double, DIM, DIM > &eigen_vec)
Definition: HenckyOps.hpp:37
constexpr double eps
Definition: HenckyOps.hpp:17
auto dd_f
Definition: HenckyOps.hpp:21
auto d_f
Definition: HenckyOps.hpp:20
auto is_eq(const double &a, const double &b)
Definition: HenckyOps.hpp:25
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:88
MatrixBoundedArray< double, 9 > MatrixDouble3by3
Definition: Types.hpp:116
MoFEMErrorCode determinantTensor2by2(T1 &t, T2 &det)
Calculate determinant 2 by 2.
Definition: Templates.hpp:1202
static auto determinantTensor3by3(T &t)
Calculate the determinant of a 3x3 matrix or a tensor of rank 2.
Definition: Templates.hpp:1161
const double r
rate factor
double scale
Definition: plastic.cpp:103
constexpr double t
plate stiffness
Definition: plate.cpp:76
FTensor::Index< 'm', 3 > m
auto getMatFirstPiolaStress()
Definition: HenckyOps.hpp:90
MatrixDouble matEigVal
Definition: HenckyOps.hpp:81
MatrixDouble matLogCdC
Definition: HenckyOps.hpp:84
MatrixDouble matEigVec
Definition: HenckyOps.hpp:82
MatrixDouble matTangent
Definition: HenckyOps.hpp:88
boost::shared_ptr< MatrixDouble > matLogCPlastic
Definition: HenckyOps.hpp:79
MatrixDouble matFirstPiolaStress
Definition: HenckyOps.hpp:85
MatrixDouble matLogC
Definition: HenckyOps.hpp:83
boost::shared_ptr< MatrixDouble > matDPtr
Definition: HenckyOps.hpp:78
MatrixDouble matSecondPiolaStress
Definition: HenckyOps.hpp:86
MatrixDouble matHenckyStress
Definition: HenckyOps.hpp:87
boost::shared_ptr< MatrixDouble > matGradPtr
Definition: HenckyOps.hpp:77
boost::shared_ptr< CommonData > commonDataPtr
Definition: HenckyOps.hpp:169
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HenckyOps.hpp:118
OpCalculateEigenVals(const std::string field_name, boost::shared_ptr< CommonData > common_data)
Definition: HenckyOps.hpp:111
OpCalculateHenckyPlasticStress(const std::string field_name, boost::shared_ptr< CommonData > common_data, boost::shared_ptr< MatrixDouble > mat_D_ptr, const double scale=1)
Definition: HenckyOps.hpp:316
boost::shared_ptr< MatrixDouble > matDPtr
Definition: HenckyOps.hpp:360
boost::shared_ptr< MatrixDouble > matLogCPlastic
Definition: HenckyOps.hpp:361
boost::shared_ptr< CommonData > commonDataPtr
Definition: HenckyOps.hpp:359
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HenckyOps.hpp:327
OpCalculateHenckyStress(const std::string field_name, boost::shared_ptr< CommonData > common_data)
Definition: HenckyOps.hpp:275
boost::shared_ptr< CommonData > commonDataPtr
Definition: HenckyOps.hpp:311
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HenckyOps.hpp:282
OpCalculateLogC_dC(const std::string field_name, boost::shared_ptr< CommonData > common_data)
Definition: HenckyOps.hpp:224
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HenckyOps.hpp:231
boost::shared_ptr< CommonData > commonDataPtr
Definition: HenckyOps.hpp:270
boost::shared_ptr< CommonData > commonDataPtr
Definition: HenckyOps.hpp:219
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HenckyOps.hpp:181
OpCalculateLogC(const std::string field_name, boost::shared_ptr< CommonData > common_data)
Definition: HenckyOps.hpp:174
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HenckyOps.hpp:374
boost::shared_ptr< CommonData > commonDataPtr
Definition: HenckyOps.hpp:417
OpCalculatePiolaStress(const std::string field_name, boost::shared_ptr< CommonData > common_data)
Definition: HenckyOps.hpp:367
boost::shared_ptr< CommonData > commonDataPtr
Definition: HenckyOps.hpp:508
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: HenckyOps.hpp:433
boost::shared_ptr< MatrixDouble > matDPtr
Definition: HenckyOps.hpp:509
OpHenckyTangent(const std::string field_name, boost::shared_ptr< CommonData > common_data, boost::shared_ptr< MatrixDouble > mat_D_ptr=nullptr)
Definition: HenckyOps.hpp:421
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
[Postprocessing]
Definition: HenckyOps.hpp:538
OpPostProcHencky(const std::string field_name, moab::Interface &post_proc_mesh, std::vector< EntityHandle > &map_gauss_pts, boost::shared_ptr< CommonData > common_data_ptr)
Definition: HenckyOps.hpp:526
std::vector< EntityHandle > & mapGaussPts
Definition: HenckyOps.hpp:521
moab::Interface & postProcMesh
Definition: HenckyOps.hpp:520
boost::shared_ptr< CommonData > commonDataPtr
Definition: HenckyOps.hpp:522
std::array< bool, MBMAXTYPE > doEntities
If true operator is executed for entity.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element