v0.15.5
Loading...
Searching...
No Matches
AdolCHuHu.cpp
Go to the documentation of this file.
1/**
2 * @file AdolCHuHu.cpp
3 * @author your name (you@domain.com)
4 * @version 0.1
5 * @date 2026-03-29
6 *
7 * @copyright Copyright (c) 2026
8 *
9 */
10
11#include <MoFEM.hpp>
12
13using namespace MoFEM;
14
15#include "AdolCOps.hpp"
16#include "AdolCHuHu.hpp"
17
18namespace AdolCOps {
19
20template <int DIM>
23
25
27 createOp(boost::shared_ptr<PhysicalEquations> physical_ptr, bool eval_stress,
28 bool eval_tangent) override;
29
30 MoFEMErrorCode getOptions(MoFEM::Interface *m_field_ptr = nullptr) override {
32 MOFEM_LOG_CHANNEL("WORLD");
33
34 PetscOptionsBegin(PETSC_COMM_WORLD, "huhu_", "", "none");
35 CHKERR PetscOptionsScalar("-K", "Bulk modulus K", "", K, &K, PETSC_NULLPTR);
36 PetscOptionsEnd();
37
38 MOFEM_TAG_AND_LOG("WORLD", Sev::inform, "Default HuHu parameters:")
39 << " K = " << K;
41
42 std::string block_name = "MAT_HUHU";
43
44 for (auto &m :
45
46 m_field_ptr->getInterface<MeshsetsManager>()->
47
48 getCubitMeshsetPtr(
49 std::regex((boost::format("%s(.*)") % block_name).str()))
50
51 ) {
52
53 std::vector<double> block_data;
54 CHKERR m->getAttributes(block_data);
55 if (block_data.size() < 1) {
56 block_data.push_back(
57 K); // Use default K if not provided in block attributes
58 }
59 auto get_block_ents = [&]() {
60 Range ents;
61 CHK_MOAB_THROW(m_field_ptr->get_moab().get_entities_by_handle(
62 m->meshset, ents, true),
63 "can not get block entities");
64 return ents;
65 };
66
67 A::paramVecByRange.push_back({get_block_ents(), {block_data[0]}});
68
69 MOFEM_TAG_AND_LOG("WORLD", Sev::inform, "MatBlock for HuHu")
70 << *m << " K = " << block_data[0];
71 }
72
74 };
75
76 MoFEMErrorCode setParams(EntityHandle ent, int gg) override {
77 set_param_vec(A::tAg, defaultMaterialParameters.size(),
79 return 0;
80 }
81
84
85 A::adolcDataPtr->insertCommonData("grad", MatrixDouble(DIM * DIM, 1));
86 A::adolcDataPtr->insertCommonData("k", MatrixDouble(1, 1));
87 A::adolcDataPtr->insertCommonData("k_dF", MatrixDouble(DIM * DIM, 1));
88
89 A::adolcDataPtr->insertActiveData("F", MatrixDouble(DIM, DIM));
90 A::adolcDataPtr->insertDependentData("k", MatrixDouble(1, 1));
91 A::adolcDataPtr->insertDependentDerivativesData("k_dF",
92 MatrixDouble(DIM, DIM));
93
94 auto t_F = getFTensor2FromPtr<DIM, DIM>(
95 A::adolcDataPtr->getActiveDataPtr("F")->data().data());
96 auto k = (*A::adolcDataPtr->getDependentDataPtr("k"))(0, 0);
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) = t_kd(i, J);
106
108 adouble ta_k;
109
110 adouble det_aF;
112
113 trace_on(A::tAg);
114 auto p_K = mkparam(K);
115 ta_F(i, J) <<= t_F(i, J);
116 ta_k = p_K;
117 ta_k >>= k;
118
119 trace_off();
120
122 }
123
124protected:
125 double K = 1e-6;
126 std::vector<double> defaultMaterialParameters = {K};
127};
128
129template <>
130boost::shared_ptr<PhysicalEquations>
132 boost::shared_ptr<ADolCData> adolc_data_ptr, int tag) {
133 return boost::make_shared<AdolCHuHu<3>>(adolc_data_ptr, tag);
134}
135
136template <>
137boost::shared_ptr<PhysicalEquations>
139 boost::shared_ptr<ADolCData> adolc_data_ptr, int tag) {
140 return boost::make_shared<AdolCHuHu<2>>(adolc_data_ptr, tag);
141}
142
143template <int DIM, typename DomainEleOp>
145
146 OpEvalAdolCHuHuImpl(boost::shared_ptr<PhysicalEquations> physical_ptr,
147 bool eval_stress, bool eval_tangent)
148 : DomainEleOp("U", DomainEleOp::OPROW), physicalPtr(physical_ptr),
149 evalStress(eval_stress), evalTangent(eval_tangent) {}
150
151 MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
152
153protected:
154 boost::shared_ptr<PhysicalEquations> physicalPtr;
155 const bool evalStress;
156 const bool evalTangent;
157};
158
159template <int DIM, typename DomainEleOp>
161 EntityType type,
162 EntData &data) {
164 int nb_dofs = data.getIndices().size();
165 if (!nb_dofs)
167 int nb_integration_pts = data.getN().size1();
168
169 constexpr int adolc_return_value = 0;
170
171 auto physics_ptr = physicalPtr;
172
173 auto get_tag = [&]() {
174 if (physics_ptr->tagVsRangePtr) {
175 for (const auto &tag_range_pair : *physics_ptr->tagVsRangePtr) {
176 if (tag_range_pair.second.find(DomainEleOp::getFEEntityHandle()) !=
177 tag_range_pair.second.end()) {
178 return tag_range_pair.first;
179 }
180 }
181 }
182 return physics_ptr->tAg; // Default tag if no range matches
183 };
184
185 auto mat_grad_ptr = physics_ptr->adolcDataPtr->getCommonDataPtr("grad");
186 auto mat_k_ptr = physics_ptr->adolcDataPtr->getCommonDataPtr("k");
187 auto mat_k_dF_ptr = physics_ptr->adolcDataPtr->getCommonDataPtr("k_dF");
188
189#ifndef NDEBUG
190 if (!mat_grad_ptr || !mat_k_ptr || !mat_k_dF_ptr) {
191 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
192 "Missing common data for ADOL-C evaluation");
193 }
194 if (mat_grad_ptr->size1() != DIM * DIM) {
195 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
196 "Inconsistent size of gradient matrix for ADOL-C evaluation");
197 }
198 if (mat_grad_ptr->size2() != nb_integration_pts) {
199 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
200 "Inconsistent size of gradient matrix data for ADOL-C evaluation");
201 }
202#endif
203
204 if (evalStress) {
205 FTENSOR_INDEX(DIM, i);
206 FTENSOR_INDEX(DIM, J);
207 FTENSOR_INDEX(DIM, K);
208 FTENSOR_INDEX(DIM, L);
209
210 auto t_grad_at_pts = getFTensor2FromMat<DIM, DIM>(*mat_grad_ptr);
211 mat_k_ptr->resize(1, nb_integration_pts, false);
212 auto t_k_at_pts = getFTensor1FromMat<1>(*mat_k_ptr);
213
214 auto next = [&]() {
215 ++t_grad_at_pts;
216 ++t_k_at_pts;
217 };
218
220
221 for (auto gg = 0; gg != nb_integration_pts; ++gg) {
222 auto t_F = getFTensor2FromPtr<DIM, DIM>(
223 physics_ptr->adolcDataPtr->getActiveDataPtr("F")->data().data());
224 t_F(i, J) = t_grad_at_pts(i, J) + t_kd(i, J);
225 CHKERR physics_ptr->setActiveContinuousVector();
226 CHKERR physics_ptr->setDependentContinuousVector();
227 CHKERR physics_ptr->setParams(this->getFEEntityHandle(), gg);
228 int r = ::function(get_tag(), physics_ptr->dependentVariables.size(),
229 physics_ptr->activeVariables.size(),
230 physics_ptr->activeVariables.data(),
231 physics_ptr->dependentVariables.data());
232 if (PetscUnlikely(r < adolc_return_value)) {
233 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
234 "ADOL-C function evaluation with error");
235 }
236 CHKERR physics_ptr->getDependentContinuousVector();
237 t_k_at_pts(0) =
238 (*physics_ptr->adolcDataPtr->getDependentDataPtr("k"))(0, 0);
239 next();
240 }
241 }
242
243 if (evalTangent) {
244 FTENSOR_INDEX(DIM, i);
245 FTENSOR_INDEX(DIM, J);
246 FTENSOR_INDEX(DIM, k);
247 FTENSOR_INDEX(DIM, L);
248
249 const int number_of_active_variables = physics_ptr->activeVariables.size();
250 const int number_of_dependent_variables =
251 physics_ptr->dependentVariables.size();
252 CHKERR physics_ptr->setDependentDerivativesContinuousVector();
253 std::vector<double *> jac_ptr(number_of_dependent_variables);
254 for (unsigned int n = 0; n != number_of_dependent_variables; ++n) {
255 jac_ptr[n] = &(
256 physics_ptr
257 ->dependentVariablesDerivatives[n * number_of_active_variables]);
258 }
259
260 auto t_grad_at_pts = getFTensor2FromMat<DIM, DIM>(*mat_grad_ptr);
261 mat_k_dF_ptr->resize(DIM * DIM, nb_integration_pts, false);
262 auto t_k_dF_at_pts = getFTensor2FromMat<DIM, DIM>(*mat_k_dF_ptr);
264
265 auto next = [&]() {
266 ++t_grad_at_pts;
267 ++t_k_dF_at_pts;
268 };
269
270
271 for (auto gg = 0; gg != nb_integration_pts; ++gg) {
272 auto t_F = getFTensor2FromPtr<DIM, DIM>(
273 physics_ptr->adolcDataPtr->getActiveDataPtr("F")->data().data());
274 t_F(i, J) = t_grad_at_pts(i, J) + t_kd(i, J);
275 CHKERR physics_ptr->setActiveContinuousVector();
276 CHKERR physics_ptr->setDependentContinuousVector();
277 CHKERR physics_ptr->setParams(this->getFEEntityHandle(), gg);
278 int r = ::jacobian(get_tag(), number_of_dependent_variables,
279 number_of_active_variables,
280 physics_ptr->activeVariables.data(), jac_ptr.data());
281 if (PetscUnlikely(r < adolc_return_value)) {
282 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
283 "ADOL-C function evaluation with error");
284 }
285 CHKERR physics_ptr->getDependentDerivativesContinuousVector();
286 auto t_k_dF = getFTensor2FromPtr<DIM, DIM>(
287 physics_ptr->adolcDataPtr->getDependentDerivativesDataPtr("k_dF")
288 ->data()
289 .data());
290 t_k_dF_at_pts(i, J) = t_k_dF(i, J);
291 next();
292 }
293 }
294
296}
297
298template <int DIM, typename EleOp>
299static EleOp *createOpImpl(boost::shared_ptr<PhysicalEquations> physical_ptr,
300 bool eval_stress, bool eval_tangent) {
301 return new OpEvalAdolCHuHuImpl<DIM, EleOp>(physical_ptr, eval_stress,
302 eval_tangent);
303}
304
305template <>
307AdolCHuHu<3>::createOp(boost::shared_ptr<PhysicalEquations> physical_ptr,
308 bool eval_stress, bool eval_tangent) {
310 return createOpImpl<3, EleOp>(physical_ptr, eval_stress, eval_tangent);
311}
312
313template <>
315AdolCHuHu<2>::createOp(boost::shared_ptr<PhysicalEquations> physical_ptr,
316 bool eval_stress, bool eval_tangent) {
318 return createOpImpl<2, EleOp>(physical_ptr, eval_stress, eval_tangent);
319}
320
321} // namespace AdolCOps
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
#define FTENSOR_INDEX(DIM, I)
Kronecker Delta class.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#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_OPERATION_UNSUCCESSFUL
Definition definitions.h:34
#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
const double n
refractive index of diffusive medium
FTensor::Index< 'J', DIM1 > J
Definition level_set.cpp:30
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
boost::shared_ptr< PhysicalEquations > createAdolCPhysicalEquationsPtr< HUHU, MODEL_2D_PLANE_STRAIN >(boost::shared_ptr< ADolCData > adolc_data_ptr, int tag)
static EleOp * createOpImpl(boost::shared_ptr< PhysicalEquations > physical_ptr, bool eval_stress, bool eval_tangent)
boost::shared_ptr< PhysicalEquations > createAdolCPhysicalEquationsPtr< HUHU, MODEL_3D >(boost::shared_ptr< ADolCData > adolc_data_ptr, int tag)
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
UBlasMatrix< double > MatrixDouble
Definition Types.hpp:77
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
constexpr IntegrationType I
FTensor::Index< 'm', 3 > m
MoFEMErrorCode recordTape() override
Definition AdolCHuHu.cpp:82
MoFEMErrorCode getOptions(MoFEM::Interface *m_field_ptr=nullptr) override
Definition AdolCHuHu.cpp:30
MoFEMErrorCode setParams(EntityHandle ent, int gg) override
Definition AdolCHuHu.cpp:76
std::vector< double > defaultMaterialParameters
PhysicalEquations()=delete
ForcesAndSourcesCore::UserDataOperator * createOp(boost::shared_ptr< PhysicalEquations > physical_ptr, bool eval_stress, bool eval_tangent) override
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
OpEvalAdolCHuHuImpl(boost::shared_ptr< PhysicalEquations > physical_ptr, bool eval_stress, bool eval_tangent)
boost::shared_ptr< PhysicalEquations > physicalPtr
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.
Data on single entity (This is passed as argument to DataOperator::doWork)
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.
Interface for managing meshsets containing materials and boundary conditions.