v0.15.5
Loading...
Searching...
No Matches
NonLinearElasticElement.cpp
Go to the documentation of this file.
1/**
2 * \brief Operators and data structures for nonlinear elastic material
3 *
4 * Implementation of nonlinear elastic element.
5 */
6
7
8
9#include <MoFEM.hpp>
10
11using namespace MoFEM;
13
14#include <adolc/adolc.h>
16
18 : VolumeElementForcesAndSourcesCore(m_field), A(PETSC_NULLPTR), F(PETSC_NULLPTR),
19 addToRule(1) {
20
21 auto create_vec = [&]() {
22 if (mField.get_comm_rank() == 0) {
23 return createVectorMPI(mField.get_comm(), 1, 1);
24 } else {
25 return createVectorMPI(mField.get_comm(), 0, 1);
26 }
27 };
28
29 V = create_vec();
30}
31
33 return 2 * (order - 1) + addToRule;
34};
35
38
39 CHKERR VolumeElementForcesAndSourcesCore::preProcess();
40
41 if (A != PETSC_NULLPTR) {
42 snes_B = A;
43 }
44
45 if (F != PETSC_NULLPTR) {
46 snes_f = F;
47 }
48
49 switch (snes_ctx) {
50 case CTX_SNESNONE:
51 CHKERR VecZeroEntries(V);
52 break;
53 default:
54 break;
55 }
56
58}
59
62
63 switch (snes_ctx) {
64 case CTX_SNESNONE:
65 CHKERR VecAssemblyBegin(V);
66 CHKERR VecAssemblyEnd(V);
67 CHKERR VecSum(V, &eNergy);
68 break;
69 default:
70 break;
71 }
72
73 CHKERR VolumeElementForcesAndSourcesCore::postProcess();
74
76}
77
79 short int tag)
80 : feRhs(m_field), feLhs(m_field), feEnergy(m_field), mField(m_field),
81 tAg(tag) {}
82
84 const std::string field_name,
85 std::vector<VectorDouble> &values_at_gauss_pts,
86 std::vector<MatrixDouble> &gardient_at_gauss_pts)
89 valuesAtGaussPts(values_at_gauss_pts),
90 gradientAtGaussPts(gardient_at_gauss_pts), zeroAtType(MBVERTEX) {}
91
93 int side, EntityType type, EntitiesFieldData::EntData &data) {
95
96 const int nb_dofs = data.getFieldData().size();
97 const int nb_base_functions = data.getN().size2();
98 if (nb_dofs == 0) {
100 }
101 const int nb_gauss_pts = data.getN().size1();
102 const int rank = data.getFieldDofs()[0]->getNbOfCoeffs();
103
104 // initialize
105 if (type == zeroAtType) {
106 valuesAtGaussPts.resize(nb_gauss_pts);
107 gradientAtGaussPts.resize(nb_gauss_pts);
108 for (int gg = 0; gg != nb_gauss_pts; gg++) {
109 valuesAtGaussPts[gg].resize(rank, false);
110 gradientAtGaussPts[gg].resize(rank, 3, false);
111 }
112 for (int gg = 0; gg != nb_gauss_pts; gg++) {
113 valuesAtGaussPts[gg].clear();
114 gradientAtGaussPts[gg].clear();
115 }
116 }
117
118 auto base_function = data.getFTensor0N();
119 auto diff_base_functions = data.getFTensor1DiffN<3>();
120 FTensor::Index<'i', 3> i;
121 FTensor::Index<'j', 3> j;
122
123 if (rank == 1) {
124
125 for (int gg = 0; gg != nb_gauss_pts; gg++) {
126 auto field_data = data.getFTensor0FieldData();
127 double &val = valuesAtGaussPts[gg][0];
128 FTensor::Tensor1<double *, 3> grad(&gradientAtGaussPts[gg](0, 0),
129 &gradientAtGaussPts[gg](0, 1),
130 &gradientAtGaussPts[gg](0, 2));
131 int bb = 0;
132 for (; bb != nb_dofs; bb++) {
133 val += base_function * field_data;
134 grad(i) += diff_base_functions(i) * field_data;
135 ++diff_base_functions;
136 ++base_function;
137 ++field_data;
138 }
139 for (; bb != nb_base_functions; bb++) {
140 ++diff_base_functions;
141 ++base_function;
142 }
143 }
144
145 } else if (rank == 3) {
146
147 for (int gg = 0; gg != nb_gauss_pts; gg++) {
148 auto field_data = data.getFTensor1FieldData<3>();
149 FTensor::Tensor1<double *, 3> values(&valuesAtGaussPts[gg][0],
150 &valuesAtGaussPts[gg][1],
151 &valuesAtGaussPts[gg][2]);
153 &gradientAtGaussPts[gg](0, 0), &gradientAtGaussPts[gg](0, 1),
154 &gradientAtGaussPts[gg](0, 2), &gradientAtGaussPts[gg](1, 0),
155 &gradientAtGaussPts[gg](1, 1), &gradientAtGaussPts[gg](1, 2),
156 &gradientAtGaussPts[gg](2, 0), &gradientAtGaussPts[gg](2, 1),
157 &gradientAtGaussPts[gg](2, 2));
158 int bb = 0;
159 for (; bb != nb_dofs / 3; bb++) {
160 values(i) += base_function * field_data(i);
161 gradient(i, j) += field_data(i) * diff_base_functions(j);
162 ++diff_base_functions;
163 ++base_function;
164 ++field_data;
165 }
166 for (; bb != nb_base_functions; bb++) {
167 ++diff_base_functions;
168 ++base_function;
169 }
170 }
171
172 } else {
173 // FIXME: THat part is inefficient
174 VectorDouble &values = data.getFieldData();
175 for (int gg = 0; gg < nb_gauss_pts; gg++) {
176 VectorAdaptor N = data.getN(gg, nb_dofs / rank);
177 MatrixAdaptor diffN = data.getDiffN(gg, nb_dofs / rank);
178 for (int dd = 0; dd < nb_dofs / rank; dd++) {
179 for (int rr1 = 0; rr1 < rank; rr1++) {
180 valuesAtGaussPts[gg][rr1] += N[dd] * values[rank * dd + rr1];
181 for (int rr2 = 0; rr2 < 3; rr2++) {
182 gradientAtGaussPts[gg](rr1, rr2) +=
183 diffN(dd, rr2) * values[rank * dd + rr1];
184 }
185 }
186 }
187 }
188 }
189
191}
192
194 const std::string field_name, CommonData &common_data)
195 : OpGetDataAtGaussPts(field_name, common_data.dataAtGaussPts[field_name],
196 common_data.gradAtGaussPts[field_name]) {}
197
200 BlockData &data, CommonData &common_data,
201 int tag, bool jacobian, bool ale,
202 bool field_disp)
205 dAta(data), commonData(common_data), tAg(tag), adlocReturnValue(0),
206 jAcobian(jacobian), fUnction(!jacobian), aLe(ale), fieldDisp(field_disp) {
207
208}
209
212 const int gg) {
214
215 CHKERR dAta.materialAdoublePtr->calculateP_PiolaKirchhoffI(
216 dAta, getNumeredEntFiniteElementPtr());
217
218 if (aLe) {
219 auto &t_P = dAta.materialAdoublePtr->t_P;
220 auto &t_invH = dAta.materialAdoublePtr->t_invH;
221 t_P(i, j) = t_P(i, k) * t_invH(j, k);
222 t_P(i, j) *= dAta.materialAdoublePtr->detH;
223 }
224
225 commonData.sTress[gg].resize(3, 3, false);
226 for (int dd1 = 0; dd1 < 3; dd1++) {
227 for (int dd2 = 0; dd2 < 3; dd2++) {
228 dAta.materialAdoublePtr->P(dd1, dd2) >>=
229 (commonData.sTress[gg])(dd1, dd2);
230 }
231 }
232
234}
235
238 const int gg) {
240
241 trace_on(tAg, 0);
242
243 dAta.materialAdoublePtr->F.resize(3, 3, false);
244
245 if (!aLe) {
246
247 nbActiveVariables = 0;
248 for (int dd1 = 0; dd1 < 3; dd1++) {
249 for (int dd2 = 0; dd2 < 3; dd2++) {
250 dAta.materialAdoublePtr->F(dd1, dd2) <<= (*ptrh)[gg](dd1, dd2);
251 if (fieldDisp) {
252 if (dd1 == dd2) {
253 dAta.materialAdoublePtr->F(dd1, dd2) += 1;
254 }
255 }
256 nbActiveVariables++;
257 }
258 }
259
260 } else {
261
262 nbActiveVariables = 0;
263
264 dAta.materialAdoublePtr->h.resize(3, 3, false);
265 for (int dd1 = 0; dd1 < 3; dd1++) {
266 for (int dd2 = 0; dd2 < 3; dd2++) {
267 dAta.materialAdoublePtr->h(dd1, dd2) <<= (*ptrh)[gg](dd1, dd2);
268 nbActiveVariables++;
269 }
270 }
271
272 dAta.materialAdoublePtr->H.resize(3, 3, false);
273 for (int dd1 = 0; dd1 < 3; dd1++) {
274 for (int dd2 = 0; dd2 < 3; dd2++) {
275 dAta.materialAdoublePtr->H(dd1, dd2) <<= (*ptrH)[gg](dd1, dd2);
276 nbActiveVariables++;
277 }
278 }
279
280 dAta.materialAdoublePtr->detH = determinantTensor3by3(dAta.materialAdoublePtr->H);
281 dAta.materialAdoublePtr->invH.resize(3, 3, false);
282 CHKERR invertTensor3by3(dAta.materialAdoublePtr->H,
283 dAta.materialAdoublePtr->detH,
284 dAta.materialAdoublePtr->invH);
285
286 auto &t_F = dAta.materialAdoublePtr->t_F;
287 auto &t_h = dAta.materialAdoublePtr->t_h;
288 auto &t_invH = dAta.materialAdoublePtr->t_invH;
289
290 t_F(i, j) = t_h(i, k) * t_invH(k, j);
291
292 }
293
294 CHKERR dAta.materialAdoublePtr->setUserActiveVariables(nbActiveVariables);
295 CHKERR calculateStress(gg);
296
297 trace_off();
298
300}
301
305
306 int r;
307
308 if (fUnction) {
309 commonData.sTress[gg].resize(3, 3, false);
310 // play recorder for values
311 r = ::function(tAg, 9, nbActiveVariables, &activeVariables[0],
312 &commonData.sTress[gg](0, 0));
313 if (r < adlocReturnValue) { // function is locally analytic
314 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
315 "ADOL-C function evaluation with error r = %d", r);
316 }
317 }
318
319 if (jAcobian) {
320 commonData.jacStress[gg].resize(9, nbActiveVariables, false);
321 double *jac_ptr[] = {
322 &(commonData.jacStress[gg](0, 0)), &(commonData.jacStress[gg](1, 0)),
323 &(commonData.jacStress[gg](2, 0)), &(commonData.jacStress[gg](3, 0)),
324 &(commonData.jacStress[gg](4, 0)), &(commonData.jacStress[gg](5, 0)),
325 &(commonData.jacStress[gg](6, 0)), &(commonData.jacStress[gg](7, 0)),
326 &(commonData.jacStress[gg](8, 0))};
327 // play recorder for jacobians
328 r = jacobian(tAg, 9, nbActiveVariables, &activeVariables[0], jac_ptr);
329 if (r < adlocReturnValue) {
330 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
331 "ADOL-C function evaluation with error");
332 }
333 }
334
336}
337
339 int row_side, EntityType row_type,
340 EntitiesFieldData::EntData &row_data) {
342
343 // do it only once, no need to repeat this for edges,faces or tets
344 if (row_type != MBVERTEX)
346
347 if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
348 dAta.tEts.end()) {
350 }
351
352 int nb_dofs = row_data.getFieldData().size();
353 if (nb_dofs == 0)
355 dAta.materialAdoublePtr->commonDataPtr = &commonData;
356 dAta.materialAdoublePtr->opPtr = this;
357
358 int nb_gauss_pts = row_data.getN().size1();
359 commonData.sTress.resize(nb_gauss_pts);
360 commonData.jacStress.resize(nb_gauss_pts);
361
363 if (aLe) {
365 }
366
367 for (int gg = 0; gg != nb_gauss_pts; gg++) {
368
369 dAta.materialAdoublePtr->gG = gg;
370
371 // Record tag and calculate stress
372 if (recordTagForIntegrationPoint(gg)) {
373 CHKERR recordTag(gg);
374 }
375
376 // Set active variables vector
377 if (jAcobian || (!recordTagForIntegrationPoint(gg))) {
378 activeVariables.resize(nbActiveVariables, false);
379 if (!aLe) {
380 for (int dd1 = 0; dd1 < 3; dd1++) {
381 for (int dd2 = 0; dd2 < 3; dd2++) {
382 activeVariables(dd1 * 3 + dd2) = (*ptrh)[gg](dd1, dd2);
383 }
384 }
385 } else {
386 for (int dd1 = 0; dd1 < 3; dd1++) {
387 for (int dd2 = 0; dd2 < 3; dd2++) {
388 activeVariables(dd1 * 3 + dd2) = (*ptrh)[gg](dd1, dd2);
389 }
390 }
391 for (int dd1 = 0; dd1 < 3; dd1++) {
392 for (int dd2 = 0; dd2 < 3; dd2++) {
393 activeVariables(9 + dd1 * 3 + dd2) = (*ptrH)[gg](dd1, dd2);
394 }
395 }
396 }
397 CHKERR dAta.materialAdoublePtr->setUserActiveVariables(activeVariables);
398
399 // Play tag and calculate stress or tangent
400 if (jAcobian || (!recordTagForIntegrationPoint(gg))) {
401 CHKERR playTag(gg);
402 }
403 }
404 }
405
407}
408
410 const std::string
411 field_name, ///< field name for spatial positions or displacements
412 BlockData &data, CommonData &common_data, int tag, bool gradient,
413 bool hessian, bool ale, bool field_disp)
416 dAta(data), commonData(common_data), tAg(tag), gRadient(gradient),
417 hEssian(hessian), aLe(ale), fieldDisp(field_disp) {}
418
422 CHKERR dAta.materialAdoublePtr->calculateElasticEnergy(
423 dAta, getNumeredEntFiniteElementPtr());
424 dAta.materialAdoublePtr->eNergy >>= commonData.eNergy[gg];
426}
427
431
432 trace_on(tAg, 0);
433
434 if (!aLe) {
435
436 nbActiveVariables = 0;
437 for (int dd1 = 0; dd1 < 3; dd1++) {
438 for (int dd2 = 0; dd2 < 3; dd2++) {
439 dAta.materialAdoublePtr->F(dd1, dd2) <<= (*ptrh)[gg](dd1, dd2);
440 if (fieldDisp) {
441 if (dd1 == dd2) {
442 dAta.materialAdoublePtr->F(dd1, dd2) += 1;
443 }
444 }
445 nbActiveVariables++;
446 }
447 }
448
449 } else {
450
451 nbActiveVariables = 0;
452
453 dAta.materialAdoublePtr->h.resize(3, 3, false);
454 for (int dd1 = 0; dd1 < 3; dd1++) {
455 for (int dd2 = 0; dd2 < 3; dd2++) {
456 dAta.materialAdoublePtr->h(dd1, dd2) <<= (*ptrh)[gg](dd1, dd2);
457 nbActiveVariables++;
458 }
459 }
460
461 dAta.materialAdoublePtr->H.resize(3, 3, false);
462 for (int dd1 = 0; dd1 < 3; dd1++) {
463 for (int dd2 = 0; dd2 < 3; dd2++) {
464 dAta.materialAdoublePtr->H(dd1, dd2) <<= (*ptrH)[gg](dd1, dd2);
465 nbActiveVariables++;
466 }
467 }
468
469 dAta.materialAdoublePtr->detH = determinantTensor3by3(dAta.materialAdoublePtr->H);
470 dAta.materialAdoublePtr->invH.resize(3, 3, false);
471 CHKERR invertTensor3by3(dAta.materialAdoublePtr->H,
472 dAta.materialAdoublePtr->detH,
473 dAta.materialAdoublePtr->invH);
474
475 auto &t_F = dAta.materialAdoublePtr->t_F;
476 auto &t_h = dAta.materialAdoublePtr->t_h;
477 auto &t_invH = dAta.materialAdoublePtr->t_invH;
478
479 t_F(i, j) = t_h(i, k) * t_invH(k, j);
480
481 }
482
483 CHKERR dAta.materialAdoublePtr->setUserActiveVariables(nbActiveVariables);
485
486 trace_off();
487
489}
490
494
495 if (gRadient) {
496 commonData.jacEnergy[gg].resize(nbActiveVariables, false);
497 int r = ::gradient(tAg, nbActiveVariables, &activeVariables[0],
498 &commonData.jacEnergy[gg][0]);
499 if (r < 0) {
500 // That means that energy function is not smooth and derivative
501 // can not be calculated,
502 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
503 "ADOL-C function evaluation with error");
504 }
505 }
506
507 if (hEssian) {
508 commonData.hessianEnergy[gg].resize(nbActiveVariables * nbActiveVariables,
509 false);
510 double *H[nbActiveVariables];
511 for (int n = 0; n != nbActiveVariables; n++) {
512 H[n] = &(commonData.hessianEnergy[gg][n * nbActiveVariables]);
513 }
514 int r = ::hessian(tAg, nbActiveVariables, &*activeVariables.begin(), H);
515 if (r < 0) {
516 // That means that energy function is not smooth and derivative
517 // can not be calculated,
518 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
519 "ADOL-C function evaluation with error");
520 }
521 }
522
524}
525
527 int row_side, EntityType row_type,
528 EntitiesFieldData::EntData &row_data) {
530
531 // do it only once, no need to repeat this for edges,faces or tets
532 if (row_type != MBVERTEX)
534
535 if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
536 dAta.tEts.end()) {
538 }
539
540 int nb_dofs = row_data.getFieldData().size();
541 if (nb_dofs == 0)
543 dAta.materialAdoublePtr->commonDataPtr = &commonData;
544 dAta.materialAdoublePtr->opPtr = this;
545
546 int nb_gauss_pts = row_data.getN().size1();
547 commonData.eNergy.resize(nb_gauss_pts);
548 commonData.jacEnergy.resize(nb_gauss_pts);
549
551 if (aLe) {
553 }
554
555 for (int gg = 0; gg != nb_gauss_pts; gg++) {
556
557 dAta.materialAdoublePtr->gG = gg;
558
559 // Record tag and calualte stress
560 if (recordTagForIntegrationPoint(gg)) {
561 CHKERR recordTag(gg);
562 }
563
564 activeVariables.resize(nbActiveVariables, false);
565 if (!aLe) {
566 for (int dd1 = 0; dd1 < 3; dd1++) {
567 for (int dd2 = 0; dd2 < 3; dd2++) {
568 activeVariables(dd1 * 3 + dd2) = (*ptrh)[gg](dd1, dd2);
569 }
570 }
571 } else {
572 for (int dd1 = 0; dd1 < 3; dd1++) {
573 for (int dd2 = 0; dd2 < 3; dd2++) {
574 activeVariables(dd1 * 3 + dd2) = (*ptrh)[gg](dd1, dd2);
575 }
576 }
577 for (int dd1 = 0; dd1 < 3; dd1++) {
578 for (int dd2 = 0; dd2 < 3; dd2++) {
579 activeVariables(9 + dd1 * 3 + dd2) = (*ptrH)[gg](dd1, dd2);
580 }
581 }
582 }
583 CHKERR dAta.materialAdoublePtr->setUserActiveVariables(activeVariables);
584
585 // Play tag and calculate stress or tangent
586 CHKERR playTag(gg);
587 }
588
590}
591
597
599 int row_side, EntityType row_type,
600 EntitiesFieldData::EntData &row_data) {
602
603 int nb_dofs = row_data.getIndices().size();
604 int *indices_ptr = &row_data.getIndices()[0];
605 if (!dAta.forcesOnlyOnEntitiesRow.empty()) {
606 iNdices.resize(nb_dofs, false);
607 noalias(iNdices) = row_data.getIndices();
608 indices_ptr = &iNdices[0];
609 VectorDofs &dofs = row_data.getFieldDofs();
610 VectorDofs::iterator dit = dofs.begin();
611 for (int ii = 0; dit != dofs.end(); dit++, ii++) {
612 if (dAta.forcesOnlyOnEntitiesRow.find((*dit)->getEnt()) ==
613 dAta.forcesOnlyOnEntitiesRow.end()) {
614 iNdices[ii] = -1;
615 }
616 }
617 }
618 CHKERR VecSetValues(getFEMethod()->snes_f, nb_dofs, indices_ptr, &nf[0],
619 ADD_VALUES);
621}
622
624 int row_side, EntityType row_type,
625 EntitiesFieldData::EntData &row_data) {
627
628 if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
629 dAta.tEts.end()) {
631 }
632
633 const int nb_dofs = row_data.getIndices().size();
634 if (nb_dofs == 0)
636 if ((unsigned int)nb_dofs > 3 * row_data.getN().size2()) {
637 SETERRQ(PETSC_COMM_SELF, 1, "data inconsistency");
638 }
639 const int nb_base_functions = row_data.getN().size2();
640 const int nb_gauss_pts = row_data.getN().size1();
641
642 nf.resize(nb_dofs, false);
643 nf.clear();
644
645 auto diff_base_functions = row_data.getFTensor1DiffN<3>();
646 FTensor::Index<'i', 3> i;
647 FTensor::Index<'j', 3> j;
648
649 for (int gg = 0; gg != nb_gauss_pts; gg++) {
650 double val = getVolume() * getGaussPts()(3, gg);
651 MatrixDouble3by3 &stress = commonData.sTress[gg];
653 &stress(0, 0), &stress(0, 1), &stress(0, 2), &stress(1, 0),
654 &stress(1, 1), &stress(1, 2), &stress(2, 0), &stress(2, 1),
655 &stress(2, 2));
656 FTensor::Tensor1<double *, 3> rhs(&nf[0], &nf[1], &nf[2], 3);
657 int bb = 0;
658 for (; bb != nb_dofs / 3; bb++) {
659 rhs(i) += val * t3(i, j) * diff_base_functions(j);
660 ++rhs;
661 ++diff_base_functions;
662 }
663 for (; bb != nb_base_functions; bb++) {
664 ++diff_base_functions;
665 }
666 }
667
668 CHKERR aSemble(row_side, row_type, row_data);
669
671}
672
674 BlockData &data,
675 CommonData &common_data,
676 SmartPetscObj<Vec> ghost_vec,
677 bool field_disp)
680 dAta(data), commonData(common_data), ghostVec(ghost_vec, true),
681 fieldDisp(field_disp) {}
682
684 int row_side, EntityType row_type,
685 EntitiesFieldData::EntData &row_data) {
687
688 if (row_type != MBVERTEX)
690 if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
691 dAta.tEts.end()) {
693 }
694
695 std::vector<MatrixDouble> &F =
697 dAta.materialDoublePtr->F.resize(3, 3, false);
698
699 double energy = 0;
700
701 for (unsigned int gg = 0; gg != row_data.getN().size1(); ++gg) {
702 double val = getVolume() * getGaussPts()(3, gg);
703 noalias(dAta.materialDoublePtr->F) = F[gg];
704 if (fieldDisp) {
705 for (int dd = 0; dd < 3; dd++) {
706 dAta.materialDoublePtr->F(dd, dd) += 1;
707 }
708 }
709 int nb_active_variables = 0;
710 CHKERR dAta.materialDoublePtr->setUserActiveVariables(nb_active_variables);
711 CHKERR dAta.materialDoublePtr->calculateElasticEnergy(
712 dAta, getNumeredEntFiniteElementPtr());
713 energy += val * dAta.materialDoublePtr->eNergy;
714 }
715
716 CHKERR VecSetValue(ghostVec, 0, energy, ADD_VALUES);
718}
719
721 const std::string vel_field, const std::string field_name, BlockData &data,
722 CommonData &common_data)
724 vel_field, field_name, UserDataOperator::OPROWCOL),
725 dAta(data), commonData(common_data), aLe(false) {}
726
727template <int S>
729 int gg, MatrixDouble &jac_stress,
730 MatrixDouble &jac) {
732 jac.clear();
733 FTensor::Index<'i', 3> i;
734 FTensor::Index<'j', 3> j;
735 FTensor::Index<'k', 3> k;
736 int nb_col = col_data.getFieldData().size();
737 double *diff_ptr =
738 const_cast<double *>(&(col_data.getDiffN(gg, nb_col / 3)(0, 0)));
739 // First two indices 'i','j' derivatives of 1st Piola-stress, third index 'k'
740 // is displacement component
742 &jac_stress(3 * 0 + 0, S + 0), &jac_stress(3 * 0 + 0, S + 1),
743 &jac_stress(3 * 0 + 0, S + 2), &jac_stress(3 * 0 + 1, S + 0),
744 &jac_stress(3 * 0 + 1, S + 1), &jac_stress(3 * 0 + 1, S + 2),
745 &jac_stress(3 * 0 + 2, S + 0), &jac_stress(3 * 0 + 2, S + 1),
746 &jac_stress(3 * 0 + 2, S + 2), &jac_stress(3 * 1 + 0, S + 0),
747 &jac_stress(3 * 1 + 0, S + 1), &jac_stress(3 * 1 + 0, S + 2),
748 &jac_stress(3 * 1 + 1, S + 0), &jac_stress(3 * 1 + 1, S + 1),
749 &jac_stress(3 * 1 + 1, S + 2), &jac_stress(3 * 1 + 2, S + 0),
750 &jac_stress(3 * 1 + 2, S + 1), &jac_stress(3 * 1 + 2, S + 2),
751 &jac_stress(3 * 2 + 0, S + 0), &jac_stress(3 * 2 + 0, S + 1),
752 &jac_stress(3 * 2 + 0, S + 2), &jac_stress(3 * 2 + 1, S + 0),
753 &jac_stress(3 * 2 + 1, S + 1), &jac_stress(3 * 2 + 1, S + 2),
754 &jac_stress(3 * 2 + 2, S + 0), &jac_stress(3 * 2 + 2, S + 1),
755 &jac_stress(3 * 2 + 2, S + 2));
757 &jac_stress(3 * 0 + 0, S + 3), &jac_stress(3 * 0 + 0, S + 4),
758 &jac_stress(3 * 0 + 0, S + 5), &jac_stress(3 * 0 + 1, S + 3),
759 &jac_stress(3 * 0 + 1, S + 4), &jac_stress(3 * 0 + 1, S + 5),
760 &jac_stress(3 * 0 + 2, S + 3), &jac_stress(3 * 0 + 2, S + 4),
761 &jac_stress(3 * 0 + 2, S + 5), &jac_stress(3 * 1 + 0, S + 3),
762 &jac_stress(3 * 1 + 0, S + 4), &jac_stress(3 * 1 + 0, S + 5),
763 &jac_stress(3 * 1 + 1, S + 3), &jac_stress(3 * 1 + 1, S + 4),
764 &jac_stress(3 * 1 + 1, S + 5), &jac_stress(3 * 1 + 2, S + 3),
765 &jac_stress(3 * 1 + 2, S + 4), &jac_stress(3 * 1 + 2, S + 5),
766 &jac_stress(3 * 2 + 0, S + 3), &jac_stress(3 * 2 + 0, S + 4),
767 &jac_stress(3 * 2 + 0, S + 5), &jac_stress(3 * 2 + 1, S + 3),
768 &jac_stress(3 * 2 + 1, S + 4), &jac_stress(3 * 2 + 1, S + 5),
769 &jac_stress(3 * 2 + 2, S + 3), &jac_stress(3 * 2 + 2, S + 4),
770 &jac_stress(3 * 2 + 2, S + 5));
772 &jac_stress(3 * 0 + 0, S + 6), &jac_stress(3 * 0 + 0, S + 7),
773 &jac_stress(3 * 0 + 0, S + 8), &jac_stress(3 * 0 + 1, S + 6),
774 &jac_stress(3 * 0 + 1, S + 7), &jac_stress(3 * 0 + 1, S + 8),
775 &jac_stress(3 * 0 + 2, S + 6), &jac_stress(3 * 0 + 2, S + 7),
776 &jac_stress(3 * 0 + 2, S + 8), &jac_stress(3 * 1 + 0, S + 6),
777 &jac_stress(3 * 1 + 0, S + 7), &jac_stress(3 * 1 + 0, S + 8),
778 &jac_stress(3 * 1 + 1, S + 6), &jac_stress(3 * 1 + 1, S + 7),
779 &jac_stress(3 * 1 + 1, S + 8), &jac_stress(3 * 1 + 2, S + 6),
780 &jac_stress(3 * 1 + 2, S + 7), &jac_stress(3 * 1 + 2, S + 8),
781 &jac_stress(3 * 2 + 0, S + 6), &jac_stress(3 * 2 + 0, S + 7),
782 &jac_stress(3 * 2 + 0, S + 8), &jac_stress(3 * 2 + 1, S + 6),
783 &jac_stress(3 * 2 + 1, S + 7), &jac_stress(3 * 2 + 1, S + 8),
784 &jac_stress(3 * 2 + 2, S + 6), &jac_stress(3 * 2 + 2, S + 7),
785 &jac_stress(3 * 2 + 2, S + 8));
786 // Derivate of 1st Piola-stress multiplied by gradient of defamation for
787 // base function (dd) and displacement component (rr)
789 &jac(0, 0), &jac(1, 0), &jac(2, 0), &jac(3, 0), &jac(4, 0), &jac(5, 0),
790 &jac(6, 0), &jac(7, 0), &jac(8, 0));
792 &jac(0, 1), &jac(1, 1), &jac(2, 1), &jac(3, 1), &jac(4, 1), &jac(5, 1),
793 &jac(6, 1), &jac(7, 1), &jac(8, 1));
795 &jac(0, 2), &jac(1, 2), &jac(2, 2), &jac(3, 2), &jac(4, 2), &jac(5, 2),
796 &jac(6, 2), &jac(7, 2), &jac(8, 2));
798 diff_ptr, &diff_ptr[1], &diff_ptr[2]);
799 for (int dd = 0; dd != nb_col / 3; ++dd) {
800 t2_1_0(i, j) += t3_1_0(i, j, k) * diff(k);
801 t2_1_1(i, j) += t3_1_1(i, j, k) * diff(k);
802 t2_1_2(i, j) += t3_1_2(i, j, k) * diff(k);
803 ++t2_1_0;
804 ++t2_1_1;
805 ++t2_1_2;
806 ++diff;
807 }
809}
810
812 EntitiesFieldData::EntData &col_data, int gg) {
813 return get_jac<0>(col_data, gg, commonData.jacStress[gg], jac);
814}
815
817 int row_side, int col_side, EntityType row_type, EntityType col_type,
819 EntitiesFieldData::EntData &col_data) {
821
822 int nb_row = row_data.getIndices().size();
823 int nb_col = col_data.getIndices().size();
824
825 int *row_indices_ptr = &row_data.getIndices()[0];
826 int *col_indices_ptr = &col_data.getIndices()[0];
827
828 if (!dAta.forcesOnlyOnEntitiesRow.empty()) {
829 rowIndices.resize(nb_row, false);
830 noalias(rowIndices) = row_data.getIndices();
831 row_indices_ptr = &rowIndices[0];
832 VectorDofs &dofs = row_data.getFieldDofs();
833 VectorDofs::iterator dit = dofs.begin();
834 for (int ii = 0; dit != dofs.end(); dit++, ii++) {
835 if (dAta.forcesOnlyOnEntitiesRow.find((*dit)->getEnt()) ==
836 dAta.forcesOnlyOnEntitiesRow.end()) {
837 rowIndices[ii] = -1;
838 }
839 }
840 }
841
842 if (!dAta.forcesOnlyOnEntitiesCol.empty()) {
843 colIndices.resize(nb_col, false);
844 noalias(colIndices) = col_data.getIndices();
845 col_indices_ptr = &colIndices[0];
846 VectorDofs &dofs = col_data.getFieldDofs();
847 VectorDofs::iterator dit = dofs.begin();
848 for (int ii = 0; dit != dofs.end(); dit++, ii++) {
849 if (dAta.forcesOnlyOnEntitiesCol.find((*dit)->getEnt()) ==
850 dAta.forcesOnlyOnEntitiesCol.end()) {
851 colIndices[ii] = -1;
852 }
853 }
854 }
855
856 CHKERR MatSetValues(getFEMethod()->snes_B, nb_row, row_indices_ptr, nb_col,
857 col_indices_ptr, &k(0, 0), ADD_VALUES);
858
859 // is symmetric
860 if (row_side != col_side || row_type != col_type) {
861
862 row_indices_ptr = &row_data.getIndices()[0];
863 col_indices_ptr = &col_data.getIndices()[0];
864
865 if (!dAta.forcesOnlyOnEntitiesCol.empty()) {
866 rowIndices.resize(nb_row, false);
867 noalias(rowIndices) = row_data.getIndices();
868 row_indices_ptr = &rowIndices[0];
869 VectorDofs &dofs = row_data.getFieldDofs();
870 VectorDofs::iterator dit = dofs.begin();
871 for (int ii = 0; dit != dofs.end(); dit++, ii++) {
872 if (dAta.forcesOnlyOnEntitiesCol.find((*dit)->getEnt()) ==
873 dAta.forcesOnlyOnEntitiesCol.end()) {
874 rowIndices[ii] = -1;
875 }
876 }
877 }
878
879 if (!dAta.forcesOnlyOnEntitiesRow.empty()) {
880 colIndices.resize(nb_col, false);
881 noalias(colIndices) = col_data.getIndices();
882 col_indices_ptr = &colIndices[0];
883 VectorDofs &dofs = col_data.getFieldDofs();
884 VectorDofs::iterator dit = dofs.begin();
885 for (int ii = 0; dit != dofs.end(); dit++, ii++) {
886 if (dAta.forcesOnlyOnEntitiesRow.find((*dit)->getEnt()) ==
887 dAta.forcesOnlyOnEntitiesRow.end()) {
888 colIndices[ii] = -1;
889 }
890 }
891 }
892
893 trans_k.resize(nb_col, nb_row, false);
894 noalias(trans_k) = trans(k);
895 CHKERR MatSetValues(getFEMethod()->snes_B, nb_col, col_indices_ptr, nb_row,
896 row_indices_ptr, &trans_k(0, 0), ADD_VALUES);
897 }
898
900}
901
903 int row_side, int col_side, EntityType row_type, EntityType col_type,
905 EntitiesFieldData::EntData &col_data) {
907
908 int nb_row = row_data.getIndices().size();
909 int nb_col = col_data.getIndices().size();
910 if (nb_row == 0)
912 if (nb_col == 0)
914
915 if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
916 dAta.tEts.end()) {
918 }
919
920 // const int nb_base_functions = row_data.getN().size2();
921 const int nb_gauss_pts = row_data.getN().size1();
922
923 FTensor::Index<'i', 3> i;
924 FTensor::Index<'j', 3> j;
925 FTensor::Index<'m', 3> m;
926
927 k.resize(nb_row, nb_col, false);
928 k.clear();
929 jac.resize(9, nb_col, false);
930
931 for (int gg = 0; gg != nb_gauss_pts; gg++) {
932 CHKERR getJac(col_data, gg);
933 double val = getVolume() * getGaussPts()(3, gg);
935 &jac(3 * 0 + 0, 0), &jac(3 * 0 + 0, 1), &jac(3 * 0 + 0, 2),
936 &jac(3 * 0 + 1, 0), &jac(3 * 0 + 1, 1), &jac(3 * 0 + 1, 2),
937 &jac(3 * 0 + 2, 0), &jac(3 * 0 + 2, 1), &jac(3 * 0 + 2, 2),
938 &jac(3 * 1 + 0, 0), &jac(3 * 1 + 0, 1), &jac(3 * 1 + 0, 2),
939 &jac(3 * 1 + 1, 0), &jac(3 * 1 + 1, 1), &jac(3 * 1 + 1, 2),
940 &jac(3 * 1 + 2, 0), &jac(3 * 1 + 2, 1), &jac(3 * 1 + 2, 2),
941 &jac(3 * 2 + 0, 0), &jac(3 * 2 + 0, 1), &jac(3 * 2 + 0, 2),
942 &jac(3 * 2 + 1, 0), &jac(3 * 2 + 1, 1), &jac(3 * 2 + 1, 2),
943 &jac(3 * 2 + 2, 0), &jac(3 * 2 + 2, 1), &jac(3 * 2 + 2, 2));
944 for (int cc = 0; cc != nb_col / 3; cc++) {
945 auto diff_base_functions = row_data.getFTensor1DiffN<3>(gg, 0);
947 &k(0, 3 * cc + 0), &k(0, 3 * cc + 1), &k(0, 3 * cc + 2),
948 &k(1, 3 * cc + 0), &k(1, 3 * cc + 1), &k(1, 3 * cc + 2),
949 &k(2, 3 * cc + 0), &k(2, 3 * cc + 1), &k(2, 3 * cc + 2), 3 * nb_col);
950 for (int rr = 0; rr != nb_row / 3; rr++) {
951 lhs(i, j) += val * t3_1(i, m, j) * diff_base_functions(m);
952 ++diff_base_functions;
953 ++lhs;
954 }
955 ++t3_1;
956 }
957 }
958
959 CHKERR aSemble(row_side, col_side, row_type, col_type, row_data, col_data);
960
962}
963
965 const std::string vel_field, const std::string field_name, BlockData &data,
966 CommonData &common_data)
967 : OpLhsPiolaKirchhoff_dx(vel_field, field_name, data, common_data) {
968 sYmm = false;
969}
970
972 EntitiesFieldData::EntData &col_data, int gg) {
973 return get_jac<9>(col_data, gg, commonData.jacStress[gg], jac);
974}
975
977 int row_side, int col_side, EntityType row_type, EntityType col_type,
979 EntitiesFieldData::EntData &col_data) {
981
982 int nb_row = row_data.getIndices().size();
983 int nb_col = col_data.getIndices().size();
984
985 int *row_indices_ptr = &row_data.getIndices()[0];
986 if (!dAta.forcesOnlyOnEntitiesRow.empty()) {
987 rowIndices.resize(nb_row, false);
988 noalias(rowIndices) = row_data.getIndices();
989 row_indices_ptr = &rowIndices[0];
990 VectorDofs &dofs = row_data.getFieldDofs();
991 VectorDofs::iterator dit = dofs.begin();
992 for (int ii = 0; dit != dofs.end(); dit++, ii++) {
993 if (dAta.forcesOnlyOnEntitiesRow.find((*dit)->getEnt()) ==
994 dAta.forcesOnlyOnEntitiesRow.end()) {
995 rowIndices[ii] = -1;
996 }
997 }
998 }
999
1000 int *col_indices_ptr = &col_data.getIndices()[0];
1001 if (!dAta.forcesOnlyOnEntitiesCol.empty()) {
1002 colIndices.resize(nb_col, false);
1003 noalias(colIndices) = col_data.getIndices();
1004 col_indices_ptr = &colIndices[0];
1005 VectorDofs &dofs = col_data.getFieldDofs();
1006 VectorDofs::iterator dit = dofs.begin();
1007 for (int ii = 0; dit != dofs.end(); dit++, ii++) {
1008 if (dAta.forcesOnlyOnEntitiesCol.find((*dit)->getEnt()) ==
1009 dAta.forcesOnlyOnEntitiesCol.end()) {
1010 colIndices[ii] = -1;
1011 }
1012 }
1013 }
1014
1015 /*for(int dd1 = 0;dd1<k.size1();dd1++) {
1016 for(int dd2 = 0;dd2<k.size2();dd2++) {
1017 if(k(dd1,dd2)!=k(dd1,dd2)) {
1018 SETERRQ(PETSC_COMM_SELF,1,"Wrong result");
1019 }
1020 }
1021 }*/
1022
1023 CHKERR MatSetValues(getFEMethod()->snes_B, nb_row, row_indices_ptr, nb_col,
1024 col_indices_ptr, &k(0, 0), ADD_VALUES);
1025
1027}
1028
1030 const std::string field_name, BlockData &data, CommonData &common_data,
1031 int tag, bool jacobian, bool ale)
1032 : OpJacobianPiolaKirchhoffStress(field_name, data, common_data, tag,
1033 jacobian, ale, false) {}
1034
1037 const int gg) {
1039
1040 CHKERR dAta.materialAdoublePtr->calculatesIGma_EshelbyStress(
1041 dAta, getNumeredEntFiniteElementPtr());
1042 if (aLe) {
1043 auto &t_sIGma = dAta.materialAdoublePtr->t_sIGma;
1044 auto &t_invH = dAta.materialAdoublePtr->t_invH;
1045 t_sIGma(i, j) = t_sIGma(i, k) * t_invH(j, k);
1046 t_sIGma(i, j) *= dAta.materialAdoublePtr->detH;
1047
1048 }
1049 commonData.sTress[gg].resize(3, 3, false);
1050 for (int dd1 = 0; dd1 < 3; dd1++) {
1051 for (int dd2 = 0; dd2 < 3; dd2++) {
1052 dAta.materialAdoublePtr->sIGma(dd1, dd2) >>=
1053 (commonData.sTress[gg])(dd1, dd2);
1054 }
1055 }
1056
1058}
1059
1063
1065 const std::string vel_field, const std::string field_name, BlockData &data,
1066 CommonData &common_data)
1067 : OpLhsPiolaKirchhoff_dX(vel_field, field_name, data, common_data) {}
1068
1070 EntitiesFieldData::EntData &col_data, int gg) {
1071 return get_jac<0>(col_data, gg, commonData.jacStress[gg], jac);
1072}
1073
1075 const std::string vel_field, const std::string field_name, BlockData &data,
1076 CommonData &common_data)
1077 : OpLhsPiolaKirchhoff_dx(vel_field, field_name, data, common_data) {}
1078
1080 EntitiesFieldData::EntData &col_data, int gg) {
1081 return get_jac<9>(col_data, gg, commonData.jacStress[gg], jac);
1082}
1083
1086 materialDoublePtr,
1088 materialAdoublePtr) {
1090
1091 if (!materialDoublePtr) {
1093 "Pointer for materialDoublePtr not allocated");
1094 }
1095 if (!materialAdoublePtr) {
1097 "Pointer for materialAdoublePtr not allocated");
1098 }
1099
1101 mField, BLOCKSET | MAT_ELASTICSET, it)) {
1102 Mat_Elastic mydata;
1103 CHKERR it->getAttributeDataStructure(mydata);
1104 int id = it->getMeshsetId();
1105 EntityHandle meshset = it->getMeshset();
1106 CHKERR mField.get_moab().get_entities_by_type(meshset, MBTET,
1107 setOfBlocks[id].tEts, true);
1108 setOfBlocks[id].iD = id;
1109 setOfBlocks[id].E = mydata.data.Young;
1110 setOfBlocks[id].PoissonRatio = mydata.data.Poisson;
1111 setOfBlocks[id].materialDoublePtr = materialDoublePtr;
1112 setOfBlocks[id].materialAdoublePtr = materialAdoublePtr;
1113 }
1114
1116}
1117
1119 const std::string element_name,
1120 const std::string spatial_position_field_name,
1121 const std::string material_position_field_name, const bool ale) {
1123
1124 CHKERR mField.add_finite_element(element_name, MF_ZERO);
1126 element_name, spatial_position_field_name);
1128 element_name, spatial_position_field_name);
1130 element_name, spatial_position_field_name);
1131 if (mField.check_field(material_position_field_name)) {
1132 if (ale) {
1134 element_name, material_position_field_name);
1136 element_name, material_position_field_name);
1137 }
1139 element_name, material_position_field_name);
1140 }
1141
1142 std::map<int, BlockData>::iterator sit = setOfBlocks.begin();
1143 for (; sit != setOfBlocks.end(); sit++) {
1144 CHKERR mField.add_ents_to_finite_element_by_type(sit->second.tEts, MBTET,
1145 element_name);
1146 }
1147
1149}
1150
1152 const std::string spatial_position_field_name,
1153 const std::string material_position_field_name, const bool ale,
1154 const bool field_disp) {
1156
1157 commonData.spatialPositions = spatial_position_field_name;
1158 commonData.meshPositions = material_position_field_name;
1159
1160 // Rhs
1161 feRhs.getOpPtrVector().push_back(
1162 new OpGetCommonDataAtGaussPts(spatial_position_field_name, commonData));
1163 if (mField.check_field(material_position_field_name)) {
1165 material_position_field_name, commonData));
1166 }
1167 std::map<int, BlockData>::iterator sit = setOfBlocks.begin();
1168 for (; sit != setOfBlocks.end(); sit++) {
1170 spatial_position_field_name, sit->second, commonData, tAg, false, ale,
1171 field_disp));
1173 spatial_position_field_name, sit->second, commonData));
1174 }
1175
1176 // Energy
1177 feEnergy.getOpPtrVector().push_back(
1178 new OpGetCommonDataAtGaussPts(spatial_position_field_name, commonData));
1179 if (mField.check_field(material_position_field_name)) {
1181 material_position_field_name, commonData));
1182 }
1183 sit = setOfBlocks.begin();
1184 for (; sit != setOfBlocks.end(); sit++) {
1185 feEnergy.getOpPtrVector().push_back(
1186 new OpEnergy(spatial_position_field_name, sit->second, commonData,
1187 feEnergy.V, field_disp));
1188 }
1189
1190 // Lhs
1191 feLhs.getOpPtrVector().push_back(
1192 new OpGetCommonDataAtGaussPts(spatial_position_field_name, commonData));
1193 if (mField.check_field(material_position_field_name)) {
1195 material_position_field_name, commonData));
1196 }
1197 sit = setOfBlocks.begin();
1198 for (; sit != setOfBlocks.end(); sit++) {
1200 spatial_position_field_name, sit->second, commonData, tAg, true, ale,
1201 field_disp));
1203 spatial_position_field_name, spatial_position_field_name, sit->second,
1204 commonData));
1205 }
1206
1208}
static MoFEMErrorCode calculateEnergy(DM dm, boost::shared_ptr< map< int, BlockData > > block_sets_ptr, const std::string x_field, const std::string X_field, const bool ale, const bool field_disp, SmartPetscObj< Vec > &v_energy_ptr)
static MoFEMErrorCode get_jac(EntitiesFieldData::EntData &col_data, int gg, MatrixDouble &jac_stress, MatrixDouble &jac)
Operators and data structures for non-linear elastic analysis.
@ MF_ZERO
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MAT_ELASTICSET
block name is "MAT_ELASTIC"
@ BLOCKSET
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition definitions.h:34
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
constexpr int order
@ F
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string name, const bool recursive=true)=0
add entities to finite element
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
virtual bool check_field(const std::string &name) const =0
check if field is in database
#define _IT_CUBITMESHSETS_BY_BCDATA_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet in a moFEM field.
const double n
refractive index of diffusive medium
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition Types.hpp:115
MatrixShallowArrayAdaptor< double > MatrixAdaptor
Matrix adaptor.
Definition Types.hpp:132
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
MoFEMErrorCode invertTensor3by3(ublas::matrix< T, L, A > &jac_data, ublas::vector< T, A > &det_data, ublas::matrix< T, L, A > &inv_jac_data)
Calculate inverse of tensor rank 2 at integration points.
auto createVectorMPI(MPI_Comm comm, PetscInt n, PetscInt N)
Create MPI Vector.
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
static auto determinantTensor3by3(T &t)
Calculate the determinant of a 3x3 matrix or a tensor of rank 2.
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
ublas::vector< FEDofEntity *, DofsAllocator > VectorDofs
constexpr AssemblyType A
constexpr auto field_name
FTensor::Index< 'm', 3 > m
const int N
Definition speed_test.cpp:3
virtual moab::Interface & get_moab()=0
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =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 & getDiffN(const FieldApproximationBase base)
get derivatives of base functions
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorDouble & getFieldData() const
Get DOF values on entity.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1FieldData()
Return FTensor of rank 1, i.e. vector from field data coefficients.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0FieldData()
Return scalar files as a FTensor of rank 0.
const VectorDofs & getFieldDofs() const
Get DOF data structures (const version)
const VectorInt & getIndices() const
Get global indices of degrees of freedom on entity.
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
Elastic material data structure.
intrusive_ptr for managing petsc objects
FTensor::Tensor2< double *, 3, 3 > & getJac()
get element Jacobian
data for calculation heat conductivity and heat capacity elements
common data used by volume elements
std::vector< MatrixDouble > jacStress
this is simply material tangent operator
std::map< std::string, std::vector< MatrixDouble > > gradAtGaussPts
Implementation of elastic (non-linear) St. Kirchhoff equation.
MoFEMErrorCode preProcess()
Pre-processing function executed at loop initialization.
MoFEMErrorCode postProcess()
Post-processing function executed at loop completion.
int getRule(int order)
it is used to calculate nb. of Gauss integration points
OpEnergy(const std::string field_name, BlockData &data, CommonData &common_data, SmartPetscObj< Vec > ghost_vec, bool field_disp)
MoFEMErrorCode doWork(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
OpGetCommonDataAtGaussPts(const std::string field_name, CommonData &common_data)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
operator calculating deformation gradient
OpGetDataAtGaussPts(const std::string field_name, std::vector< VectorDouble > &values_at_gauss_pts, std::vector< MatrixDouble > &gradient_at_gauss_pts)
virtual MoFEMErrorCode playTag(const int gg)
Play ADOL-C tape.
OpJacobianEnergy(const std::string field_name, BlockData &data, CommonData &common_data, int tag, bool gradient, bool hessian, bool ale, bool field_disp)
virtual MoFEMErrorCode recordTag(const int gg)
Record ADOL-C tape.
MoFEMErrorCode doWork(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
virtual MoFEMErrorCode calculateEnergy(const int gg)
Calculate Paola-Kirchhoff I stress.
MoFEMErrorCode calculateStress(const int gg)
Calculate Paola-Kirchhoff I stress.
OpJacobianEshelbyStress(const std::string field_name, BlockData &data, CommonData &common_data, int tag, bool jacobian, bool ale)
virtual MoFEMErrorCode recordTag(const int gg)
Record ADOL-C tape.
virtual MoFEMErrorCode playTag(const int gg)
Play ADOL-C tape.
OpJacobianPiolaKirchhoffStress(const std::string field_name, BlockData &data, CommonData &common_data, int tag, bool jacobian, bool ale, bool field_disp)
Construct operator to calculate Piola-Kirchhoff stress or its derivatives over gradient deformation.
MoFEMErrorCode doWork(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
Calculate stress or jacobian at gauss points.
virtual MoFEMErrorCode calculateStress(const int gg)
Calculate Paola-Kirchhoff I stress.
OpLhsEshelby_dX(const std::string vel_field, const std::string field_name, BlockData &data, CommonData &common_data)
OpLhsEshelby_dx(const std::string vel_field, const std::string field_name, BlockData &data, CommonData &common_data)
OpLhsPiolaKirchhoff_dX(const std::string vel_field, const std::string field_name, BlockData &data, CommonData &common_data)
MoFEMErrorCode aSemble(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
virtual MoFEMErrorCode aSemble(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
OpLhsPiolaKirchhoff_dx(const std::string vel_field, const std::string field_name, BlockData &data, CommonData &common_data)
OpRhsEshelbyStress(const std::string field_name, BlockData &data, CommonData &common_data)
virtual MoFEMErrorCode aSemble(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
MoFEMErrorCode doWork(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
OpRhsPiolaKirchhoff(const std::string field_name, BlockData &data, CommonData &common_data)
MyVolumeFE feRhs
calculate right hand side for tetrahedral elements
NonlinearElasticElement(MoFEM::Interface &m_field, short int tag)
MoFEMErrorCode addElement(const std::string element_name, const std::string spatial_position_field_name, const std::string material_position_field_name="MESH_NODE_POSITIONS", const bool ale=false)
std::map< int, BlockData > setOfBlocks
maps block set id with appropriate BlockData
MoFEMErrorCode setBlocks(boost::shared_ptr< FunctionsToCalculatePiolaKirchhoffI< double > > materialDoublePtr, boost::shared_ptr< FunctionsToCalculatePiolaKirchhoffI< adouble > > materialAdoublePtr)
MyVolumeFE feEnergy
calculate elastic energy
MoFEMErrorCode setOperators(const std::string spatial_position_field_name, const std::string material_position_field_name="MESH_NODE_POSITIONS", const bool ale=false, const bool field_disp=false)
Set operators to calculate left hand tangent matrix and right hand residual.
double H
Hardening.
Definition plastic.cpp:128