v0.14.0
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 
11 namespace HenckyOps {
12 
13 constexpr double eps = std::numeric_limits<float>::epsilon();
14 
15 auto f = [](double v) { return 0.5 * std::log(v); };
16 auto d_f = [](double v) { return 0.5 / v; };
17 auto dd_f = [](double v) { return -0.5 / (v * v); };
18 
19 struct 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 
26 double isEq::absMax = 1;
27 
28 inline auto is_eq(const double &a, const double &b) {
29  return isEq::check(a, b);
30 };
31 
32 template <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 
40 template <int DIM>
43 
44  isEq::absMax =
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  {
75  eig(i) = eig_c(i);
76  eigen_vec(i, j) = eigen_vec_c(i, j);
77  }
78 };
79 
80 struct 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 
93 
94  inline auto getMatFirstPiolaStress() {
95  return boost::shared_ptr<MatrixDouble>(shared_from_this(),
97  }
98 
99  inline auto getMatHenckyStress() {
100  return boost::shared_ptr<MatrixDouble>(shared_from_this(),
101  &matHenckyStress);
102  }
103 
104  inline auto getMatLogC() {
105  return boost::shared_ptr<MatrixDouble>(shared_from_this(), &matLogC);
106  }
107 
108  inline auto getMatTangent() {
109  return boost::shared_ptr<MatrixDouble>(shared_from_this(), &matTangent);
110  }
111 };
112 
113 template <int DIM, IntegrationType I, typename DomainEleOp>
115 
116 template <int DIM, IntegrationType I, typename DomainEleOp>
118 
119 template <int DIM, IntegrationType I, typename DomainEleOp>
121 
122 template <int DIM, IntegrationType I, typename DomainEleOp, int S>
124 
125 template <int DIM, IntegrationType I, typename DomainEleOp, int S>
127 
128 template <int DIM, IntegrationType I, typename DomainEleOp, int S>
130 
131 template <int DIM, IntegrationType I, typename DomainEleOp, int S>
133 
134 template <int DIM, typename DomainEleOp>
136 
138  boost::shared_ptr<CommonData> common_data)
140  commonDataPtr(common_data) {
141  std::fill(&DomainEleOp::doEntities[MBEDGE],
142  &DomainEleOp::doEntities[MBMAXTYPE], false);
143  }
144 
145  MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
147 
151 
152  constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
153  // const size_t nb_gauss_pts = matGradPtr->size2();
154  const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
155 
156  auto t_grad = getFTensor2FromMat<DIM, DIM>(*(commonDataPtr->matGradPtr));
157 
158  commonDataPtr->matEigVal.resize(DIM, nb_gauss_pts, false);
159  commonDataPtr->matEigVec.resize(DIM * DIM, nb_gauss_pts, false);
160  auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
161  auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
162 
163  for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
164 
169 
170  F(i, j) = t_grad(i, j) + t_kd(i, j);
171  C(i, j) = F(k, i) ^ F(k, j);
172 
173  for (int ii = 0; ii != DIM; ii++)
174  for (int jj = 0; jj != DIM; jj++)
175  eigen_vec(ii, jj) = C(ii, jj);
176 
177  CHKERR computeEigenValuesSymmetric(eigen_vec, eig);
178  for (auto ii = 0; ii != DIM; ++ii)
179  eig(ii) = std::max(std::numeric_limits<double>::epsilon(), eig(ii));
180 
181  // rare case when two eigen values are equal
182  auto nb_uniq = get_uniq_nb<DIM>(&eig(0));
183  if constexpr (DIM == 3) {
184  if (nb_uniq == 2) {
185  sort_eigen_vals<DIM>(eig, eigen_vec);
186  }
187  }
188 
189  t_eig_val(i) = eig(i);
190  t_eig_vec(i, j) = eigen_vec(i, j);
191 
192  ++t_grad;
193  ++t_eig_val;
194  ++t_eig_vec;
195  }
196 
198  }
199 
200 private:
201  boost::shared_ptr<CommonData> commonDataPtr;
202 };
203 
204 template <int DIM, typename DomainEleOp>
206 
207  OpCalculateLogCImpl(const std::string field_name,
208  boost::shared_ptr<CommonData> common_data)
210  commonDataPtr(common_data) {
211  std::fill(&DomainEleOp::doEntities[MBEDGE],
212  &DomainEleOp::doEntities[MBMAXTYPE], false);
213  }
214 
215  MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
217 
220 
221  // const size_t nb_gauss_pts = matGradPtr->size2();
222  const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
223  constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
224  commonDataPtr->matLogC.resize(size_symm, nb_gauss_pts, false);
225 
226  auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
227  auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
228 
229  auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
230 
231  for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
232  t_logC(i, j) = EigenMatrix::getMat(t_eig_val, t_eig_vec, f)(i, j);
233  ++t_eig_val;
234  ++t_eig_vec;
235  ++t_logC;
236  }
237 
239  }
240 
241 private:
242  boost::shared_ptr<CommonData> commonDataPtr;
243 };
244 
245 template <int DIM, typename DomainEleOp>
247 
249  boost::shared_ptr<CommonData> common_data)
251  commonDataPtr(common_data) {
252  std::fill(&DomainEleOp::doEntities[MBEDGE],
253  &DomainEleOp::doEntities[MBMAXTYPE], false);
254  }
255 
256  MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
258 
263 
264  const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
265  constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
266  commonDataPtr->matLogCdC.resize(size_symm * size_symm, nb_gauss_pts, false);
267  auto t_logC_dC = getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->matLogCdC);
268  auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
269  auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
270 
271  for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
272  // rare case when two eigen values are equal
273  auto nb_uniq = get_uniq_nb<DIM>(&t_eig_val(0));
274  t_logC_dC(i, j, k, l) =
275  2 * EigenMatrix::getDiffMat(t_eig_val, t_eig_vec, f, d_f,
276  nb_uniq)(i, j, k, l);
277 
278  ++t_logC_dC;
279  ++t_eig_val;
280  ++t_eig_vec;
281  }
282 
284  }
285 
286 private:
287  boost::shared_ptr<CommonData> commonDataPtr;
288 };
289 
290 template <int DIM, typename DomainEleOp, int S>
292  : public DomainEleOp {
293 
295  boost::shared_ptr<CommonData> common_data)
297  commonDataPtr(common_data) {
298  std::fill(&DomainEleOp::doEntities[MBEDGE],
299  &DomainEleOp::doEntities[MBMAXTYPE], false);
300  }
301 
302  MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
304 
309 
310  // const size_t nb_gauss_pts = matGradPtr->size2();
311  const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
312  auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*commonDataPtr->matDPtr);
313  auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
314  constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
315  commonDataPtr->matHenckyStress.resize(size_symm, nb_gauss_pts, false);
316  auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
317 
318  for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
319  t_T(i, j) = t_D(i, j, k, l) * t_logC(k, l);
320  ++t_logC;
321  ++t_T;
322  ++t_D;
323  }
324 
326  }
327 
328 private:
329  boost::shared_ptr<CommonData> commonDataPtr;
330 };
331 
332 template <int DIM, typename DomainEleOp, int S>
334  : public DomainEleOp {
335 
337  boost::shared_ptr<CommonData> common_data,
338  boost::shared_ptr<MatrixDouble> mat_D_ptr,
339  const double scale = 1)
340  : DomainEleOp(field_name, DomainEleOp::OPROW), commonDataPtr(common_data),
341  scaleStress(scale), matDPtr(mat_D_ptr) {
342  std::fill(&DomainEleOp::doEntities[MBEDGE],
343  &DomainEleOp::doEntities[MBMAXTYPE], false);
344 
345  matLogCPlastic = commonDataPtr->matLogCPlastic;
346  }
347 
348  MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
350 
355 
356  // const size_t nb_gauss_pts = matGradPtr->size2();
357  const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
358  auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*matDPtr);
359  auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
360  auto t_logCPlastic = getFTensor2SymmetricFromMat<DIM>(*matLogCPlastic);
361  constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
362  commonDataPtr->matHenckyStress.resize(size_symm, nb_gauss_pts, false);
363  auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
364 
365  for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
366  t_T(i, j) = t_D(i, j, k, l) * (t_logC(k, l) - t_logCPlastic(k, l));
367  t_T(i, j) /= scaleStress;
368  ++t_logC;
369  ++t_T;
370  ++t_D;
371  ++t_logCPlastic;
372  }
373 
375  }
376 
377 private:
378  boost::shared_ptr<CommonData> commonDataPtr;
379  boost::shared_ptr<MatrixDouble> matDPtr;
380  boost::shared_ptr<MatrixDouble> matLogCPlastic;
381  const double scaleStress;
382 };
383 
384 template <int DIM, typename DomainEleOp, int S>
386  : public DomainEleOp {
387 
389  boost::shared_ptr<CommonData> common_data)
391  commonDataPtr(common_data) {
392  std::fill(&DomainEleOp::doEntities[MBEDGE],
393  &DomainEleOp::doEntities[MBMAXTYPE], false);
394  }
395 
396  MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
398 
403 
404  constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
405 
406  // const size_t nb_gauss_pts = matGradPtr->size2();
407  const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
408  #ifdef HENCKY_SMALL_STRAIN
409  auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*commonDataPtr->matDPtr);
410  #endif
411  auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
412  auto t_logC_dC = getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->matLogCdC);
413  constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
414  commonDataPtr->matFirstPiolaStress.resize(DIM * DIM, nb_gauss_pts, false);
415  commonDataPtr->matSecondPiolaStress.resize(size_symm, nb_gauss_pts, false);
416  auto t_P = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matFirstPiolaStress);
417  auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
418  auto t_S =
419  getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matSecondPiolaStress);
420  auto t_grad = getFTensor2FromMat<DIM, DIM>(*(commonDataPtr->matGradPtr));
421 
422  for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
423 
424 #ifdef HENCKY_SMALL_STRAIN
425  t_P(i, j) = t_D(i, j, k, l) * t_grad(k, l);
426 #else
428  t_F(i, j) = t_grad(i, j) + t_kd(i, j);
429  t_S(k, l) = t_T(i, j) * t_logC_dC(i, j, k, l);
430  t_P(i, l) = t_F(i, k) * t_S(k, l);
431 #endif
432 
433  ++t_grad;
434  ++t_logC;
435  ++t_logC_dC;
436  ++t_P;
437  ++t_T;
438  ++t_S;
439 #ifdef HENCKY_SMALL_STRAIN
440  ++t_D;
441 #endif
442  }
443 
445  }
446 
447 private:
448  boost::shared_ptr<CommonData> commonDataPtr;
449 };
450 
451 template <int DIM, typename DomainEleOp, int S>
453  OpHenckyTangentImpl(const std::string field_name,
454  boost::shared_ptr<CommonData> common_data,
455  boost::shared_ptr<MatrixDouble> mat_D_ptr = nullptr)
457  commonDataPtr(common_data) {
458  std::fill(&DomainEleOp::doEntities[MBEDGE],
459  &DomainEleOp::doEntities[MBMAXTYPE], false);
460  if (mat_D_ptr)
461  matDPtr = mat_D_ptr;
462  else
463  matDPtr = commonDataPtr->matDPtr;
464  }
465 
466  MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
468 
477 
478  constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
479  // const size_t nb_gauss_pts = matGradPtr->size2();
480  const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
481  commonDataPtr->matTangent.resize(DIM * DIM * DIM * DIM, nb_gauss_pts);
482  auto dP_dF =
483  getFTensor4FromMat<DIM, DIM, DIM, DIM, 1>(commonDataPtr->matTangent);
484 
485  auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*matDPtr);
486  auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
487  auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
488  auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
489  auto t_S =
490  getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matSecondPiolaStress);
491  auto t_grad = getFTensor2FromMat<DIM, DIM>(*(commonDataPtr->matGradPtr));
492  auto t_logC_dC = getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->matLogCdC);
493 
494  for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
495 
496 #ifdef HENCKY_SMALL_STRAIN
497  dP_dF(i, j, k, l) = t_D(i, j, k, l);
498 #else
499 
501  t_F(i, j) = t_grad(i, j) + t_kd(i, j);
502 
503  // rare case when two eigen values are equal
504  auto nb_uniq = get_uniq_nb<DIM>(&t_eig_val(0));
506  dC_dF(i, j, k, l) = (t_kd(i, l) * t_F(k, j)) + (t_kd(j, l) * t_F(k, i));
507 
508  auto TL = EigenMatrix::getDiffDiffMat(t_eig_val, t_eig_vec, f, d_f, dd_f,
509  t_T, nb_uniq);
510  TL(i, j, k, l) *= 4;
511  FTensor::Ddg<double, DIM, DIM> P_D_P_plus_TL;
512  P_D_P_plus_TL(i, j, k, l) =
513  TL(i, j, k, l) +
514  (t_logC_dC(i, j, o, p) * t_D(o, p, m, n)) * t_logC_dC(m, n, k, l);
515  P_D_P_plus_TL(i, j, k, l) *= 0.5;
516  dP_dF(i, j, m, n) = t_kd(i, m) * (t_kd(k, n) * t_S(k, j));
517  dP_dF(i, j, m, n) +=
518  t_F(i, k) * (P_D_P_plus_TL(k, j, o, p) * dC_dF(o, p, m, n));
519 
520 #endif
521 
522  ++dP_dF;
523 
524  ++t_grad;
525  ++t_eig_val;
526  ++t_eig_vec;
527  ++t_logC_dC;
528  ++t_S;
529  ++t_T;
530  ++t_D;
531  }
532 
534  }
535 
536 private:
537  boost::shared_ptr<CommonData> commonDataPtr;
538  boost::shared_ptr<MatrixDouble> matDPtr;
539 };
540 
541 template <typename DomainEleOp> struct HenckyIntegrators {
542  template <int DIM, IntegrationType I>
544 
545  template <int DIM, IntegrationType I>
547 
548  template <int DIM, IntegrationType I>
550 
551  template <int DIM, IntegrationType I, int S>
554 
555  template <int DIM, IntegrationType I, int S = 0>
558 
559  template <int DIM, IntegrationType I, int S>
560  using OpCalculatePiolaStress =
562 
563  template <int DIM, IntegrationType I, int S>
565 };
566 
567 template <int DIM>
570  boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
571  std::string block_name,
572  boost::shared_ptr<MatrixDouble> mat_D_Ptr, Sev sev,
573  double scale = 1) {
575 
576  struct OpMatBlocks : public DomainEleOp {
577  OpMatBlocks(boost::shared_ptr<MatrixDouble> m, double bulk_modulus_K,
578  double shear_modulus_G, MoFEM::Interface &m_field, Sev sev,
579  std::vector<const CubitMeshSets *> meshset_vec_ptr,
580  double scale)
581  : DomainEleOp(NOSPACE, DomainEleOp::OPSPACE), matDPtr(m),
582  bulkModulusKDefault(bulk_modulus_K),
583  shearModulusGDefault(shear_modulus_G), scaleYoungModulus(scale) {
584  CHK_THROW_MESSAGE(extractBlockData(m_field, meshset_vec_ptr, sev),
585  "Can not get data from block");
586  }
587 
588  MoFEMErrorCode doWork(int side, EntityType type,
591 
592  for (auto &b : blockData) {
593 
594  if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
595  CHKERR getMatDPtr(matDPtr, b.bulkModulusK * scaleYoungModulus,
596  b.shearModulusG * scaleYoungModulus);
598  }
599  }
600 
601  CHKERR getMatDPtr(matDPtr, bulkModulusKDefault * scaleYoungModulus,
602  shearModulusGDefault * scaleYoungModulus);
604  }
605 
606  private:
607  boost::shared_ptr<MatrixDouble> matDPtr;
608  const double scaleYoungModulus;
609 
610  struct BlockData {
611  double bulkModulusK;
612  double shearModulusG;
613  Range blockEnts;
614  };
615 
616  double bulkModulusKDefault;
617  double shearModulusGDefault;
618  std::vector<BlockData> blockData;
619 
621  extractBlockData(MoFEM::Interface &m_field,
622  std::vector<const CubitMeshSets *> meshset_vec_ptr,
623  Sev sev) {
625 
626  for (auto m : meshset_vec_ptr) {
627  MOFEM_TAG_AND_LOG("WORLD", sev, "MatBlock") << *m;
628  std::vector<double> block_data;
629  CHKERR m->getAttributes(block_data);
630  if (block_data.size() != 2) {
631  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
632  "Expected that block has two attribute");
633  }
634  auto get_block_ents = [&]() {
635  Range ents;
636  CHKERR
637  m_field.get_moab().get_entities_by_handle(m->meshset, ents, true);
638  return ents;
639  };
640 
641  double young_modulus = block_data[0];
642  double poisson_ratio = block_data[1];
643  double bulk_modulus_K = young_modulus / (3 * (1 - 2 * poisson_ratio));
644  double shear_modulus_G = young_modulus / (2 * (1 + poisson_ratio));
645 
646  MOFEM_TAG_AND_LOG("WORLD", sev, "MatBlock")
647  << "E = " << young_modulus << " nu = " << poisson_ratio;
648 
649  blockData.push_back(
650  {bulk_modulus_K, shear_modulus_G, get_block_ents()});
651  }
652  MOFEM_LOG_CHANNEL("WORLD");
654  }
655 
656  MoFEMErrorCode getMatDPtr(boost::shared_ptr<MatrixDouble> mat_D_ptr,
657  double bulk_modulus_K, double shear_modulus_G) {
659  //! [Calculate elasticity tensor]
660  auto set_material_stiffness = [&]() {
666  double A = (DIM == 2)
667  ? 2 * shear_modulus_G /
668  (bulk_modulus_K + (4. / 3.) * shear_modulus_G)
669  : 1;
670  auto t_D = getFTensor4DdgFromMat<DIM, DIM, 0>(*mat_D_ptr);
671  t_D(i, j, k, l) =
672  2 * shear_modulus_G * ((t_kd(i, k) ^ t_kd(j, l)) / 4.) +
673  A * (bulk_modulus_K - (2. / 3.) * shear_modulus_G) * t_kd(i, j) *
674  t_kd(k, l);
675  };
676  //! [Calculate elasticity tensor]
677  constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
678  mat_D_ptr->resize(size_symm * size_symm, 1);
679  set_material_stiffness();
681  }
682  };
683 
684  double E = 1.0;
685  double nu = 0.3;
686 
687  CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "", "", "none");
688  CHKERR PetscOptionsScalar("-young_modulus", "Young modulus", "", E, &E,
689  PETSC_NULL);
690  CHKERR PetscOptionsScalar("-poisson_ratio", "poisson ratio", "", nu, &nu,
691  PETSC_NULL);
692  ierr = PetscOptionsEnd();
693 
694  double bulk_modulus_K = E / (3 * (1 - 2 * nu));
695  double shear_modulus_G = E / (2 * (1 + nu));
696  pip.push_back(new OpMatBlocks(
697  mat_D_Ptr, bulk_modulus_K, shear_modulus_G, m_field, sev,
698 
699  // Get blockset using regular expression
700  m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
701 
702  (boost::format("%s(.*)") % block_name).str()
703 
704  )),
705  scale
706 
707  ));
708 
710 }
711 
712 template <int DIM, IntegrationType I, typename DomainEleOp>
714  MoFEM::Interface &m_field,
715  boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
716  std::string field_name, std::string block_name, Sev sev, double scale = 1) {
717 
718  auto common_ptr = boost::make_shared<HenckyOps::CommonData>();
719  common_ptr->matDPtr = boost::make_shared<MatrixDouble>();
720  common_ptr->matGradPtr = boost::make_shared<MatrixDouble>();
721 
722  CHK_THROW_MESSAGE(addMatBlockOps<DIM>(m_field, pip, block_name,
723  common_ptr->matDPtr, sev, scale),
724  "addMatBlockOps");
725 
727 
729  field_name, common_ptr->matGradPtr));
730  pip.push_back(new typename H::template OpCalculateEigenVals<DIM, I>(
731  field_name, common_ptr));
732  pip.push_back(
733  new typename H::template OpCalculateLogC<DIM, I>(field_name, common_ptr));
734  pip.push_back(new typename H::template OpCalculateLogC_dC<DIM, I>(
735  field_name, common_ptr));
736  // Assumes constant D matrix per entity
737  pip.push_back(new typename H::template OpCalculateHenckyStress<DIM, I, 0>(
738  field_name, common_ptr));
739  pip.push_back(new typename H::template OpCalculatePiolaStress<DIM, I, 0>(
740  field_name, common_ptr));
741 
742  return common_ptr;
743 }
744 
745 template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
747  MoFEM::Interface &m_field,
748  boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
749  std::string field_name, boost::shared_ptr<HenckyOps::CommonData> common_ptr,
750  Sev sev) {
752 
753  using B = typename FormsIntegrators<DomainEleOp>::template Assembly<
754  A>::template LinearForm<I>;
755  using OpInternalForcePiola =
756  typename B::template OpGradTimesTensor<1, DIM, DIM>;
757  pip.push_back(
758  new OpInternalForcePiola("U", common_ptr->getMatFirstPiolaStress()));
759 
761 }
762 
763 template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
765  MoFEM::Interface &m_field,
766  boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
767  std::string field_name, std::string block_name, Sev sev, double scale = 1) {
769 
770  auto common_ptr = commonDataFactory<DIM, I, DomainEleOp>(
771  m_field, pip, field_name, block_name, sev, scale);
772  CHKERR opFactoryDomainRhs<DIM, A, I, DomainEleOp>(m_field, pip, field_name,
773  common_ptr, sev);
774 
776 }
777 
778 template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
780  MoFEM::Interface &m_field,
781  boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
782  std::string field_name, boost::shared_ptr<HenckyOps::CommonData> common_ptr,
783  Sev sev) {
785 
786  using B = typename FormsIntegrators<DomainEleOp>::template Assembly<
787  A>::template BiLinearForm<I>;
788  using OpKPiola = typename B::template OpGradTensorGrad<1, DIM, DIM, 1>;
789 
791  // Assumes constant D matrix per entity
792  pip.push_back(
793  new typename H::template OpHenckyTangent<DIM, I, 0>(field_name, common_ptr));
794  pip.push_back(
795  new OpKPiola(field_name, field_name, common_ptr->getMatTangent()));
796 
798 }
799 
800 template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
802  MoFEM::Interface &m_field,
803  boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
804  std::string field_name, std::string block_name, Sev sev, double scale = 1) {
806 
807  auto common_ptr = commonDataFactory<DIM, I, DomainEleOp>(
808  m_field, pip, field_name, block_name, sev, scale);
809  CHKERR opFactoryDomainLhs<DIM, A, I, DomainEleOp>(m_field, pip, field_name,
810  common_ptr, sev);
811 
813 }
814 } // namespace HenckyOps
815 
816 #endif // __HENCKY_OPS_HPP__
NOSPACE
@ NOSPACE
Definition: definitions.h:83
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:128
HenckyOps::CommonData::matLogCPlastic
boost::shared_ptr< MatrixDouble > matLogCPlastic
Definition: HenckyOps.hpp:83
HenckyOps::OpCalculatePiolaStressImpl< DIM, GAUSS, DomainEleOp, S >::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Definition: HenckyOps.hpp:396
HenckyOps::CommonData::getMatFirstPiolaStress
auto getMatFirstPiolaStress()
Definition: HenckyOps.hpp:94
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
young_modulus
double young_modulus
Young modulus.
Definition: plastic.cpp:125
HenckyOps::OpCalculateHenckyStressImpl
Definition: HenckyOps.hpp:123
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:609
HenckyOps::OpCalculatePiolaStressImpl< DIM, GAUSS, DomainEleOp, S >::commonDataPtr
boost::shared_ptr< CommonData > commonDataPtr
Definition: HenckyOps.hpp:448
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
E
HenckyOps::isEq::absMax
static double absMax
Definition: HenckyOps.hpp:23
HenckyOps::CommonData::matEigVec
MatrixDouble matEigVec
Definition: HenckyOps.hpp:86
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
HenckyOps::OpCalculateLogC_dCImpl< DIM, GAUSS, DomainEleOp >::OpCalculateLogC_dCImpl
OpCalculateLogC_dCImpl(const std::string field_name, boost::shared_ptr< CommonData > common_data)
Definition: HenckyOps.hpp:248
HenckyOps::d_f
auto d_f
Definition: HenckyOps.hpp:16
HenckyOps::CommonData::getMatLogC
auto getMatLogC()
Definition: HenckyOps.hpp:104
HenckyOps::opFactoryDomainLhs
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)
Definition: HenckyOps.hpp:779
FTensor::Kronecker_Delta
Kronecker Delta class.
Definition: Kronecker_Delta.hpp:15
OpKPiola
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradTensorGrad< 1, SPACE_DIM, SPACE_DIM, 1 > OpKPiola
[Only used for dynamics]
Definition: seepage.cpp:64
HenckyOps::OpCalculateLogCImpl< DIM, GAUSS, DomainEleOp >::OpCalculateLogCImpl
OpCalculateLogCImpl(const std::string field_name, boost::shared_ptr< CommonData > common_data)
Definition: HenckyOps.hpp:207
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
bulk_modulus_K
double bulk_modulus_K
Definition: dynamic_first_order_con_law.cpp:95
HenckyOps::OpCalculateEigenValsImpl< DIM, GAUSS, DomainEleOp >::OpCalculateEigenValsImpl
OpCalculateEigenValsImpl(const std::string field_name, boost::shared_ptr< CommonData > common_data)
Definition: HenckyOps.hpp:137
HenckyOps::CommonData
Definition: HenckyOps.hpp:80
HenckyOps::CommonData::matEigVal
MatrixDouble matEigVal
Definition: HenckyOps.hpp:85
HenckyOps
Definition: HenckyOps.hpp:11
FTensor::Tensor2_symmetric
Definition: Tensor2_symmetric_value.hpp:13
scale
double scale
Definition: plastic.cpp:123
HenckyOps::CommonData::matLogC
MatrixDouble matLogC
Definition: HenckyOps.hpp:87
HenckyOps::OpCalculatePiolaStressImpl< DIM, GAUSS, DomainEleOp, S >
Definition: HenckyOps.hpp:385
FTensor::Tensor2
Definition: Tensor2_value.hpp:16
HenckyOps::OpCalculateHenckyPlasticStressImpl< DIM, GAUSS, DomainEleOp, S >::scaleStress
const double scaleStress
Definition: HenckyOps.hpp:381
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
HenckyOps::OpHenckyTangentImpl< DIM, GAUSS, DomainEleOp, S >::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Definition: HenckyOps.hpp:466
HenckyOps::OpCalculateHenckyPlasticStressImpl
Definition: HenckyOps.hpp:126
HenckyOps::OpCalculatePiolaStressImpl
Definition: HenckyOps.hpp:129
HenckyOps::isEq::check
static auto check(const double &a, const double &b)
Definition: HenckyOps.hpp:20
HenckyOps::is_eq
auto is_eq(const double &a, const double &b)
Definition: HenckyOps.hpp:28
HenckyOps::OpCalculateLogC_dCImpl< DIM, GAUSS, DomainEleOp >::commonDataPtr
boost::shared_ptr< CommonData > commonDataPtr
Definition: HenckyOps.hpp:287
EigenMatrix::getMat
auto getMat(A &&t_val, B &&t_vec, Fun< double > f)
Get the Mat object.
Definition: MatrixFunction.hpp:114
EigenMatrix::getDiffMat
auto getDiffMat(A &&t_val, B &&t_vec, Fun< double > f, Fun< double > d_f, const int nb)
Get the Diff Mat object.
Definition: MatrixFunction.hpp:166
MoFEM::Sev
MoFEM::LogManager::SeverityLevel Sev
Definition: CoreTemplates.hpp:17
HenckyOps::OpCalculateLogC_dCImpl
Definition: HenckyOps.hpp:120
HenckyOps::OpHenckyTangentImpl
Definition: HenckyOps.hpp:132
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
DomainEleOp
DomainEle::UserDataOperator DomainEleOp
Finire element operator type.
Definition: child_and_parent.cpp:36
OpInternalForcePiola
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesTensor< 1, SPACE_DIM, SPACE_DIM > OpInternalForcePiola
Definition: seepage.cpp:66
BlockData
Definition: ElasticityMixedFormulation.hpp:12
a
constexpr double a
Definition: approx_sphere.cpp:30
HenckyOps::OpCalculatePiolaStressImpl< DIM, GAUSS, DomainEleOp, S >::OpCalculatePiolaStressImpl
OpCalculatePiolaStressImpl(const std::string field_name, boost::shared_ptr< CommonData > common_data)
Definition: HenckyOps.hpp:388
HenckyOps::OpHenckyTangentImpl< DIM, GAUSS, DomainEleOp, S >
Definition: HenckyOps.hpp:452
convert.type
type
Definition: convert.py:64
HenckyOps::OpCalculateLogCImpl< DIM, GAUSS, DomainEleOp >::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Definition: HenckyOps.hpp:215
HenckyOps::OpCalculateLogCImpl
Definition: HenckyOps.hpp:117
HenckyOps::eps
constexpr double eps
Definition: HenckyOps.hpp:13
poisson_ratio
double poisson_ratio
Poisson ratio.
Definition: plastic.cpp:126
HenckyOps::CommonData::getMatHenckyStress
auto getMatHenckyStress()
Definition: HenckyOps.hpp:99
MoFEM::GAUSS
@ GAUSS
Definition: FormsIntegrators.hpp:136
HenckyOps::CommonData::matFirstPiolaStress
MatrixDouble matFirstPiolaStress
Definition: HenckyOps.hpp:89
HenckyOps::OpCalculateHenckyPlasticStressImpl< DIM, GAUSS, DomainEleOp, S >::matDPtr
boost::shared_ptr< MatrixDouble > matDPtr
Definition: HenckyOps.hpp:379
HenckyOps::CommonData::matSecondPiolaStress
MatrixDouble matSecondPiolaStress
Definition: HenckyOps.hpp:90
HenckyOps::dd_f
auto dd_f
Definition: HenckyOps.hpp:17
FTensor::Tensor4
Definition: Tensor4_value.hpp:18
size_symm
constexpr auto size_symm
Definition: plastic.cpp:42
HenckyOps::OpHenckyTangentImpl< DIM, GAUSS, DomainEleOp, S >::matDPtr
boost::shared_ptr< MatrixDouble > matDPtr
Definition: HenckyOps.hpp:538
HenckyOps::OpCalculateEigenValsImpl< DIM, GAUSS, DomainEleOp >::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Definition: HenckyOps.hpp:145
HenckyOps::OpCalculateEigenValsImpl
Definition: HenckyOps.hpp:114
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
t_kd
constexpr auto t_kd
Definition: free_surface.cpp:139
BiLinearForm
HenckyOps::isEq
Definition: HenckyOps.hpp:19
OpGradTimesTensor
FormsIntegrators< DomainEleOp >::Assembly< A >::LinearForm< I >::OpGradTimesTensor< 1, FIELD_DIM, SPACE_DIM > OpGradTimesTensor
Definition: operators_tests.cpp:48
HenckyOps::OpCalculateHenckyStressImpl< DIM, GAUSS, DomainEleOp, S >::commonDataPtr
boost::shared_ptr< CommonData > commonDataPtr
Definition: HenckyOps.hpp:329
EntData
EntitiesFieldData::EntData EntData
Definition: child_and_parent.cpp:37
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
FTensor::Index< 'i', DIM >
convert.n
n
Definition: convert.py:82
HenckyOps::CommonData::matLogCdC
MatrixDouble matLogCdC
Definition: HenckyOps.hpp:88
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
Range
DomainEleOp
HenckyOps::sort_eigen_vals
auto sort_eigen_vals(FTensor::Tensor1< double, DIM > &eig, FTensor::Tensor2< double, DIM, DIM > &eigen_vec)
Definition: HenckyOps.hpp:41
MOFEM_TAG_AND_LOG
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
Definition: LogManager.hpp:362
HenckyOps::CommonData::matGradPtr
boost::shared_ptr< MatrixDouble > matGradPtr
Definition: HenckyOps.hpp:81
HenckyOps::CommonData::matDPtr
boost::shared_ptr< MatrixDouble > matDPtr
Definition: HenckyOps.hpp:82
shear_modulus_G
double shear_modulus_G
Definition: dynamic_first_order_con_law.cpp:96
HenckyOps::CommonData::getMatTangent
auto getMatTangent()
Definition: HenckyOps.hpp:108
HenckyOps::OpCalculateLogCImpl< DIM, GAUSS, DomainEleOp >::commonDataPtr
boost::shared_ptr< CommonData > commonDataPtr
Definition: HenckyOps.hpp:242
HenckyOps::OpCalculateHenckyPlasticStressImpl< DIM, GAUSS, DomainEleOp, S >::OpCalculateHenckyPlasticStressImpl
OpCalculateHenckyPlasticStressImpl(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:336
HenckyOps::opFactoryDomainRhs
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)
Definition: HenckyOps.hpp:746
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
FTensor::Ddg
Definition: Ddg_value.hpp:7
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
HenckyOps::OpCalculateEigenValsImpl< DIM, GAUSS, DomainEleOp >::commonDataPtr
boost::shared_ptr< CommonData > commonDataPtr
Definition: HenckyOps.hpp:201
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
HenckyOps::OpHenckyTangentImpl< DIM, GAUSS, DomainEleOp, S >::commonDataPtr
boost::shared_ptr< CommonData > commonDataPtr
Definition: HenckyOps.hpp:537
HenckyOps::OpCalculateHenckyStressImpl< DIM, GAUSS, DomainEleOp, S >::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Definition: HenckyOps.hpp:302
HenckyOps::OpCalculateHenckyPlasticStressImpl< DIM, GAUSS, DomainEleOp, S >::matLogCPlastic
boost::shared_ptr< MatrixDouble > matLogCPlastic
Definition: HenckyOps.hpp:380
HenckyOps::OpHenckyTangentImpl< DIM, GAUSS, DomainEleOp, S >::OpHenckyTangentImpl
OpHenckyTangentImpl(const std::string field_name, boost::shared_ptr< CommonData > common_data, boost::shared_ptr< MatrixDouble > mat_D_ptr=nullptr)
Definition: HenckyOps.hpp:453
HenckyOps::OpCalculateHenckyPlasticStressImpl< DIM, GAUSS, DomainEleOp, S >::commonDataPtr
boost::shared_ptr< CommonData > commonDataPtr
Definition: HenckyOps.hpp:378
HenckyOps::OpCalculateLogC_dCImpl< DIM, GAUSS, DomainEleOp >::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Definition: HenckyOps.hpp:256
HenckyOps::commonDataFactory
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)
Definition: HenckyOps.hpp:713
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
FTensor::Kronecker_Delta_symmetric
Kronecker Delta class symmetric.
Definition: Kronecker_Delta.hpp:49
HenckyOps::OpCalculateHenckyPlasticStressImpl< DIM, GAUSS, DomainEleOp, S >::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Definition: HenckyOps.hpp:348
HenckyOps::OpCalculateHenckyStressImpl< DIM, GAUSS, DomainEleOp, S >::OpCalculateHenckyStressImpl
OpCalculateHenckyStressImpl(const std::string field_name, boost::shared_ptr< CommonData > common_data)
Definition: HenckyOps.hpp:294
HenckyOps::addMatBlockOps
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)
Definition: HenckyOps.hpp:569
HenckyOps::CommonData::matTangent
MatrixDouble matTangent
Definition: HenckyOps.hpp:92
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
HenckyOps::HenckyIntegrators
Definition: HenckyOps.hpp:541
HenckyOps::CommonData::matHenckyStress
MatrixDouble matHenckyStress
Definition: HenckyOps.hpp:91
HenckyOps::get_uniq_nb
auto get_uniq_nb(double *ptr)
Definition: HenckyOps.hpp:32
H
double H
Hardening.
Definition: plastic.cpp:128
OpCalculateVectorFieldGradient
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
l
FTensor::Index< 'l', 3 > l
Definition: matrix_function.cpp:21
MoFEM::computeEigenValuesSymmetric
MoFEMErrorCode computeEigenValuesSymmetric(const MatrixDouble &mat, VectorDouble &eig, MatrixDouble &eigen_vec)
compute eigenvalues of a symmetric matrix using lapack dsyev
Definition: Templates.hpp:1452
F
@ F
Definition: free_surface.cpp:396
EigenMatrix::getDiffDiffMat
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.
Definition: MatrixFunction.hpp:279