v0.15.5
Loading...
Searching...
No Matches
AdolCNeohookean.cpp
Go to the documentation of this file.
1/**
2 * @file AdolCNeohookean.cpp
3 * @author your name (you@domain.com)
4 * @brief
5 * @version 0.1
6 * @date 2026-03-26
7 *
8 * @copyright Copyright (c) 2026
9 *
10 */
11
12namespace AdolCOps {
13
14template <int DIM> struct AdolCNeohookean : public AdolCElasticImpl<DIM> {
15 using AdolCElasticImpl<DIM>::AdolCElasticImpl;
16
18
19 MoFEMErrorCode getOptions(MoFEM::Interface *m_field_ptr = nullptr) override {
21
22 MOFEM_LOG_CHANNEL("WORLD");
23
24 PetscOptionsBegin(PETSC_COMM_WORLD, "neo_hookean_", "", "none");
25 CHKERR PetscOptionsScalar("-c10", "C10", "", C10, &C10, PETSC_NULLPTR);
26 CHKERR PetscOptionsScalar("-K", "Bulk modulus K", "", K, &K, PETSC_NULLPTR);
27 PetscOptionsEnd();
28
29 MOFEM_TAG_AND_LOG("WORLD", Sev::inform, "Default Neohookean parameters:")
30 << " C10 = " << C10 << " K = " << K;
32
33 std::string block_name = "MAT_NEOHOOKEAN";
34
35 for (auto &m :
36
37 m_field_ptr->getInterface<MeshsetsManager>()->
38
39 getCubitMeshsetPtr(
40 std::regex((boost::format("%s(.*)") % block_name).str()))
41
42 ) {
43
44 std::vector<double> block_data;
45 CHKERR m->getAttributes(block_data);
46 if (block_data.size() < 2) {
48 "Expected that block has two attribute");
49 }
50 auto get_block_ents = [&]() {
51 Range ents;
52 CHK_MOAB_THROW(m_field_ptr->get_moab().get_entities_by_handle(
53 m->meshset, ents, true),
54 "can not get block entities");
55 return ents;
56 };
57
58 A::paramVecByRange.push_back(
59 {get_block_ents(), {block_data[0], block_data[1]}});
60
61 MOFEM_TAG_AND_LOG("WORLD", Sev::inform, "MatBlock for Neohookean")
62 << *m << " C10 = " << block_data[0] << " K = " << block_data[1];
63 }
64
66 };
67
68 MoFEMErrorCode setParams(EntityHandle ent, int gg) override {
69 for (auto &[range, param_vec] : A::paramVecByRange) {
70 if (std::find(range.begin(), range.end(), ent) != range.end()) {
71 set_param_vec(A::tAg, param_vec.size(), param_vec.data());
72 return 0;
73 }
74 }
75 set_param_vec(A::tAg, defaultMaterialParameters.size(),
77 return 0;
78 }
79
80 MoFEMErrorCode recordTape() override {
82
83 A::adolcDataPtr->insertCommonData("grad", MatrixDouble(DIM * DIM, 1));
84 A::adolcDataPtr->insertCommonData("P", MatrixDouble(DIM * DIM, 1));
85 A::adolcDataPtr->insertCommonData("P_dF",
86 MatrixDouble(DIM * DIM, DIM * DIM));
87
88 A::adolcDataPtr->insertActiveData("F", MatrixDouble(DIM, DIM));
89 A::adolcDataPtr->insertDependentData("P", MatrixDouble(DIM, DIM));
90 A::adolcDataPtr->insertDependentDerivativesData(
91 "P_dF", MatrixDouble(DIM * DIM, DIM * DIM));
92
93 auto t_F = getFTensor2FromPtr<DIM, DIM>(
94 A::adolcDataPtr->getActiveDataPtr("F")->data().data());
95 auto t_P = getFTensor2FromPtr<DIM, DIM>(
96 A::adolcDataPtr->getDependentDataPtr("P")->data().data());
97
99
100 FTENSOR_INDEX(DIM, i);
101 FTENSOR_INDEX(DIM, j);
102 FTENSOR_INDEX(DIM, I);
103 FTENSOR_INDEX(DIM, J);
104
105 t_F(i, J) = 0;
106
109
110 adouble det_aF;
112
113 trace_on(A::tAg);
114 auto p_c10 = mkparam(C10);
115 auto p_K = mkparam(K);
116
117 ta_F(i, J) <<= t_F(i, J);
118 // assume that gradient from approximated displacement, that why we add diagonal
119 ta_F(i, J) += t_kd(i, J);
120
121 det_aF = determinantTensor(ta_F);
122 CHKERR invertTensor(ta_F, det_aF, ta_invF);
123
124 ta_P(i, I) = 2. * p_c10 * (ta_F(i, I) - ta_invF(i, I)) +
125 p_K * log(det_aF) * ta_invF(i, I);
126
127 ta_P(i, I) >>= t_P(i, I);
128
129 trace_off();
130
132 }
133
134protected:
135 double C10 = 1;
136 double K = 1;
137 std::vector<double> defaultMaterialParameters = {C10, K};
138};
139
140template <>
141boost::shared_ptr<PhysicalEquations>
143 boost::shared_ptr<ADolCData> adolc_data_ptr, int tag) {
144 return boost::make_shared<AdolCNeohookean<3>>(adolc_data_ptr, tag);
145}
146
147template <>
148boost::shared_ptr<PhysicalEquations>
150 boost::shared_ptr<ADolCData> adolc_data_ptr, int tag) {
151 return boost::make_shared<AdolCNeohookean<2>>(adolc_data_ptr, tag);
152}
153} // namespace AdolCOps
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
#define FTENSOR_INDEX(DIM, I)
Kronecker Delta class.
#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 ...
#define CHK_MOAB_THROW(err, msg)
Check error code of MoAB function and throw MoFEM exception.
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
constexpr auto t_kd
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'J', DIM1 > J
Definition level_set.cpp:30
FTensor::Index< 'j', 3 > j
boost::shared_ptr< PhysicalEquations > createAdolCPhysicalEquationsPtr< ELASTICITY::NEOHOOKEAN, MODEL_3D >(boost::shared_ptr< ADolCData > adolc_data_ptr, int tag)
boost::shared_ptr< PhysicalEquations > createAdolCPhysicalEquationsPtr< ELASTICITY::NEOHOOKEAN, MODEL_2D_PLANE_STRAIN >(boost::shared_ptr< ADolCData > adolc_data_ptr, int tag)
constexpr IntegrationType I
FTensor::Index< 'm', 3 > m
MoFEMErrorCode setParams(EntityHandle ent, int gg) override
std::vector< double > defaultMaterialParameters
MoFEMErrorCode getOptions(MoFEM::Interface *m_field_ptr=nullptr) override
MoFEMErrorCode recordTape() override
std::vector< std::pair< Range, std::vector< double > > > paramVecByRange
Definition AdolCOps.hpp:165
boost::shared_ptr< ADolCData > adolcDataPtr
Definition AdolCOps.hpp:167
Deprecated interface functions.