v0.14.0
Loading...
Searching...
No Matches
EshelbianPlasticity.hpp
Go to the documentation of this file.
1/**
2 * \file EshelbianPlasticity.hpp
3 * \brief Eshelbian plasticity interface
4 *
5 * \brief Problem implementation for mix element for large-strain elasticity
6 *
7 * \todo Implementation of plasticity
8 */
9
10constexpr int SPACE_DIM = 3;
11constexpr AssemblyType A = AssemblyType::SCHUR;
12
13#ifndef __ESHELBIAN_PLASTICITY_HPP__
14#define __ESHELBIAN_PLASTICITY_HPP__
15
16namespace EshelbianPlasticity {
17
19inline static enum RotSelector rotSelector = LARGE_ROT;
20inline static double exponentBase = exp(1);
22
23inline double f_log(const double v) { return pow(exponentBase, v); }
24inline double d_f_log(const double v) {
25 return pow(exponentBase, v) * log(exponentBase);
26}
27inline double dd_f_log(const double v) {
28 return pow(exponentBase, v) * log(exponentBase) * log(exponentBase);
29}
30
31inline double f_linear(const double v) { return 1 + v; }
32inline double d_f_linear(const double v) { return 1; }
33inline double dd_f_linear(const double v) { return 0; }
34
35static inline boost::function<double(const double)> f = f_log;
36static inline boost::function<double(const double)> d_f = d_f_log;
37static inline boost::function<double(const double)> dd_f = dd_f_log;
38
39using MatrixPtr = boost::shared_ptr<MatrixDouble>;
40using VectorPtr = boost::shared_ptr<VectorDouble>;
41
42using EntData = EntitiesFieldData::EntData;
43using UserDataOperator = ForcesAndSourcesCore::UserDataOperator;
44using VolUserDataOperator = VolumeElementForcesAndSourcesCore::UserDataOperator;
45using FaceUserDataOperator = FaceElementForcesAndSourcesCore::UserDataOperator;
46
49 : public boost::enable_shared_from_this<DataAtIntegrationPts> {
50
51 MatrixDouble approxPAtPts;
52 MatrixDouble approxSigmaAtPts;
53 MatrixDouble divPAtPts;
54 MatrixDouble divSigmaAtPts;
55 MatrixDouble wAtPts;
56 MatrixDouble wDotAtPts;
57 MatrixDouble wDotDotAtPts;
59 MatrixDouble streachTensorAtPts;
60
65 MatrixDouble rotAxisAtPts;
66 MatrixDouble rotAxisDotAtPts;
67 MatrixDouble WAtPts;
68 MatrixDouble W0AtPts;
69 MatrixDouble GAtPts;
70 MatrixDouble G0AtPts;
71 MatrixDouble UAtPts;
72
73 MatrixDouble hAtPts;
74 MatrixDouble rotMatAtPts;
75 MatrixDouble diffRotMatAtPts;
76 MatrixDouble diff2RotMatAtPts;
77 MatrixDouble PAtPts;
78 MatrixDouble SigmaAtPts;
79 VectorDouble phiAtPts;
80 MatrixDouble flowAtPts;
81
82 MatrixDouble P_dh0;
83 MatrixDouble P_dh1;
84 MatrixDouble P_dh2;
85 MatrixDouble P_dH0;
86 MatrixDouble P_dH1;
87 MatrixDouble P_dH2;
88 MatrixDouble Sigma_dh0;
89 MatrixDouble Sigma_dh1;
90 MatrixDouble Sigma_dh2;
91 MatrixDouble Sigma_dH0;
92 MatrixDouble Sigma_dH1;
93 MatrixDouble Sigma_dH2;
94 MatrixDouble phi_dh;
95 MatrixDouble phi_dH;
96 MatrixDouble Flow_dh0;
97 MatrixDouble Flow_dh1;
98 MatrixDouble Flow_dh2;
99 MatrixDouble Flow_dH0;
100 MatrixDouble Flow_dH1;
101 MatrixDouble Flow_dH2;
102
103 MatrixDouble eigenVals;
104 MatrixDouble eigenVecs;
105 VectorInt nbUniq;
106
107 MatrixDouble matD;
108
110 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
112 }
114 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &approxPAtPts);
115 }
116
118 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &divPAtPts);
119 }
120
122 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &divSigmaAtPts);
123 }
124
126 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &wAtPts);
127 }
128
130 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &wDotAtPts);
131 }
132
134 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &wDotDotAtPts);
135 }
136
138 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
140 }
141
143 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
145 }
146
148 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
150 }
151
153 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &rotAxisAtPts);
154 }
155
157 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
159 }
160
162 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &GAtPts);
163 }
164
166 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &GAtPts);
167 }
168
170 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &matD);
171 }
172
174 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &UAtPts);
175 }
176
177 boost::shared_ptr<PhysicalEquations> physicsPtr;
178};
179
180struct OpJacobian;
181
183
190
192
195
197 PhysicalEquations(const int size_active, const int size_dependent)
198 : activeVariables(size_active, 0),
199 dependentVariablesPiola(size_dependent, 0),
200 dependentVariablesPiolaDirevatives(size_dependent * size_active, 0) {}
201 virtual ~PhysicalEquations() = default;
202
203 virtual MoFEMErrorCode recordTape(const int tag, DTensor2Ptr *t_h) = 0;
204
205 virtual OpJacobian *
206 returnOpJacobian(const int tag, const bool eval_rhs, const bool eval_lhs,
207 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
208 boost::shared_ptr<PhysicalEquations> &physics_ptr) = 0;
209
210 std::vector<double> activeVariables;
211 std::vector<double> dependentVariablesPiola;
213
214 /** \name Active variables */
215
216 /**@{*/
217
218 template <int S>
219 inline static DTensor2Ptr get_VecTensor2(std::vector<double> &v) {
220 return DTensor2Ptr(&v[S + 0], &v[S + 1], &v[S + 2], &v[S + 3], &v[S + 4],
221 &v[S + 5], &v[S + 6], &v[S + 7], &v[S + 8]);
222 }
223
224 template <int S>
225 inline static DTensor0Ptr get_VecTensor0(std::vector<double> &v) {
226 return DTensor0Ptr(&v[S + 0]);
227 }
228
229 template <int S0>
230 inline static DTensor3Ptr get_vecTensor3(std::vector<double> &v,
231 const int nba) {
232
233 const int A00 = nba * 0 + S0;
234 const int A01 = nba * 1 + S0;
235 const int A02 = nba * 2 + S0;
236 const int A10 = nba * 3 + S0;
237 const int A11 = nba * 4 + S0;
238 const int A12 = nba * 5 + S0;
239 const int A20 = nba * 6 + S0;
240 const int A21 = nba * 7 + S0;
241 const int A22 = nba * 8 + S0;
242
243 return DTensor3Ptr(
244
245 &v[A00 + 0], &v[A00 + 1], &v[A00 + 2], &v[A01 + 0], &v[A01 + 1],
246 &v[A01 + 2], &v[A02 + 0], &v[A02 + 1], &v[A02 + 2],
247
248 &v[A10 + 0], &v[A10 + 1], &v[A10 + 2], &v[A11 + 0], &v[A11 + 1],
249 &v[A11 + 2], &v[A12 + 0], &v[A12 + 1], &v[A12 + 2],
250
251 &v[A20 + 0], &v[A20 + 1], &v[A20 + 2], &v[A21 + 0], &v[A21 + 1],
252 &v[A21 + 2], &v[A22 + 0], &v[A22 + 1], &v[A22 + 2]
253
254 );
255 }
256
257 /**@}*/
258
259 /** \name Active variables */
260
261 /**@{*/
262
263 inline DTensor2Ptr get_h() { return get_VecTensor2<0>(activeVariables); }
264
265 /**@}*/
266
267 /** \name Dependent variables */
268
269 /**@{*/
270
272 return get_VecTensor2<0>(dependentVariablesPiola);
273 }
274
275 /**@}*/
276
277 /** \name Derivatives of dependent variables */
278
279 /**@{*/
280
282 return get_vecTensor3<0>(dependentVariablesPiolaDirevatives,
283 activeVariables.size());
284 }
286 return get_vecTensor3<3>(dependentVariablesPiolaDirevatives,
287 activeVariables.size());
288 }
290 return get_vecTensor3<6>(dependentVariablesPiolaDirevatives,
291 activeVariables.size());
292 }
293
294 /**@}*/
295};
296
297struct BcDisp {
298 BcDisp(std::string name, std::vector<double> &attr, Range &faces);
299 std::string blockName;
301 VectorDouble3 vals;
302 VectorInt3 flags;
303};
304typedef std::vector<BcDisp> BcDispVec;
305
306struct BcRot {
307 BcRot(std::string name, std::vector<double> &attr, Range &faces);
308 std::string blockName;
310 VectorDouble3 vals;
311 double theta;
312};
313typedef std::vector<BcRot> BcRotVec;
314
315typedef std::vector<Range> TractionFreeBc;
316
318 TractionBc(std::string name, std::vector<double> &attr, Range &faces);
319 std::string blockName;
321 VectorDouble3 vals;
322 VectorInt3 flags;
323};
324typedef std::vector<TractionBc> TractionBcVec;
325
327
328 OpJacobian(const int tag, const bool eval_rhs, const bool eval_lhs,
329 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
330 boost::shared_ptr<PhysicalEquations> &physics_ptr)
331 : UserDataOperator(NOSPACE, OPSPACE), tAg(tag), evalRhs(eval_rhs),
332 evalLhs(eval_lhs), dataAtPts(data_ptr), physicsPtr(physics_ptr) {}
333
334 virtual MoFEMErrorCode evaluateRhs(EntData &data) = 0;
335 virtual MoFEMErrorCode evaluateLhs(EntData &data) = 0;
336
337 MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
338
339protected:
341
342 int tAg = -1; ///< adol-c tape
343 bool evalRhs = false;
344 bool evalLhs = false;
345
346 boost::shared_ptr<DataAtIntegrationPts>
347 dataAtPts; ///< data at integration pts
348 boost::shared_ptr<PhysicalEquations>
349 physicsPtr; ///< material physical equations
350};
351
352template <typename T> struct OpAssembleBasic : public T {
353
354 using ScaleOff = boost::function<double()>;
356
357 boost::shared_ptr<DataAtIntegrationPts>
358 dataAtPts; ///< data at integration pts
359
360 OpAssembleBasic(const std::string &field_name,
361 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
362 const char type)
363 : T(field_name, type), dataAtPts(data_ptr), assembleSymmetry(false) {}
364
366 const std::string &row_field, const std::string &col_field,
367 boost::shared_ptr<DataAtIntegrationPts> data_ptr, const char type,
368 const bool assemble_symmetry, ScaleOff scale_off = []() { return 1; })
369 : T(row_field, col_field, type, false), dataAtPts(data_ptr),
370 assembleSymmetry(assemble_symmetry), scaleOff(scale_off) {}
371
372 VectorDouble nF; ///< local right hand side vector
373 MatrixDouble K; ///< local tangent matrix
374 MatrixDouble transposeK;
375
377
378 virtual MoFEMErrorCode integrate(EntData &data) {
380 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
382 }
383
384 virtual MoFEMErrorCode integrate(int row_side, EntityType row_type,
385 EntData &data) {
387 CHKERR integrate(data);
389 }
390
391 virtual MoFEMErrorCode assemble(EntData &data) {
393 double *vec_ptr = &*nF.begin();
394 int nb_dofs = data.getIndices().size();
395 int *ind_ptr = &*data.getIndices().begin();
396 CHKERR VecSetValues(this->getTSf(), nb_dofs, ind_ptr, vec_ptr, ADD_VALUES);
398 }
399
400 virtual MoFEMErrorCode assemble(int row_side, EntityType row_type,
401 EntData &data) {
403 CHKERR assemble(data);
405 }
406
407 virtual MoFEMErrorCode integrate(EntData &row_data, EntData &col_data) {
409 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
411 }
412
413 virtual MoFEMErrorCode assemble(int row_side, int col_side,
414 EntityType row_type, EntityType col_type,
415 EntData &row_data, EntData &col_data) {
417
418 if (assembleSymmetry) {
419 const auto row_nb_dofs = row_data.getIndices().size();
420 const auto col_nb_dofs = col_data.getIndices().size();
421 transposeK.resize(col_nb_dofs, row_nb_dofs, false);
422 noalias(transposeK) = trans(K);
423 transposeK *= scaleOff();
424 }
425
426 CHKERR MatSetValues<SchurL2Mats>(this->getTSA(), row_data, col_data, K,
427 ADD_VALUES);
428 if (assembleSymmetry) {
429 CHKERR MatSetValues<SchurL2Mats>(this->getTSA(), col_data, row_data,
430 transposeK, ADD_VALUES);
431 }
432
434 }
435
436 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
438 if (data.getIndices().empty())
440 nF.resize(data.getIndices().size(), false);
441 nF.clear();
442 CHKERR integrate(side, type, data);
443 CHKERR assemble(side, type, data);
445 }
446
447 MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
448 EntityType col_type, EntData &row_data,
449 EntData &col_data) {
451 if (row_data.getIndices().empty())
453 if (col_data.getIndices().empty())
455 K.resize(row_data.getIndices().size(), col_data.getIndices().size(), false);
456 K.clear();
457 CHKERR integrate(row_data, col_data);
458 CHKERR assemble(row_side, col_side, row_type, col_type, row_data, col_data);
460 }
461};
462
463struct OpAssembleVolume : public OpAssembleBasic<VolUserDataOperator> {
464 OpAssembleVolume(const std::string &field,
465 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
466 const char type)
467 : OpAssembleBasic<VolUserDataOperator>(field, data_ptr, type) {}
468
469 OpAssembleVolume(const std::string &row_field, const std::string &col_field,
470 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
471 const char type, const bool assemble_symmetry)
472 : OpAssembleBasic<VolUserDataOperator>(row_field, col_field, data_ptr,
473 type, assemble_symmetry) {}
474};
475
476struct OpAssembleFace : public OpAssembleBasic<FaceUserDataOperator> {
477 OpAssembleFace(const std::string &field,
478 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
479 const char type)
480 : OpAssembleBasic<FaceUserDataOperator>(field, data_ptr, type) {}
481
482 OpAssembleFace(const std::string &row_field, const std::string &col_field,
483 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
484 const char type, const bool assemble_symmetry)
485 : OpAssembleBasic<FaceUserDataOperator>(row_field, col_field, data_ptr,
486 type, assemble_symmetry) {}
487};
488
490 boost::shared_ptr<DataAtIntegrationPts>
491 dataAtPts; ///< data at integration pts
493 boost::shared_ptr<DataAtIntegrationPts> &data_ptr)
494 : VolUserDataOperator(NOSPACE, OPSPACE), dataAtPts(data_ptr) {}
495 MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
496};
497
499 const double alphaW;
500 const double alphaRho;
502 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
503 const double alpha, const double rho)
504 : OpAssembleVolume(field_name, data_ptr, OPROW), alphaW(alpha),
505 alphaRho(rho) {}
506 MoFEMErrorCode integrate(EntData &data);
507};
508
510 OpSpatialRotation(const std::string &field_name,
511 boost::shared_ptr<DataAtIntegrationPts> &data_ptr)
512 : OpAssembleVolume(field_name, data_ptr, OPROW) {}
513 MoFEMErrorCode integrate(EntData &data);
514};
515
517
518 OpSpatialPhysical(const std::string &field_name,
519 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
520 const double alpha_u)
521 : OpAssembleVolume(field_name, data_ptr, OPROW), alphaU(alpha_u) {}
522
523 MoFEMErrorCode integrate(EntData &data);
524
525 MoFEMErrorCode integratePiola(EntData &data);
526 MoFEMErrorCode integrateHencky(EntData &data);
527
528private:
529 const double alphaU;
530};
531
534 boost::shared_ptr<DataAtIntegrationPts> &data_ptr)
535 : OpAssembleVolume(field_name, data_ptr, OPROW) {}
536 MoFEMErrorCode integrate(EntData &data);
537
538private:
539};
540
543 boost::shared_ptr<DataAtIntegrationPts> &data_ptr)
544 : OpAssembleVolume(field_name, data_ptr, OPROW) {}
545 MoFEMErrorCode integrate(EntData &data);
546
547private:
548};
549
552 boost::shared_ptr<DataAtIntegrationPts> &data_ptr)
553 : OpAssembleVolume(field_name, data_ptr, OPROW) {}
554 MoFEMErrorCode integrate(EntData &data);
555};
556
557struct OpDispBc : public OpAssembleFace {
558 boost::shared_ptr<BcDispVec> bcDispPtr;
559 OpDispBc(const std::string &field_name,
560 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
561 boost::shared_ptr<BcDispVec> &bc_disp_ptr)
562 : OpAssembleFace(field_name, data_ptr, OPROW), bcDispPtr(bc_disp_ptr) {}
563 MoFEMErrorCode integrate(EntData &data);
564};
565
567 boost::shared_ptr<BcDispVec> bcDispPtr;
568 OpDispBc_dx(const std::string &row_field_name,
569 const std::string &col_field_name,
570 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
571 boost::shared_ptr<BcDispVec> &bc_disp_ptr)
572 : OpAssembleFace(row_field_name, col_field_name, data_ptr, OPROWCOL,
573 false),
574 bcDispPtr(bc_disp_ptr) {}
575 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
576};
577
579 boost::shared_ptr<BcRotVec> bcRotPtr;
580 OpRotationBc(const std::string &field_name,
581 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
582 boost::shared_ptr<BcRotVec> &bc_rot_ptr)
583 : OpAssembleFace(field_name, data_ptr, OPROW), bcRotPtr(bc_rot_ptr) {}
584 MoFEMErrorCode integrate(EntData &data);
585};
586
588 boost::shared_ptr<BcRotVec> bcRotPtr;
589 OpRotationBc_dx(const std::string &row_field_name,
590 const std::string &col_field_name,
591 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
592 boost::shared_ptr<BcRotVec> &bc_rot_ptr)
593 : OpAssembleFace(row_field_name, col_field_name, data_ptr, OPROWCOL,
594 false),
595 bcRotPtr(bc_rot_ptr) {}
596 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
597};
598
599struct FeTractionBc;
600
602 OpTractionBc(std::string row_field, FeTractionBc &bc_fe)
603 : FaceUserDataOperator(row_field, FaceUserDataOperator::OPROW),
604 bcFe(bc_fe) {}
605 MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
606
607protected:
609 MatrixDouble matPP;
610 MatrixDouble vecPv;
611};
612
613struct FeTractionBc : public FEMethod {
615 FeTractionBc(MoFEM::Interface &m_field, const std::string field_name,
616 boost::shared_ptr<TractionBcVec> &bc)
617 : FEMethod(), mField(m_field), bcData(bc), fieldName(field_name) {}
618 MoFEMErrorCode preProcess();
619 MoFEMErrorCode postProcess() { return 0; }
620
621 friend MoFEMErrorCode OpTractionBc::doWork(int side, EntityType type,
622 EntData &data);
623
624protected:
626 boost::shared_ptr<TractionBcVec> bcData;
627 std::string fieldName;
628};
629
631 OpSpatialEquilibrium_dw_dP(const std::string &row_field,
632 const std::string &col_field,
633 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
634 const bool assemble_off = false)
635 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL,
636 assemble_off) {
637 sYmm = false;
638 }
639 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
640};
641
643 const double alphaW;
644 const double alphaRho;
645 OpSpatialEquilibrium_dw_dw(const std::string &row_field,
646 const std::string &col_field,
647 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
648 const double alpha, const double rho)
649 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL, false),
650 alphaW(alpha), alphaRho(rho) {
651 sYmm = true;
652 }
653 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
654};
655
657 const double alphaU;
658 OpSpatialPhysical_du_du(const std::string &row_field,
659 const std::string &col_field,
660 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
661 const double alpha)
662 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL, false),
663 alphaU(alpha) {
664 sYmm = false;
665 }
666 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
667 MoFEMErrorCode integratePiola(EntData &row_data, EntData &col_data);
668 MoFEMErrorCode integrateHencky(EntData &row_data, EntData &col_data);
669};
670
672 OpSpatialPhysical_du_dP(const std::string &row_field,
673 const std::string &col_field,
674 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
675 const bool assemble_off = false)
676 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL,
677 assemble_off) {
678 sYmm = false;
679 }
680
681 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
682};
683
686 const std::string &row_field, const std::string &col_field,
687 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
688 const bool assemble_off = false)
689 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL,
690 assemble_off) {
691 sYmm = false;
692 }
693
694 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
695};
696
698 OpSpatialPhysical_du_domega(const std::string &row_field,
699 const std::string &col_field,
700 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
701 const bool assemble_off = false)
702 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL,
703 assemble_off) {
704 sYmm = false;
705 }
706 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
707};
708
710 OpSpatialPhysical_du_dx(const std::string &row_field,
711 const std::string &col_field,
712 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
713 const bool assemble_off = false)
714 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL,
715 assemble_off) {
716 sYmm = false;
717 }
718 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
719};
720
722 OpSpatialRotation_domega_dP(const std::string &row_field,
723 const std::string &col_field,
724 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
725 const bool assemble_off)
726 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL,
727 assemble_off) {
728 sYmm = false;
729 }
730 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
731};
732
735 const std::string &row_field, const std::string &col_field,
736 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
737 const bool assemble_off)
738 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL,
739 assemble_off) {
740 sYmm = false;
741 }
742 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
743};
744
747 const std::string &row_field, const std::string &col_field,
748 boost::shared_ptr<DataAtIntegrationPts> &data_ptr)
749 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL, false) {
750 sYmm = false;
751 }
752 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
753};
754
756 OpSpatialRotation_domega_du(const std::string &row_field,
757 const std::string &col_field,
758 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
759 const bool assemble_off = false)
760 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL,
761 assemble_off) {
762 sYmm = false;
763 }
764 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
765};
766
769 const std::string &row_field, const std::string &col_field,
770 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
771 const bool assemble_off)
772 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL,
773 assemble_off) {
774 sYmm = false;
775 }
776 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
777
778private:
779};
780
783 const std::string &row_field, const std::string &col_field,
784 boost::shared_ptr<DataAtIntegrationPts> &data_ptr,
785 const bool assemble_off)
786 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL,
787 assemble_off) {
788 sYmm = false;
789 }
790 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
791
792private:
793};
794
796
797 moab::Interface &postProcMesh;
798 std::vector<EntityHandle> &mapGaussPts;
799 boost::shared_ptr<DataAtIntegrationPts> dataAtPts;
800
801 OpPostProcDataStructure(moab::Interface &post_proc_mesh,
802 std::vector<EntityHandle> &map_gauss_pts,
803 const std::string field_name,
804 boost::shared_ptr<DataAtIntegrationPts> &data_ptr)
806 postProcMesh(post_proc_mesh), mapGaussPts(map_gauss_pts),
807 dataAtPts(data_ptr) {}
808
809 MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
810};
811
813 OpSpatialPrj(const std::string &row_field,
814 boost::shared_ptr<DataAtIntegrationPts> &data_ptr)
815 : OpAssembleVolume(row_field, data_ptr, OPROW) {}
816 MoFEMErrorCode integrate(EntData &row_data);
817};
818
820 OpSpatialPrj_dx_dx(const std::string &row_field, const std::string &col_field,
821 boost::shared_ptr<DataAtIntegrationPts> &data_ptr)
822 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL, false) {
823 // FIXME: That is symmetric
824 sYmm = false;
825 }
826 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
827};
828
830 OpSpatialPrj_dx_dw(const std::string &row_field, const std::string &col_field,
831 boost::shared_ptr<DataAtIntegrationPts> &data_ptr)
832 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL, false) {
833 sYmm = false;
834 }
835 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
836};
837
839
841 boost::shared_ptr<DataAtIntegrationPts> data_at_pts,
842 boost::shared_ptr<double> &e)
844 dataAtPts(data_at_pts), energy(e) {}
845
846 MoFEMErrorCode doWork(int side, EntityType type,
847 EntitiesFieldData::EntData &data);
848
849private:
850 boost::shared_ptr<DataAtIntegrationPts> dataAtPts;
851 boost::shared_ptr<double> energy;
852};
853
855
856 /**
857 * \brief Getting interface of core database
858 * @param uuid unique ID of interface
859 * @param iface returned pointer to interface
860 * @return error code
861 */
862 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
863 UnknownInterface **iface) const;
864
866
867 boost::shared_ptr<DataAtIntegrationPts> dataAtPts;
868 boost::shared_ptr<PhysicalEquations> physicalEquations;
869
870 boost::shared_ptr<VolumeElementForcesAndSourcesCore> elasticFeRhs;
871 boost::shared_ptr<VolumeElementForcesAndSourcesCore> elasticFeLhs;
872 boost::shared_ptr<FaceElementForcesAndSourcesCore> elasticBcLhs;
873 boost::shared_ptr<FaceElementForcesAndSourcesCore> elasticBcRhs;
874 boost::shared_ptr<FaceElementForcesAndSourcesCore> contactRhs;
875
876 SmartPetscObj<DM> dM; ///< Coupled problem all fields
877 SmartPetscObj<DM> dmElastic; ///< Elastic problem
878 SmartPetscObj<DM> dmPrjSpatial; ///< Projection spatial displacement
879
880 const std::string piolaStress;
881 const std::string eshelbyStress;
882 const std::string spatialL2Disp;
883 const std::string spatialH1Disp;
884 const std::string materialL2Disp;
885 const std::string streachTensor;
886 const std::string rotAxis;
887 const std::string materialGradient;
888 const std::string tauField;
889 const std::string lambdaField;
890 const std::string bubbleField;
891
892 const std::string elementVolumeName;
893 const std::string naturalBcElement;
894 const std::string essentialBcElement;
895 const std::string skinElement;
896 const std::string contactElement;
897
899 virtual ~EshelbianCore();
900
903 double alphaU;
904 double alphaW;
905 double alphaRho;
906 double precEps;
908 double precEpsW;
909
910 MoFEMErrorCode getOptions();
911
912 boost::shared_ptr<BcDispVec> bcSpatialDispVecPtr;
913 boost::shared_ptr<BcRotVec> bcSpatialRotationVecPtr;
914 boost::shared_ptr<TractionBcVec> bcSpatialTraction;
915 boost::shared_ptr<TractionFreeBc> bcSpatialFreeTraction;
916
917 template <typename BC>
918 MoFEMErrorCode getBc(boost::shared_ptr<BC> &bc_vec_ptr,
919 const std::string block_name, const int nb_attributes) {
921 for (auto it :
922 mField.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
923
924 (boost::format("%s(.*)") % block_name).str()
925
926 ))
927
928 ) {
929 std::vector<double> block_attributes;
930 CHKERR it->getAttributes(block_attributes);
931 if (block_attributes.size() != nb_attributes) {
932 SETERRQ3(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
933 "In block %s expected %d attributes, but given %d",
934 it->getName().c_str(), nb_attributes, block_attributes.size());
935 }
936 Range faces;
937 CHKERR it->getMeshsetIdEntitiesByDimension(mField.get_moab(), 2, faces,
938 true);
939 bc_vec_ptr->emplace_back(it->getName(), block_attributes, faces);
940 }
942 }
943
944 MoFEMErrorCode getSpatialDispBc();
945
946 inline MoFEMErrorCode getSpatialRotationBc() {
947 bcSpatialRotationVecPtr = boost::make_shared<BcRotVec>();
948 return getBc(bcSpatialRotationVecPtr, "SPATIAL_ROTATION_BC", 4);
949 }
950
951 MoFEMErrorCode getSpatialTractionBc();
952
953 /**
954 * @brief Remove all, but entities where kinematic constrains are applied.
955 *
956 * @param meshset
957 * @param bc_ptr
958 * @param disp_block_set_name
959 * @param rot_block_set_name
960 * @param contact_set_name
961 * @return MoFEMErrorCode
962 */
963 MoFEMErrorCode getTractionFreeBc(const EntityHandle meshset,
964 boost::shared_ptr<TractionFreeBc> &bc_ptr,
965 const std::string contact_set_name);
966
967 inline MoFEMErrorCode
970 boost::shared_ptr<TractionFreeBc>(new TractionFreeBc());
971 return getTractionFreeBc(meshset, bcSpatialFreeTraction, "CONTACT");
972 }
973
974 MoFEMErrorCode addFields(const EntityHandle meshset = 0);
975 MoFEMErrorCode addVolumeFiniteElement(const EntityHandle meshset = 0);
976 MoFEMErrorCode addBoundaryFiniteElement(const EntityHandle meshset = 0);
977 MoFEMErrorCode addDMs(const BitRefLevel bit = BitRefLevel().set(0));
978
979 MoFEMErrorCode addMaterial_HMHHStVenantKirchhoff(const int tape,
980 const double lambda,
981 const double mu,
982 const double sigma_y);
983
984 MoFEMErrorCode addMaterial_HMHMooneyRivlin(const int tape, const double alpha,
985 const double beta,
986 const double lambda,
987 const double sigma_y);
988
989 MoFEMErrorCode addMaterial_Hencky(double E, double nu);
990
991 MoFEMErrorCode setBaseVolumeElementOps(
992 const int tag, const bool do_rhs, const bool do_lhs,
993 boost::shared_ptr<VolumeElementForcesAndSourcesCore> &fe);
994
995 MoFEMErrorCode setVolumeElementOps(
996 const int tag, const bool add_elastic, const bool add_material,
997 boost::shared_ptr<VolumeElementForcesAndSourcesCore> &fe_rhs,
998 boost::shared_ptr<VolumeElementForcesAndSourcesCore> &fe_lhs);
999
1000 MoFEMErrorCode
1001 setFaceElementOps(const bool add_elastic, const bool add_material,
1002 boost::shared_ptr<FaceElementForcesAndSourcesCore> &fe_rhs,
1003 boost::shared_ptr<FaceElementForcesAndSourcesCore> &fe_lhs);
1004
1005 MoFEMErrorCode setContactElementOps(
1006 boost::shared_ptr<FaceElementForcesAndSourcesCore> &fe_rhs,
1007 boost::shared_ptr<FaceElementForcesAndSourcesCore> &fe_lhs);
1008
1009 MoFEMErrorCode setElasticElementOps(const int tag);
1010 MoFEMErrorCode setElasticElementToTs(DM dm);
1011
1012 MoFEMErrorCode solveElastic(TS ts, Mat m, Vec f, Vec x);
1013
1014 MoFEMErrorCode postProcessResults(const int tag, const std::string file);
1015 MoFEMErrorCode gettingNorms();
1016
1017 struct SetUpSchur {
1018 static boost::shared_ptr<SetUpSchur> createSetUpSchur(
1019
1020 MoFEM::Interface &m_field, SmartPetscObj<Mat> m, SmartPetscObj<Mat> p,
1021 EshelbianCore *ep_core_ptr
1022
1023 );
1024 virtual MoFEMErrorCode setUp(KSP solver) = 0;
1025
1026 virtual MoFEMErrorCode preProc() = 0;
1027 virtual MoFEMErrorCode postProc() = 0;
1028
1029 protected:
1030 SetUpSchur() = default;
1031 };
1032};
1033
1034} // namespace EshelbianPlasticity
1035
1036#endif //__ESHELBIAN_PLASTICITY_HPP__
static Index< 'p', 3 > p
constexpr int SPACE_DIM
constexpr AssemblyType A
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
@ NOSPACE
Definition: definitions.h:83
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ 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()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
FTensor::Index< 'm', SPACE_DIM > m
auto bit
set bit
static double lambda
const double v
phase velocity of light in medium (cm/ns)
const double T
std::vector< BcDisp > BcDispVec
double dd_f_linear(const double v)
boost::shared_ptr< MatrixDouble > MatrixPtr
double f_log(const double v)
static boost::function< double(const double)> d_f
double d_f_linear(const double v)
VolumeElementForcesAndSourcesCore::UserDataOperator VolUserDataOperator
static boost::function< double(const double)> f
double dd_f_log(const double v)
double f_linear(const double v)
static boost::function< double(const double)> dd_f
std::vector< Range > TractionFreeBc
ForcesAndSourcesCore::UserDataOperator UserDataOperator
double d_f_log(const double v)
static enum RotSelector rotSelector
std::vector< BcRot > BcRotVec
std::vector< TractionBc > TractionBcVec
FaceElementForcesAndSourcesCore::UserDataOperator FaceUserDataOperator
boost::shared_ptr< VectorDouble > VectorPtr
double rho
Definition: plastic.cpp:182
constexpr auto field_name
boost::shared_ptr< PhysicalEquations > physicsPtr
static boost::shared_ptr< SetUpSchur > createSetUpSchur(MoFEM::Interface &m_field, SmartPetscObj< Mat > m, SmartPetscObj< Mat > p, EshelbianCore *ep_core_ptr)
virtual MoFEMErrorCode setUp(KSP solver)=0
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Getting interface of core database.
boost::shared_ptr< FaceElementForcesAndSourcesCore > elasticBcLhs
boost::shared_ptr< VolumeElementForcesAndSourcesCore > elasticFeLhs
MoFEMErrorCode addMaterial_Hencky(double E, double nu)
MoFEMErrorCode postProcessResults(const int tag, const std::string file)
MoFEMErrorCode addBoundaryFiniteElement(const EntityHandle meshset=0)
MoFEMErrorCode addDMs(const BitRefLevel bit=BitRefLevel().set(0))
MoFEMErrorCode addMaterial_HMHHStVenantKirchhoff(const int tape, const double lambda, const double mu, const double sigma_y)
boost::shared_ptr< VolumeElementForcesAndSourcesCore > elasticFeRhs
MoFEMErrorCode setElasticElementOps(const int tag)
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
MoFEMErrorCode solveElastic(TS ts, Mat m, Vec f, Vec x)
MoFEMErrorCode getSpatialTractionFreeBc(const EntityHandle meshset=0)
MoFEMErrorCode getBc(boost::shared_ptr< BC > &bc_vec_ptr, const std::string block_name, const int nb_attributes)
MoFEMErrorCode setBaseVolumeElementOps(const int tag, const bool do_rhs, const bool do_lhs, boost::shared_ptr< VolumeElementForcesAndSourcesCore > &fe)
boost::shared_ptr< BcRotVec > bcSpatialRotationVecPtr
MoFEMErrorCode setContactElementOps(boost::shared_ptr< FaceElementForcesAndSourcesCore > &fe_rhs, boost::shared_ptr< FaceElementForcesAndSourcesCore > &fe_lhs)
MoFEMErrorCode setFaceElementOps(const bool add_elastic, const bool add_material, boost::shared_ptr< FaceElementForcesAndSourcesCore > &fe_rhs, boost::shared_ptr< FaceElementForcesAndSourcesCore > &fe_lhs)
MoFEMErrorCode addMaterial_HMHMooneyRivlin(const int tape, const double alpha, const double beta, const double lambda, const double sigma_y)
SmartPetscObj< DM > dM
Coupled problem all fields.
MoFEMErrorCode addVolumeFiniteElement(const EntityHandle meshset=0)
SmartPetscObj< DM > dmElastic
Elastic problem.
boost::shared_ptr< BcDispVec > bcSpatialDispVecPtr
boost::shared_ptr< PhysicalEquations > physicalEquations
MoFEMErrorCode setVolumeElementOps(const int tag, const bool add_elastic, const bool add_material, boost::shared_ptr< VolumeElementForcesAndSourcesCore > &fe_rhs, boost::shared_ptr< VolumeElementForcesAndSourcesCore > &fe_lhs)
boost::shared_ptr< FaceElementForcesAndSourcesCore > contactRhs
MoFEMErrorCode gettingNorms()
[Getting norms]
boost::shared_ptr< TractionFreeBc > bcSpatialFreeTraction
boost::shared_ptr< TractionBcVec > bcSpatialTraction
MoFEMErrorCode addFields(const EntityHandle meshset=0)
boost::shared_ptr< FaceElementForcesAndSourcesCore > elasticBcRhs
SmartPetscObj< DM > dmPrjSpatial
Projection spatial displacement.
MoFEMErrorCode getTractionFreeBc(const EntityHandle meshset, boost::shared_ptr< TractionFreeBc > &bc_ptr, const std::string contact_set_name)
Remove all, but entities where kinematic constrains are applied.
boost::shared_ptr< TractionBcVec > bcData
FeTractionBc(MoFEM::Interface &m_field, const std::string field_name, boost::shared_ptr< TractionBcVec > &bc)
virtual MoFEMErrorCode assemble(EntData &data)
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
virtual MoFEMErrorCode integrate(int row_side, EntityType row_type, EntData &data)
virtual MoFEMErrorCode integrate(EntData &data)
OpAssembleBasic(const std::string &row_field, const std::string &col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const char type, const bool assemble_symmetry, ScaleOff scale_off=[]() { return 1;})
virtual MoFEMErrorCode assemble(int row_side, EntityType row_type, EntData &data)
MatrixDouble K
local tangent matrix
VectorDouble nF
local right hand side vector
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
data at integration pts
virtual MoFEMErrorCode assemble(int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
virtual MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
OpAssembleBasic(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const char type)
OpAssembleFace(const std::string &row_field, const std::string &col_field, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, const char type, const bool assemble_symmetry)
OpAssembleFace(const std::string &field, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, const char type)
OpAssembleVolume(const std::string &field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const char type)
OpAssembleVolume(const std::string &row_field, const std::string &col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const char type, const bool assemble_symmetry)
OpCalculateRotationAndSpatialGradient(boost::shared_ptr< DataAtIntegrationPts > &data_ptr)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
data at integration pts
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
OpCalculateStrainEnergy(const std::string field_name, boost::shared_ptr< DataAtIntegrationPts > data_at_pts, boost::shared_ptr< double > &e)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpDispBc_dx(const std::string &row_field_name, const std::string &col_field_name, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, boost::shared_ptr< BcDispVec > &bc_disp_ptr)
boost::shared_ptr< BcDispVec > bcDispPtr
boost::shared_ptr< BcDispVec > bcDispPtr
OpDispBc(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, boost::shared_ptr< BcDispVec > &bc_disp_ptr)
MoFEMErrorCode integrate(EntData &data)
boost::shared_ptr< PhysicalEquations > physicsPtr
material physical equations
virtual MoFEMErrorCode evaluateRhs(EntData &data)=0
OpJacobian(const int tag, const bool eval_rhs, const bool eval_lhs, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, boost::shared_ptr< PhysicalEquations > &physics_ptr)
virtual MoFEMErrorCode evaluateLhs(EntData &data)=0
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
data at integration pts
moab::Interface & postProcMesh
OpPostProcDataStructure(moab::Interface &post_proc_mesh, std::vector< EntityHandle > &map_gauss_pts, const std::string field_name, boost::shared_ptr< DataAtIntegrationPts > &data_ptr)
std::vector< EntityHandle > & mapGaussPts
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
OpRotationBc_dx(const std::string &row_field_name, const std::string &col_field_name, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, boost::shared_ptr< BcRotVec > &bc_rot_ptr)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
boost::shared_ptr< BcRotVec > bcRotPtr
MoFEMErrorCode integrate(EntData &data)
OpRotationBc(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, boost::shared_ptr< BcRotVec > &bc_rot_ptr)
boost::shared_ptr< BcRotVec > bcRotPtr
OpSpatialConsistency_dBubble_domega(const std::string &row_field, const std::string &col_field, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, const bool assemble_off)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialConsistency_dP_domega(const std::string &row_field, const std::string &col_field, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, const bool assemble_off)
OpSpatialConsistencyBubble(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > &data_ptr)
OpSpatialConsistencyDivTerm(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > &data_ptr)
OpSpatialConsistencyP(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > &data_ptr)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialEquilibrium_dw_dP(const std::string &row_field, const std::string &col_field, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, const bool assemble_off=false)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialEquilibrium_dw_dw(const std::string &row_field, const std::string &col_field, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, const double alpha, const double rho)
OpSpatialEquilibrium(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, const double alpha, const double rho)
OpSpatialPhysical_du_dBubble(const std::string &row_field, const std::string &col_field, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, const bool assemble_off=false)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialPhysical_du_dP(const std::string &row_field, const std::string &col_field, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, const bool assemble_off=false)
OpSpatialPhysical_du_domega(const std::string &row_field, const std::string &col_field, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, const bool assemble_off=false)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
MoFEMErrorCode integrateHencky(EntData &row_data, EntData &col_data)
OpSpatialPhysical_du_du(const std::string &row_field, const std::string &col_field, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, const double alpha)
MoFEMErrorCode integratePiola(EntData &row_data, EntData &col_data)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialPhysical_du_dx(const std::string &row_field, const std::string &col_field, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, const bool assemble_off=false)
MoFEMErrorCode integratePiola(EntData &data)
MoFEMErrorCode integrate(EntData &data)
MoFEMErrorCode integrateHencky(EntData &data)
OpSpatialPhysical(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, const double alpha_u)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialPrj_dx_dw(const std::string &row_field, const std::string &col_field, boost::shared_ptr< DataAtIntegrationPts > &data_ptr)
OpSpatialPrj_dx_dx(const std::string &row_field, const std::string &col_field, boost::shared_ptr< DataAtIntegrationPts > &data_ptr)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
MoFEMErrorCode integrate(EntData &row_data)
OpSpatialPrj(const std::string &row_field, boost::shared_ptr< DataAtIntegrationPts > &data_ptr)
OpSpatialRotation_domega_dBubble(const std::string &row_field, const std::string &col_field, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, const bool assemble_off)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialRotation_domega_dP(const std::string &row_field, const std::string &col_field, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, const bool assemble_off)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialRotation_domega_domega(const std::string &row_field, const std::string &col_field, boost::shared_ptr< DataAtIntegrationPts > &data_ptr)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialRotation_domega_du(const std::string &row_field, const std::string &col_field, boost::shared_ptr< DataAtIntegrationPts > &data_ptr, const bool assemble_off=false)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
MoFEMErrorCode integrate(EntData &data)
OpSpatialRotation(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > &data_ptr)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
OpTractionBc(std::string row_field, FeTractionBc &bc_fe)
FTensor::Tensor3< double, 3, 3, 3 > DTensor3
FTensor::Tensor3< FTensor::PackPtr< double *, 1 >, 3, 3, 3 > DTensor3Ptr
FTensor::Tensor2< FTensor::PackPtr< double *, 1 >, 3, 3 > DTensor2Ptr
FTensor::Tensor1< double, 3 > DTensor1
static DTensor3Ptr get_vecTensor3(std::vector< double > &v, const int nba)
FTensor::Tensor2< adouble, 3, 3 > ATensor2
FTensor::Tensor2< double, 3, 3 > DTensor2
PhysicalEquations(const int size_active, const int size_dependent)
std::vector< double > dependentVariablesPiolaDirevatives
static DTensor2Ptr get_VecTensor2(std::vector< double > &v)
FTensor::Tensor1< adouble, 3 > ATensor1
FTensor::Tensor3< adouble, 3, 3, 3 > ATensor3
virtual MoFEMErrorCode recordTape(const int tag, DTensor2Ptr *t_h)=0
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > DTensor0Ptr
static DTensor0Ptr get_VecTensor0(std::vector< double > &v)
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)=0
virtual moab::Interface & get_moab()=0
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
const VectorInt & getIndices() const
Get global indices of dofs on entity.
base class for all interface classes
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.