v0.15.0
Loading...
Searching...
No Matches
HenckyOps.hpp
Go to the documentation of this file.
1/**
2 * \file HenckyOps.hpp
3 * \example HenckyOps.hpp
4 *
5 * @copyright Copyright (c) 2023
6 */
7
8#ifndef __HENCKY_OPS_HPP__
9#define __HENCKY_OPS_HPP__
10
11namespace HenckyOps {
12
13constexpr double eps = std::numeric_limits<float>::epsilon();
14
15auto f = [](double v) { return 0.5 * std::log(v); };
16auto d_f = [](double v) { return 0.5 / v; };
17auto dd_f = [](double v) { return -0.5 / (v * v); };
18
19struct isEq {
20 static inline auto check(const double &a, const double &b) {
21 return std::abs(a - b) / absMax < eps;
22 }
23 static double absMax;
24};
25
26double isEq::absMax = 1;
27
28inline auto is_eq(const double &a, const double &b) {
29 return isEq::check(a, b);
30};
31
32template <int DIM> inline auto get_uniq_nb(double *ptr) {
33 std::array<double, DIM> tmp;
34 std::copy(ptr, &ptr[DIM], tmp.begin());
35 std::sort(tmp.begin(), tmp.end());
36 isEq::absMax = std::max(std::abs(tmp[0]), std::abs(tmp[DIM - 1]));
37 return std::distance(tmp.begin(), std::unique(tmp.begin(), tmp.end(), is_eq));
38};
39
40template <int DIM>
43
45 std::max(std::max(std::abs(eig(0)), std::abs(eig(1))), std::abs(eig(2)));
46
47 int i = 0, j = 1, k = 2;
48
49 if (is_eq(eig(0), eig(1))) {
50 i = 0;
51 j = 2;
52 k = 1;
53 } else if (is_eq(eig(0), eig(2))) {
54 i = 0;
55 j = 1;
56 k = 2;
57 } else if (is_eq(eig(1), eig(2))) {
58 i = 1;
59 j = 0;
60 k = 2;
61 }
62
64 eigen_vec(i, 0), eigen_vec(i, 1), eigen_vec(i, 2),
65
66 eigen_vec(j, 0), eigen_vec(j, 1), eigen_vec(j, 2),
67
68 eigen_vec(k, 0), eigen_vec(k, 1), eigen_vec(k, 2)};
69
70 FTensor::Tensor1<double, 3> eig_c{eig(i), eig(j), eig(k)};
71
72 {
73 FTensor::Index<'i', DIM> i;
74 FTensor::Index<'j', DIM> j;
75 eig(i) = eig_c(i);
76 eigen_vec(i, j) = eigen_vec_c(i, j);
77 }
78};
79
80struct CommonData : public boost::enable_shared_from_this<CommonData> {
81 boost::shared_ptr<MatrixDouble> matGradPtr;
82 boost::shared_ptr<MatrixDouble> matDPtr;
83 boost::shared_ptr<MatrixDouble> matLogCPlastic;
84
85 MatrixDouble matEigVal;
86 MatrixDouble matEigVec;
87 MatrixDouble matLogC;
88 MatrixDouble matLogCdC;
89 MatrixDouble matCdF;
90 MatrixDouble matFirstPiolaStress;
92 MatrixDouble matHenckyStress;
93 MatrixDouble matTangent;
95 inline auto getMatFirstPiolaStress() {
96 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
98 }
100 inline auto getMatHenckyStress() {
101 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
103 }
105 inline auto getMatLogC() {
106 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &matLogC);
107 }
109 inline auto getMatTangent() {
110 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &matTangent);
111 }
112};
113
114template <int DIM, IntegrationType I, typename DomainEleOp>
116
117template <int DIM, IntegrationType I, typename DomainEleOp>
119
120template <int DIM, IntegrationType I, typename DomainEleOp>
122
123template <int DIM, IntegrationType I, typename DomainEleOp, int S>
125
126template <int DIM, IntegrationType I, typename DomainEleOp, int S>
128
129template <int DIM, IntegrationType I, typename DomainEleOp, int S>
131
132template <int DIM, IntegrationType I, typename DomainEleOp, int S>
134
135template <int DIM, IntegrationType I, typename DomainEleOp, int S>
137
138template <int DIM, IntegrationType I, typename DomainEleOp, int S>
140
141template <int DIM, IntegrationType I, typename DomainEleOp, int S>
143
144template <int DIM, typename DomainEleOp>
145struct OpCalculateEigenValsImpl<DIM, GAUSS, DomainEleOp> : public DomainEleOp {
148 boost::shared_ptr<CommonData> common_data)
150 commonDataPtr(common_data) {
151 std::fill(&DomainEleOp::doEntities[MBEDGE],
152 &DomainEleOp::doEntities[MBMAXTYPE], false);
153 }
155 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
157
158 FTensor::Index<'i', DIM> i;
159 FTensor::Index<'j', DIM> j;
160 FTensor::Index<'k', DIM> k;
161
162 constexpr auto t_kd = FTensor::Kronecker_Delta<double>();
163 // const size_t nb_gauss_pts = matGradPtr->size2();
164 const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
165
166 auto t_grad = getFTensor2FromMat<DIM, DIM>(*(commonDataPtr->matGradPtr));
167
168 commonDataPtr->matEigVal.resize(DIM, nb_gauss_pts, false);
169 commonDataPtr->matEigVec.resize(DIM * DIM, nb_gauss_pts, false);
170 auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
171 auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
172
173 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
174
179
180 F(i, j) = t_grad(i, j) + t_kd(i, j);
181 C(i, j) = F(k, i) ^ F(k, j);
182
183 for (int ii = 0; ii != DIM; ii++)
184 for (int jj = 0; jj != DIM; jj++)
185 eigen_vec(ii, jj) = C(ii, jj);
186
187 CHKERR computeEigenValuesSymmetric(eigen_vec, eig);
188 for (auto ii = 0; ii != DIM; ++ii)
189 eig(ii) = std::max(std::numeric_limits<double>::epsilon(), eig(ii));
190
191 // rare case when two eigen values are equal
192 auto nb_uniq = get_uniq_nb<DIM>(&eig(0));
193 if constexpr (DIM == 3) {
194 if (nb_uniq == 2) {
195 sort_eigen_vals<DIM>(eig, eigen_vec);
196 }
197 }
198
199 t_eig_val(i) = eig(i);
200 t_eig_vec(i, j) = eigen_vec(i, j);
201
202 ++t_grad;
203 ++t_eig_val;
204 ++t_eig_vec;
205 }
206
208 }
209
210private:
211 boost::shared_ptr<CommonData> commonDataPtr;
212};
213
214template <int DIM, typename DomainEleOp>
215struct OpCalculateLogCImpl<DIM, GAUSS, DomainEleOp> : public DomainEleOp {
218 boost::shared_ptr<CommonData> common_data)
220 commonDataPtr(common_data) {
221 std::fill(&DomainEleOp::doEntities[MBEDGE],
222 &DomainEleOp::doEntities[MBMAXTYPE], false);
223 }
225 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
227
228 FTensor::Index<'i', DIM> i;
229 FTensor::Index<'j', DIM> j;
230
231 // const size_t nb_gauss_pts = matGradPtr->size2();
232 const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
233 constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
234 commonDataPtr->matLogC.resize(size_symm, nb_gauss_pts, false);
235
236 auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
237 auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
238
239 auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
240
241 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
242 t_logC(i, j) = EigenMatrix::getMat(t_eig_val, t_eig_vec, f)(i, j);
243 ++t_eig_val;
244 ++t_eig_vec;
245 ++t_logC;
246 }
247
249 }
250
251private:
252 boost::shared_ptr<CommonData> commonDataPtr;
253};
254
255template <int DIM, typename DomainEleOp>
256struct OpCalculateLogC_dCImpl<DIM, GAUSS, DomainEleOp> : public DomainEleOp {
259 boost::shared_ptr<CommonData> common_data)
261 commonDataPtr(common_data) {
262 std::fill(&DomainEleOp::doEntities[MBEDGE],
263 &DomainEleOp::doEntities[MBMAXTYPE], false);
264 }
266 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
268
269 FTensor::Index<'i', DIM> i;
270 FTensor::Index<'j', DIM> j;
271 FTensor::Index<'k', DIM> k;
272 FTensor::Index<'l', DIM> l;
273
274 const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
275 constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
276 commonDataPtr->matLogCdC.resize(size_symm * size_symm, nb_gauss_pts, false);
277 auto t_logC_dC = getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->matLogCdC);
278 auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
279 auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
280
281 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
282 // rare case when two eigen values are equal
283 auto nb_uniq = get_uniq_nb<DIM>(&t_eig_val(0));
284 t_logC_dC(i, j, k, l) =
285 2 * EigenMatrix::getDiffMat(t_eig_val, t_eig_vec, f, d_f,
286 nb_uniq)(i, j, k, l);
287
288 ++t_logC_dC;
289 ++t_eig_val;
290 ++t_eig_vec;
291 }
292
294 }
295
296private:
297 boost::shared_ptr<CommonData> commonDataPtr;
298};
299
300template <int DIM, typename DomainEleOp, int S>
301struct OpCalculateHenckyStressImpl<DIM, GAUSS, DomainEleOp, S>
302 : public DomainEleOp {
305 boost::shared_ptr<CommonData> common_data)
307 commonDataPtr(common_data) {
308 std::fill(&DomainEleOp::doEntities[MBEDGE],
309 &DomainEleOp::doEntities[MBMAXTYPE], false);
310 }
312 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
314
315 FTensor::Index<'i', DIM> i;
316 FTensor::Index<'j', DIM> j;
317 FTensor::Index<'k', DIM> k;
318 FTensor::Index<'l', DIM> l;
319
320 // const size_t nb_gauss_pts = matGradPtr->size2();
321 const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
322 auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*commonDataPtr->matDPtr);
323 auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
324 constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
325 commonDataPtr->matHenckyStress.resize(size_symm, nb_gauss_pts, false);
326 auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
327
328 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
329 t_T(i, j) = t_D(i, j, k, l) * t_logC(k, l);
330 ++t_logC;
331 ++t_T;
332 ++t_D;
333 }
334
336 }
337
338private:
339 boost::shared_ptr<CommonData> commonDataPtr;
340};
341
342template <int DIM, typename DomainEleOp, int S>
344 : public DomainEleOp {
347 const std::string field_name, boost::shared_ptr<VectorDouble> temperature,
348 boost::shared_ptr<CommonData> common_data,
349 boost::shared_ptr<VectorDouble> coeff_expansion_ptr,
350 boost::shared_ptr<double> ref_temp_ptr)
351 : DomainEleOp(field_name, DomainEleOp::OPROW), tempPtr(temperature),
352 commonDataPtr(common_data), coeffExpansionPtr(coeff_expansion_ptr),
353 refTempPtr(ref_temp_ptr) {
354 std::fill(&DomainEleOp::doEntities[MBEDGE],
355 &DomainEleOp::doEntities[MBMAXTYPE], false);
356 }
358 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
360
361 FTensor::Index<'i', DIM> i;
362 FTensor::Index<'j', DIM> j;
363 FTensor::Index<'k', DIM> k;
364 FTensor::Index<'l', DIM> l;
365
366 constexpr auto t_kd = FTensor::Kronecker_Delta<double>();
367
368 // const size_t nb_gauss_pts = matGradPtr->size2();
369 const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
370 auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*commonDataPtr->matDPtr);
371 auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
372 auto t_logC_dC = getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->matLogCdC);
373 constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
374 commonDataPtr->matHenckyStress.resize(size_symm, nb_gauss_pts, false);
375 commonDataPtr->matFirstPiolaStress.resize(DIM * DIM, nb_gauss_pts, false);
376 commonDataPtr->matSecondPiolaStress.resize(size_symm, nb_gauss_pts, false);
377 auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
378 auto t_P = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matFirstPiolaStress);
379 auto t_S =
380 getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matSecondPiolaStress);
381 auto t_grad = getFTensor2FromMat<DIM, DIM>(*(commonDataPtr->matGradPtr));
382 auto t_temp = getFTensor0FromVec(*tempPtr);
383
385 t_coeff_exp(i, j) = 0;
386 for (auto d = 0; d != SPACE_DIM; ++d) {
387 t_coeff_exp(d, d) = (*coeffExpansionPtr)[d];
388 }
389
390 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
391#ifdef HENCKY_SMALL_STRAIN
392 t_P(i, j) = t_D(i, j, k, l) *
393 (t_grad(k, l) - t_coeff_exp(k, l) * (t_temp - (*refTempPtr)));
394#else
395 t_T(i, j) = t_D(i, j, k, l) *
396 (t_logC(k, l) - t_coeff_exp(k, l) * (t_temp - (*refTempPtr)));
398 t_F(i, j) = t_grad(i, j) + t_kd(i, j);
399 t_S(k, l) = t_T(i, j) * t_logC_dC(i, j, k, l);
400 t_P(i, l) = t_F(i, k) * t_S(k, l);
401#endif
402 ++t_grad;
403 ++t_logC;
404 ++t_logC_dC;
405 ++t_P;
406 ++t_T;
407 ++t_S;
408 ++t_D;
409 ++t_temp;
410 }
411
413 }
414
415private:
416 boost::shared_ptr<CommonData> commonDataPtr;
417 boost::shared_ptr<VectorDouble> tempPtr;
418 boost::shared_ptr<VectorDouble> coeffExpansionPtr;
419 boost::shared_ptr<double> refTempPtr;
420};
421
422template <int DIM, typename DomainEleOp, int S>
424 : public DomainEleOp {
427 boost::shared_ptr<CommonData> common_data,
428 boost::shared_ptr<MatrixDouble> mat_D_ptr,
429 const double scale = 1)
430 : DomainEleOp(field_name, DomainEleOp::OPROW), commonDataPtr(common_data),
431 scaleStress(scale), matDPtr(mat_D_ptr) {
432 std::fill(&DomainEleOp::doEntities[MBEDGE],
433 &DomainEleOp::doEntities[MBMAXTYPE], false);
434
435 matLogCPlastic = commonDataPtr->matLogCPlastic;
436 }
438 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
440
441 FTensor::Index<'i', DIM> i;
442 FTensor::Index<'j', DIM> j;
443 FTensor::Index<'k', DIM> k;
444 FTensor::Index<'l', DIM> l;
445
446 // const size_t nb_gauss_pts = matGradPtr->size2();
447 const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
448 auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*matDPtr);
449 auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
450 auto t_logCPlastic = getFTensor2SymmetricFromMat<DIM>(*matLogCPlastic);
451 constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
452 commonDataPtr->matHenckyStress.resize(size_symm, nb_gauss_pts, false);
453 auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
454
455 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
456 t_T(i, j) = t_D(i, j, k, l) * (t_logC(k, l) - t_logCPlastic(k, l));
457 t_T(i, j) /= scaleStress;
458 ++t_logC;
459 ++t_T;
460 ++t_D;
461 ++t_logCPlastic;
462 }
463
465 }
466
467private:
468 boost::shared_ptr<CommonData> commonDataPtr;
469 boost::shared_ptr<MatrixDouble> matDPtr;
470 boost::shared_ptr<MatrixDouble> matLogCPlastic;
471 const double scaleStress;
472};
473
474template <int DIM, typename DomainEleOp, int S>
476 : public DomainEleOp {
479 const std::string field_name, boost::shared_ptr<VectorDouble> temperature,
480 boost::shared_ptr<CommonData> common_data,
481 boost::shared_ptr<VectorDouble> coeff_expansion_ptr,
482 boost::shared_ptr<double> ref_temp_ptr)
483 : DomainEleOp(field_name, DomainEleOp::OPROW), tempPtr(temperature),
484 commonDataPtr(common_data), coeffExpansionPtr(coeff_expansion_ptr),
485 refTempPtr(ref_temp_ptr) {
486 std::fill(&DomainEleOp::doEntities[MBEDGE],
487 &DomainEleOp::doEntities[MBMAXTYPE], false);
488
489 matLogCPlastic = commonDataPtr->matLogCPlastic;
490 }
491
492 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
494
495 FTensor::Index<'i', DIM> i;
496 FTensor::Index<'j', DIM> j;
497 FTensor::Index<'k', DIM> k;
498 FTensor::Index<'l', DIM> l;
499
500 constexpr auto t_kd = FTensor::Kronecker_Delta<double>();
501
502 // const size_t nb_gauss_pts = matGradPtr->size2();
503 const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
504 auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*commonDataPtr->matDPtr);
505 auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
506 auto t_logC_dC = getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->matLogCdC);
507 auto t_logCPlastic = getFTensor2SymmetricFromMat<DIM>(*matLogCPlastic);
508 constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
509 commonDataPtr->matHenckyStress.resize(size_symm, nb_gauss_pts, false);
510 commonDataPtr->matFirstPiolaStress.resize(DIM * DIM, nb_gauss_pts, false);
511 commonDataPtr->matSecondPiolaStress.resize(size_symm, nb_gauss_pts, false);
512 auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
513 auto t_P = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matFirstPiolaStress);
514 auto t_S =
515 getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matSecondPiolaStress);
516 auto t_grad = getFTensor2FromMat<DIM, DIM>(*(commonDataPtr->matGradPtr));
517 auto t_temp = getFTensor0FromVec(*tempPtr);
518
520 t_coeff_exp(i, j) = 0;
521 for (auto d = 0; d != SPACE_DIM; ++d) {
522 t_coeff_exp(d, d) = (*coeffExpansionPtr)[d];
523 }
524
525 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
526#ifdef HENCKY_SMALL_STRAIN
527 t_P(i, j) =
528 t_D(i, j, k, l) * (t_grad(k, l) - t_logCPlastic(k, l) -
529 t_coeff_exp(k, l) * (t_temp - (*refTempPtr)));
530#else
531 t_T(i, j) =
532 t_D(i, j, k, l) * (t_logC(k, l) - t_logCPlastic(k, l) -
533 t_coeff_exp(k, l) * (t_temp - (*refTempPtr)));
535 t_F(i, j) = t_grad(i, j) + t_kd(i, j);
536 t_S(k, l) = t_T(i, j) * t_logC_dC(i, j, k, l);
537 t_P(i, l) = t_F(i, k) * t_S(k, l);
538#endif
539 ++t_grad;
540 ++t_logC;
541 ++t_P;
542 ++t_T;
543 ++t_S;
544 ++t_D;
545 ++t_temp;
546 ++t_logCPlastic;
547 }
548
550 }
551
552private:
553 boost::shared_ptr<CommonData> commonDataPtr;
554 boost::shared_ptr<VectorDouble> tempPtr;
555 boost::shared_ptr<MatrixDouble> matLogCPlastic;
556 boost::shared_ptr<VectorDouble> coeffExpansionPtr;
557 boost::shared_ptr<double> refTempPtr;
558};
559
560template <int DIM, typename DomainEleOp, int S>
561struct OpCalculatePiolaStressImpl<DIM, GAUSS, DomainEleOp, S>
562 : public DomainEleOp {
563
565 boost::shared_ptr<CommonData> common_data)
567 commonDataPtr(common_data) {
568 std::fill(&DomainEleOp::doEntities[MBEDGE],
569 &DomainEleOp::doEntities[MBMAXTYPE], false);
570 }
571
572 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
574
575 FTensor::Index<'i', DIM> i;
576 FTensor::Index<'j', DIM> j;
577 FTensor::Index<'k', DIM> k;
578 FTensor::Index<'l', DIM> l;
579
580 constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
581
582 // const size_t nb_gauss_pts = matGradPtr->size2();
583 const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
584#ifdef HENCKY_SMALL_STRAIN
585 auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*commonDataPtr->matDPtr);
586#endif
587 auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
588 auto t_logC_dC = getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->matLogCdC);
589 constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
590 commonDataPtr->matFirstPiolaStress.resize(DIM * DIM, nb_gauss_pts, false);
591 commonDataPtr->matSecondPiolaStress.resize(size_symm, nb_gauss_pts, false);
592 auto t_P = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matFirstPiolaStress);
593 auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
594 auto t_S =
595 getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matSecondPiolaStress);
596 auto t_grad = getFTensor2FromMat<DIM, DIM>(*(commonDataPtr->matGradPtr));
597
598 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
599
600#ifdef HENCKY_SMALL_STRAIN
601 t_P(i, j) = t_D(i, j, k, l) * t_grad(k, l);
602#else
604 t_F(i, j) = t_grad(i, j) + t_kd(i, j);
605 t_S(k, l) = t_T(i, j) * t_logC_dC(i, j, k, l);
606 t_P(i, l) = t_F(i, k) * t_S(k, l);
607#endif
608
609 ++t_grad;
610 ++t_logC;
611 ++t_logC_dC;
612 ++t_P;
613 ++t_T;
614 ++t_S;
615#ifdef HENCKY_SMALL_STRAIN
616 ++t_D;
617#endif
618 }
619
621 }
622
623private:
624 boost::shared_ptr<CommonData> commonDataPtr;
625};
627template <int DIM, typename DomainEleOp, int S>
628struct OpHenckyTangentImpl<DIM, GAUSS, DomainEleOp, S> : public DomainEleOp {
630 boost::shared_ptr<CommonData> common_data,
631 boost::shared_ptr<MatrixDouble> mat_D_ptr = nullptr)
633 commonDataPtr(common_data) {
634 std::fill(&DomainEleOp::doEntities[MBEDGE],
635 &DomainEleOp::doEntities[MBMAXTYPE], false);
636 if (mat_D_ptr)
637 matDPtr = mat_D_ptr;
638 else
639 matDPtr = commonDataPtr->matDPtr;
640 }
641
642 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
644
645 FTensor::Index<'i', DIM> i;
646 FTensor::Index<'j', DIM> j;
647 FTensor::Index<'k', DIM> k;
648 FTensor::Index<'l', DIM> l;
649 FTensor::Index<'m', DIM> m;
650 FTensor::Index<'n', DIM> n;
651 FTensor::Index<'o', DIM> o;
652 FTensor::Index<'p', DIM> p;
653
654 constexpr auto t_kd = FTensor::Kronecker_Delta<double>();
655 // const size_t nb_gauss_pts = matGradPtr->size2();
656 const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
657 commonDataPtr->matTangent.resize(DIM * DIM * DIM * DIM, nb_gauss_pts);
658 auto dP_dF =
659 getFTensor4FromMat<DIM, DIM, DIM, DIM, 1>(commonDataPtr->matTangent);
660
661 auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*matDPtr);
662 auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
663 auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
664 auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
665 auto t_S =
666 getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matSecondPiolaStress);
667 auto t_grad = getFTensor2FromMat<DIM, DIM>(*(commonDataPtr->matGradPtr));
668 auto t_logC_dC = getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->matLogCdC);
669 commonDataPtr->matCdF.resize(DIM * DIM * DIM * DIM, nb_gauss_pts);
670 auto dC_dF =
671 getFTensor4FromMat<DIM, DIM, DIM, DIM, 1>(commonDataPtr->matCdF);
672
673 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
674
675#ifdef HENCKY_SMALL_STRAIN
676 dP_dF(i, j, k, l) = t_D(i, j, k, l);
677#else
678
680 t_F(i, j) = t_grad(i, j) + t_kd(i, j);
681
682 // rare case when two eigen values are equal
683 auto nb_uniq = get_uniq_nb<DIM>(&t_eig_val(0));
684
685 dC_dF(i, j, k, l) = (t_kd(i, l) * t_F(k, j)) + (t_kd(j, l) * t_F(k, i));
686
687 auto TL = EigenMatrix::getDiffDiffMat(t_eig_val, t_eig_vec, f, d_f, dd_f,
688 t_T, nb_uniq);
689 TL(i, j, k, l) *= 4;
690 FTensor::Ddg<double, DIM, DIM> P_D_P_plus_TL;
691 P_D_P_plus_TL(i, j, k, l) =
692 TL(i, j, k, l) +
693 (t_logC_dC(i, j, o, p) * t_D(o, p, m, n)) * t_logC_dC(m, n, k, l);
694 P_D_P_plus_TL(i, j, k, l) *= 0.5;
695 dP_dF(i, j, m, n) = t_kd(i, m) * (t_kd(k, n) * t_S(k, j));
696 dP_dF(i, j, m, n) +=
697 t_F(i, k) * (P_D_P_plus_TL(k, j, o, p) * dC_dF(o, p, m, n));
698
699#endif
700
701 ++dP_dF;
702
703 ++t_grad;
704 ++t_eig_val;
705 ++t_eig_vec;
706 ++t_logC_dC;
707 ++t_S;
708 ++t_T;
709 ++t_D;
710 ++dC_dF;
711 }
712
714 }
715
716private:
717 boost::shared_ptr<CommonData> commonDataPtr;
718 boost::shared_ptr<MatrixDouble> matDPtr;
719};
720
721template <int DIM, typename AssemblyDomainEleOp, int S>
722struct OpCalculateHenckyThermalStressdTImpl<DIM, GAUSS, AssemblyDomainEleOp, S>
723 : public AssemblyDomainEleOp {
725 const std::string row_field_name, const std::string col_field_name,
726 boost::shared_ptr<CommonData> elastic_common_data_ptr,
727 boost::shared_ptr<VectorDouble> coeff_expansion_ptr);
728
729 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data,
730 EntitiesFieldData::EntData &col_data);
731
732private:
733 boost::shared_ptr<CommonData> elasticCommonDataPtr;
734 boost::shared_ptr<VectorDouble> coeffExpansionPtr;
735};
737template <int DIM, typename AssemblyDomainEleOp, int S>
740 const std::string row_field_name, const std::string col_field_name,
741 boost::shared_ptr<CommonData> elastic_common_data_ptr,
742 boost::shared_ptr<VectorDouble> coeff_expansion_ptr)
743 : AssemblyDomainEleOp(row_field_name, col_field_name,
745 elasticCommonDataPtr(elastic_common_data_ptr),
746 coeffExpansionPtr(coeff_expansion_ptr) {
747 this->sYmm = false;
748}
749
750template <int DIM, typename AssemblyDomainEleOp, int S>
751MoFEMErrorCode
753 iNtegrate(EntitiesFieldData::EntData &row_data,
754 EntitiesFieldData::EntData &col_data) {
756
757 auto &locMat = AssemblyDomainEleOp::locMat;
758
759 const auto nb_integration_pts = row_data.getN().size1();
760 const auto nb_row_base_functions = row_data.getN().size2();
761 auto t_w = this->getFTensor0IntegrationWeight();
762
764 auto t_row_diff_base = row_data.getFTensor1DiffN<DIM>();
765 auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*elasticCommonDataPtr->matDPtr);
766 auto t_grad =
767 getFTensor2FromMat<DIM, DIM>(*(elasticCommonDataPtr->matGradPtr));
768 auto t_logC_dC =
769 getFTensor4DdgFromMat<DIM, DIM>(elasticCommonDataPtr->matLogCdC);
773 FTensor::Index<'i', DIM> i;
774 FTensor::Index<'j', DIM> j;
775 FTensor::Index<'k', DIM> k;
776 FTensor::Index<'l', DIM> l;
777 FTensor::Index<'m', DIM> m;
778 FTensor::Index<'n', DIM> n;
779 FTensor::Index<'o', DIM> o;
780
782 t_coeff_exp(i, j) = 0;
783 for (auto d = 0; d != SPACE_DIM; ++d) {
784 t_coeff_exp(d, d) = (*coeffExpansionPtr)[d];
785 }
786
787 t_eigen_strain(i, j) = t_D(i, j, k, l) * t_coeff_exp(k, l);
788
789 for (auto gg = 0; gg != nb_integration_pts; ++gg) {
790
791 t_F(i, j) = t_grad(i, j) + t_kd(i, j);
792
793 double alpha = this->getMeasure() * t_w;
794 auto rr = 0;
795 for (; rr != AssemblyDomainEleOp::nbRows / DIM; ++rr) {
796 auto t_mat = getFTensor1FromMat<DIM, 1>(locMat, rr * DIM);
797 auto t_col_base = col_data.getFTensor0N(gg, 0);
798 for (auto cc = 0; cc != AssemblyDomainEleOp::nbCols; cc++) {
799#ifdef HENCKY_SMALL_STRAIN
800 t_mat(i) -=
801 (t_row_diff_base(j) * t_eigen_strain(i, j)) * (t_col_base * alpha);
802#else
803 t_mat(i) -= (t_row_diff_base(j) *
804 (t_F(i, o) * ((t_D(m, n, k, l) * t_coeff_exp(k, l)) *
805 t_logC_dC(m, n, o, j)))) *
806 (t_col_base * alpha);
807#endif
808
809 ++t_mat;
810 ++t_col_base;
811 }
812
813 ++t_row_diff_base;
814 }
815 for (; rr != nb_row_base_functions; ++rr)
816 ++t_row_diff_base;
817
818 ++t_w;
819 ++t_grad;
820 ++t_logC_dC;
821 ++t_D;
822 }
823
825}
826
827template <typename DomainEleOp> struct HenckyIntegrators {
828 template <int DIM, IntegrationType I>
830
831 template <int DIM, IntegrationType I>
833
834 template <int DIM, IntegrationType I>
836
837 template <int DIM, IntegrationType I, int S>
840
841 template <int DIM, IntegrationType I, int S>
844
845 template <int DIM, IntegrationType I, int S>
848
849 template <int DIM, IntegrationType I, int S>
852
853 template <int DIM, IntegrationType I, int S>
856
857 template <int DIM, IntegrationType I, int S>
859
860 template <int DIM, IntegrationType I, typename AssemblyDomainEleOp, int S>
863};
864
865template <int DIM>
866MoFEMErrorCode
868 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
869 std::string block_name,
870 boost::shared_ptr<MatrixDouble> mat_D_Ptr, Sev sev,
871 double scale = 1) {
873
874 PetscBool plane_strain_flag = PETSC_FALSE;
875 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-plane_strain",
876 &plane_strain_flag, PETSC_NULLPTR);
877
878 struct OpMatBlocks : public DomainEleOp {
879 OpMatBlocks(boost::shared_ptr<MatrixDouble> m, double bulk_modulus_K,
880 double shear_modulus_G, MoFEM::Interface &m_field, Sev sev,
881 std::vector<const CubitMeshSets *> meshset_vec_ptr,
882 double scale, PetscBool plane_strain_flag)
883 : DomainEleOp(NOSPACE, DomainEleOp::OPSPACE), matDPtr(m),
884 bulkModulusKDefault(bulk_modulus_K),
885 shearModulusGDefault(shear_modulus_G), scaleYoungModulus(scale),
886 planeStrainFlag(plane_strain_flag) {
887 CHK_THROW_MESSAGE(extractBlockData(m_field, meshset_vec_ptr, sev),
888 "Can not get data from block");
889 }
890
891 MoFEMErrorCode doWork(int side, EntityType type,
892 EntitiesFieldData::EntData &data) {
894
895 for (auto &b : blockData) {
896
897 if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
898 CHKERR getMatDPtr(matDPtr, b.bulkModulusK * scaleYoungModulus,
899 b.shearModulusG * scaleYoungModulus,
900 planeStrainFlag);
902 }
903 }
904
905 CHKERR getMatDPtr(matDPtr, bulkModulusKDefault * scaleYoungModulus,
906 shearModulusGDefault * scaleYoungModulus,
907 planeStrainFlag);
909 }
910
911 private:
912 boost::shared_ptr<MatrixDouble> matDPtr;
913 const double scaleYoungModulus;
914 const PetscBool planeStrainFlag;
915
916 struct BlockData {
917 double bulkModulusK;
918 double shearModulusG;
919 Range blockEnts;
920 };
921
922 double bulkModulusKDefault;
923 double shearModulusGDefault;
924 std::vector<BlockData> blockData;
926 MoFEMErrorCode
927 extractBlockData(MoFEM::Interface &m_field,
928 std::vector<const CubitMeshSets *> meshset_vec_ptr,
929 Sev sev) {
931
932 for (auto m : meshset_vec_ptr) {
933 MOFEM_TAG_AND_LOG("WORLD", sev, "MatBlock") << *m;
934 std::vector<double> block_data;
935 CHKERR m->getAttributes(block_data);
936 if (block_data.size() < 2) {
937 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
938 "Expected that block has two attribute");
939 }
940 auto get_block_ents = [&]() {
941 Range ents;
942 CHKERR
943 m_field.get_moab().get_entities_by_handle(m->meshset, ents, true);
944 return ents;
945 };
946
947 double young_modulus = block_data[0];
948 double poisson_ratio = block_data[1];
949 double bulk_modulus_K = young_modulus / (3 * (1 - 2 * poisson_ratio));
950 double shear_modulus_G = young_modulus / (2 * (1 + poisson_ratio));
951
952 MOFEM_TAG_AND_LOG("WORLD", sev, "MatBlock")
953 << "E = " << young_modulus << " nu = " << poisson_ratio;
954
955 blockData.push_back(
956 {bulk_modulus_K, shear_modulus_G, get_block_ents()});
957 }
960 }
961
962 MoFEMErrorCode getMatDPtr(boost::shared_ptr<MatrixDouble> mat_D_ptr,
963 double bulk_modulus_K, double shear_modulus_G,
964 PetscBool is_plane_strain) {
966 //! [Calculate elasticity tensor]
967 auto set_material_stiffness = [&]() {
968 FTensor::Index<'i', DIM> i;
969 FTensor::Index<'j', DIM> j;
970 FTensor::Index<'k', DIM> k;
971 FTensor::Index<'l', DIM> l;
973 double A = (SPACE_DIM == 2 && !is_plane_strain)
974 ? 2 * shear_modulus_G /
975 (bulk_modulus_K + (4. / 3.) * shear_modulus_G)
976 : 1;
977 auto t_D = getFTensor4DdgFromMat<DIM, DIM, 0>(*mat_D_ptr);
978 t_D(i, j, k, l) =
979 2 * shear_modulus_G * ((t_kd(i, k) ^ t_kd(j, l)) / 4.) +
980 A * (bulk_modulus_K - (2. / 3.) * shear_modulus_G) * t_kd(i, j) *
981 t_kd(k, l);
982 };
983 //! [Calculate elasticity tensor]
984 constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
985 mat_D_ptr->resize(size_symm * size_symm, 1);
986 set_material_stiffness();
988 }
989 };
990
991 double E = young_modulus;
992 double nu = poisson_ratio;
993
994 PetscOptionsBegin(PETSC_COMM_WORLD, "", "", "none");
995 CHKERR PetscOptionsScalar("-young_modulus", "Young modulus", "", E, &E,
996 PETSC_NULLPTR);
997 CHKERR PetscOptionsScalar("-poisson_ratio", "poisson ratio", "", nu, &nu,
998 PETSC_NULLPTR);
999 PetscOptionsEnd();
1000
1001 double bulk_modulus_K = E / (3 * (1 - 2 * nu));
1002 double shear_modulus_G = E / (2 * (1 + nu));
1003 pip.push_back(new OpMatBlocks(
1004 mat_D_Ptr, bulk_modulus_K, shear_modulus_G, m_field, sev,
1005
1006 // Get blockset using regular expression
1007 m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
1008
1009 (boost::format("%s(.*)") % block_name).str()
1010
1011 )),
1012 scale, plane_strain_flag
1014 ));
1015
1017}
1018
1019template <int DIM, IntegrationType I, typename DomainEleOp>
1021 MoFEM::Interface &m_field,
1022 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
1023 std::string field_name, std::string block_name, Sev sev, double scale = 1) {
1024
1025 auto common_ptr = boost::make_shared<HenckyOps::CommonData>();
1026 common_ptr->matDPtr = boost::make_shared<MatrixDouble>();
1027 common_ptr->matGradPtr = boost::make_shared<MatrixDouble>();
1028
1029 CHK_THROW_MESSAGE(addMatBlockOps<DIM>(m_field, pip, block_name,
1030 common_ptr->matDPtr, sev, scale),
1031 "addMatBlockOps");
1032
1033 using H = HenckyIntegrators<DomainEleOp>;
1034
1035 pip.push_back(new OpCalculateVectorFieldGradient<DIM, DIM>(
1036 field_name, common_ptr->matGradPtr));
1037 pip.push_back(new typename H::template OpCalculateEigenVals<DIM, I>(
1038 field_name, common_ptr));
1039 pip.push_back(
1040 new typename H::template OpCalculateLogC<DIM, I>(field_name, common_ptr));
1041 pip.push_back(new typename H::template OpCalculateLogC_dC<DIM, I>(
1042 field_name, common_ptr));
1043 // Assumes constant D matrix per entity
1044 pip.push_back(new typename H::template OpCalculateHenckyStress<DIM, I, 0>(
1045 field_name, common_ptr));
1046 pip.push_back(new typename H::template OpCalculatePiolaStress<DIM, I, 0>(
1047 field_name, common_ptr));
1048
1049 return common_ptr;
1050}
1051
1052template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
1054 MoFEM::Interface &m_field,
1055 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
1056 std::string field_name, boost::shared_ptr<HenckyOps::CommonData> common_ptr,
1057 Sev sev) {
1059
1060 using B = typename FormsIntegrators<DomainEleOp>::template Assembly<
1061 A>::template LinearForm<I>;
1062 using OpInternalForcePiola =
1063 typename B::template OpGradTimesTensor<1, DIM, DIM>;
1064 pip.push_back(
1065 new OpInternalForcePiola("U", common_ptr->getMatFirstPiolaStress()));
1066
1068}
1069
1070template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
1072 MoFEM::Interface &m_field,
1073 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
1074 std::string field_name, std::string block_name, Sev sev, double scale = 1) {
1076
1077 auto common_ptr = commonDataFactory<DIM, I, DomainEleOp>(
1078 m_field, pip, field_name, block_name, sev, scale);
1079 CHKERR opFactoryDomainRhs<DIM, A, I, DomainEleOp>(m_field, pip, field_name,
1080 common_ptr, sev);
1081
1083}
1084
1085template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
1087 MoFEM::Interface &m_field,
1088 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
1089 std::string field_name, boost::shared_ptr<HenckyOps::CommonData> common_ptr,
1090 Sev sev) {
1092
1093 using B = typename FormsIntegrators<DomainEleOp>::template Assembly<
1094 A>::template BiLinearForm<I>;
1095 using OpKPiola = typename B::template OpGradTensorGrad<1, DIM, DIM, 1>;
1096
1097 using H = HenckyIntegrators<DomainEleOp>;
1098 // Assumes constant D matrix per entity
1099 pip.push_back(new typename H::template OpHenckyTangent<DIM, I, 0>(
1100 field_name, common_ptr));
1101 pip.push_back(
1102 new OpKPiola(field_name, field_name, common_ptr->getMatTangent()));
1103
1105}
1106
1107template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
1109 MoFEM::Interface &m_field,
1110 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
1111 std::string field_name, std::string block_name, Sev sev, double scale = 1) {
1113
1114 auto common_ptr = commonDataFactory<DIM, I, DomainEleOp>(
1115 m_field, pip, field_name, block_name, sev, scale);
1116 CHKERR opFactoryDomainLhs<DIM, A, I, DomainEleOp>(m_field, pip, field_name,
1117 common_ptr, sev);
1118
1120}
1121} // namespace HenckyOps
1122
1123#endif // __HENCKY_OPS_HPP__
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
PetscBool is_plane_strain
Definition adjoint.cpp:42
constexpr double a
constexpr int SPACE_DIM
[Define dimension]
DomainEle::UserDataOperator DomainEleOp
Finire element operator type.
Kronecker Delta class symmetric.
Kronecker Delta class.
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ NOSPACE
Definition definitions.h:83
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ 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.
double bulk_modulus_K
double shear_modulus_G
@ F
constexpr auto t_kd
@ GAUSS
Gaussian quadrature integration.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
FTensor::Index< 'i', SPACE_DIM > i
const double v
phase velocity of light in medium (cm/ns)
const double n
refractive index of diffusive medium
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.
auto get_uniq_nb(double *ptr)
Definition HenckyOps.hpp:32
auto sort_eigen_vals(FTensor::Tensor1< double, DIM > &eig, FTensor::Tensor2< double, DIM, DIM > &eigen_vec)
Definition HenckyOps.hpp:41
MoFEMErrorCode opFactoryDomainLhs(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string field_name, boost::shared_ptr< HenckyOps::CommonData > common_ptr, Sev sev)
constexpr double eps
Definition HenckyOps.hpp:13
MoFEMErrorCode opFactoryDomainRhs(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string field_name, boost::shared_ptr< HenckyOps::CommonData > common_ptr, Sev sev)
MoFEMErrorCode addMatBlockOps(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string block_name, boost::shared_ptr< MatrixDouble > mat_D_Ptr, Sev sev, double scale=1)
auto commonDataFactory(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string field_name, std::string block_name, Sev sev, double scale=1)
auto is_eq(const double &a, const double &b)
Definition HenckyOps.hpp:28
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
FormsIntegrators< DomainEleOp >::Assembly< A >::LinearForm< I >::OpGradTimesTensor< 1, FIELD_DIM, SPACE_DIM > OpGradTimesTensor
constexpr AssemblyType A
[Define dimension]
constexpr auto field_name
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesTensor< 1, SPACE_DIM, SPACE_DIM > OpInternalForcePiola
Definition seepage.cpp:65
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradTensorGrad< 1, SPACE_DIM, SPACE_DIM, 1 > OpKPiola
[Only used for dynamics]
Definition seepage.cpp:63
FTensor::Index< 'm', 3 > m
MatrixDouble matEigVec
Definition HenckyOps.hpp:86
boost::shared_ptr< MatrixDouble > matLogCPlastic
Definition HenckyOps.hpp:83
MatrixDouble matCdF
Definition HenckyOps.hpp:89
MatrixDouble matHenckyStress
Definition HenckyOps.hpp:91
MatrixDouble matLogCdC
Definition HenckyOps.hpp:88
MatrixDouble matEigVal
Definition HenckyOps.hpp:85
MatrixDouble matFirstPiolaStress
Definition HenckyOps.hpp:89
boost::shared_ptr< MatrixDouble > matDPtr
Definition HenckyOps.hpp:82
MatrixDouble matSecondPiolaStress
Definition HenckyOps.hpp:90
MatrixDouble matLogC
Definition HenckyOps.hpp:87
MatrixDouble matTangent
Definition HenckyOps.hpp:92
boost::shared_ptr< MatrixDouble > matGradPtr
Definition HenckyOps.hpp:81
OpCalculateEigenValsImpl(const std::string field_name, boost::shared_ptr< CommonData > common_data)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
OpCalculateHenckyPlasticStressImpl(const std::string field_name, boost::shared_ptr< CommonData > common_data, boost::shared_ptr< MatrixDouble > mat_D_ptr, const double scale=1)
OpCalculateHenckyStressImpl(const std::string field_name, boost::shared_ptr< CommonData > common_data)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
OpCalculateHenckyThermalStressImpl(const std::string field_name, boost::shared_ptr< VectorDouble > temperature, boost::shared_ptr< CommonData > common_data, boost::shared_ptr< VectorDouble > coeff_expansion_ptr, boost::shared_ptr< double > ref_temp_ptr)
OpCalculateHenckyThermalStressdTImpl(const std::string row_field_name, const std::string col_field_name, boost::shared_ptr< CommonData > elastic_common_data_ptr, boost::shared_ptr< VectorDouble > coeff_expansion_ptr)
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
OpCalculateHenckyThermoPlasticStressImpl(const std::string field_name, boost::shared_ptr< VectorDouble > temperature, boost::shared_ptr< CommonData > common_data, boost::shared_ptr< VectorDouble > coeff_expansion_ptr, boost::shared_ptr< double > ref_temp_ptr)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
OpCalculateLogCImpl(const std::string field_name, boost::shared_ptr< CommonData > common_data)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
OpCalculateLogC_dCImpl(const std::string field_name, boost::shared_ptr< CommonData > common_data)
OpCalculatePiolaStressImpl(const std::string field_name, boost::shared_ptr< CommonData > common_data)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
OpHenckyTangentImpl(const std::string field_name, boost::shared_ptr< CommonData > common_data, boost::shared_ptr< MatrixDouble > mat_D_ptr=nullptr)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
static double absMax
Definition HenckyOps.hpp:23
static auto check(const double &a, const double &b)
Definition HenckyOps.hpp:20
virtual moab::Interface & get_moab()=0
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
double young_modulus
Young modulus.
Definition plastic.cpp:125
double poisson_ratio
Poisson ratio.
Definition plastic.cpp:126
double scale
Definition plastic.cpp:123
constexpr auto size_symm
Definition plastic.cpp:42
double H
Hardening.
Definition plastic.cpp:128