v0.14.0
Hencky.cpp
Go to the documentation of this file.
1 /**
2  * @file Hencky.cpp
3  * @brief Implementation of Hencky material
4  * @date 2024-08-31
5  *
6  * @copyright Copyright (c) 2024
7  *
8  */
9 
10 namespace EshelbianPlasticity {
11 
12 struct HMHHencky : public PhysicalEquations {
13 
14  HMHHencky(MoFEM::Interface &m_field, const double E, const double nu)
15  : PhysicalEquations(0, 0), mField(m_field), E(E), nu(nu) {}
16 
17  MoFEMErrorCode recordTape(const int tag, DTensor2Ptr *t_h) { return 0; }
18 
19  struct OpHenckyJacobian : public OpJacobian {
20  OpHenckyJacobian(boost::shared_ptr<DataAtIntegrationPts> data_ptr,
21  boost::shared_ptr<HMHHencky> hencky_ptr)
22  : OpJacobian(), dataAtGaussPts(data_ptr), henckyPtr(hencky_ptr) {
24  "getOptions failed");
25  CHK_THROW_MESSAGE(henckyPtr->extractBlockData(Sev::inform),
26  "Can not get data from block");
27  }
28 
29  MoFEMErrorCode doWork(int side, EntityType type,
32 
33  for (auto &b : henckyPtr->blockData) {
34 
35  if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
36  dataAtGaussPts->mu = b.bulkModulusK;
37  dataAtGaussPts->lambda = b.bulkModulusK;
38  CHKERR henckyPtr->getMatDPtr(dataAtGaussPts->getMatDPtr(),
39  dataAtGaussPts->getMatAxiatorDPtr(),
40  dataAtGaussPts->getMatDeviatorDPtr(),
41  b.bulkModulusK, b.shearModulusG);
43  }
44  }
45 
46  const auto E = henckyPtr->E;
47  const auto nu = henckyPtr->nu;
48 
49  double bulk_modulus_K = E / (3 * (1 - 2 * nu));
50  double shear_modulus_G = E / (2 * (1 + nu));
51 
54 
55  CHKERR henckyPtr->getMatDPtr(dataAtGaussPts->getMatDPtr(),
56  dataAtGaussPts->getMatAxiatorDPtr(),
57  dataAtGaussPts->getMatDeviatorDPtr(),
59 
61  }
62 
63  MoFEMErrorCode evaluateRhs(EntData &data) { return 0; }
64  MoFEMErrorCode evaluateLhs(EntData &data) { return 0; }
65 
66  private:
67  boost::shared_ptr<DataAtIntegrationPts> dataAtGaussPts;
68  boost::shared_ptr<HMHHencky> henckyPtr;
69  };
70 
71  virtual OpJacobian *
72  returnOpJacobian(const int tag, const bool eval_rhs, const bool eval_lhs,
73  boost::shared_ptr<DataAtIntegrationPts> data_ptr,
74  boost::shared_ptr<PhysicalEquations> physics_ptr) {
75  return (new OpHenckyJacobian(
76  data_ptr, boost::dynamic_pointer_cast<HMHHencky>(physics_ptr)));
77  }
78 
80 
81  OpSpatialPhysical(const std::string &field_name,
82  boost::shared_ptr<DataAtIntegrationPts> data_ptr,
83  const double alpha_u);
84 
86 
88 
90 
91  private:
92  const double alphaU;
93  PetscBool polyConvex = PETSC_FALSE;
94  };
95 
96  virtual VolUserDataOperator *
97  returnOpSpatialPhysical(const std::string &field_name,
98  boost::shared_ptr<DataAtIntegrationPts> data_ptr,
99  const double alpha_u) {
100  return new OpSpatialPhysical(field_name, data_ptr, alpha_u);
101  }
102 
104  const double alphaU;
105  OpSpatialPhysical_du_du(std::string row_field, std::string col_field,
106  boost::shared_ptr<DataAtIntegrationPts> data_ptr,
107  const double alpha);
108  MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
109  MoFEMErrorCode integrateHencky(EntData &row_data, EntData &col_data);
111  EntData &col_data);
112 
113  private:
114  PetscBool polyConvex = PETSC_FALSE;
115 
117  };
118 
120  std::string row_field, std::string col_field,
121  boost::shared_ptr<DataAtIntegrationPts> data_ptr, const double alpha) {
122  return new OpSpatialPhysical_du_du(row_field, col_field, data_ptr, alpha);
123  }
124 
125  /**
126  * @brief Calculate energy density for Hencky material model
127  *
128  *
129  * \f[
130  *
131  * \Psi(\log{\mathbf{U}}) = \frac{1}{2} U_{IJ} D_{IJKL} U_{KL} = \frac{1}{2}
132  * U_{IJ} T_{IJ}
133  *
134  * \f]
135  * where \f$T_{IJ} = D_{IJKL} U_{KL}\f$ is a a Hencky stress.
136  *
137  */
139 
140  OpCalculateEnergy(boost::shared_ptr<DataAtIntegrationPts> data_ptr,
141  boost::shared_ptr<double> total_energy_ptr);
142  MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
143 
144  private:
145  boost::shared_ptr<DataAtIntegrationPts> dataAtPts;
146  boost::shared_ptr<double> totalEnergyPtr;
147  };
148 
150  returnOpCalculateEnergy(boost::shared_ptr<DataAtIntegrationPts> data_ptr,
151  boost::shared_ptr<double> total_energy_ptr) {
152  return new OpCalculateEnergy(data_ptr, total_energy_ptr);
153  }
154 
157  boost::shared_ptr<DataAtIntegrationPts> data_ptr,
158  boost::shared_ptr<HMHHencky> hencky_ptr);
159  MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
160 
161  private:
162  boost::shared_ptr<DataAtIntegrationPts>
163  dataAtPts; ///< data at integration pts
164  boost::shared_ptr<HMHHencky> henckyPtr;
165  };
166 
168  boost::shared_ptr<DataAtIntegrationPts> data_ptr,
169  boost::shared_ptr<PhysicalEquations> physics_ptr) {
170  return new OpCalculateStretchFromStress(
171  data_ptr, boost::dynamic_pointer_cast<HMHHencky>(physics_ptr));
172  }
173 
174  MoFEMErrorCode getOptions(boost::shared_ptr<DataAtIntegrationPts> data_ptr) {
176  CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "hencky_", "", "none");
177 
178  CHKERR PetscOptionsScalar("-young_modulus", "Young modulus", "", E, &E,
179  PETSC_NULL);
180  CHKERR PetscOptionsScalar("-poisson_ratio", "poisson ratio", "", nu, &nu,
181  PETSC_NULL);
182 
183  ierr = PetscOptionsEnd();
184  CHKERRG(ierr);
185 
187  }
188 
190  return extractBlockData(
191 
192  mField.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
193 
194  (boost::format("%s(.*)") % "MAT_ELASTIC").str()
195 
196  )),
197 
198  sev);
199  }
200 
202  extractBlockData(std::vector<const CubitMeshSets *> meshset_vec_ptr,
203  Sev sev) {
205 
206  for (auto m : meshset_vec_ptr) {
207  MOFEM_TAG_AND_LOG("WORLD", sev, "MatBlock") << *m;
208  std::vector<double> block_data;
209  CHKERR m->getAttributes(block_data);
210  if (block_data.size() != 2) {
211  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
212  "Expected that block has two attribute");
213  }
214  auto get_block_ents = [&]() {
215  Range ents;
216  CHKERR mField.get_moab().get_entities_by_handle(m->meshset, ents, true);
217  return ents;
218  };
219 
220  double young_modulus = block_data[0];
221  double poisson_ratio = block_data[1];
222  double bulk_modulus_K = young_modulus / (3 * (1 - 2 * poisson_ratio));
223  double shear_modulus_G = young_modulus / (2 * (1 + poisson_ratio));
224 
225  MOFEM_TAG_AND_LOG("WORLD", sev, "MatBlock")
226  << "E = " << young_modulus << " nu = " << poisson_ratio;
227 
228  blockData.push_back({bulk_modulus_K, shear_modulus_G, get_block_ents()});
229  }
230  MOFEM_LOG_CHANNEL("WORLD");
232  }
233 
234  MoFEMErrorCode getMatDPtr(boost::shared_ptr<MatrixDouble> mat_D_ptr,
235  boost::shared_ptr<MatrixDouble> mat_axiator_D_ptr,
236  boost::shared_ptr<MatrixDouble> mat_deviator_D_ptr,
237  double bulk_modulus_K, double shear_modulus_G) {
239  //! [Calculate elasticity tensor]
240  auto set_material_stiffness = [&]() {
246  double A = (SPACE_DIM == 2)
247  ? 2 * shear_modulus_G /
248  (bulk_modulus_K + (4. / 3.) * shear_modulus_G)
249  : 1;
250 
251  auto t_D = getFTensor4DdgFromPtr<SPACE_DIM, SPACE_DIM, 0>(
252  &*(mat_D_ptr->data().begin()));
253  auto t_axiator_D = getFTensor4DdgFromPtr<SPACE_DIM, SPACE_DIM, 0>(
254  &*mat_axiator_D_ptr->data().begin());
255  auto t_deviator_D = getFTensor4DdgFromPtr<SPACE_DIM, SPACE_DIM, 0>(
256  &*mat_deviator_D_ptr->data().begin());
257  t_axiator_D(i, j, k, l) = A *
258  (bulk_modulus_K - (2. / 3.) * shear_modulus_G) *
259  t_kd(i, j) * t_kd(k, l);
260  t_deviator_D(i, j, k, l) =
261  2 * shear_modulus_G * ((t_kd(i, k) ^ t_kd(j, l)) / 4.);
262  t_D(i, j, k, l) = t_axiator_D(i, j, k, l) + t_deviator_D(i, j, k, l);
263  };
264  //! [Calculate elasticity tensor]
265  constexpr auto size_symm = (SPACE_DIM * (SPACE_DIM + 1)) / 2;
266  mat_D_ptr->resize(size_symm, size_symm, false);
267  mat_axiator_D_ptr->resize(size_symm, size_symm, false);
268  mat_deviator_D_ptr->resize(size_symm, size_symm, false);
269  set_material_stiffness();
271  }
272 
273 private:
275 
276  struct BlockData {
277  double bulkModulusK;
280  };
281  std::vector<BlockData> blockData;
282 
283  double E;
284  double nu;
285 };
286 
288  const std::string &field_name,
289  boost::shared_ptr<DataAtIntegrationPts> data_ptr, const double alpha_u)
290  : OpAssembleVolume(field_name, data_ptr, OPROW), alphaU(alpha_u) {
291 
292  CHK_MOAB_THROW(PetscOptionsGetBool(PETSC_NULL, "", "-poly_convex",
293  &polyConvex, PETSC_NULL),
294  "get polyconvex option failed");
295 }
296 
299  if (polyConvex) {
300  CHKERR integratePolyconvexHencky(data);
301  } else {
302  CHKERR integrateHencky(data);
303  }
305 }
306 
309 
311  auto t_L = symm_L_tensor();
312 
313  int nb_dofs = data.getIndices().size();
314  int nb_integration_pts = data.getN().size1();
315  auto v = getVolume();
316  auto t_w = getFTensor0IntegrationWeight();
317  auto t_approx_P_adjont_log_du =
318  getFTensor1FromMat<size_symm>(dataAtPts->adjointPdUAtPts);
319  auto t_log_stretch_h1 =
320  getFTensor2SymmetricFromMat<3>(dataAtPts->logStretchTotalTensorAtPts);
321  auto t_dot_log_u =
322  getFTensor2SymmetricFromMat<3>(dataAtPts->logStretchDotTensorAtPts);
323 
324  auto t_D = getFTensor4DdgFromPtr<3, 3, 0>(&*dataAtPts->matD.data().begin());
325 
330  auto get_ftensor2 = [](auto &v) {
332  &v[0], &v[1], &v[2], &v[3], &v[4], &v[5]);
333  };
334 
335  int nb_base_functions = data.getN().size2();
336  auto t_row_base_fun = data.getFTensor0N();
337  for (int gg = 0; gg != nb_integration_pts; ++gg) {
338  double a = v * t_w;
339  auto t_nf = get_ftensor2(nF);
340 
342  t_T(i, j) =
343  t_D(i, j, k, l) * (t_log_stretch_h1(k, l) + alphaU * t_dot_log_u(k, l));
345  t_residual(L) =
346  a * (t_approx_P_adjont_log_du(L) - t_L(i, j, L) * t_T(i, j));
347 
348  int bb = 0;
349  for (; bb != nb_dofs / 6; ++bb) {
350  t_nf(L) += t_row_base_fun * t_residual(L);
351  ++t_nf;
352  ++t_row_base_fun;
353  }
354  for (; bb != nb_base_functions; ++bb)
355  ++t_row_base_fun;
356 
357  ++t_w;
358  ++t_approx_P_adjont_log_du;
359  ++t_dot_log_u;
360  ++t_log_stretch_h1;
361  }
363 }
364 
368 
370  auto t_L = symm_L_tensor();
371 
372  int nb_dofs = data.getIndices().size();
373  int nb_integration_pts = data.getN().size1();
374  auto v = getVolume();
375  auto t_w = getFTensor0IntegrationWeight();
376  auto t_approx_P_adjont_log_du =
377  getFTensor1FromMat<size_symm>(dataAtPts->adjointPdUAtPts);
378  auto t_log_stretch_h1 =
379  getFTensor2SymmetricFromMat<3>(dataAtPts->logStretchTotalTensorAtPts);
380  auto t_dot_log_u =
381  getFTensor2SymmetricFromMat<3>(dataAtPts->logStretchDotTensorAtPts);
382 
383  auto t_D = getFTensor4DdgFromPtr<3, 3, 0>(&*dataAtPts->matD.data().begin());
384 
389  auto get_ftensor2 = [](auto &v) {
391  &v[0], &v[1], &v[2], &v[3], &v[4], &v[5]);
392  };
393 
394  constexpr double nohat_k = 1. / 4;
395  constexpr double hat_k = 1. / 8;
396  double mu = dataAtPts->mu;
397  double lambda = dataAtPts->lambda;
398 
399  constexpr double third = boost::math::constants::third<double>();
401  auto t_diff_deviator = diff_deviator(diff_tensor());
402 
403  int nb_base_functions = data.getN().size2();
404  auto t_row_base_fun = data.getFTensor0N();
405  for (int gg = 0; gg != nb_integration_pts; ++gg) {
406  double a = v * t_w;
407  auto t_nf = get_ftensor2(nF);
408 
409  double log_det = t_log_stretch_h1(i, i);
410  double log_det2 = log_det * log_det;
412  t_dev(i, j) = t_log_stretch_h1(i, j) - t_kd(i, j) * (third * log_det);
413  double dev_norm2 = t_dev(i, j) * t_dev(i, j);
414 
416  auto A = 2 * mu * std::exp(nohat_k * dev_norm2);
417  auto B = lambda * std::exp(hat_k * log_det2) * log_det;
418  t_T(i, j) =
419 
420  A * (t_dev(k, l) * t_diff_deviator(k, l, i, j))
421 
422  +
423 
424  B * t_kd(i, j)
425 
426  +
427 
428  alphaU * t_D(i, j, k, l) * t_dot_log_u(k, l);
429 
431  t_residual(L) =
432  a * (t_approx_P_adjont_log_du(L) - t_L(i, j, L) * t_T(i, j));
433 
434  int bb = 0;
435  for (; bb != nb_dofs / size_symm; ++bb) {
436  t_nf(L) += t_row_base_fun * t_residual(L);
437  ++t_nf;
438  ++t_row_base_fun;
439  }
440  for (; bb != nb_base_functions; ++bb)
441  ++t_row_base_fun;
442 
443  ++t_w;
444  ++t_approx_P_adjont_log_du;
445  ++t_dot_log_u;
446  ++t_log_stretch_h1;
447  }
449 }
450 
452  std::string row_field, std::string col_field,
453  boost::shared_ptr<DataAtIntegrationPts> data_ptr, const double alpha)
454  : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL, false),
455  alphaU(alpha) {
456  sYmm = false;
457 
458  CHK_MOAB_THROW(PetscOptionsGetBool(PETSC_NULL, "", "-poly_convex",
459  &polyConvex, PETSC_NULL),
460  "get polyconvex option failed");
461 }
462 
465  EntData &col_data) {
467 
468  if (polyConvex) {
469  CHKERR integratePolyconvexHencky(row_data, col_data);
470  } else {
471  CHKERR integrateHencky(row_data, col_data);
472  }
474 }
475 
478  EntData &col_data) {
480 
483  auto t_L = symm_L_tensor();
484  auto t_diff = diff_tensor();
485 
486  int nb_integration_pts = row_data.getN().size1();
487  int row_nb_dofs = row_data.getIndices().size();
488  int col_nb_dofs = col_data.getIndices().size();
489 
490  auto get_ftensor2 = [](MatrixDouble &m, const int r, const int c) {
492  size_symm>(
493 
494  &m(r + 0, c + 0), &m(r + 0, c + 1), &m(r + 0, c + 2), &m(r + 0, c + 3),
495  &m(r + 0, c + 4), &m(r + 0, c + 5),
496 
497  &m(r + 1, c + 0), &m(r + 1, c + 1), &m(r + 1, c + 2), &m(r + 1, c + 3),
498  &m(r + 1, c + 4), &m(r + 1, c + 5),
499 
500  &m(r + 2, c + 0), &m(r + 2, c + 1), &m(r + 2, c + 2), &m(r + 2, c + 3),
501  &m(r + 2, c + 4), &m(r + 2, c + 5),
502 
503  &m(r + 3, c + 0), &m(r + 3, c + 1), &m(r + 3, c + 2), &m(r + 3, c + 3),
504  &m(r + 3, c + 4), &m(r + 3, c + 5),
505 
506  &m(r + 4, c + 0), &m(r + 4, c + 1), &m(r + 4, c + 2), &m(r + 4, c + 3),
507  &m(r + 4, c + 4), &m(r + 4, c + 5),
508 
509  &m(r + 5, c + 0), &m(r + 5, c + 1), &m(r + 5, c + 2), &m(r + 5, c + 3),
510  &m(r + 5, c + 4), &m(r + 5, c + 5)
511 
512  );
513  };
520 
521  auto v = getVolume();
522  auto t_w = getFTensor0IntegrationWeight();
523 
524  auto t_approx_P_adjont_dstretch =
525  getFTensor2FromMat<3, 3>(dataAtPts->adjointPdstretchAtPts);
526  auto t_eigen_vals = getFTensor1FromMat<3>(dataAtPts->eigenVals);
527  auto t_eigen_vecs = getFTensor2FromMat<3, 3>(dataAtPts->eigenVecs);
528  auto &nbUniq = dataAtPts->nbUniq;
529 
530  int row_nb_base_functions = row_data.getN().size2();
531  auto t_row_base_fun = row_data.getFTensor0N();
532 
533  auto get_dP = [&]() {
534  dP.resize(size_symm * size_symm, nb_integration_pts, false);
535  auto ts_a = getTSa();
536 
537  auto t_D = getFTensor4DdgFromPtr<3, 3, 0>(&*dataAtPts->matD.data().begin());
539  t_dP_tmp(L, J) = -(1 + alphaU * ts_a) *
540  (t_L(i, j, L) *
541  ((t_D(i, j, m, n) * t_diff(m, n, k, l)) * t_L(k, l, J)));
542 
544  auto t_approx_P_adjont_dstretch =
545  getFTensor2FromMat<3, 3>(dataAtPts->adjointPdstretchAtPts);
546  auto t_eigen_vals = getFTensor1FromMat<3>(dataAtPts->eigenVals);
547  auto t_eigen_vecs = getFTensor2FromMat<3, 3>(dataAtPts->eigenVecs);
548  auto &nbUniq = dataAtPts->nbUniq;
549 
550  auto t_dP = getFTensor2FromMat<size_symm, size_symm>(dP);
551  for (auto gg = 0; gg != nb_integration_pts; ++gg) {
552 
553  // Work of symmetric tensor on undefined tensor is equal to the work
554  // of the symmetric part of it
556  t_sym(i, j) = (t_approx_P_adjont_dstretch(i, j) ||
557  t_approx_P_adjont_dstretch(j, i));
558  t_sym(i, j) /= 2.0;
559  auto t_diff2_uP2 = EigenMatrix::getDiffDiffMat(
560  t_eigen_vals, t_eigen_vecs, EshelbianCore::f, EshelbianCore::d_f,
561  EshelbianCore::dd_f, t_sym, nbUniq[gg]);
562  t_dP(L, J) = t_L(i, j, L) *
563  ((t_diff2_uP2(i, j, k, l) + t_diff2_uP2(k, l, i, j)) *
564  t_L(k, l, J)) /
565  2. +
566  t_dP_tmp(L, J);
567 
568  ++t_dP;
569  ++t_approx_P_adjont_dstretch;
570  ++t_eigen_vals;
571  ++t_eigen_vecs;
572  }
573  } else {
574  auto t_dP = getFTensor2FromMat<size_symm, size_symm>(dP);
575  for (auto gg = 0; gg != nb_integration_pts; ++gg) {
576  t_dP(L, J) = t_dP_tmp(L, J);
577  ++t_dP;
578  }
579  }
580 
581  return getFTensor2FromMat<size_symm, size_symm>(dP);
582  };
583 
584  auto t_dP = get_dP();
585  for (int gg = 0; gg != nb_integration_pts; ++gg) {
586  double a = v * t_w;
587 
588  int rr = 0;
589  for (; rr != row_nb_dofs / 6; ++rr) {
590  auto t_col_base_fun = col_data.getFTensor0N(gg, 0);
591  auto t_m = get_ftensor2(K, 6 * rr, 0);
592  for (int cc = 0; cc != col_nb_dofs / 6; ++cc) {
593  const double b = a * t_row_base_fun * t_col_base_fun;
594  t_m(L, J) += b * t_dP(L, J);
595  ++t_m;
596  ++t_col_base_fun;
597  }
598  ++t_row_base_fun;
599  }
600 
601  for (; rr != row_nb_base_functions; ++rr) {
602  ++t_row_base_fun;
603  }
604 
605  ++t_w;
606  ++t_dP;
607  }
609 }
610 
612  EntData &row_data, EntData &col_data) {
614 
617  auto t_L = symm_L_tensor();
618  auto t_diff = diff_tensor();
619 
620  int nb_integration_pts = row_data.getN().size1();
621  int row_nb_dofs = row_data.getIndices().size();
622  int col_nb_dofs = col_data.getIndices().size();
623 
624  auto get_ftensor2 = [](MatrixDouble &m, const int r, const int c) {
626  size_symm>(
627 
628  &m(r + 0, c + 0), &m(r + 0, c + 1), &m(r + 0, c + 2), &m(r + 0, c + 3),
629  &m(r + 0, c + 4), &m(r + 0, c + 5),
630 
631  &m(r + 1, c + 0), &m(r + 1, c + 1), &m(r + 1, c + 2), &m(r + 1, c + 3),
632  &m(r + 1, c + 4), &m(r + 1, c + 5),
633 
634  &m(r + 2, c + 0), &m(r + 2, c + 1), &m(r + 2, c + 2), &m(r + 2, c + 3),
635  &m(r + 2, c + 4), &m(r + 2, c + 5),
636 
637  &m(r + 3, c + 0), &m(r + 3, c + 1), &m(r + 3, c + 2), &m(r + 3, c + 3),
638  &m(r + 3, c + 4), &m(r + 3, c + 5),
639 
640  &m(r + 4, c + 0), &m(r + 4, c + 1), &m(r + 4, c + 2), &m(r + 4, c + 3),
641  &m(r + 4, c + 4), &m(r + 4, c + 5),
642 
643  &m(r + 5, c + 0), &m(r + 5, c + 1), &m(r + 5, c + 2), &m(r + 5, c + 3),
644  &m(r + 5, c + 4), &m(r + 5, c + 5)
645 
646  );
647  };
654 
655  auto v = getVolume();
656  auto t_w = getFTensor0IntegrationWeight();
657 
658  int row_nb_base_functions = row_data.getN().size2();
659  auto t_row_base_fun = row_data.getFTensor0N();
660 
661  auto get_dP = [&]() {
662  dP.resize(size_symm * size_symm, nb_integration_pts, false);
663  auto ts_a = getTSa();
664 
665  auto t_D = getFTensor4DdgFromPtr<3, 3, 0>(&*dataAtPts->matD.data().begin());
666 
667  constexpr double nohat_k = 1. / 4;
668  constexpr double hat_k = 1. / 8;
669  double mu = dataAtPts->mu;
670  double lambda = dataAtPts->lambda;
671 
672  constexpr double third = boost::math::constants::third<double>();
674  auto t_diff_deviator = diff_deviator(diff_tensor());
675 
676  auto t_approx_P_adjont_dstretch =
677  getFTensor2FromMat<3, 3>(dataAtPts->adjointPdstretchAtPts);
678  auto t_log_stretch_h1 =
679  getFTensor2SymmetricFromMat<3>(dataAtPts->logStretchTotalTensorAtPts);
680  auto t_eigen_vals = getFTensor1FromMat<3>(dataAtPts->eigenVals);
681  auto t_eigen_vecs = getFTensor2FromMat<3, 3>(dataAtPts->eigenVecs);
682  auto &nbUniq = dataAtPts->nbUniq;
683 
684  auto t_dP = getFTensor2FromMat<size_symm, size_symm>(dP);
685  for (auto gg = 0; gg != nb_integration_pts; ++gg) {
686 
687  double log_det = t_log_stretch_h1(i, i);
688  double log_det2 = log_det * log_det;
690  t_dev(i, j) = t_log_stretch_h1(i, j) - t_kd(i, j) * (third * log_det);
691  double dev_norm2 = t_dev(i, j) * t_dev(i, j);
692 
693  auto A = 2 * mu * std::exp(nohat_k * dev_norm2);
694  auto B = lambda * std::exp(hat_k * log_det2) * log_det;
695 
696  FTensor::Tensor2_symmetric<double, 3> t_A_diff, t_B_diff;
697  t_A_diff(i, j) =
698  (A * 2 * nohat_k) * (t_dev(k, l) * t_diff_deviator(k, l, i, j));
699  t_B_diff(i, j) = (B * 2 * hat_k) * log_det * t_kd(i, j) +
700  lambda * std::exp(hat_k * log_det2) * t_kd(i, j);
702  t_dT(i, j, k, l) =
703  t_A_diff(i, j) * (t_dev(m, n) * t_diff_deviator(m, n, k, l))
704 
705  +
706 
707  A * t_diff_deviator(m, n, i, j) * t_diff_deviator(m, n, k, l)
708 
709  +
710 
711  t_B_diff(i, j) * t_kd(k, l);
712 
713  t_dP(L, J) = -t_L(i, j, L) *
714  ((
715 
716  t_dT(i, j, k, l)
717 
718  +
719 
720  (alphaU * ts_a) * (t_D(i, j, m, n) * t_diff(m, n, k, l)
721 
722  )) *
723  t_L(k, l, J));
724 
725  // Work of symmetric tensor on undefined tensor is equal to the work
726  // of the symmetric part of it
729  t_sym(i, j) = (t_approx_P_adjont_dstretch(i, j) ||
730  t_approx_P_adjont_dstretch(j, i));
731  t_sym(i, j) /= 2.0;
732  auto t_diff2_uP2 = EigenMatrix::getDiffDiffMat(
733  t_eigen_vals, t_eigen_vecs, EshelbianCore::f, EshelbianCore::d_f,
734  EshelbianCore::dd_f, t_sym, nbUniq[gg]);
735  t_dP(L, J) += t_L(i, j, L) *
736  ((t_diff2_uP2(i, j, k, l) + t_diff2_uP2(k, l, i, j)) *
737  t_L(k, l, J)) /
738  2.;
739  }
740 
741  ++t_dP;
742  ++t_approx_P_adjont_dstretch;
743  ++t_log_stretch_h1;
744  ++t_eigen_vals;
745  ++t_eigen_vecs;
746  }
747 
748  return getFTensor2FromMat<size_symm, size_symm>(dP);
749  };
750 
751  auto t_dP = get_dP();
752  for (int gg = 0; gg != nb_integration_pts; ++gg) {
753  double a = v * t_w;
754 
755  int rr = 0;
756  for (; rr != row_nb_dofs / 6; ++rr) {
757  auto t_col_base_fun = col_data.getFTensor0N(gg, 0);
758  auto t_m = get_ftensor2(K, 6 * rr, 0);
759  for (int cc = 0; cc != col_nb_dofs / 6; ++cc) {
760  const double b = a * t_row_base_fun * t_col_base_fun;
761  t_m(L, J) += b * t_dP(L, J);
762  ++t_m;
763  ++t_col_base_fun;
764  }
765  ++t_row_base_fun;
766  }
767 
768  for (; rr != row_nb_base_functions; ++rr) {
769  ++t_row_base_fun;
770  }
771 
772  ++t_w;
773  ++t_dP;
774  }
776 }
777 
779  boost::shared_ptr<DataAtIntegrationPts> data_ptr,
780  boost::shared_ptr<double> total_energy_ptr)
781  : VolUserDataOperator(NOSPACE, OPSPACE), dataAtPts(data_ptr),
782  totalEnergyPtr(total_energy_ptr) {}
783 
785  EntData &data) {
787 
792 
793  int nb_integration_pts = getGaussPts().size2();
794  auto t_log_u =
795  getFTensor2SymmetricFromMat<3>(dataAtPts->logStretchTotalTensorAtPts);
796  auto t_D = getFTensor4DdgFromPtr<3, 3, 0>(&*dataAtPts->matD.data().begin());
797 
798  dataAtPts->energyAtPts.resize(nb_integration_pts, false);
799  auto t_energy = getFTensor0FromVec(dataAtPts->energyAtPts);
800 
801  for (auto gg = 0; gg != nb_integration_pts; ++gg) {
802 
803  t_energy = 0.5 * (t_log_u(i, j) * (t_D(i, j, k, l) * t_log_u(k, l)));
804 
805  ++t_log_u;
806  ++t_energy;
807  }
808 
809  if (totalEnergyPtr) {
810  auto t_w = getFTensor0IntegrationWeight();
811  auto t_energy = getFTensor0FromVec(dataAtPts->energyAtPts);
812  double loc_energy = 0;
813  for (auto gg = 0; gg != nb_integration_pts; ++gg) {
814  loc_energy += t_energy * t_w;
815  ++t_w;
816  ++t_energy;
817  }
818  *totalEnergyPtr += getMeasure() * loc_energy;
819  }
820 
822 }
823 
825  boost::shared_ptr<DataAtIntegrationPts> data_ptr,
826  boost::shared_ptr<HMHHencky> hencky_ptr)
827  : VolUserDataOperator(NOSPACE, OPSPACE), dataAtPts(data_ptr),
828  henckyPtr(hencky_ptr) {}
829 
831  EntityType type,
832  EntData &data) {
834 
841 
842  auto nb_integration_pts = dataAtPts->approxPAtPts.size2();
843 #ifdef NDEBUG
844  if (nb_integration_pts != getGaussPts().size2()) {
845  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
846  "inconsistent number of integration points");
847  }
848 #endif // NDEBUG
849 
850  auto get_D = [&]() {
852  for (auto &b : henckyPtr->blockData) {
853 
854  if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
855  dataAtPts->mu = b.bulkModulusK;
856  dataAtPts->lambda = b.bulkModulusK;
857  CHKERR henckyPtr->getMatDPtr(
858  dataAtPts->getMatDPtr(), dataAtPts->getMatAxiatorDPtr(),
859  dataAtPts->getMatDeviatorDPtr(), b.bulkModulusK, b.shearModulusG);
861  }
862  }
863 
864  const auto E = henckyPtr->E;
865  const auto nu = henckyPtr->nu;
866 
867  double bulk_modulus_K = E / (3 * (1 - 2 * nu));
868  double shear_modulus_G = E / (2 * (1 + nu));
869 
870  dataAtPts->mu = shear_modulus_G;
871  dataAtPts->lambda = bulk_modulus_K;
872 
873  CHKERR henckyPtr->getMatDPtr(
874  dataAtPts->getMatDPtr(), dataAtPts->getMatAxiatorDPtr(),
875  dataAtPts->getMatDeviatorDPtr(), bulk_modulus_K, shear_modulus_G);
877  };
878 
879  auto get_invert_D = [&]() {
881  dataAtPts->matInvD.resize(size_symm, size_symm, false);
882  noalias(dataAtPts->matInvD) = dataAtPts->matD;
883  CHKERR computeMatrixInverse(dataAtPts->matInvD);
885  };
886 
887  CHKERR get_D();
888  CHKERR get_invert_D();
889 
890  dataAtPts->logStretchTensorAtPts.resize(size_symm, nb_integration_pts, false);
891 
892  auto t_log_u =
893  getFTensor2SymmetricFromMat<3>(dataAtPts->logStretchTensorAtPts);
894  auto t_approx_P =
895  getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(dataAtPts->approxPAtPts);
896  auto t_inv_D = getFTensor4DdgFromPtr<SPACE_DIM, SPACE_DIM, 0>(
897  &*dataAtPts->matInvD.data().begin());
898 
899  // note: add rotation, so we can extract rigid body motion, work then with
900  // symmetric part.
901  for (auto gg = 0; gg != nb_integration_pts; ++gg) {
902  t_log_u(i, j) = t_inv_D(i, j, k, l) * t_approx_P(k, l);
903  ++t_log_u;
904  ++t_approx_P;
905  ++t_inv_D;
906  }
907 
909 }
910 
911 } // namespace EshelbianPlasticity
NOSPACE
@ NOSPACE
Definition: definitions.h:83
EshelbianPlasticity::HMHHencky::OpSpatialPhysical_du_du::integratePolyconvexHencky
MoFEMErrorCode integratePolyconvexHencky(EntData &row_data, EntData &col_data)
Definition: Hencky.cpp:611
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
CHK_MOAB_THROW
#define CHK_MOAB_THROW(err, msg)
Check error code of MoAB function and throw MoFEM exception.
Definition: definitions.h:589
get_D
auto get_D
Definition: free_surface.cpp:252
MoFEM::K
VectorDouble K
Definition: Projection10NodeCoordsOnField.cpp:125
EshelbianPlasticity::HMHHencky::OpHenckyJacobian::evaluateRhs
MoFEMErrorCode evaluateRhs(EntData &data)
Definition: Hencky.cpp:63
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:128
EshelbianPlasticity::HMHHencky::getOptions
MoFEMErrorCode getOptions(boost::shared_ptr< DataAtIntegrationPts > data_ptr)
Definition: Hencky.cpp:174
EshelbianPlasticity::HMHHencky::OpCalculateEnergy::totalEnergyPtr
boost::shared_ptr< double > totalEnergyPtr
Definition: Hencky.cpp:146
EshelbianPlasticity::LINEAR
@ LINEAR
Definition: EshelbianPlasticity.hpp:44
EshelbianPlasticity::size_symm
constexpr static auto size_symm
Definition: EshelbianAux.hpp:39
EshelbianPlasticity::EshelbianCore::dd_f
static boost::function< double(const double)> dd_f
Definition: EshelbianPlasticity.hpp:905
EshelbianPlasticity::HMHHencky::recordTape
MoFEMErrorCode recordTape(const int tag, DTensor2Ptr *t_h)
Definition: Hencky.cpp:17
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:121
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:609
EshelbianPlasticity::HMHHencky::OpSpatialPhysical::OpSpatialPhysical
OpSpatialPhysical(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha_u)
Definition: Hencky.cpp:287
EshelbianPlasticity::HMHHencky::OpHenckyJacobian
Definition: Hencky.cpp:19
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
E
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
EshelbianPlasticity::HMHHencky::OpHenckyJacobian::henckyPtr
boost::shared_ptr< HMHHencky > henckyPtr
Definition: Hencky.cpp:68
EshelbianPlasticity
Definition: CGGTonsorialBubbleBase.hpp:11
EshelbianPlasticity::diff_tensor
auto diff_tensor()
Definition: EshelbianAux.hpp:41
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
bulk_modulus_K
double bulk_modulus_K
Definition: dynamic_first_order_con_law.cpp:96
EshelbianPlasticity::HMHHencky::blockData
std::vector< BlockData > blockData
Definition: Hencky.cpp:281
EshelbianPlasticity::PhysicalEquations
Definition: EshelbianPlasticity.hpp:219
EshelbianPlasticity::HMHHencky::OpSpatialPhysical::integrate
MoFEMErrorCode integrate(EntData &data)
Definition: Hencky.cpp:297
J
FTensor::Index< 'J', DIM1 > J
Definition: level_set.cpp:30
EshelbianPlasticity::HMHHencky::BlockData
Definition: Hencky.cpp:276
EshelbianPlasticity::HMHHencky::OpCalculateEnergy
Calculate energy density for Hencky material model.
Definition: Hencky.cpp:138
FTensor::Tensor2_symmetric
Definition: Tensor2_symmetric_value.hpp:13
EshelbianPlasticity::EshelbianCore::d_f
static boost::function< double(const double)> d_f
Definition: EshelbianPlasticity.hpp:904
FTENSOR_INDEX
#define FTENSOR_INDEX(DIM, I)
Definition: Templates.hpp:2011
EshelbianPlasticity::HMHHencky::OpCalculateStretchFromStress::henckyPtr
boost::shared_ptr< HMHHencky > henckyPtr
Definition: Hencky.cpp:164
EshelbianPlasticity::HMHHencky::nu
double nu
Definition: Hencky.cpp:284
EshelbianPlasticity::HMHHencky::returnOpSpatialPhysical
virtual VolUserDataOperator * returnOpSpatialPhysical(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha_u)
Definition: Hencky.cpp:97
sdf.r
int r
Definition: sdf.py:8
FTensor::Tensor2
Definition: Tensor2_value.hpp:16
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
EshelbianPlasticity::HMHHencky::OpSpatialPhysical_du_du::alphaU
const double alphaU
Definition: Hencky.cpp:104
EshelbianPlasticity::HMHHencky::extractBlockData
MoFEMErrorCode extractBlockData(std::vector< const CubitMeshSets * > meshset_vec_ptr, Sev sev)
Definition: Hencky.cpp:202
MoFEM::EntitiesFieldData::EntData::getFTensor0N
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
Definition: EntitiesFieldData.hpp:1502
EshelbianPlasticity::HMHHencky::OpHenckyJacobian::evaluateLhs
MoFEMErrorCode evaluateLhs(EntData &data)
Definition: Hencky.cpp:64
EshelbianPlasticity::HMHHencky::OpSpatialPhysical::integrateHencky
MoFEMErrorCode integrateHencky(EntData &data)
Definition: Hencky.cpp:307
c
const double c
speed of light (cm/ns)
Definition: initial_diffusion.cpp:39
EshelbianPlasticity::HMHHencky::OpHenckyJacobian::OpHenckyJacobian
OpHenckyJacobian(boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< HMHHencky > hencky_ptr)
Definition: Hencky.cpp:20
MoFEM::Sev
MoFEM::LogManager::SeverityLevel Sev
Definition: CoreTemplates.hpp:17
EshelbianPlasticity::HMHHencky::OpCalculateStretchFromStress::dataAtPts
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
data at integration pts
Definition: Hencky.cpp:163
EshelbianPlasticity::HMHHencky::OpSpatialPhysical_du_du::integrateHencky
MoFEMErrorCode integrateHencky(EntData &row_data, EntData &col_data)
Definition: Hencky.cpp:477
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
EshelbianPlasticity::HMHHencky::OpSpatialPhysical
Definition: Hencky.cpp:79
SPACE_DIM
constexpr int SPACE_DIM
Definition: child_and_parent.cpp:16
a
constexpr double a
Definition: approx_sphere.cpp:30
EshelbianPlasticity::HMHHencky::BlockData::bulkModulusK
double bulkModulusK
Definition: Hencky.cpp:277
EshelbianPlasticity::EshelbianCore::stretchSelector
static enum StretchSelector stretchSelector
Definition: EshelbianPlasticity.hpp:898
convert.type
type
Definition: convert.py:64
MoFEM::getFTensor0FromVec
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
EshelbianPlasticity::OpAssembleVolume
Definition: EshelbianPlasticity.hpp:517
EshelbianPlasticity::HMHHencky::OpCalculateEnergy::dataAtPts
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
Definition: Hencky.cpp:145
EshelbianPlasticity::HMHHencky
Definition: Hencky.cpp:12
EshelbianPlasticity::EshelbianCore::f
static boost::function< double(const double)> f
Definition: EshelbianPlasticity.hpp:903
MoFEM::EntitiesFieldData::EntData::getIndices
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Definition: EntitiesFieldData.hpp:1214
EshelbianPlasticity::HMHHencky::getMatDPtr
MoFEMErrorCode getMatDPtr(boost::shared_ptr< MatrixDouble > mat_D_ptr, boost::shared_ptr< MatrixDouble > mat_axiator_D_ptr, boost::shared_ptr< MatrixDouble > mat_deviator_D_ptr, double bulk_modulus_K, double shear_modulus_G)
Definition: Hencky.cpp:234
poisson_ratio
double poisson_ratio
Poisson ratio.
Definition: plastic.cpp:122
EshelbianPlasticity::HMHHencky::OpCalculateStretchFromStress
Definition: Hencky.cpp:155
EshelbianPlasticity::HMHHencky::E
double E
Definition: Hencky.cpp:283
EshelbianPlasticity::HMHHencky::OpCalculateEnergy::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Definition: Hencky.cpp:784
EshelbianPlasticity::HMHHencky::OpHenckyJacobian::dataAtGaussPts
boost::shared_ptr< DataAtIntegrationPts > dataAtGaussPts
Definition: Hencky.cpp:67
EshelbianPlasticity::HMHHencky::returnOpCalculateStretchFromStress
VolUserDataOperator * returnOpCalculateStretchFromStress(boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
Definition: Hencky.cpp:167
MoFEM::L
VectorDouble L
Definition: Projection10NodeCoordsOnField.cpp:124
EshelbianPlasticity::HMHHencky::mField
MoFEM::Interface & mField
Definition: Hencky.cpp:274
EshelbianPlasticity::diff_deviator
auto diff_deviator(FTensor::Ddg< double, 3, 3 > &&t_diff_stress)
Definition: EshelbianAux.hpp:12
EshelbianPlasticity::HMHHencky::returnOpCalculateEnergy
VolUserDataOperator * returnOpCalculateEnergy(boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< double > total_energy_ptr)
Definition: Hencky.cpp:150
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
t_kd
constexpr auto t_kd
Definition: free_surface.cpp:137
EshelbianPlasticity::symm_L_tensor
auto symm_L_tensor()
Definition: EshelbianAux.hpp:52
EshelbianPlasticity::HMHHencky::extractBlockData
MoFEMErrorCode extractBlockData(Sev sev)
Definition: Hencky.cpp:189
EshelbianPlasticity::HMHHencky::HMHHencky
HMHHencky(MoFEM::Interface &m_field, const double E, const double nu)
Definition: Hencky.cpp:14
MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
Definition: VolumeElementForcesAndSourcesCore.hpp:108
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', SPACE_DIM >
EshelbianPlasticity::HMHHencky::OpSpatialPhysical_du_du::polyConvex
PetscBool polyConvex
Definition: Hencky.cpp:114
convert.n
n
Definition: convert.py:82
EshelbianPlasticity::HMHHencky::OpCalculateStretchFromStress::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Definition: Hencky.cpp:830
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
Range
MOFEM_TAG_AND_LOG
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
Definition: LogManager.hpp:362
shear_modulus_G
double shear_modulus_G
Definition: dynamic_first_order_con_law.cpp:97
MoFEM::EntitiesFieldData::EntData::getN
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
Definition: EntitiesFieldData.hpp:1318
EshelbianPlasticity::HMHHencky::BlockData::shearModulusG
double shearModulusG
Definition: Hencky.cpp:278
EshelbianPlasticity::OpJacobian
Definition: EshelbianPlasticity.hpp:380
EshelbianPlasticity::HMHHencky::OpSpatialPhysical::integratePolyconvexHencky
MoFEMErrorCode integratePolyconvexHencky(EntData &data)
Definition: Hencky.cpp:366
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
FTensor::Ddg< double, 3, 3 >
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
mu
double mu
Definition: dynamic_first_order_con_law.cpp:98
lambda
static double lambda
Definition: incompressible_elasticity.cpp:199
EshelbianPlasticity::HMHHencky::returnOpSpatialPhysical_du_du
VolUserDataOperator * returnOpSpatialPhysical_du_du(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha)
Definition: Hencky.cpp:119
EshelbianPlasticity::HMHHencky::OpSpatialPhysical::alphaU
const double alphaU
Definition: Hencky.cpp:92
EigenMatrix::getDiffDiffMat
FTensor::Ddg< double, 3, 3 > getDiffDiffMat(Val< double, 3 > &t_val, Vec< double, 3 > &t_vec, Fun< double > f, Fun< double > d_f, Fun< double > dd_f, FTensor::Tensor2< double, 3, 3 > &t_S, const int nb)
Definition: MatrixFunction.cpp:78
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
EshelbianPlasticity::HMHHencky::returnOpJacobian
virtual OpJacobian * returnOpJacobian(const int tag, const bool eval_rhs, const bool eval_lhs, boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
Definition: Hencky.cpp:72
MoFEM::computeMatrixInverse
MoFEMErrorCode computeMatrixInverse(MatrixDouble &mat)
compute matrix inverse with lapack dgetri
Definition: Templates.hpp:1365
EshelbianPlasticity::HMHHencky::OpSpatialPhysical_du_du::dP
MatrixDouble dP
Definition: Hencky.cpp:116
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
EshelbianPlasticity::HMHHencky::OpCalculateStretchFromStress::OpCalculateStretchFromStress
OpCalculateStretchFromStress(boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< HMHHencky > hencky_ptr)
Definition: Hencky.cpp:824
FTensor::Kronecker_Delta_symmetric
Kronecker Delta class symmetric.
Definition: Kronecker_Delta.hpp:49
third
constexpr double third
Definition: EshelbianADOL-C.cpp:17
EshelbianPlasticity::HMHHencky::OpSpatialPhysical::polyConvex
PetscBool polyConvex
Definition: Hencky.cpp:93
EshelbianPlasticity::HMHHencky::BlockData::blockEnts
Range blockEnts
Definition: Hencky.cpp:279
EshelbianPlasticity::HMHHencky::OpSpatialPhysical_du_du::OpSpatialPhysical_du_du
OpSpatialPhysical_du_du(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha)
Definition: Hencky.cpp:451
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
EshelbianPlasticity::HMHHencky::OpCalculateEnergy::OpCalculateEnergy
OpCalculateEnergy(boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< double > total_energy_ptr)
Definition: Hencky.cpp:778
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
EshelbianPlasticity::HMHHencky::OpSpatialPhysical_du_du
Definition: Hencky.cpp:103
CHKERRG
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:496
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
EshelbianPlasticity::HMHHencky::OpSpatialPhysical_du_du::integrate
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
Definition: Hencky.cpp:464
EshelbianPlasticity::HMHHencky::OpHenckyJacobian::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Definition: Hencky.cpp:29
MoFEM::PetscOptionsGetBool
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:182