v0.15.5
Loading...
Searching...
No Matches
HMHNeohookean.cpp
Go to the documentation of this file.
1/**
2 * @file HMHNeohookean.cpp
3 * @brief Implementation of NeoHookean material
4 * @date 2024-08-31
5 *
6 * @copyright Copyright (c) 2024
7 *
8 */
9
10namespace EshelbianPlasticity {
11
12// #define NEOHOOKEAN_SCALING
14
15 /**
16 * @brief Definition of Neo-hookean function
17 *
18 * P(i, I) = 2. * c10_default * (tF(i, I) - tInvF(i, I)) +
19 * K_default * log(detF) * tInvF(i, I)
20 * tCof(i, I) = detF * tInvF(I, i)
21 * Psi = c10* (trace(tCof) + 2ln(det(F)) - 3) + (K/2.) * (log(detF))^2
22 *
23 * @param c10
24 * @param v
25 * @return double
26 */
27 static inline double fun_neohookean(double c10, double v) {
28 return (2.0 * c10) * (-1. + EshelbianCore::f(2 * v));
29 }
30
31 /**
32 * @brief Definition of derivative of Neo-hookean function
33 *
34 * @param c10
35 * @param v
36 * @return * double
37 */
38 static inline double fun_d_neohookean(double c10, double v) {
39 return
40
41 (4.0 * c10) * EshelbianCore::d_f(2 * v);
42 }
43
44 /**
45 * @brief Definition of axiator of Neo-hookean function
46 *
47 * Psi = (K/2.) * (log(J))^2 = (K/2.) * (log(exp tr H))^2 = (K/2.) * (tr H)^2
48 *
49 * @param c10
50 * @param v
51 * @return double
52 */
53 static inline double fun_neohookean_bulk(double K, double tr) {
54 return K * tr;
55 }
56
57 /**
58 * @brief Definition of derivative of axiator of Neo-hookean function
59 *
60 * @param K
61 * @param tr
62 * @return double
63 */
64 static inline double fun_diff_neohookean_bulk(double K, double tr) {
65 return K;//K * exp(tr) * (1. + tr); // K*2*ln(J)/J
66 }
67
68 static constexpr int numberOfActiveVariables = 9;
69 static constexpr int numberOfDependentVariables = 9;
70
71 HMHNeohookean(MoFEM::Interface &m_field, const double c10, const double K)
73 mField(m_field), c10_default(c10), K_default(K) {
74
75 CHK_THROW_MESSAGE(getOptions(), "get options failed");
77 "extract block data failed");
78
79#ifdef NEOHOOKEAN_SCALING
80 if (blockData.size()) {
81 double ave_K = 0;
82 for (auto b : blockData) {
83 ave_K += b.K;
84 }
85 eqScaling = ave_K / blockData.size();
86 }
87#endif
88
89 MOFEM_LOG("EP", Sev::inform) << "Neo-Hookean scaling = " << eqScaling;
90
91 if (adolCOn == PETSC_FALSE) {
94 "Stretch selector is not equal to LOG");
95 } else {
96 if (EshelbianCore::exponentBase != exp(1)) {
98 "Exponent base is not equal to exp(1)");
99 }
100 }
101 }
102 }
103
105 for (auto &b : blockData) {
106 if (b.blockEnts.find(ent) != b.blockEnts.end()) {
107 return std::make_pair(b.c10, b.K);
108 }
109 }
110 if (blockData.size() != 0)
112 "Block not found for entity handle. If you mat set "
113 "block, set it to all elements");
114 return std::make_pair(c10_default, K_default);
115 }
116
118 OpGetScale(boost::shared_ptr<double> scale_ptr,
119 boost::shared_ptr<PhysicalEquations> physics_ptr)
121 scalePtr(scale_ptr), physicsPtr(physics_ptr) {}
122
123 MoFEMErrorCode doWork(int side, EntityType type,
124 EntitiesFieldData::EntData &data) {
126 auto neohookean_ptr =
127 boost::dynamic_pointer_cast<HMHNeohookean>(physicsPtr);
128 *(scalePtr) = neohookean_ptr->eqScaling;
130 }
131
132 private:
133 boost::shared_ptr<double> scalePtr;
134 boost::shared_ptr<PhysicalEquations> physicsPtr;
135 };
136
138 returnOpSetScale(boost::shared_ptr<double> scale_ptr,
139 boost::shared_ptr<PhysicalEquations> physics_ptr) {
140
141 return new OpGetScale(scale_ptr, physics_ptr);
142 };
143
144 struct OpJacobian : public EshelbianPlasticity::OpJacobian {
145 using EshelbianPlasticity::OpJacobian::OpJacobian;
146 MoFEMErrorCode evaluateRhs(EntData &data) { return 0; }
147 MoFEMErrorCode evaluateLhs(EntData &data) { return 0; }
148 };
149
151 returnOpJacobian(const int tag, const bool eval_rhs, const bool eval_lhs,
152 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
153 boost::shared_ptr<PhysicalEquations> physics_ptr) {
154 if (adolCOn) {
155 return PhysicalEquations::returnOpJacobian(tag, eval_rhs, eval_lhs,
156 data_ptr, physics_ptr);
157 } else {
158 return (new OpJacobian(tag, eval_rhs, eval_lhs, data_ptr, physics_ptr));
159 }
160 }
161
162 MoFEMErrorCode getOptions() {
164 PetscOptionsBegin(PETSC_COMM_WORLD, "neo_hookean_", "", "none");
165
166 c10_default = 1;
167 CHKERR PetscOptionsScalar("-c10", "C10", "", c10_default, &c10_default,
168 PETSC_NULLPTR);
169 K_default = 1;
170 CHKERR PetscOptionsScalar("-K", "Bulk modulus K", "", K_default, &K_default,
171 PETSC_NULLPTR);
172
173 alphaGradU = 0;
174 CHKERR PetscOptionsScalar("-viscosity_alpha_grad_u", "viscosity", "",
175 alphaGradU, &alphaGradU, PETSC_NULLPTR);
176
177 CHKERR PetscOptionsBool("-adolc_on", "adolc ON", "", adolCOn, &adolCOn,
178 PETSC_NULLPTR);
179
180 PetscOptionsEnd();
181
182 MOFEM_LOG_CHANNEL("WORLD");
183 MOFEM_TAG_AND_LOG("WORLD", Sev::inform, "MatBlock Neo-Hookean (default)")
184 << "c10 = " << c10_default << " K = " << K_default
185 << " grad alpha u = " << alphaGradU;
186 MOFEM_TAG_AND_LOG("WORLD", Sev::inform, "MatBlock Neo-Hookean (default)")
187 << "ADOL-C On " << (adolCOn)
188 ? "Yes"
189 : "No";
191 }
192
193 MoFEMErrorCode extractBlockData(Sev sev) {
194 return extractBlockData(
195
196 mField.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
197
198 (boost::format("%s(.*)") % "MAT_NEOHOOKEAN").str()
199
200 )),
201
202 sev);
203 }
204
205 MoFEMErrorCode
206 extractBlockData(std::vector<const CubitMeshSets *> meshset_vec_ptr,
207 Sev sev) {
209
210 for (auto m : meshset_vec_ptr) {
211 MOFEM_LOG("EP", sev) << *m;
212 std::vector<double> block_data;
213 CHKERR m->getAttributes(block_data);
214 if (block_data.size() < 2) {
215 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
216 "Expected that block has atleast two attributes");
217 }
218 auto get_block_ents = [&]() {
219 Range ents;
220 CHKERR mField.get_moab().get_entities_by_handle(m->meshset, ents, true);
221 return ents;
222 };
223
224 double c10 = block_data[0];
225 double K = block_data[1];
226
227 blockData.push_back({c10, K, get_block_ents()});
228
229 MOFEM_LOG("EP", sev) << "MatBlock Neo-Hookean c10 = "
230 << blockData.back().c10
231 << " K = " << blockData.back().K << " nb ents. = "
232 << blockData.back().blockEnts.size();
233 }
235 }
236
237 MoFEMErrorCode recordTape(const int tape, DTensor2Ptr *t_h_ptr) {
239
240 FTENSOR_INDEX(3, i);
241 FTENSOR_INDEX(3, j);
242 FTENSOR_INDEX(3, I);
243 FTENSOR_INDEX(3, J);
244
245 auto large = [&]() {
247
248 constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
249
250 auto ih = get_h();
251 if (t_h_ptr)
252 ih(i, j) = (*t_h_ptr)(i, j);
253 else {
254 ih(i, j) = t_kd(i, j);
255 }
256
257 auto r_P = get_P();
258
259 enableMinMaxUsingAbs();
260 trace_on(tape);
261
262 // Set active variables to ADOL-C
263 th(i, j) <<= get_h()(i, j);
264
265 tF(i, I) = th(i, I);
266 CHKERR determinantTensor3by3(tF, detF);
267 CHKERR invertTensor3by3(tF, detF, tInvF);
268
269 // tCof(i, I) = detF * tInvF(I, i);
270
271 // Psi = c10* (trace(tCof) + 2ln(det(F)) - 3) + (K/2.) * (log(detF))^2;
272
273 tP(i, I) = 2. * c10_default * (tF(i, I) - tInvF(i, I)) +
274 K_default * log(detF) * tInvF(i, I);
275
276 // Set dependent variables to ADOL-C
277 tP(i, j) >>= r_P(i, j);
278
279 trace_off();
280
282 };
283
286 case LARGE_ROT:
287 case MODERATE_ROT:
288 CHKERR large();
289 break;
290 default:
291 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
292 "gradApproximator not handled");
293 };
294
296 }
297
299
300 OpSpatialPhysical(const std::string &field_name,
301 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
302 const double alpha_u);
303
304 MoFEMErrorCode integrate(EntData &data);
305
306 private:
307 const double alphaU;
308 };
309
310 virtual VolUserDataOperator *
312 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
313 const double alpha_u) {
314 if (adolCOn) {
316 alpha_u);
317 } else {
318 return new OpSpatialPhysical(field_name, data_ptr, alpha_u);
319 }
320 }
321
323
325 const std::string &field_name,
326 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
327 boost::shared_ptr<ExternalStrainVec> external_strain_vec_ptr,
328 std::map<std::string, boost::shared_ptr<ScalingMethod>> smv);
329
330 MoFEMErrorCode integrate(EntData &data);
331
332 private:
333 boost::shared_ptr<ExternalStrainVec> externalStrainVecPtr;
334 std::map<std::string, boost::shared_ptr<ScalingMethod>> scalingMethodsMap;
335 };
336
338 const std::string &field_name,
339 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
340 boost::shared_ptr<ExternalStrainVec> external_strain_vec_ptr,
341 std::map<std::string, boost::shared_ptr<ScalingMethod>> smv) {
342 return new OpSpatialPhysicalExternalStrain(field_name, data_ptr,
343 external_strain_vec_ptr, smv);
344 }
345
347 OpSpatialPhysical_du_du(std::string row_field, std::string col_field,
348 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
349 const double alpha);
350 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
351
352 private:
353 const double alphaU;
354 };
355
357 std::string row_field, std::string col_field,
358 boost::shared_ptr<DataAtIntegrationPts> data_ptr, const double alpha) {
359 if (adolCOn) {
361 row_field, col_field, data_ptr, alpha);
362 } else {
363 return new OpSpatialPhysical_du_du(row_field, col_field, data_ptr, alpha);
364 }
365 }
366
367private:
369
371 double K_default;
373 struct BlockData {
374 double c10;
375 double K;
377 };
378 std::vector<BlockData> blockData;
379
380 double eqScaling = 1.;
381
382 // ADOL-C
383 PetscBool adolCOn = PETSC_FALSE;
384
391};
392
394 const std::string &field_name,
395 boost::shared_ptr<DataAtIntegrationPts> data_ptr, const double alpha_u)
396 : OpAssembleVolume(field_name, data_ptr, OPROW), alphaU(alpha_u) {}
397
398extern "C" {
399void tetcircumcenter_tp(double a[3], double b[3], double c[3], double d[3],
400 double circumcenter[3], double *xi, double *eta,
401 double *zeta);
402}
403
406
407 auto neohookean_ptr =
408 boost::dynamic_pointer_cast<HMHNeohookean>(dataAtPts->physicsPtr);
409 if (!neohookean_ptr) {
410 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
411 "Pointer to HMHNeohookean is null");
412 }
413 auto [def_c10, def_K] =
414 neohookean_ptr->getMaterialParameters(getFEEntityHandle());
415
416 double c10 = def_c10 / neohookean_ptr->eqScaling;
417 double alpha_u = alphaU / neohookean_ptr->eqScaling;
418 double K = def_K / neohookean_ptr->eqScaling;
419
420 double alpha_grad_u =
421 neohookean_ptr->alphaGradU / neohookean_ptr->eqScaling;
422
424 auto t_L = symm_L_tensor();
425
426 constexpr auto t_kd_sym = FTensor::Kronecker_Delta_symmetric<int>();
427 constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
428
429 int nb_dofs = data.getIndices().size();
430 int nb_integration_pts = data.getN().size1();
431 auto v = getVolume();
432 auto t_w = getFTensor0IntegrationWeight();
433 auto t_approx_P_adjoint_log_du =
434 getFTensor1FromMat<size_symm>(dataAtPts->adjointPdUAtPts);
435 auto t_u = getFTensor2SymmetricFromMat<3>(dataAtPts->stretchTensorAtPts);
436 auto t_grad_h1 = getFTensor2FromMat<3, 3>(dataAtPts->wGradH1AtPts);
437 auto t_dot_log_u =
438 getFTensor2SymmetricFromMat<3>(dataAtPts->logStretchDotTensorAtPts);
439 auto t_diff_u =
440 getFTensor4DdgFromMat<3, 3, 1>(dataAtPts->diffStretchTensorAtPts);
441 auto t_log_u =
442 getFTensor2SymmetricFromMat<3>(dataAtPts->logStretchTensorAtPts);
443 auto t_grad_log_u =
444 getFTensor2FromMat<6, 3>(dataAtPts->gradLogStretchDotTensorAtPts);
445 auto t_log_u2_h1 =
446 getFTensor2SymmetricFromMat<3>(dataAtPts->logStretch2H1AtPts);
447
448 auto t_eigen_vals = getFTensor1FromMat<3>(dataAtPts->eigenVals);
449 auto t_eigen_vecs = getFTensor2FromMat<3, 3>(dataAtPts->eigenVecs);
450 auto &nbUniq = dataAtPts->nbUniq;
451 auto t_nb_uniq =
452 FTensor::Tensor0<FTensor::PackPtr<int *, 1>>(nbUniq.data().data());
453
454 auto t_diff = diff_tensor();
455
462
463 auto get_ftensor2 = [](auto &v) {
465 &v[0], &v[1], &v[2], &v[3], &v[4], &v[5]);
466 };
467
468 int nb_base_functions = data.getN().size2();
469 auto t_row_base_fun = data.getFTensor0N();
470 auto t_grad_base_fun = data.getFTensor1DiffN<3>();
471
472 auto no_h1 = [&]() {
474
475 for (int gg = 0; gg != nb_integration_pts; ++gg) {
476 double a = v * t_w;
477 ++t_w;
478
479 auto neohookean = [c10](auto v) { return fun_neohookean(c10, v); };
480 auto t_neohookean_hencky =
481 EigenMatrix::getMat(t_eigen_vals, t_eigen_vecs, neohookean);
482 ++t_eigen_vals;
483 ++t_eigen_vecs;
484 ++t_nb_uniq;
485
486 const double tr = t_log_u(i, i);
488 t_P(L) = t_L(i, j, L) * (t_neohookean_hencky(i, j) +
489 t_kd(i, j) * fun_neohookean_bulk(K, tr));
491 t_viscous_P(L) = alpha_u * (t_L(i, j, L) * t_dot_log_u(i, j));
492
494 t_residual(L) = t_approx_P_adjoint_log_du(L) - t_P(L) - t_viscous_P(L);
495 t_residual(L) *= a;
496
498 t_grad_residual(L, i) = alpha_grad_u * t_grad_log_u(L, i);
499 t_grad_residual(L, i) *= a;
500
501 ++t_approx_P_adjoint_log_du;
502 ++t_log_u;
503 ++t_dot_log_u;
504 ++t_grad_log_u;
505
506 auto t_nf = getFTensor1FromPtr<size_symm>(&*nF.data().begin());
507 int bb = 0;
508 for (; bb != nb_dofs / size_symm; ++bb) {
509 t_nf(L) -= t_row_base_fun * t_residual(L);
510 t_nf(L) += t_grad_base_fun(i) * t_grad_residual(L, i);
511 ++t_nf;
512 ++t_row_base_fun;
513 ++t_grad_base_fun;
514 }
515 for (; bb != nb_base_functions; ++bb) {
516 ++t_row_base_fun;
517 ++t_grad_base_fun;
518 }
519 }
520
522 };
523
524 auto large = [&]() {
526 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
527 "Not implemented for Neo-Hookean (used ADOL-C)");
529 };
530
533 CHKERR no_h1();
534 break;
535 case LARGE_ROT:
536 case MODERATE_ROT:
537 CHKERR large();
538 break;
539 default:
540 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
541 "gradApproximator not handled");
542 };
543
545}
546
548 const std::string &field_name,
549 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
550 boost::shared_ptr<ExternalStrainVec> external_strain_vec_ptr,
551 std::map<std::string, boost::shared_ptr<ScalingMethod>> smv)
552 : OpAssembleVolume(field_name, data_ptr, OPROW),
553 externalStrainVecPtr(external_strain_vec_ptr), scalingMethodsMap{smv} {}
554
555MoFEMErrorCode
558
559 auto neohookean_ptr =
560 boost::dynamic_pointer_cast<HMHNeohookean>(dataAtPts->physicsPtr);
561 if (!neohookean_ptr) {
562 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
563 "Pointer to HMHNeohookean is null");
564 }
565
566 double time = OpAssembleVolume::getFEMethod()->ts_t;
569 }
570 // get entity of tet
572 // iterate over all block data
573
574 for (auto &ext_strain_block : (*externalStrainVecPtr)) {
575 auto block_name = "(.*)ANALYTICAL_EXTERNALSTRAIN(.*)";
576 std::regex reg_name(block_name);
577 if (std::regex_match(ext_strain_block.blockName, reg_name)) {
578 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
579 "Analytical external strain not implemented for Neo-Hookean "
580 "material.");
581 }
582 // check if finite element entity is part of the EXTERNALSTRAIN block
583 if (ext_strain_block.ents.find(fe_ent) != ext_strain_block.ents.end()) {
584 double scale = 1;
585 if (scalingMethodsMap.find(ext_strain_block.blockName) !=
586 scalingMethodsMap.end()) {
587 scale *=
588 scalingMethodsMap.at(ext_strain_block.blockName)->getScale(time);
589 } else {
590 MOFEM_LOG("SELF", Sev::warning)
591 << "No scaling method found for " << ext_strain_block.blockName;
592 }
593
594 // get ExternalStrain block data
595 double external_strain_val = scale * ext_strain_block.val;
596 double K = ext_strain_block.bulkModulusK;
597
599 auto t_L = symm_L_tensor();
600 constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
601
602 int nb_dofs = data.getIndices().size();
603 int nb_integration_pts = data.getN().size1();
604 auto v = getVolume();
605 auto t_w = getFTensor0IntegrationWeight();
606
613
614 int nb_base_functions = data.getN().size2();
615 auto t_row_base_fun = data.getFTensor0N();
616
617 const double tr = 3.0 * external_strain_val;
618 const double sigma_J = K * tr;
619
620 for (int gg = 0; gg != nb_integration_pts; ++gg) {
621 double a = v * t_w;
622 ++t_w;
623
625 t_residual(L) = 0.0;
626 t_residual(L) += (t_L(i, j, L) * t_kd(i, j)) * sigma_J;
627 t_residual(L) *= a;
628
629 auto t_nf = getFTensor1FromPtr<size_symm>(&*nF.data().begin());
630 int bb = 0;
631 for (; bb != nb_dofs / size_symm; ++bb) {
632 t_nf(L) -= t_row_base_fun * t_residual(L);
633 ++t_nf;
634 ++t_row_base_fun;
635 }
636 for (; bb != nb_base_functions; ++bb) {
637 ++t_row_base_fun;
638 }
639 }
640 }
641 }
643}
644
646 std::string row_field, std::string col_field,
647 boost::shared_ptr<DataAtIntegrationPts> data_ptr, const double alpha)
648 : OpAssembleVolumePositiveDefine(row_field, col_field, data_ptr, OPROWCOL,
649 false),
650 alphaU(alpha) {
651 sYmm = false;
652}
653
654MoFEMErrorCode
656 EntData &col_data) {
658
659 auto neohookean_ptr =
660 boost::dynamic_pointer_cast<HMHNeohookean>(dataAtPts->physicsPtr);
661 if (!neohookean_ptr) {
662 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
663 "Pointer to HMHNeohookean is null");
664 }
665 auto [def_c10, def_K] =
666 neohookean_ptr->getMaterialParameters(getFEEntityHandle());
667
668 double c10 = def_c10 / neohookean_ptr->eqScaling;
669 double alpha_u = alphaU / neohookean_ptr->eqScaling;
670 double lambda = def_K / neohookean_ptr->eqScaling;
671
672 double alpha_grad_u =
673 neohookean_ptr->alphaGradU / neohookean_ptr->eqScaling;
674
677
678 constexpr auto t_kd_sym = FTensor::Kronecker_Delta_symmetric<int>();
679 constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
680
681 auto t_L = symm_L_tensor();
682 auto t_diff = diff_tensor();
683
684 int nb_integration_pts = row_data.getN().size1();
685 int row_nb_dofs = row_data.getIndices().size();
686 int col_nb_dofs = col_data.getIndices().size();
687
688 auto get_ftensor2 = [](MatrixDouble &m, const int r, const int c) {
690 size_symm>(
691
692 &m(r + 0, c + 0), &m(r + 0, c + 1), &m(r + 0, c + 2), &m(r + 0, c + 3),
693 &m(r + 0, c + 4), &m(r + 0, c + 5),
694
695 &m(r + 1, c + 0), &m(r + 1, c + 1), &m(r + 1, c + 2), &m(r + 1, c + 3),
696 &m(r + 1, c + 4), &m(r + 1, c + 5),
697
698 &m(r + 2, c + 0), &m(r + 2, c + 1), &m(r + 2, c + 2), &m(r + 2, c + 3),
699 &m(r + 2, c + 4), &m(r + 2, c + 5),
700
701 &m(r + 3, c + 0), &m(r + 3, c + 1), &m(r + 3, c + 2), &m(r + 3, c + 3),
702 &m(r + 3, c + 4), &m(r + 3, c + 5),
703
704 &m(r + 4, c + 0), &m(r + 4, c + 1), &m(r + 4, c + 2), &m(r + 4, c + 3),
705 &m(r + 4, c + 4), &m(r + 4, c + 5),
706
707 &m(r + 5, c + 0), &m(r + 5, c + 1), &m(r + 5, c + 2), &m(r + 5, c + 3),
708 &m(r + 5, c + 4), &m(r + 5, c + 5)
709
710 );
711 };
712
719
720 auto v = getVolume();
721 auto ts_a = getTSa();
722 auto t_w = getFTensor0IntegrationWeight();
723
724 int row_nb_base_functions = row_data.getN().size2();
725 auto t_row_base_fun = row_data.getFTensor0N();
726 auto t_row_grad_fun = row_data.getFTensor1DiffN<3>();
727
728 auto t_grad_h1 = getFTensor2FromMat<3, 3>(dataAtPts->wGradH1AtPts);
729 auto t_diff_u =
730 getFTensor4DdgFromMat<3, 3, 1>(dataAtPts->diffStretchTensorAtPts);
731 auto t_log_u =
732 getFTensor2SymmetricFromMat<3>(dataAtPts->logStretchTensorAtPts);
733 auto t_log_u2_h1 =
734 getFTensor2SymmetricFromMat<3>(dataAtPts->logStretch2H1AtPts);
735 auto t_u = getFTensor2SymmetricFromMat<3>(dataAtPts->stretchTensorAtPts);
736 auto t_approx_P_adjoint__dstretch =
737 getFTensor2FromMat<3, 3>(dataAtPts->adjointPdstretchAtPts);
738 auto t_eigen_vals = getFTensor1FromMat<3>(dataAtPts->eigenVals);
739 auto t_eigen_vecs = getFTensor2FromMat<3, 3>(dataAtPts->eigenVecs);
740 auto &nbUniq = dataAtPts->nbUniq;
741 auto t_nb_uniq =
742 FTensor::Tensor0<FTensor::PackPtr<int *, 1>>(nbUniq.data().data());
743
744 auto no_h1 = [&]() {
746
747 for (int gg = 0; gg != nb_integration_pts; ++gg) {
748 double a = v * t_w;
749 ++t_w;
750
751 auto neohookean = [c10](auto v) { return fun_neohookean(c10, v); };
752 auto d_neohookean = [c10, lambda](auto v) {
753 return fun_d_neohookean(c10, v);
754 };
755
756 auto t_diff_neohookean = EigenMatrix::getDiffMat(
757 t_eigen_vals, t_eigen_vecs, neohookean, d_neohookean, t_nb_uniq);
758
759 const auto tr = t_log_u(i, i);
761 t_dP(L, J) = -t_L(i, j, L) * ((t_diff_neohookean(i, j, k, l) +
763 t_kd_sym(i, j) * t_kd_sym(k, l)) *
764 t_L(k, l, J));
765 t_dP(L, J) -= (alpha_u * ts_a) *
766 (t_L(i, j, L) * (t_diff(i, j, k, l) * t_L(k, l, J)));
767
768 if constexpr (1) {
770 t_deltaP(i, j) = (t_approx_P_adjoint__dstretch(i, j) ||
771 t_approx_P_adjoint__dstretch(j, i)) /
772 2.;
773 auto t_diff2_uP = EigenMatrix::getDiffDiffMat(
774 t_eigen_vals, t_eigen_vecs, EshelbianCore::f, EshelbianCore::d_f,
775 EshelbianCore::dd_f, t_deltaP, t_nb_uniq);
776 t_dP(L, J) += t_L(i, j, L) * (t_diff2_uP(i, j, k, l) * t_L(k, l, J));
777 }
778 ++t_approx_P_adjoint__dstretch;
779 ++t_log_u;
780 ++t_eigen_vals;
781 ++t_eigen_vecs;
782 ++t_nb_uniq;
783
784 int rr = 0;
785 for (; rr != row_nb_dofs / size_symm; ++rr) {
786 auto t_col_base_fun = col_data.getFTensor0N(gg, 0);
787 auto t_col_grad_fun = col_data.getFTensor1DiffN<3>(gg, 0);
788
789 auto t_m = get_ftensor2(K, 6 * rr, 0);
790 for (int cc = 0; cc != col_nb_dofs / size_symm; ++cc) {
791 double b = a * t_row_base_fun * t_col_base_fun;
792 double c = (a * alpha_grad_u * ts_a) *
793 (t_row_grad_fun(i) * t_col_grad_fun(i));
794 t_m(L, J) -= b * t_dP(L, J);
795 t_m(L, J) += c * t_kd_sym(L, J);
796
797 ++t_m;
798 ++t_col_base_fun;
799 ++t_col_grad_fun;
800 }
801 ++t_row_base_fun;
802 ++t_row_grad_fun;
803 }
804
805 for (; rr != row_nb_base_functions; ++rr) {
806 ++t_row_base_fun;
807 ++t_row_grad_fun;
808 }
809
810 }
812 };
813
814 auto large = [&]() {
816 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
817 "Not implemented for Neo-Hookean (used ADOL-C)");
819 };
820
823 CHKERR no_h1();
824 break;
825 case LARGE_ROT:
826 case MODERATE_ROT:
827 CHKERR large();
828 break;
829 default:
830 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
831 "gradApproximator not handled");
832 };
833
835}
836
837} // namespace EshelbianPlasticity
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
#define FTENSOR_INDEX(DIM, I)
constexpr double a
constexpr int SPACE_DIM
Kronecker Delta class symmetric.
Kronecker Delta class.
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
@ 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
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
constexpr auto t_kd
double eta
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
FTensor::Index< 'i', SPACE_DIM > i
static double lambda
const double c
speed of light (cm/ns)
const double v
phase velocity of light in medium (cm/ns)
const double n
refractive index of diffusive medium
FTensor::Index< 'J', DIM1 > J
Definition level_set.cpp:30
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.
void tetcircumcenter_tp(double a[3], double b[3], double c[3], double d[3], double circumcenter[3], double *xi, double *eta, double *zeta)
ForcesAndSourcesCore::UserDataOperator UserDataOperator
static constexpr auto size_symm
constexpr IntegrationType I
constexpr auto field_name
FTensor::Index< 'm', 3 > m
static enum StretchSelector stretchSelector
static double dynamicTime
static enum RotSelector gradApproximator
static PetscBool dynamicRelaxation
static double exponentBase
static boost::function< double(const double)> f
static boost::function< double(const double)> dd_f
static boost::function< double(const double)> d_f
OpGetScale(boost::shared_ptr< double > scale_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
boost::shared_ptr< PhysicalEquations > physicsPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
std::map< std::string, boost::shared_ptr< ScalingMethod > > scalingMethodsMap
OpSpatialPhysicalExternalStrain(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< ExternalStrainVec > external_strain_vec_ptr, std::map< std::string, boost::shared_ptr< ScalingMethod > > smv)
OpSpatialPhysical_du_du(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialPhysical(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha_u)
static double fun_diff_neohookean_bulk(double K, double tr)
Definition of derivative of axiator of Neo-hookean function.
static double fun_d_neohookean(double c10, double v)
Definition of derivative of Neo-hookean function.
static double fun_neohookean_bulk(double K, double tr)
Definition of axiator of Neo-hookean function.
virtual VolUserDataOperator * returnOpSpatialPhysical(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha_u)
std::vector< BlockData > blockData
MoFEMErrorCode extractBlockData(std::vector< const CubitMeshSets * > meshset_vec_ptr, Sev sev)
MoFEMErrorCode extractBlockData(Sev sev)
static constexpr int numberOfDependentVariables
VolUserDataOperator * returnOpSetScale(boost::shared_ptr< double > scale_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
static double fun_neohookean(double c10, double v)
Definition of Neo-hookean function.
auto getMaterialParameters(EntityHandle ent)
static constexpr int numberOfActiveVariables
virtual VolUserDataOperator * returnOpSpatialPhysicalExternalStrain(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< ExternalStrainVec > external_strain_vec_ptr, std::map< std::string, boost::shared_ptr< ScalingMethod > > smv)
HMHNeohookean(MoFEM::Interface &m_field, const double c10, const double K)
VolUserDataOperator * returnOpSpatialPhysical_du_du(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha)
MoFEMErrorCode recordTape(const int tape, DTensor2Ptr *t_h_ptr)
UserDataOperator * returnOpJacobian(const int tag, const bool eval_rhs, const bool eval_lhs, boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
virtual VolUserDataOperator * returnOpSpatialPhysical(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha_u)
virtual VolUserDataOperator * returnOpSpatialPhysical_du_du(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha)
virtual UserDataOperator * returnOpJacobian(const int tag, const bool eval_rhs, const bool eval_lhs, boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< PhysicalEquations > physics_ptr)
virtual moab::Interface & get_moab()=0
bool sYmm
If true assume that matrix is symmetric structure.
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorInt & getIndices() const
Get global indices of degrees of freedom on entity.
EntityHandle getFEEntityHandle() const
Return finite element entity handle.
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
PetscReal ts_t
Current time value.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
double scale
Definition plastic.cpp:123
double zeta
Viscous hardening.
Definition plastic.cpp:130