v0.15.5
Loading...
Searching...
No Matches
EshelbianOperators.hpp
Go to the documentation of this file.
1/**
2 * @file EshelbianOperators.hpp
3 * @author your name (you@domain.com)
4 * @brief
5 * @version 0.1
6 * @date 2024-12-31
7 *
8 * @copyright Copyright (c) 2024
9 *
10 */
11
13
14 OpJacobian(const int tag, const bool eval_rhs, const bool eval_lhs,
15 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
16 boost::shared_ptr<PhysicalEquations> physics_ptr)
17 : UserDataOperator(NOSPACE, OPSPACE), tAg(tag), evalRhs(eval_rhs),
18 evalLhs(eval_lhs), dataAtPts(data_ptr), physicsPtr(physics_ptr) {}
19
20 virtual MoFEMErrorCode evaluateRhs(EntData &data) = 0;
21 virtual MoFEMErrorCode evaluateLhs(EntData &data) = 0;
22
23 MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
24
25protected:
27
28 int tAg = -1; ///< adol-c tape
29 bool evalRhs = false;
30 bool evalLhs = false;
31
32 boost::shared_ptr<DataAtIntegrationPts>
33 dataAtPts; ///< data at integration pts
34 boost::shared_ptr<PhysicalEquations>
35 physicsPtr; ///< material physical equations
36};
37
38template <typename T> struct OpAssembleBasic : public T {
39
40 using ScaleOff = boost::function<double()>;
41 const bool assembleSymmetry;
42
43 boost::shared_ptr<DataAtIntegrationPts>
44 dataAtPts; ///< data at integration pts
45
46 OpAssembleBasic(const std::string &field_name,
47 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
48 const char type)
49 : T(field_name, type), assembleSymmetry(false), dataAtPts(data_ptr) {}
50
52 std::string row_field, std::string col_field,
53 boost::shared_ptr<DataAtIntegrationPts> data_ptr, const char type,
54 const bool assemble_symmetry, ScaleOff scale_off = []() { return 1; })
55 : T(row_field, col_field, type, false),
56 assembleSymmetry(assemble_symmetry), dataAtPts(data_ptr),
57 scaleOff(scale_off) {}
58
60 : T(space, T::OPSPACE), assembleSymmetry(false) {}
61
62 VectorDouble nF; ///< local right hand side vector
63 MatrixDouble K; ///< local tangent matrix
64 MatrixDouble transposeK;
65
67
68 virtual MoFEMErrorCode integrate(EntData &data) {
70 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
72 }
73
74 virtual MoFEMErrorCode integrate(int row_side, EntityType row_type,
75 EntData &data) {
77 CHKERR integrate(data);
79 }
80
81 virtual MoFEMErrorCode assemble(EntData &data) {
83 double *vec_ptr = nF.data().data();
84 int nb_dofs = data.getIndices().size();
85 int *ind_ptr = data.getIndices().data().data();
86 CHKERR VecSetValues(this->getTSf(), nb_dofs, ind_ptr, vec_ptr, ADD_VALUES);
88 }
89
90 virtual MoFEMErrorCode assemble(int row_side, EntityType row_type,
91 EntData &data) {
93 CHKERR assemble(data);
95 }
96
97 virtual MoFEMErrorCode integrate(EntData &row_data, EntData &col_data) {
99 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
101 }
102
103 virtual MoFEMErrorCode assemble(int row_side, int col_side,
104 EntityType row_type, EntityType col_type,
105 EntData &row_data, EntData &col_data) {
107
108 if (assembleSymmetry) {
109 const auto row_nb_dofs = row_data.getIndices().size();
110 const auto col_nb_dofs = col_data.getIndices().size();
111 transposeK.resize(col_nb_dofs, row_nb_dofs, false);
112 noalias(transposeK) = trans(K);
113 transposeK *= scaleOff();
114 }
115
116 CHKERR MatSetValues<AssemblyTypeSelector<A>>(this->getTSB(), row_data,
117 col_data, K, ADD_VALUES);
118 if (assembleSymmetry) {
119 CHKERR MatSetValues<AssemblyTypeSelector<A>>(
120 this->getTSB(), col_data, row_data, transposeK, ADD_VALUES);
121 }
122
124 }
125
126 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
128 if (data.getIndices().empty())
130 nF.resize(data.getIndices().size(), false);
131 nF.clear();
132 CHKERR integrate(side, type, data);
133 CHKERR assemble(side, type, data);
135 }
136
137 MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
138 EntityType col_type, EntData &row_data,
139 EntData &col_data) {
141 if (row_data.getIndices().empty())
143 if (col_data.getIndices().empty())
145 K.resize(row_data.getIndices().size(), col_data.getIndices().size(), false);
146 K.clear();
147 CHKERR integrate(row_data, col_data);
148 CHKERR assemble(row_side, col_side, row_type, col_type, row_data, col_data);
150 }
151};
152
153struct OpAssembleVolume : public OpAssembleBasic<VolUserDataOperator> {
155 using ScaleOff = typename OP::ScaleOff;
156
158
159 MoFEMErrorCode assemble(int row_side, int col_side, EntityType row_type,
160 EntityType col_type, EntData &row_data,
161 EntData &col_data) {
163 CHKERR OP::assemble(row_side, col_side, row_type, col_type, row_data,
164 col_data);
165 constexpr bool store_matrices = false;
166 if constexpr (store_matrices) {
168 m.resize(K.size1(), K.size2(), false);
169 noalias(m) = K;
170 if (assembleSymmetry) {
172 m.resize(K.size2(), K.size1(), false);
173 noalias(m) = transposeK;
174 }
175 }
177 };
178
179protected:
180 static inline std::map<std::pair<std::string, std::string>, MatrixDouble>
182};
183
184struct OpAssembleFace : public OpAssembleBasic<FaceUserDataOperator> {
185 using OpAssembleBasic<FaceUserDataOperator>::OpAssembleBasic;
186};
187
189 using OpAssembleVolume::OpAssembleVolume;
190 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
192 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
193 "OpAssembleVolumePositiveDefine not implemented");
195 }
196
197 MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
198 EntityType col_type, EntData &row_data,
199 EntData &col_data) {
201 if (row_data.getIndices().empty())
203 if (col_data.getIndices().empty())
205 K.resize(row_data.getIndices().size(), col_data.getIndices().size(), false);
206 K.clear();
207 CHKERR integrate(row_data, col_data);
208
209 constexpr bool run_psd = true;
210
211 if (run_psd) {
212
213 auto symmetrize = [](MatrixDouble &mat) {
214 MatrixDouble r = (mat + trans(mat)) / 2.0;
215 return r;
216 };
217 auto check_positive_definite = [](MatrixDouble &mat) {
218 int n = mat.size1();
219 return (lapack_dpotrf('U', n, mat.data().data(), n) == 0);
220 };
221 auto frobenius_norm = [](MatrixDouble &mat) {
222 double sum = 0.0;
223 for (std::size_t i = 0; i < mat.size1(); ++i)
224 for (std::size_t j = 0; j < mat.size2(); ++j)
225 sum += mat(i, j) * mat(i, j);
226 return std::sqrt(sum);
227 };
228
229 auto higham_method_one_pass = [&](MatrixDouble &mat) {
230 int n = mat.size1();
231 const int lda = n;
232 const int size = (n + 2) * n;
233 int lwork = size;
234 double *work = new double[size];
235
236 VectorDouble eig(n, 0.0);
237 MatrixDouble eigen_vec = prod(trans(mat), mat);
238 if (lapack_dsyev('V', 'U', n, eigen_vec.data().data(), lda,
239 eig.data().data(), work, lwork) > 0)
241 "The algorithm failed to compute eigenvalues.");
242 delete[] work;
243
244 for (auto &e : eig) {
245 constexpr double eps = std::numeric_limits<double>::epsilon();
246 e = std::max(std::sqrt(std::abs(e)), eps);
247 }
248
249 auto diagonal_matrix =
250 ublas::diagonal_matrix<double, ublas::row_major,
251 VecAllocator<double>>(n, eig.data());
252 MatrixDouble symm_mat;
253 symm_mat = prod(trans(eigen_vec), diagonal_matrix);
254 symm_mat = prod(symm_mat, eigen_vec);
255
256 MatrixDouble hat_m = (mat + symm_mat) / 2.0;
257
258 return hat_m;
259 };
260
261 auto hat_k = higham_method_one_pass(K);
262 K.swap(hat_k);
263 }
264
265 CHKERR assemble(row_side, col_side, row_type, col_type, row_data, col_data);
267 }
268
269private:
270};
271
272// Add operator to calculate Eshelby tensor
273
275 : public ForcesAndSourcesCore::UserDataOperator {
276 OpCalculateEshelbyStress(boost::shared_ptr<DataAtIntegrationPts> data_ptr)
278 dataAtPts(data_ptr) {}
279 MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
280
281private:
282 boost::shared_ptr<DataAtIntegrationPts>
283 dataAtPts; ///< data at integration pts
284};
285
287
289 boost::shared_ptr<DataAtIntegrationPts> data_ptr, double alpha_omega = 0)
291 alphaOmega(alpha_omega) {}
292 MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
293
294private:
295 boost::shared_ptr<DataAtIntegrationPts>
296 dataAtPts; ///< data at integration pts
298};
299
302 boost::shared_ptr<DataAtIntegrationPts> data_ptr)
303 : SideEleOp(NOSPACE, SideEleOp::OPSPACE), dataAtPts(data_ptr) {
304 std::fill(&doEntities[MBVERTEX], &doEntities[MBMAXTYPE], false);
305 for (auto t = moab::CN::TypeDimensionMap[SPACE_DIM].first;
306 t <= moab::CN::TypeDimensionMap[SPACE_DIM].second; ++t)
307 doEntities[t] = true;
308 sYmm = false;
309 }
310 MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
311
312private:
313 boost::shared_ptr<DataAtIntegrationPts>
314 dataAtPts; ///< data at integration pts
315};
316
317struct OpCalculateReactionForces : public FaceUserDataOperator {
318 OpCalculateReactionForces(boost::shared_ptr<DataAtIntegrationPts> data_ptr,
319 std::string block_name, Range block_entities,
320 std::array<double, 6> &reaction_vec)
321 : FaceUserDataOperator(NOSPACE, OPSPACE), dataAtPts(data_ptr),
322 blockName(block_name), blockEntities(block_entities),
323 reactionVec(reaction_vec) {}
324 MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
325
326private:
327 boost::shared_ptr<DataAtIntegrationPts>
328 dataAtPts; ///< data at integration pts
330 std::string blockName;
331 std::array<double, 6> &reactionVec;
332};
333
336 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
337 const double alpha, const double rho)
338 : OpAssembleVolume(field_name, data_ptr, OPROW), alphaW(alpha),
339 alphaRho(rho) {}
340 MoFEMErrorCode integrate(EntData &data);
341
342private:
343 const double alphaW;
344 const double alphaRho;
345};
346
348 OpSpatialRotation(const std::string &field_name,
349 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
350 double alpha_omega)
351 : OpAssembleVolume(field_name, data_ptr, OPROW), alphaOmega(alpha_omega) {
352 }
353 MoFEMErrorCode integrate(EntData &data);
354
355private:
357};
358
361 boost::shared_ptr<DataAtIntegrationPts> data_ptr)
362 : OpAssembleVolume(field_name, data_ptr, OPROW) {}
363 MoFEMErrorCode integrate(EntData &data);
364
365private:
366};
367
370 boost::shared_ptr<DataAtIntegrationPts> data_ptr)
371 : OpAssembleVolume(field_name, data_ptr, OPROW) {}
372 MoFEMErrorCode integrate(EntData &data);
373
374private:
375};
376
379 boost::shared_ptr<DataAtIntegrationPts> data_ptr)
380 : OpAssembleVolume(field_name, data_ptr, OPROW) {}
381 MoFEMErrorCode integrate(EntData &data);
382};
383
384template <int INTERP_ORDER>
386 OpGetInternalStress(boost::shared_ptr<DataAtIntegrationPts> data_ptr,
387 std::string tag_name = "")
388 : UserDataOperator(H1, OPLAST), dataAtPts(data_ptr), tagName(tag_name) {
389 std::fill(&doEntities[MBVERTEX], &doEntities[MBMAXTYPE], false);
390 doEntities[MBVERTEX] = true;
391 }
392 MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
393
394private:
395 boost::shared_ptr<DataAtIntegrationPts> dataAtPts;
396 std::string tagName;
397};
398
399template <bool VOIGT>
402 const std::string &field_name,
403 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
404 boost::shared_ptr<ScalingMethod> scaling_method_ptr)
405 : OpAssembleVolume(field_name, data_ptr, OPROW),
406 scalingMethodPtr(scaling_method_ptr) {}
407 MoFEMErrorCode integrate(EntData &row_data);
408
409protected:
410 boost::shared_ptr<ScalingMethod> scalingMethodPtr;
411};
412
413template <AssemblyType A, IntegrationType I> struct OpDispBcImpl;
414
415template <AssemblyType A>
416struct OpDispBcImpl<A, GAUSS>
417 : public FormsIntegrators<FaceUserDataOperator>::Assembly<A>::OpBrokenBase {
418
419 using OP = typename FormsIntegrators<FaceUserDataOperator>::Assembly<
421
423 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
424 boost::shared_ptr<BcDispVec> &bc_disp_ptr,
425 std::map<std::string, boost::shared_ptr<ScalingMethod>> smv,
426 boost::shared_ptr<Range> ents_ptr = nullptr)
427 : OP(broken_base_side_data, ents_ptr), bcDispPtr(bc_disp_ptr),
428 scalingMethodsMap(smv) {}
429
430protected:
431 MoFEMErrorCode iNtegrate(EntData &data);
432 boost::shared_ptr<BcDispVec> bcDispPtr;
433 std::map<std::string, boost::shared_ptr<ScalingMethod>> scalingMethodsMap;
434};
435
436struct OpDispBc : public OpDispBcImpl<PETSC, GAUSS> {
438 using OP::OpDispBcImpl;
439
440protected:
441 MoFEMErrorCode iNtegrate(EntData &data);
442};
443
444template <AssemblyType A, IntegrationType I> struct OpRotationBcImpl;
445
446template <AssemblyType A>
447struct OpRotationBcImpl<A, GAUSS>
448 : public FormsIntegrators<FaceUserDataOperator>::Assembly<A>::OpBrokenBase {
449
450 using OP = typename FormsIntegrators<FaceUserDataOperator>::Assembly<
452
454 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
455 boost::shared_ptr<BcRotVec> &bc_rot_ptr,
456 std::map<std::string, boost::shared_ptr<ScalingMethod>> smv,
457 boost::shared_ptr<Range> ents_ptr = nullptr)
458 : OP(broken_base_side_data, ents_ptr), bcRotPtr(bc_rot_ptr),
459 scalingMethodsMap(smv) {}
460
461protected:
462 MoFEMErrorCode iNtegrate(EntData &data);
463 boost::shared_ptr<BcRotVec> bcRotPtr;
464 std::map<std::string, boost::shared_ptr<ScalingMethod>> scalingMethodsMap;
465};
466
467/**
468 * @brief Apply rotation boundary condition
469 *
470 * The boundary condition is applied getting data from the SPATIAL_ROTATION_BC
471 * blockset.
472 *
473 * Centre of rotation is given by the first three attributes of the blockset.
474 * The fort parameter is angle of rotation
475 *
476 */
477struct OpRotationBc : public OpRotationBcImpl<PETSC, GAUSS> {
479 using OP::OpRotationBcImpl;
480
481protected:
482 MoFEMErrorCode iNtegrate(EntData &data);
483};
484
486 : public FormsIntegrators<FaceUserDataOperator>::Assembly<PETSC>::OpBase {
488 const std::string row_field, boost::shared_ptr<TractionBcVec> bc_data,
489 boost::shared_ptr<double> piola_scale_ptr,
490 std::map<std::string, boost::shared_ptr<ScalingMethod>> smv)
491 : FormsIntegrators<FaceUserDataOperator>::Assembly<PETSC>::OpBase(
492 row_field, row_field, FaceUserDataOperator::OPROW),
493 bcData(bc_data), piolaScalePtr(), scalingMethodsMap(smv) {}
494
495 MoFEMErrorCode iNtegrate(EntData &data);
496
497protected:
498 boost::shared_ptr<TractionBcVec> bcData;
499 boost::shared_ptr<double> piolaScalePtr;
500 std::map<std::string, boost::shared_ptr<ScalingMethod>> scalingMethodsMap;
501};
502
504 : public FormsIntegrators<FaceUserDataOperator>::Assembly<PETSC>::OpBase {
506 const std::string row_field, boost::shared_ptr<PressureBcVec> bc_data,
507 boost::shared_ptr<double> piola_scale_ptr,
508 boost::shared_ptr<MatrixDouble> hybrid_grad_disp,
509 std::map<std::string, boost::shared_ptr<ScalingMethod>> smv)
510 : FormsIntegrators<FaceUserDataOperator>::Assembly<PETSC>::OpBase(
511 row_field, row_field, FaceUserDataOperator::OPROW),
512 bcData(bc_data), piolaScalePtr(), hybridGradDispPtr(hybrid_grad_disp),
513 scalingMethodsMap(smv) {}
514
515 MoFEMErrorCode iNtegrate(EntData &data);
516
517protected:
518 boost::shared_ptr<PressureBcVec> bcData;
519 boost::shared_ptr<double> piolaScalePtr;
520 boost::shared_ptr<MatrixDouble> hybridGradDispPtr;
521 std::map<std::string, boost::shared_ptr<ScalingMethod>> scalingMethodsMap;
522};
523
524template <AssemblyType A, IntegrationType I>
526
527template <AssemblyType A>
529 : public FormsIntegrators<FaceUserDataOperator>::Assembly<A>::OpBase {
530
531 using OP =
532 typename FormsIntegrators<FaceUserDataOperator>::Assembly<A>::OpBase;
533
535 std::string row_field, boost::shared_ptr<PressureBcVec> bc_data,
536 boost::shared_ptr<MatrixDouble> hybrid_grad_disp,
537 std::map<std::string, boost::shared_ptr<ScalingMethod>> smv)
538 : OP(row_field, row_field, OP::OPROWCOL), bcData(bc_data),
539 hybridGradDispPtr(hybrid_grad_disp), scalingMethodsMap(smv) {}
540
541protected:
542 MoFEMErrorCode iNtegrate(EntData &row_data, EntData &col_data);
543 boost::shared_ptr<PressureBcVec> bcData;
544 boost::shared_ptr<MatrixDouble> hybridGradDispPtr;
545 std::map<std::string, boost::shared_ptr<ScalingMethod>> scalingMethodsMap;
546};
547
549 : public OpBrokenPressureBcLhsImpl_dU<A, GAUSS> {
552
553protected:
554 MoFEMErrorCode iNtegrate(EntData &row_data, EntData &col_data);
555};
556
558 : public FormsIntegrators<FaceUserDataOperator>::Assembly<PETSC>::OpBase {
560 const std::string row_field,
561 boost::shared_ptr<AnalyticalTractionBcVec> bc_data,
562 boost::shared_ptr<double> piola_scale_ptr,
563 std::map<std::string, boost::shared_ptr<ScalingMethod>> smv)
564 : FormsIntegrators<FaceUserDataOperator>::Assembly<PETSC>::OpBase(
565 row_field, row_field, FaceUserDataOperator::OPROW),
566 bcData(bc_data), piolaScalePtr(), scalingMethodsMap(smv) {}
567
568 MoFEMErrorCode iNtegrate(EntData &data);
569
570protected:
571 boost::shared_ptr<AnalyticalTractionBcVec> bcData;
572 boost::shared_ptr<double> piolaScalePtr;
573 std::map<std::string, boost::shared_ptr<ScalingMethod>> scalingMethodsMap;
574};
575
576template <AssemblyType A, IntegrationType I> struct OpNormalDispBcRhsImpl;
577template <AssemblyType A, IntegrationType I> struct OpNormalDispBcLhsImpl_dU;
578template <AssemblyType A, IntegrationType I> struct OpNormalDispBcLhsImpl_dP;
579
580template <AssemblyType A>
582 : public FormsIntegrators<FaceUserDataOperator>::Assembly<PETSC>::OpBase {
583
584 using OP =
585 typename FormsIntegrators<FaceUserDataOperator>::Assembly<PETSC>::OpBase;
586
588 std::string row_field, boost::shared_ptr<MatrixDouble> hybrid_disp_ptr,
589 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_side_data_ptr,
590 // boost::shared_ptr<MatrixDouble> piola_stress_ptr,
591 boost::shared_ptr<NormalDisplacementBcVec> &bc_disp_ptr,
592 std::map<std::string, boost::shared_ptr<ScalingMethod>> smv,
593 boost::shared_ptr<Range> ents_ptr = nullptr)
594 : FormsIntegrators<FaceUserDataOperator>::Assembly<PETSC>::OpBase(
595 row_field, row_field, OP::OPROW),
596 hybridDispPtr(hybrid_disp_ptr),
597 brokenBaseSideDataPtr(
598 broken_side_data_ptr), /*piolaStressPtr(piola_stress_ptr),*/
599 bcDispPtr(bc_disp_ptr), scalingMethodsMap(smv) {}
600
601protected:
602 MoFEMErrorCode iNtegrate(EntData &data);
603 boost::shared_ptr<MatrixDouble> hybridDispPtr;
604 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenBaseSideDataPtr;
605 // boost::shared_ptr<MatrixDouble> piolaStressPtr;
606 boost::shared_ptr<NormalDisplacementBcVec> bcDispPtr;
607 std::map<std::string, boost::shared_ptr<ScalingMethod>> scalingMethodsMap;
608};
609
610template <AssemblyType A>
612 : public FormsIntegrators<FaceUserDataOperator>::Assembly<A>::OpBase {
613
614 using OP =
615 typename FormsIntegrators<FaceUserDataOperator>::Assembly<A>::OpBase;
616
618 std::string row_field,
619 boost::shared_ptr<NormalDisplacementBcVec> &bc_disp_ptr,
620 std::map<std::string, boost::shared_ptr<ScalingMethod>> smv,
621 boost::shared_ptr<Range> ents_ptr = nullptr)
622 : OP(row_field, row_field, OP::OPROWCOL), bcDispPtr(bc_disp_ptr),
623 scalingMethodsMap(smv) {}
624
625protected:
626 MoFEMErrorCode iNtegrate(EntData &row_data, EntData &col_data);
627 boost::shared_ptr<NormalDisplacementBcVec> bcDispPtr;
628 std::map<std::string, boost::shared_ptr<ScalingMethod>> scalingMethodsMap;
629};
630
631template <AssemblyType A>
633 : public FormsIntegrators<FaceUserDataOperator>::Assembly<A>::OpBrokenBase {
634
635 using OP = typename FormsIntegrators<FaceUserDataOperator>::Assembly<
637
639 std::string row_field,
640 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
641 boost::shared_ptr<NormalDisplacementBcVec> &bc_disp_ptr,
642 std::map<std::string, boost::shared_ptr<ScalingMethod>> smv,
643 boost::shared_ptr<Range> ents_ptr = nullptr)
644 : OP(row_field, broken_base_side_data, false, false, ents_ptr),
645 bcDispPtr(bc_disp_ptr), scalingMethodsMap(smv) {}
646
647protected:
648 MoFEMErrorCode iNtegrate(EntData &row_data, EntData &col_data);
649 boost::shared_ptr<NormalDisplacementBcVec> bcDispPtr;
650 std::map<std::string, boost::shared_ptr<ScalingMethod>> scalingMethodsMap;
651};
652
653struct OpNormalDispRhsBc : public OpNormalDispBcRhsImpl<PETSC, GAUSS> {
655 using OP::OpNormalDispBcRhsImpl;
656
657protected:
658 MoFEMErrorCode iNtegrate(EntData &data);
659};
660
664
665protected:
666 MoFEMErrorCode iNtegrate(EntData &row_data, EntData &col_data);
667};
668
672
673protected:
674 MoFEMErrorCode iNtegrate(EntData &row_data, EntData &col_data);
675};
676
677template <AssemblyType A, IntegrationType I> struct OpAnalyticalDispBcImpl;
678
679template <AssemblyType A>
681 : public FormsIntegrators<FaceUserDataOperator>::Assembly<A>::OpBrokenBase {
682
683 using OP = typename FormsIntegrators<FaceUserDataOperator>::Assembly<
685
687 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
688 boost::shared_ptr<AnalyticalDisplacementBcVec> &bc_disp_ptr,
689 std::map<std::string, boost::shared_ptr<ScalingMethod>> smv,
690 boost::shared_ptr<Range> ents_ptr = nullptr)
691 : OP(broken_base_side_data, ents_ptr), bcDispPtr(bc_disp_ptr),
692 scalingMethodsMap(smv) {}
693
694protected:
695 MoFEMErrorCode iNtegrate(EntData &data);
696 boost::shared_ptr<AnalyticalDisplacementBcVec> bcDispPtr;
697 std::map<std::string, boost::shared_ptr<ScalingMethod>> scalingMethodsMap;
698};
699
700struct OpAnalyticalDispBc : public OpAnalyticalDispBcImpl<PETSC, GAUSS> {
702 using OP::OpAnalyticalDispBcImpl;
703
704protected:
705 MoFEMErrorCode iNtegrate(EntData &data);
706};
707
709 OpSpatialEquilibrium_dw_dP(std::string row_field, std::string col_field,
710 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
711 const bool assemble_off = false)
712 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL,
713 assemble_off) {
714 sYmm = false;
715 }
716 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
717};
718
720 const double alphaW;
721 const double alphaRho;
722 OpSpatialEquilibrium_dw_dw(std::string row_field, std::string col_field,
723 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
724 const double alpha, const double rho)
725 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL, false),
726 alphaW(alpha), alphaRho(rho) {
727 sYmm = true;
728 }
729 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
730};
731
733 OpSpatialPhysical_du_dP(std::string row_field, std::string col_field,
734 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
735 const bool assemble_off = false)
736 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL,
737 assemble_off) {
738 sYmm = false;
739 }
740
741 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
742};
743
745 OpSpatialPhysical_du_dBubble(std::string row_field, std::string col_field,
746 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
747 const bool assemble_off = false)
748 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL,
749 assemble_off) {
750 sYmm = false;
751 }
752
753 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
754};
755
757 using OpAssembleVolume::OpAssembleVolume;
758 OpSpatialPhysical_du_domega(std::string row_field, std::string col_field,
759 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
760 const bool assemble_off)
761 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL,
762 assemble_off) {
763 sYmm = false;
764 }
765 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
766};
767
769 OpSpatialRotation_domega_du(std::string row_field, std::string col_field,
770 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
771 double alpha_omega)
772 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL, false) {
773 sYmm = false;
774 }
775 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
776
777private:
778};
779
781 OpSpatialRotation_domega_dP(std::string row_field, std::string col_field,
782 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
783 const bool assemble_off)
784 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL,
785 assemble_off) {
786 sYmm = false;
787 }
788 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
789};
790
793 std::string row_field, std::string col_field,
794 boost::shared_ptr<DataAtIntegrationPts> data_ptr, const bool assemble_off)
795 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL,
796 assemble_off) {
797 sYmm = false;
798 }
799 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
800};
801
804 std::string row_field, std::string col_field,
805 boost::shared_ptr<DataAtIntegrationPts> data_ptr, double alpha_omega)
806 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL, false),
807 alphaOmega(alpha_omega) {
808 sYmm = false;
809 }
810 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
811
812private:
814};
815
817 OpSpatialConsistency_dP_dP(std::string row_field, std::string col_field,
818 boost::shared_ptr<DataAtIntegrationPts> data_ptr)
819 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL, false) {
820 sYmm = false;
821 }
822 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
823
824private:
825 template <int S>
826 MoFEMErrorCode integrateImpl(EntData &row_data, EntData &col_data);
827};
828
831 std::string row_field, std::string col_field,
832 boost::shared_ptr<DataAtIntegrationPts> data_ptr)
833 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL, false) {
834 sYmm = false;
835 }
836 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
837
838private:
839 template <int S>
840 MoFEMErrorCode integrateImpl(EntData &row_data, EntData &col_data);
841};
842
845 std::string row_field, std::string col_field,
846 boost::shared_ptr<DataAtIntegrationPts> data_ptr)
847 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL, true) {
848 sYmm = false;
849 }
850 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
851
852private:
853 template <int S>
854 MoFEMErrorCode integrateImpl(EntData &row_data, EntData &col_data);
855};
856
859 std::string row_field, std::string col_field,
860 boost::shared_ptr<DataAtIntegrationPts> data_ptr, const bool assemble_off)
861 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL,
862 assemble_off) {
863 sYmm = false;
864 }
865 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
866
867private:
868};
869
872 std::string row_field, std::string col_field,
873 boost::shared_ptr<DataAtIntegrationPts> data_ptr, const bool assemble_off)
874 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL,
875 assemble_off) {
876 sYmm = false;
877 }
878 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
879
880private:
881};
882
884 : public VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator {
885
886 using OP = VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator;
887
888 moab::Interface &postProcMesh;
889 std::vector<EntityHandle> &mapGaussPts;
890 boost::shared_ptr<DataAtIntegrationPts> dataAtPts;
891
892 OpPostProcDataStructure(moab::Interface &post_proc_mesh,
893 std::vector<EntityHandle> &map_gauss_pts,
894 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
895 int sense)
896 : OP(NOSPACE, UserDataOperator::OPSPACE), postProcMesh(post_proc_mesh),
897 mapGaussPts(map_gauss_pts), dataAtPts(data_ptr), tagSense(sense) {}
898
899 MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
900
901private:
903};
904
906 OpSpatialPrj(std::string row_field,
907 boost::shared_ptr<DataAtIntegrationPts> data_ptr)
908 : OpAssembleVolume(row_field, data_ptr, OPROW) {}
909 MoFEMErrorCode integrate(EntData &row_data);
910};
911
913 OpSpatialPrj_dx_dx(std::string row_field, std::string col_field,
914 boost::shared_ptr<DataAtIntegrationPts> data_ptr)
915 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL, false) {
916 // FIXME: That is symmetric
917 sYmm = false;
918 }
919 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
920};
921
923 OpSpatialPrj_dx_dw(std::string row_field, std::string col_field,
924 boost::shared_ptr<DataAtIntegrationPts> data_ptr)
925 : OpAssembleVolume(row_field, col_field, data_ptr, OPROWCOL, false) {
926 sYmm = false;
927 }
928 MoFEMErrorCode integrate(EntData &row_data, EntData &col_data);
929};
930
932 : public VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator {
933
934 // @note This is not on the side, since, integrate energy in the volume
935 using OP = VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator;
936
937 OpFaceSideMaterialForce(boost::shared_ptr<DataAtIntegrationPts> data_ptr)
938 : OP(NOSPACE, OPSPACE), dataAtPts(data_ptr) {}
939
940 MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
941
942private:
943 boost::shared_ptr<DataAtIntegrationPts>
944 dataAtPts; ///< data at integration pts
945};
946
948 : public FaceElementForcesAndSourcesCore::UserDataOperator {
949
950 // @note This is not on the side, since, integrate energy in the volume
951 using OP = FaceElementForcesAndSourcesCore::UserDataOperator;
952
953 OpFaceMaterialForce(boost::shared_ptr<DataAtIntegrationPts> data_ptr)
954 : OP(NOSPACE, OPSPACE), dataAtPts(data_ptr) {}
955
956 MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
957
958private:
959 boost::shared_ptr<DataAtIntegrationPts>
960 dataAtPts; ///< data at integration pts
961};
962
963template <int FE_DIM, int PROBLEM_DIM, int SPACE_DIM> struct AddHOOps;
964
965template <> struct AddHOOps<2, 3, 3> {
966 AddHOOps() = delete;
967 static MoFEMErrorCode
968 add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
969 std::vector<FieldSpace> space, std::string geom_field_name,
970 boost::shared_ptr<Range> crack_front_edges_ptr);
971};
972
973template <> struct AddHOOps<2, 2, 3> {
974 AddHOOps() = delete;
975 static MoFEMErrorCode
976 add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
977 std::vector<FieldSpace> space, std::string geom_field_name,
978 boost::shared_ptr<Range> crack_front_edges_ptr);
979};
980
981template <> struct AddHOOps<3, 3, 3> {
982 AddHOOps() = delete;
983 static MoFEMErrorCode
984 add(boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
985 std::vector<FieldSpace> space, std::string geom_field_name,
986 boost::shared_ptr<Range> crack_front_edges_ptr,
987 boost::shared_ptr<MatrixDouble> jac = nullptr,
988 boost::shared_ptr<VectorDouble> det = nullptr,
989 boost::shared_ptr<MatrixDouble> inv_jac = nullptr);
990};
991
993 FormsIntegrators<FaceElementForcesAndSourcesCore::UserDataOperator>::
994 Assembly<A>::LinearForm<GAUSS>::OpBaseTimesVector<1, SPACE_DIM, 1>;
995
996template <typename OP>
997struct OpStabBrokenBaseImpl : public OpBrokenBaseImpl<OP> {
998
999 using BASE = OpBrokenBaseImpl<OP>;
1000 using BASE::BASE;
1001
1002 MoFEMErrorCode doWork(int row_side, EntityType row_type,
1003 EntitiesFieldData::EntData &row_data) {
1005
1006 if (OP::entsPtr) {
1007 if (OP::entsPtr->find(this->getFEEntityHandle()) == OP::entsPtr->end())
1009 }
1010
1011#ifndef NDEBUG
1012 if (!BASE::brokenBaseSideData) {
1013 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "space not set");
1014 }
1015#endif // NDEBUG
1016
1017 auto do_work_rhs = [this](int row_side, EntityType row_type,
1018 EntitiesFieldData::EntData &row_data) {
1020 // get number of dofs on row
1021 OP::nbRows = row_data.getIndices().size();
1022 if (!OP::nbRows)
1024 // get number of integration points
1025 OP::nbIntegrationPts = OP::getGaussPts().size2();
1026 // get row base functions
1027 OP::nbRowBaseFunctions = OP::getNbOfBaseFunctions(row_data);
1028 // resize and clear the right hand side vector
1029 OP::locF.resize(OP::nbRows, false);
1030 OP::locF.clear();
1031 // integrate local vector
1032 CHKERR this->iNtegrate(row_data);
1033 // assemble local vector
1034 CHKERR this->aSsemble(row_data);
1036 };
1037
1038 auto do_work_lhs = [this](int row_side, int col_side, EntityType row_type,
1039 EntityType col_type,
1040 EntitiesFieldData::EntData &row_data,
1041 EntitiesFieldData::EntData &col_data) {
1043
1044 auto check_if_assemble_transpose = [&] {
1045 if (this->sYmm) {
1046 if (OP::rowSide != OP::colSide || OP::rowType != OP::colType)
1047 return true;
1048 else
1049 return false;
1050 } else if (OP::assembleTranspose) {
1051 return true;
1052 }
1053 return false;
1054 };
1055
1056 OP::rowSide = row_side;
1057 OP::rowType = row_type;
1058 OP::colSide = col_side;
1059 OP::colType = col_type;
1060 OP::nbCols = col_data.getIndices().size();
1061 OP::locMat.resize(OP::nbRows, OP::nbCols, false);
1062 OP::locMat.clear();
1063 CHKERR this->iNtegrate(row_data, col_data);
1064 CHKERR this->aSsemble(row_data, col_data, check_if_assemble_transpose());
1066 };
1067
1068 switch (OP::opType) {
1069 case OP::OPROW:
1070
1071 OP::nbRows = row_data.getIndices().size();
1072 if (!OP::nbRows)
1074 OP::nbIntegrationPts = OP::getGaussPts().size2();
1075 OP::nbRowBaseFunctions = OP::getNbOfBaseFunctions(row_data);
1076
1077 if (!OP::nbRows)
1079
1080 for (auto &bd : *BASE::brokenBaseSideData) {
1081
1082#ifndef NDEBUG
1083 if (!bd.getData().getNSharedPtr(bd.getData().getBase())) {
1084 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1085 "base functions not set");
1086 }
1087#endif
1088
1089 CHKERR do_work_lhs(
1090
1091 // side
1092 row_side, bd.getSide(),
1093
1094 // type
1095 row_type, bd.getType(),
1096
1097 // row_data
1098 row_data, bd.getData()
1099
1100 );
1101 }
1102
1103 break;
1104 case OP::OPSPACE:
1105 for (auto &bd : *BASE::brokenBaseSideData) {
1106 CHKERR do_work_rhs(bd.getSide(), bd.getType(), bd.getData());
1107 }
1108 break;
1109 default:
1111 (std::string("wrong op type ") +
1112 OpBaseDerivativesBase::OpTypeNames[OP::opType])
1113 .c_str());
1114 }
1115
1117 }
1118};
1119
1121 : public OpStabBrokenBaseImpl<OpBaseTimesVectorFace> {
1122
1125 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
1126 boost::shared_ptr<MatrixDouble> vec,
1127 ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
1128 boost::shared_ptr<Range> ents_ptr = nullptr);
1129};
1130
1132 : public OpStabBrokenBaseImpl<OpBaseTimesVectorFace> {
1133
1136 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
1137 ScalarFun beta_coeff = [](double, double, double) constexpr { return 1; },
1138 boost::shared_ptr<Range> ents_ptr = nullptr);
1139
1140protected:
1141 MoFEMErrorCode doWork(int row_side, EntityType row_type,
1142 EntitiesFieldData::EntData &row_data);
1143};
1144
1146
1149 const std::string row_field,
1150 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
1151 ScalarFun beta_coeff, boost::shared_ptr<Range> ents_ptr = nullptr);
1152
1153protected:
1154 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &data);
1155 boost::shared_ptr<std::vector<BrokenBaseSideData>> brokenBaseSideDataPtr;
1156};
1157
1159 FormsIntegrators<FaceElementForcesAndSourcesCore::UserDataOperator>::
1160 Assembly<A>::BiLinearForm<GAUSS>::OpMass<1, SPACE_DIM>;
1161
1162struct OpHyrbridBaseBrokenBase : public OpStabBrokenBaseImpl<OpMassVectorFace> {
1163
1166 std::string row_field,
1167 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
1168 ScalarFun beta, const bool assmb_transpose, const bool only_transpose,
1169 boost::shared_ptr<Range> ents_ptr = nullptr);
1170};
1171
1172struct OpBrokenBaseBrokenBase : public OpStabBrokenBaseImpl<OpMassVectorFace> {
1173
1176 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
1177 ScalarFun beta, boost::shared_ptr<Range> ents_ptr = nullptr);
1178
1179 MoFEMErrorCode doWork(int row_side, EntityType row_type,
1180 EntitiesFieldData::EntData &row_data);
1181};
FormsIntegrators< FaceElementForcesAndSourcesCore::UserDataOperator >::Assembly< A >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpMassVectorFace
FormsIntegrators< FaceElementForcesAndSourcesCore::UserDataOperator >::Assembly< A >::LinearForm< GAUSS >::OpBaseTimesVector< 1, SPACE_DIM, 1 > OpBaseTimesVectorFace
ForcesAndSourcesCore::UserDataOperator UserDataOperator
static const double eps
constexpr int SPACE_DIM
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
FieldSpace
approximation spaces
Definition definitions.h:82
@ H1
continuous field
Definition definitions.h:85
@ NOSPACE
Definition definitions.h:83
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define CHK_MOAB_THROW(err, msg)
Check error code of MoAB function and throw MoFEM exception.
@ MOFEM_IMPOSSIBLE_CASE
Definition definitions.h:35
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
@ MOFEM_INVALID_DATA
Definition definitions.h:36
@ 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.
FTensor::Index< 'i', SPACE_DIM > i
const double n
refractive index of diffusive medium
static __CLPK_integer lapack_dpotrf(char uplo, __CLPK_integer n, __CLPK_doublereal *a, __CLPK_integer lda)
static __CLPK_integer lapack_dsyev(char jobz, char uplo, __CLPK_integer n, __CLPK_doublereal *a, __CLPK_integer lda, __CLPK_doublereal *w, __CLPK_doublereal *work, __CLPK_integer lwork)
FTensor::Index< 'j', 3 > j
constexpr AssemblyType A
constexpr double t
plate stiffness
Definition plate.cpp:58
constexpr auto field_name
OpBaseImpl< PETSC, EdgeEleOp > OpBase
Definition radiation.cpp:29
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpMass
[Only used with Hooke equation (linear material model)]
Definition seepage.cpp:56
FTensor::Index< 'm', 3 > m
static MoFEMErrorCode add(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, std::vector< FieldSpace > space, std::string geom_field_name, boost::shared_ptr< Range > crack_front_edges_ptr)
static MoFEMErrorCode add(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, std::vector< FieldSpace > space, std::string geom_field_name, boost::shared_ptr< Range > crack_front_edges_ptr)
static MoFEMErrorCode add(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, std::vector< FieldSpace > space, std::string geom_field_name, boost::shared_ptr< Range > crack_front_edges_ptr, boost::shared_ptr< MatrixDouble > jac=nullptr, boost::shared_ptr< VectorDouble > det=nullptr, boost::shared_ptr< MatrixDouble > inv_jac=nullptr)
bool sYmm
If true assume that matrix is symmetric structure.
Data on single entity (This is passed as argument to DataOperator::doWork)
const VectorInt & getIndices() const
Get global indices of degrees of freedom on entity.
@ OPROW
operator doWork function is executed on FE rows
@ OPROWCOL
operator doWork is executed on FE rows &columns
@ OPSPACE
operator do Work is execute on space data
boost::shared_ptr< AnalyticalDisplacementBcVec > bcDispPtr
std::map< std::string, boost::shared_ptr< ScalingMethod > > scalingMethodsMap
MoFEMErrorCode iNtegrate(EntData &data)
OpAnalyticalDispBcImpl(boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, boost::shared_ptr< AnalyticalDisplacementBcVec > &bc_disp_ptr, std::map< std::string, boost::shared_ptr< ScalingMethod > > smv, boost::shared_ptr< Range > ents_ptr=nullptr)
MoFEMErrorCode iNtegrate(EntData &data)
MatrixDouble K
local tangent matrix
virtual MoFEMErrorCode integrate(int row_side, EntityType row_type, EntData &data)
VectorDouble nF
local right hand side vector
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
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 &data)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
virtual MoFEMErrorCode assemble(EntData &data)
OpAssembleBasic(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const char type)
virtual MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
boost::function< double()> ScaleOff
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
data at integration pts
OpAssembleBasic(const FieldSpace space)
virtual MoFEMErrorCode assemble(int row_side, EntityType row_type, EntData &data)
OpAssembleBasic(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const char type, const bool assemble_symmetry, ScaleOff scale_off=[]() { return 1;})
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
Operator for bi-linear form, usually to calculate values on left hand side.
typename OP::ScaleOff ScaleOff
static std::map< std::pair< std::string, std::string >, MatrixDouble > mapMatrix
MoFEMErrorCode assemble(int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
boost::shared_ptr< AnalyticalTractionBcVec > bcData
OpBrokenAnalyticalTractionBc(const std::string row_field, boost::shared_ptr< AnalyticalTractionBcVec > bc_data, boost::shared_ptr< double > piola_scale_ptr, std::map< std::string, boost::shared_ptr< ScalingMethod > > smv)
boost::shared_ptr< double > piolaScalePtr
MoFEMErrorCode iNtegrate(EntData &data)
std::map< std::string, boost::shared_ptr< ScalingMethod > > scalingMethodsMap
MoFEMErrorCode doWork(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
MoFEMErrorCode doWork(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
MoFEMErrorCode iNtegrate(EntData &row_data, EntData &col_data)
boost::shared_ptr< MatrixDouble > hybridGradDispPtr
OpBrokenPressureBcLhsImpl_dU(std::string row_field, boost::shared_ptr< PressureBcVec > bc_data, boost::shared_ptr< MatrixDouble > hybrid_grad_disp, std::map< std::string, boost::shared_ptr< ScalingMethod > > smv)
std::map< std::string, boost::shared_ptr< ScalingMethod > > scalingMethodsMap
boost::shared_ptr< PressureBcVec > bcData
MoFEMErrorCode iNtegrate(EntData &row_data, EntData &col_data)
boost::shared_ptr< MatrixDouble > hybridGradDispPtr
boost::shared_ptr< double > piolaScalePtr
MoFEMErrorCode iNtegrate(EntData &data)
OpBrokenPressureBc(const std::string row_field, boost::shared_ptr< PressureBcVec > bc_data, boost::shared_ptr< double > piola_scale_ptr, boost::shared_ptr< MatrixDouble > hybrid_grad_disp, std::map< std::string, boost::shared_ptr< ScalingMethod > > smv)
std::map< std::string, boost::shared_ptr< ScalingMethod > > scalingMethodsMap
boost::shared_ptr< PressureBcVec > bcData
std::map< std::string, boost::shared_ptr< ScalingMethod > > scalingMethodsMap
OpBrokenTractionBc(const std::string row_field, boost::shared_ptr< TractionBcVec > bc_data, boost::shared_ptr< double > piola_scale_ptr, std::map< std::string, boost::shared_ptr< ScalingMethod > > smv)
boost::shared_ptr< double > piolaScalePtr
boost::shared_ptr< TractionBcVec > bcData
MoFEMErrorCode iNtegrate(EntData &data)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
data at integration pts
OpCalculateEshelbyStress(boost::shared_ptr< DataAtIntegrationPts > data_ptr)
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
data at integration pts
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
std::array< double, 6 > & reactionVec
OpCalculateReactionForces(boost::shared_ptr< DataAtIntegrationPts > data_ptr, std::string block_name, Range block_entities, std::array< double, 6 > &reaction_vec)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
data at integration pts
OpCalculateRotationAndSpatialGradient(boost::shared_ptr< DataAtIntegrationPts > data_ptr, double alpha_omega=0)
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
data at integration pts
OpCalculateTractionFromSideEl(boost::shared_ptr< DataAtIntegrationPts > data_ptr)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
std::map< std::string, boost::shared_ptr< ScalingMethod > > scalingMethodsMap
OpDispBcImpl(boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, boost::shared_ptr< BcDispVec > &bc_disp_ptr, std::map< std::string, boost::shared_ptr< ScalingMethod > > smv, boost::shared_ptr< Range > ents_ptr=nullptr)
MoFEMErrorCode iNtegrate(EntData &data)
boost::shared_ptr< BcDispVec > bcDispPtr
MoFEMErrorCode iNtegrate(EntData &data)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
data at integration pts
FaceElementForcesAndSourcesCore::UserDataOperator OP
OpFaceMaterialForce(boost::shared_ptr< DataAtIntegrationPts > data_ptr)
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
Caluclate face material force and normal pressure at gauss points.
OpFaceSideMaterialForce(boost::shared_ptr< DataAtIntegrationPts > data_ptr)
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
data at integration pts
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
OpGetInternalStress(boost::shared_ptr< DataAtIntegrationPts > data_ptr, std::string tag_name="")
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &data)
boost::shared_ptr< std::vector< BrokenBaseSideData > > brokenBaseSideDataPtr
virtual MoFEMErrorCode evaluateLhs(EntData &data)=0
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)
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
data at integration pts
int tAg
adol-c tape
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
boost::shared_ptr< NormalDisplacementBcVec > bcDispPtr
OpNormalDispBcLhsImpl_dP(std::string row_field, boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, boost::shared_ptr< NormalDisplacementBcVec > &bc_disp_ptr, std::map< std::string, boost::shared_ptr< ScalingMethod > > smv, boost::shared_ptr< Range > ents_ptr=nullptr)
std::map< std::string, boost::shared_ptr< ScalingMethod > > scalingMethodsMap
MoFEMErrorCode iNtegrate(EntData &row_data, EntData &col_data)
MoFEMErrorCode iNtegrate(EntData &row_data, EntData &col_data)
std::map< std::string, boost::shared_ptr< ScalingMethod > > scalingMethodsMap
OpNormalDispBcLhsImpl_dU(std::string row_field, boost::shared_ptr< NormalDisplacementBcVec > &bc_disp_ptr, std::map< std::string, boost::shared_ptr< ScalingMethod > > smv, boost::shared_ptr< Range > ents_ptr=nullptr)
boost::shared_ptr< NormalDisplacementBcVec > bcDispPtr
boost::shared_ptr< MatrixDouble > hybridDispPtr
OpNormalDispBcRhsImpl(std::string row_field, boost::shared_ptr< MatrixDouble > hybrid_disp_ptr, boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_side_data_ptr, boost::shared_ptr< NormalDisplacementBcVec > &bc_disp_ptr, std::map< std::string, boost::shared_ptr< ScalingMethod > > smv, boost::shared_ptr< Range > ents_ptr=nullptr)
std::map< std::string, boost::shared_ptr< ScalingMethod > > scalingMethodsMap
MoFEMErrorCode iNtegrate(EntData &data)
boost::shared_ptr< NormalDisplacementBcVec > bcDispPtr
boost::shared_ptr< std::vector< BrokenBaseSideData > > brokenBaseSideDataPtr
MoFEMErrorCode iNtegrate(EntData &row_data, EntData &col_data)
MoFEMErrorCode iNtegrate(EntData &row_data, EntData &col_data)
MoFEMErrorCode iNtegrate(EntData &data)
std::vector< EntityHandle > & mapGaussPts
OpPostProcDataStructure(moab::Interface &post_proc_mesh, std::vector< EntityHandle > &map_gauss_pts, boost::shared_ptr< DataAtIntegrationPts > data_ptr, int sense)
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
MoFEMErrorCode doWork(int side, EntityType type, EntData &data)
VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator OP
OpRotationBcImpl(boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, boost::shared_ptr< BcRotVec > &bc_rot_ptr, std::map< std::string, boost::shared_ptr< ScalingMethod > > smv, boost::shared_ptr< Range > ents_ptr=nullptr)
MoFEMErrorCode iNtegrate(EntData &data)
std::map< std::string, boost::shared_ptr< ScalingMethod > > scalingMethodsMap
boost::shared_ptr< BcRotVec > bcRotPtr
Apply rotation boundary condition.
MoFEMErrorCode iNtegrate(EntData &data)
MoFEMErrorCode integrate(EntData &data)
OpSpatialConsistencyBubble(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr)
MoFEMErrorCode integrate(EntData &data)
OpSpatialConsistencyDivTerm(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr)
MoFEMErrorCode integrate(EntData &data)
OpSpatialConsistencyP(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialConsistency_dBubble_dBubble(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr)
MoFEMErrorCode integrateImpl(EntData &row_data, EntData &col_data)
OpSpatialConsistency_dBubble_dP(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr)
MoFEMErrorCode integrateImpl(EntData &row_data, EntData &col_data)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialConsistency_dBubble_domega(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const bool assemble_off)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
MoFEMErrorCode integrateImpl(EntData &row_data, EntData &col_data)
OpSpatialConsistency_dP_dP(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialConsistency_dP_domega(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const bool assemble_off)
OpSpatialEquilibrium_dw_dP(std::string row_field, 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)
OpSpatialEquilibrium_dw_dw(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha, const double rho)
MoFEMErrorCode integrate(EntData &data)
OpSpatialEquilibrium(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha, const double rho)
MoFEMErrorCode integrate(EntData &row_data)
OpSpatialPhysicalInternalStress(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, boost::shared_ptr< ScalingMethod > scaling_method_ptr)
boost::shared_ptr< ScalingMethod > scalingMethodPtr
OpSpatialPhysical_du_dBubble(std::string row_field, 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(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const bool assemble_off=false)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialPhysical_du_domega(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const bool assemble_off)
OpSpatialPrj_dx_dw(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialPrj_dx_dx(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr)
OpSpatialPrj(std::string row_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr)
MoFEMErrorCode integrate(EntData &row_data)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialRotation_domega_dBubble(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const bool assemble_off)
OpSpatialRotation_domega_dP(std::string row_field, 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(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, double alpha_omega)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialRotation_domega_du(std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, double alpha_omega)
MoFEMErrorCode integrate(EntData &row_data, EntData &col_data)
OpSpatialRotation(const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, double alpha_omega)
MoFEMErrorCode integrate(EntData &data)
MoFEMErrorCode doWork(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
OpBrokenBaseImpl< OP > BASE
double rho
Definition plastic.cpp:144