v0.15.5
Loading...
Searching...
No Matches
AdolCHuHu.hpp
Go to the documentation of this file.
1/**
2 * @file AdolCHuHu.hpp
3 * @author your name (you@domain.com)
4 * @brief
5 * @version 0.1
6 * @date 2026-03-29
7 *
8 * @copyright Copyright (c) 2026
9 *
10 */
11
12#ifndef __ADOLC_HUHU_HPP__
13#define __ADOLC_HUHU_HPP__
14
15namespace AdolCOps {
16
17struct HUHU {
18 HUHU() = delete;
19};
20
21template <int FIELD_DIM, int SPACE_DIM, AssemblyType A, IntegrationType I,
22 typename OpBase>
23struct OpRhsHu;
24
25template <int FIELD_DIM, int SPACE_DIM, AssemblyType A, IntegrationType I,
26 typename OpBase>
27struct OpLhsHuHu;
28
29template <int FIELD_DIM, int SPACE_DIM, AssemblyType A, IntegrationType I,
30 typename OpBase>
32
33
34template <int MODEL_TYPE> struct OpMaterialFactory<HUHU, MODEL_TYPE> {
36
37 template <AssemblyType A, IntegrationType I, typename DomainEleOp>
38 static MoFEMErrorCode
40 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
41 std::string fe_name, std::string field_name,
42 boost::shared_ptr<PhysicalEquations> physical_equations_ptr,
43 Sev sev = Sev::noisy) {
45 constexpr int DIM = (MODEL_TYPE == MODEL_3D) ? 3 : 2;
46 auto op_this = getPipThis<DomainEleOp>(m_field, pip, fe_name, field_name,
47 physical_equations_ptr, sev);
48 auto m_grad_grad = boost::make_shared<MatrixDouble>();
49 op_this->getOpPtrVector().push_back(
50 new OpCalculateVectorFieldHessian<DIM, DIM>(field_name, m_grad_grad));
51 op_this->getOpPtrVector().push_back(
52 physical_equations_ptr->createOp(physical_equations_ptr, true, false));
53 auto m_k = physical_equations_ptr->adolcDataPtr->getCommonDataPtr("k");
54 op_this->getOpPtrVector().push_back(
57 }
58
59 template <AssemblyType A, IntegrationType I, typename DomainEleOp>
60 static MoFEMErrorCode
62 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
63 std::string fe_name, std::string field_name,
64 boost::shared_ptr<PhysicalEquations> physical_equations_ptr,
65 Sev sev = Sev::noisy) {
67 constexpr int DIM = (MODEL_TYPE == MODEL_3D) ? 3 : 2;
68 auto op_this = getPipThis<DomainEleOp>(m_field, pip, fe_name, field_name,
69 physical_equations_ptr, sev);
70 auto m_grad_grad = boost::make_shared<MatrixDouble>();
71 op_this->getOpPtrVector().push_back(
72 new OpCalculateVectorFieldHessian<DIM, DIM>(field_name, m_grad_grad));
73 op_this->getOpPtrVector().push_back(
74 physical_equations_ptr->createOp(physical_equations_ptr, true, true));
75 auto m_k = physical_equations_ptr->adolcDataPtr->getCommonDataPtr("k");
76 auto m_diff_k =
77 physical_equations_ptr->adolcDataPtr->getCommonDataPtr("k_dF");
78 op_this->getOpPtrVector().push_back(
80 op_this->getOpPtrVector().push_back(
82 m_diff_k));
84 }
85
86private:
87 template <typename DomainEleOp>
88 static auto
90 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
91 std::string fe_name, std::string field_name,
92 boost::shared_ptr<PhysicalEquations> physical_equations_ptr,
93 Sev sev) {
94
95 auto &param_vec_by_range = physical_equations_ptr->paramVecByRange;
96 // range only if parameters are set via blockset
97 auto r = boost::make_shared<Range>();
98 for (auto &p : param_vec_by_range) {
99 r->merge(p.first);
100 }
101 MOFEM_LOG("WORLD", Sev::inform) << "HuHu number of entities " << r->size();
102
103 constexpr int DIM = (MODEL_TYPE == MODEL_3D) ? 3 : 2;
104
105 auto field_structure = m_field.get_field_structure(field_name);
106 auto base = field_structure->getApproxBase();
107 auto space = field_structure->getSpace();
108
109 using DomainEle = typename parent_of<DomainEleOp>::type;
110 auto op_this = new OpLoopThis<DomainEle>(m_field, fe_name, r, sev);
111 pip.push_back(op_this);
112 auto this_fe_ptr = op_this->getThisFEPtr();
113 auto &this_pip = op_this->getOpPtrVector();
114 this_fe_ptr->getRuleHook = [](int order_row, int order_col,
115 int order_data) {
116 return 2 * (order_data - 1);
117 };
118
119 auto base_mass = boost::make_shared<MatrixDouble>();
120 auto data_l2 = boost::make_shared<EntitiesFieldData>(MBENTITYSET);
121 auto jac_ptr = boost::make_shared<MatrixDouble>();
122 auto det_ptr = boost::make_shared<VectorDouble>();
123 auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
124
125 // calculate jacobian at integration points
126 this_pip.push_back(new OpCalculateHOJac<DIM>(jac_ptr));
127 // calculate jacobian at integration points
128 this_pip.push_back(new OpInvertMatrix<DIM>(jac_ptr, det_ptr, inv_jac_ptr));
129 // calculate mass matrix to project derivatives
130 this_pip.push_back(
131 new OpBaseDerivativesMass<1>(base_mass, data_l2, base, L2));
132 // calculate second derivative of base functions, i.e. hessian
133 this_pip.push_back(new OpBaseDerivativesNext<1>(
134 BaseDerivatives::SecondDerivative, base_mass, data_l2, base, space));
135
136 switch (space) {
137 case H1:
138 // push first base derivatives tp physical element shape
139 this_pip.push_back(
140 new OpSetHOInvJacToScalarBases<DIM, 1>(space, inv_jac_ptr));
141 // push second base derivatives tp physical element shape
142 this_pip.push_back(
143 new OpSetHOInvJacToScalarBases<DIM, 2>(space, inv_jac_ptr));
144 break;
145 default:
146 MOFEM_LOG("WORLD", Sev::error)
147 << "Unsupported space: " << space << " for field: " << field_name;
149 }
150
151 auto m_grad =
152 physical_equations_ptr->adolcDataPtr->getCommonDataPtr("grad");
153 this_pip.push_back(
154 new OpCalculateVectorFieldGradient<DIM, DIM>(field_name, m_grad));
155
156 return op_this;
157 }
158};
159
160template <>
161boost::shared_ptr<PhysicalEquations>
163 boost::shared_ptr<ADolCData> adolc_data_ptr, int tag);
164
165template <>
166boost::shared_ptr<PhysicalEquations>
168 boost::shared_ptr<ADolCData> adolc_data_ptr, int tag);
169
170template <int FIELD_DIM, int SPACE_DIM, AssemblyType A, typename OpBase>
172 : public FormsIntegrators<OpBase>::template Assembly<A>::OpBase {
173
174 using OP = typename FormsIntegrators<OpBase>::template Assembly<A>::OpBase;
175
176 OpRhsHu(const std::string field_name,
177 boost::shared_ptr<MatrixDouble> mat_vals,
178 boost::shared_ptr<MatrixDouble> mat_K,
179 boost::shared_ptr<Range> ents_ptr = nullptr)
180 : OP(field_name, field_name, OP::OPROW), matVals(mat_vals), matK(mat_K) {}
181
182protected:
183 boost::shared_ptr<MatrixDouble> matVals;
184 boost::shared_ptr<MatrixDouble> matK;
185 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &data);
186};
187
188template <int FIELD_DIM, int SPACE_DIM, AssemblyType A, typename OpBase>
190 : public FormsIntegrators<OpBase>::template Assembly<A>::OpBase {
191
192 using OP = typename FormsIntegrators<OpBase>::template Assembly<A>::OpBase;
193 OpLhsHuHu(const std::string field_name, boost::shared_ptr<MatrixDouble> mat_K,
194 boost::shared_ptr<Range> ents_ptr = nullptr)
195 : OP(field_name, field_name, OP::OPROWCOL), matK(mat_K) {}
196
197protected:
198 boost::shared_ptr<MatrixDouble> matK;
199 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data,
200 EntitiesFieldData::EntData &col_data);
201};
202
203template <int FIELD_DIM, int SPACE_DIM, AssemblyType A, typename OpBase>
205 : public FormsIntegrators<OpBase>::template Assembly<A>::OpBase {
206 using OP = typename FormsIntegrators<OpBase>::template Assembly<A>::OpBase;
207 OpLhsHuGrad(const std::string field_name,
208 boost::shared_ptr<MatrixDouble> mat_vals,
209 boost::shared_ptr<MatrixDouble> mat_diff_K,
210 boost::shared_ptr<Range> ents_ptr = nullptr)
211 : OP(field_name, field_name, OP::OPROWCOL), matVals(mat_vals),
212 matDiffK(mat_diff_K) {
213 OP::sYmm = false;
214 }
215
216protected:
217 boost::shared_ptr<MatrixDouble> matVals;
218 boost::shared_ptr<MatrixDouble> matDiffK;
219 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data,
220 EntitiesFieldData::EntData &col_data);
221};
222
223template <int FIELD_DIM, int SPACE_DIM, AssemblyType A, typename OpBase>
225 EntitiesFieldData::EntData &row_data) {
227
231
232 // get element volume
233 const double vol = OP::getMeasure();
234 // get integration weights
235 auto t_w = OP::getFTensor0IntegrationWeight();
236 // get base function gradient on rows
237 auto t_row_grad = row_data.getFTensor2DiffN2<SPACE_DIM>();
238 // get field gradient values
239 auto t_val_grad = getFTensor3FromMat<FIELD_DIM, SPACE_DIM, SPACE_DIM>(
240 *(matVals)); // tensor of second derivatives
241 // get K values
242 auto t_K = getFTensor1FromMat<1>(*(matK)); // tensor of second derivatives
243 // get coordinate at integration points
244 auto t_coords = OP::getFTensor1CoordsAtGaussPts();
245 // loop over integration points
246 for (int gg = 0; gg != OP::nbIntegrationPts; gg++) {
247 // take into account Jacobian
248 const double alpha = t_w * vol * t_K(0);
249 // calculate rhs
250 auto t_nf = getFTensor1FromPtr<FIELD_DIM>(OP::locF.data().data());
251 // loop over rows base functions
252 int rr = 0;
253 for (; rr != OP::nbRows / FIELD_DIM; rr++) {
254 // calculate element of local rhs vector
255 t_nf(i) += alpha * (t_row_grad(J, K) * t_val_grad(i, J, K));
256 ++t_row_grad; // move to another element of gradient of base
257 // function on row
258 }
259 for (; rr < OP::nbRowBaseFunctions; ++rr)
260 ++t_row_grad;
261
262 ++t_coords;
263 ++t_val_grad;
264 ++t_K;
265 ++t_w; // move to another integration weight
266 }
268}
269
270template <int FIELD_DIM, int SPACE_DIM, AssemblyType A, typename OpBase>
272 EntitiesFieldData::EntData &row_data,
273 EntitiesFieldData::EntData &col_data) {
275
280
281 constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
282
283 // get element volume
284 const double vol = OP::getMeasure();
285 // get integration weights
286 auto t_w = OP::getFTensor0IntegrationWeight();
287 // get base function gradient on rows
288 auto t_row_grad = row_data.getFTensor2DiffN2<SPACE_DIM>();
289 // get k values
290 auto t_K = getFTensor1FromMat<1>(*(matK)); // tensor of second derivatives
291 // get coordinate at integration points
292 auto t_coords = OP::getFTensor1CoordsAtGaussPts();
293 // loop over integration points
294 for (int gg = 0; gg != OP::nbIntegrationPts; gg++) {
295 // take into account Jacobian
296 const double alpha = t_w * vol * t_K(0);
297 // loop over rows base functions
298 int rr = 0;
299 for (; rr != OP::nbRows / FIELD_DIM; rr++) {
300
301 auto t_mat = getFTensor2FromPtr<FIELD_DIM, FIELD_DIM, FIELD_DIM>(
302 &OP::locMat(rr * FIELD_DIM, 0));
303 auto t_col_grad = col_data.getFTensor2DiffN2<SPACE_DIM>(gg, 0);
304
305 // calculate element of local matrix
306 for (int cc = 0; cc != OP::nbCols / FIELD_DIM; cc++) {
307 t_mat(i, j) +=
308 alpha * (t_row_grad(J, K) * t_col_grad(J, K)) * t_kd(i, j);
309 ++t_col_grad; // move to another element of gradient of base
310 // function on column
311 ++t_mat; // move to another element of local matrix
312 }
313 ++t_row_grad; // move to another element of gradient of base
314 // function on row
315 }
316 for (; rr < OP::nbRowBaseFunctions; ++rr)
317 ++t_row_grad;
318
319 ++t_coords;
320 ++t_K;
321 ++t_w; // move to another integration weight
322 }
324}
325
326template <int FIELD_DIM, int SPACE_DIM, AssemblyType A, typename OpBase>
328 EntitiesFieldData::EntData &row_data,
329 EntitiesFieldData::EntData &col_data) {
331
337
338 // get element volume
339 const double vol = OP::getMeasure();
340 // get integration weights
341 auto t_w = OP::getFTensor0IntegrationWeight();
342 // get base function gradient on rows
343 auto t_row_grad = row_data.getFTensor2DiffN2<SPACE_DIM>();
344 // get field gradient values
345 auto t_val_grad = getFTensor3FromMat<FIELD_DIM, SPACE_DIM, SPACE_DIM>(
346 *(matVals)); // tensor of second derivatives
347 // get K values
348 auto t_diff_K = getFTensor2FromMat<FIELD_DIM, SPACE_DIM>(
349 *(matDiffK)); // tensor of second derivatives
350 // get coordinate at integration points
351 auto t_coords = OP::getFTensor1CoordsAtGaussPts();
352 // loop over integration points
353 for (int gg = 0; gg != OP::nbIntegrationPts; gg++) {
354 // take into account Jacobian
355 const double alpha = t_w * vol;
356 // loop over rows base functions
357 int rr = 0;
358 for (; rr != OP::nbRows / FIELD_DIM; rr++) {
359
360 auto t_mat = getFTensor2FromPtr<FIELD_DIM, FIELD_DIM, FIELD_DIM>(
361 &OP::locMat(rr * FIELD_DIM, 0));
362 auto t_col_grad = col_data.getFTensor1DiffN<SPACE_DIM>(gg, 0);
363 for (int bb = 0; bb != OP::nbCols / FIELD_DIM; ++bb) {
364 t_mat(i, j) += alpha * (t_row_grad(J, K) * t_val_grad(i, J, K)) *
365 (t_diff_K(j, M) * t_col_grad(M));
366
367 ++t_mat; // move to another element of local matrix
368 ++t_col_grad; // move to another element of gradient of base
369 // function on column
370 }
371
372 ++t_row_grad; // move to another element of gradient of base
373 // function on row
374 }
375 for (; rr < OP::nbRowBaseFunctions; ++rr)
376 ++t_row_grad;
377
378 ++t_coords;
379 ++t_val_grad;
380 ++t_diff_K;
381 ++t_w; // move to another integration weight
382 }
384}
385
386} // namespace AdolCOps
387
388#endif // __ADOLC_HUHU_HPP__
#define FTENSOR_INDEX(DIM, I)
constexpr int SPACE_DIM
constexpr int FIELD_DIM
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
Kronecker Delta class.
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
@ L2
field with C-1 continuity
Definition definitions.h:88
@ H1
continuous field
Definition definitions.h:85
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition definitions.h:34
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
constexpr auto t_kd
virtual const Field * get_field_structure(const std::string &name, enum MoFEMTypes bh=MF_EXIST) const =0
get field structure
#define MOFEM_LOG(channel, severity)
Log.
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< HUHU, MODEL_2D_PLANE_STRAIN >(boost::shared_ptr< ADolCData > adolc_data_ptr, int tag)
boost::shared_ptr< PhysicalEquations > createAdolCPhysicalEquationsPtr< HUHU, MODEL_3D >(boost::shared_ptr< ADolCData > adolc_data_ptr, int tag)
constexpr IntegrationType I
constexpr AssemblyType A
constexpr auto field_name
OpLhsHuGrad(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals, boost::shared_ptr< MatrixDouble > mat_diff_K, boost::shared_ptr< Range > ents_ptr=nullptr)
OpLhsHuHu(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_K, boost::shared_ptr< Range > ents_ptr=nullptr)
static MoFEMErrorCode opRhsFactory(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string fe_name, std::string field_name, boost::shared_ptr< PhysicalEquations > physical_equations_ptr, Sev sev=Sev::noisy)
Definition AdolCHuHu.hpp:39
static auto getPipThis(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string fe_name, std::string field_name, boost::shared_ptr< PhysicalEquations > physical_equations_ptr, Sev sev)
Definition AdolCHuHu.hpp:89
static MoFEMErrorCode opLhsFactory(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string fe_name, std::string field_name, boost::shared_ptr< PhysicalEquations > physical_equations_ptr, Sev sev=Sev::noisy)
Definition AdolCHuHu.hpp:61
OpRhsHu(const std::string field_name, boost::shared_ptr< MatrixDouble > mat_vals, boost::shared_ptr< MatrixDouble > mat_K, boost::shared_ptr< Range > ents_ptr=nullptr)
Deprecated interface functions.
FieldApproximationBase getApproxBase() const
Get approximation basis type.