v0.15.5
Loading...
Searching...
No Matches
dynamic_first_order_con_law.cpp
Go to the documentation of this file.
1/**
2 * \file dynamic_first_order_con_law.cpp
3 * \example mofem/tutorials/adv-4_dynamic_first_order_con_law/dynamic_first_order_con_law.cpp
4 *
5 * Explicit first order conservation laws for solid dynamics
6 *
7 */
8
9#include <MoFEM.hpp>
10#include <MatrixFunction.hpp>
11
12using namespace MoFEM;
13
14template <typename T> inline double trace(FTensor::Tensor2<T, 2, 2> &t_stress) {
15 constexpr double third = boost::math::constants::third<double>();
16 return (t_stress(0, 0) + t_stress(1, 1));
17};
18
19template <typename T> inline double trace(FTensor::Tensor2<T, 3, 3> &t_stress) {
20 return (t_stress(0, 0) + t_stress(1, 1) + t_stress(2, 2));
21};
22
23constexpr int SPACE_DIM =
24 EXECUTABLE_DIMENSION; //< Space dimension of problem, mesh
25
28using DomainEleOp = DomainEle::UserDataOperator;
32
35using BoundaryEleOp = BoundaryEle::UserDataOperator;
37
38template <int DIM> struct PostProcEleByDim;
39
45
51
55
60
63
65 GAUSS>::OpBaseTimesVector<1, SPACE_DIM, 0>;
66
71 SPACE_DIM>;
72
76
79 IntegrationType::GAUSS>::OpBaseTimesVector<1, SPACE_DIM * SPACE_DIM, 1>;
80
84
88
89/** \brief Save field DOFS on vertices/tags
90 */
91
92constexpr double omega = 1.;
93constexpr double young_modulus = 1.;
94constexpr double poisson_ratio = 0.;
95double bulk_modulus_K = young_modulus / (3. * (1. - 2. * poisson_ratio));
97double mu = young_modulus / (2. * (1. + poisson_ratio));
99 ((1. + poisson_ratio) * (1. - 2. * poisson_ratio));
100
101// Operator to Calculate F
102template <int DIM_0, int DIM_1>
104 OpCalculateFStab(boost::shared_ptr<MatrixDouble> def_grad_ptr,
105 boost::shared_ptr<MatrixDouble> def_grad_stab_ptr,
106 boost::shared_ptr<MatrixDouble> def_grad_dot_ptr,
107 double tau_F_ptr, double xi_F_ptr,
108 boost::shared_ptr<MatrixDouble> grad_x_ptr,
109 boost::shared_ptr<MatrixDouble> grad_vel_ptr)
111 defGradPtr(def_grad_ptr), defGradStabPtr(def_grad_stab_ptr),
112 defGradDotPtr(def_grad_dot_ptr), tauFPtr(tau_F_ptr), xiF(xi_F_ptr),
113 gradxPtr(grad_x_ptr), gradVelPtr(grad_vel_ptr) {}
114
115 MoFEMErrorCode doWork(int side, EntityType type,
116 DataForcesAndSourcesCore::EntData &data) {
118 // Define Indicies
121
122 // Number of Gauss points
123 const size_t nb_gauss_pts = getGaussPts().size2();
124
125 defGradStabPtr->resize(DIM_0 * DIM_1, nb_gauss_pts, false);
126 defGradStabPtr->clear();
127
128 // Extract matrix from data matrix
129 auto t_F = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*defGradPtr);
130 auto t_Fstab = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*defGradStabPtr);
131 auto t_F_dot = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*defGradDotPtr);
132
133 // tau_F = alpha deltaT
134 auto tau_F = tauFPtr;
135 double xi_F = xiF;
136 auto t_gradx = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*gradxPtr);
137 auto t_gradVel = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*gradVelPtr);
138
139 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
140 // Stabilised Deformation Gradient
141 t_Fstab(i, j) = t_F(i, j) + tau_F * (t_gradVel(i, j) - t_F_dot(i, j)) +
142 xi_F * (t_gradx(i, j) - t_F(i, j));
143
144 ++t_F;
145 ++t_Fstab;
146 ++t_gradVel;
147 ++t_F_dot;
148
149 ++t_gradx;
150 }
151
153 }
154
155private:
156 double tauFPtr;
157 double xiF;
158 boost::shared_ptr<MatrixDouble> defGradPtr;
159 boost::shared_ptr<MatrixDouble> defGradStabPtr;
160 boost::shared_ptr<MatrixDouble> defGradDotPtr;
161 boost::shared_ptr<MatrixDouble> gradxPtr;
162 boost::shared_ptr<MatrixDouble> gradVelPtr;
163};
164
165// Operator to Calculate P
166template <int DIM_0, int DIM_1>
168 OpCalculatePiola(double shear_modulus, double bulk_modulus, double m_u,
169 double lambda_lamme,
170 boost::shared_ptr<MatrixDouble> first_piola_ptr,
171 boost::shared_ptr<MatrixDouble> def_grad_ptr)
173 shearModulus(shear_modulus), bulkModulus(bulk_modulus), mU(m_u),
174 lammeLambda(lambda_lamme), firstPiolaPtr(first_piola_ptr),
175 defGradPtr(def_grad_ptr) {}
176
177 MoFEMErrorCode doWork(int side, EntityType type,
178 DataForcesAndSourcesCore::EntData &data) {
180 // Define Indicies
184
185 // Define Kronecker Delta
186 constexpr auto t_kd = FTensor::Kronecker_Delta<double>();
187
188 // Number of Gauss points
189 const size_t nb_gauss_pts = getGaussPts().size2();
190
191 // Resize Piola
192 firstPiolaPtr->resize(DIM_0 * DIM_1, nb_gauss_pts, false); // ignatios check
193 firstPiolaPtr->clear();
194
195 // Extract matrix from data matrix
196 auto t_P = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*firstPiolaPtr);
197 auto t_F = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*defGradPtr);
198 const double two_o_three = 2. / 3.;
199 const double trace_t_dk = DIM_0;
200 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
201
202 t_P(i, j) = shearModulus * (t_F(i, j) + t_F(j, i) - 2. * t_kd(i, j) -
203 two_o_three * trace(t_F) * t_kd(i, j) +
204 two_o_three * trace_t_dk * t_kd(i, j)) +
205 bulkModulus * trace(t_F) * t_kd(i, j) -
206 bulkModulus * trace_t_dk * t_kd(i, j);
207
208 ++t_F;
209 ++t_P;
210 }
211
213 }
214
215private:
218 double mU;
220 boost::shared_ptr<MatrixDouble> firstPiolaPtr;
221 boost::shared_ptr<MatrixDouble> defGradPtr;
222};
223
224template <int DIM>
226 OpCalculateDisplacement(boost::shared_ptr<MatrixDouble> spatial_pos_ptr,
227 boost::shared_ptr<MatrixDouble> reference_pos_ptr,
228 boost::shared_ptr<MatrixDouble> u_ptr)
230 xPtr(spatial_pos_ptr), XPtr(reference_pos_ptr), uPtr(u_ptr) {}
231
232 MoFEMErrorCode doWork(int side, EntityType type,
233 DataForcesAndSourcesCore::EntData &data) {
235 // Define Indicies
236 FTensor::Index<'i', DIM> i;
237
238 // Number of Gauss points
239 const size_t nb_gauss_pts = getGaussPts().size2();
240
241 uPtr->resize(DIM, nb_gauss_pts, false); // ignatios check
242 uPtr->clear();
243
244 // Extract matrix from data matrix
245 auto t_x = getFTensor1FromMat<DIM>(*xPtr);
246 auto t_X = getFTensor1FromMat<DIM>(*XPtr);
247 auto t_u = getFTensor1FromMat<DIM>(*uPtr);
248 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
249
250 t_u(i) = t_x(i) - t_X(i);
251 ++t_u;
252 ++t_x;
253 ++t_X;
254 }
255
257 }
258
259private:
260 boost::shared_ptr<MatrixDouble> xPtr;
261 boost::shared_ptr<MatrixDouble> XPtr;
262 boost::shared_ptr<MatrixDouble> uPtr;
263};
264
265template <int DIM_0, int DIM_1>
269 double shear_modulus, double bulk_modulus, double m_u,
270 double lambda_lamme, boost::shared_ptr<MatrixDouble> first_piola_ptr,
271 boost::shared_ptr<MatrixDouble> def_grad_ptr,
272 boost::shared_ptr<MatrixDouble> inv_def_grad_ptr,
273 boost::shared_ptr<VectorDouble> det)
275 shearModulus(shear_modulus), bulkModulus(bulk_modulus), mU(m_u),
276 lammeLambda(lambda_lamme), firstPiolaPtr(first_piola_ptr),
277 defGradPtr(def_grad_ptr), invDefGradPtr(inv_def_grad_ptr), dEt(det) {}
278
279 MoFEMErrorCode doWork(int side, EntityType type,
280 DataForcesAndSourcesCore::EntData &data) {
282 // Define Indicies
287
288 // Define Kronecker Delta
289 constexpr auto t_kd = FTensor::Kronecker_Delta<double>();
290
291 // Number of Gauss points
292 const size_t nb_gauss_pts = getGaussPts().size2();
293
294 // Resize Piola
295 firstPiolaPtr->resize(DIM_0 * DIM_1, nb_gauss_pts, false); // ignatios check
296 firstPiolaPtr->clear();
297
298 // Extract matrix from data matrix
299 auto t_P = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*firstPiolaPtr);
300 auto t_F = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*defGradPtr);
301 auto t_inv_F = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*invDefGradPtr);
302 auto t_det = getFTensor0FromVec<1>(*dEt);
303 const double two_o_three = 2. / 3.;
304 const double one_o_three = 1. / 3.;
305 const double bulk_mod = bulkModulus;
306 const double shear_mod = shearModulus;
307 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
308
309 // Nearly incompressible NH
310 // volumetric part
311 t_P(i, j) = bulk_mod * (t_det - 1.) * t_det * t_inv_F(j, i);
312 // deviatoric part
313 t_P(i, j) +=
314 shear_mod * pow(t_det, two_o_three) *
315 (t_F(i, j) - one_o_three * (t_F(l, k) * t_F(l, k)) * t_inv_F(j, i));
316
317 ++t_F;
318 ++t_P;
319 ++t_inv_F;
320 ++t_det;
321 }
322
324 }
325
326private:
329 double mU;
331 boost::shared_ptr<MatrixDouble> firstPiolaPtr;
332 boost::shared_ptr<MatrixDouble> defGradPtr;
333 boost::shared_ptr<MatrixDouble> invDefGradPtr;
334 boost::shared_ptr<VectorDouble> dEt;
335};
336
337template <int DIM_0, int DIM_1>
341 boost::shared_ptr<MatrixDouble> def_grad_ptr,
342 boost::shared_ptr<MatrixDouble> grad_tensor_ptr)
344 defGradPtr(def_grad_ptr), gradTensorPtr(grad_tensor_ptr) {}
345
346 MoFEMErrorCode doWork(int side, EntityType type,
347 DataForcesAndSourcesCore::EntData &data) {
349 // Define Indicies
352
353 // Define Kronecker Delta
354 constexpr auto t_kd = FTensor::Kronecker_Delta<double>();
355
356 // Number of Gauss points
357 const size_t nb_gauss_pts = getGaussPts().size2();
358
359 // Resize Piola
360 defGradPtr->resize(DIM_0 * DIM_1, nb_gauss_pts, false);
361 defGradPtr->clear();
362
363 // Extract matrix from data matrix
364 auto t_F = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*defGradPtr);
365 auto t_H = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*gradTensorPtr);
366 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
367
368 t_F(i, j) = t_H(i, j) + t_kd(i, j);
369
370 ++t_F;
371 ++t_H;
372 }
373
375 }
376
377private:
378 boost::shared_ptr<MatrixDouble> gradTensorPtr;
379 boost::shared_ptr<MatrixDouble> defGradPtr;
380};
381
382struct Example;
384
385 TSPrePostProc() = default;
386 virtual ~TSPrePostProc() = default;
387
388 /**
389 * @brief Used to setup TS solver
390 *
391 * @param ts
392 * @return MoFEMErrorCode
393 */
394 MoFEMErrorCode tsSetUp(TS ts);
395
396 // SmartPetscObj<VecScatter> getScatter(Vec x, Vec y, enum FR fr);
398 static MoFEMErrorCode tsPostStage(TS ts, PetscReal stagetime,
399 PetscInt stageindex, Vec *Y);
400 static MoFEMErrorCode tsPostStep(TS ts);
401 static MoFEMErrorCode tsPreStep(TS ts);
402};
403
404static boost::weak_ptr<TSPrePostProc> tsPrePostProc;
405
408 double getScale(const double time) {
409 return sin(2. * M_PI * MoFEM::TimeScale::getScale(time));
410 };
411};
412
413struct CommonData {
414 SmartPetscObj<Mat> M; ///< Mass matrix
415 SmartPetscObj<KSP> ksp; ///< Linear solver
416};
417
418struct Example {
419
420 Example(MoFEM::Interface &m_field) : mField(m_field) {}
421
423
424private:
426
434 friend struct TSPrePostProc;
435
438 double getScale(const double time) { return 0.001 * sin(0.1 * time); };
439 };
440
443 double getScale(const double time) { return 0.001; };
444 };
445};
446
447//! [Run problem]
458}
459//! [Run problem]
460
461//! [Read mesh]
470//! [Read mesh]
471
472//! [Set up problem]
476 enum bases { AINSWORTH, DEMKOWICZ, LASBASETOPT };
477 const char *list_bases[LASBASETOPT] = {"ainsworth", "demkowicz"};
478 PetscInt choice_base_value = AINSWORTH;
479 CHKERR PetscOptionsGetEList(PETSC_NULLPTR, NULL, "-base", list_bases,
480 LASBASETOPT, &choice_base_value, PETSC_NULLPTR);
481
483 switch (choice_base_value) {
484 case AINSWORTH:
486 MOFEM_LOG("WORLD", Sev::inform)
487 << "Set AINSWORTH_LEGENDRE_BASE for displacements";
488 break;
489 case DEMKOWICZ:
491 MOFEM_LOG("WORLD", Sev::inform)
492 << "Set DEMKOWICZ_JACOBI_BASE for displacements";
493 break;
494 default:
495 base = LASTBASE;
496 break;
497 }
498 // Add field
499 CHKERR simple->addDomainField("V", H1, base, SPACE_DIM);
500 CHKERR simple->addBoundaryField("V", H1, base, SPACE_DIM);
501 CHKERR simple->addDomainField("F", H1, base, SPACE_DIM * SPACE_DIM);
502 CHKERR simple->addDataField("x_1", H1, base, SPACE_DIM);
503 CHKERR simple->addDataField("x_2", H1, base, SPACE_DIM);
504 CHKERR simple->addDataField("F_0", H1, base, SPACE_DIM * SPACE_DIM);
505 CHKERR simple->addDataField("F_dot", H1, base, SPACE_DIM * SPACE_DIM);
506
507 CHKERR simple->addDataField("GEOMETRY", H1, base, SPACE_DIM);
508 int order = 2;
509 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-order", &order, PETSC_NULLPTR);
510 CHKERR simple->setFieldOrder("V", order);
511 CHKERR simple->setFieldOrder("F", order);
512 CHKERR simple->setFieldOrder("F_0", order);
513 CHKERR simple->setFieldOrder("F_dot", order);
514 CHKERR simple->setFieldOrder("x_1", order);
515 CHKERR simple->setFieldOrder("x_2", order);
516 CHKERR simple->setFieldOrder("GEOMETRY", order);
517 CHKERR simple->setUp();
518
519 auto project_ho_geometry = [&]() {
520 Projection10NodeCoordsOnField ent_method_x(mField, "x_1");
521 CHKERR mField.loop_dofs("x_1", ent_method_x);
522 Projection10NodeCoordsOnField ent_method_x_2(mField, "x_2");
523 CHKERR mField.loop_dofs("x_2", ent_method_x_2);
524
525 Projection10NodeCoordsOnField ent_method(mField, "GEOMETRY");
526 return mField.loop_dofs("GEOMETRY", ent_method);
527 };
528 CHKERR project_ho_geometry();
529
531}
532//! [Set up problem]
533
534//! [Boundary condition]
537
539 auto bc_mng = mField.getInterface<BcManager>();
540 auto *pipeline_mng = mField.getInterface<PipelineManager>();
541 auto time_scale = boost::make_shared<TimeScale>();
542
543 PetscBool sin_time_function = PETSC_FALSE;
544 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-sin_time_function",
545 &sin_time_function, PETSC_NULLPTR);
546
547 if (sin_time_function)
548 time_scale = boost::make_shared<DynamicFirstOrderConsSinusTimeScale>();
549 else
550 time_scale = boost::make_shared<DynamicFirstOrderConsConstantTimeScale>();
551
552 pipeline_mng->getBoundaryExplicitRhsFE().reset();
554 pipeline_mng->getOpBoundaryExplicitRhsPipeline(), {NOSPACE}, "GEOMETRY");
555
557 pipeline_mng->getOpBoundaryExplicitRhsPipeline(), mField, "V",
558 {time_scale}, "FORCE", "PRESSURE", Sev::inform);
559
560 auto integration_rule = [](int, int, int approx_order) {
561 return 2 * approx_order;
562 };
563
564 CHKERR pipeline_mng->setBoundaryExplicitRhsIntegrationRule(integration_rule);
565 CHKERR pipeline_mng->setDomainExplicitRhsIntegrationRule(integration_rule);
566
567 CHKERR bc_mng->removeBlockDOFsOnEntities<DisplacementCubitBcData>(
568 simple->getProblemName(), "V");
569
570 auto get_pre_proc_hook = [&]() {
572 mField, pipeline_mng->getDomainExplicitRhsFE(), {time_scale});
573 };
574 pipeline_mng->getDomainExplicitRhsFE()->preProcessHook = get_pre_proc_hook();
575
577}
578//! [Boundary condition]
579
581 PetscInt stageindex, Vec *Y) {
583 // cerr << "tsPostStage " <<"\n";
584 if (auto ptr = tsPrePostProc.lock()) {
585 auto &m_field = ptr->fsRawPtr->mField;
586
587 auto fb = m_field.getInterface<FieldBlas>();
588 double dt;
589 CHKERR TSGetTimeStep(ts, &dt);
590 double time;
591 CHKERR TSGetTime(ts, &time);
592 PetscInt num_stages;
593 Vec *stage_solutions;
594
595 CHKERR TSGetStages(ts, &num_stages, &stage_solutions);
596 PetscPrintf(PETSC_COMM_WORLD, "Check timestep %d time %e dt %e\n",
597 num_stages, time, dt);
598
599 const double inv_num_step = (double)num_stages;
600 CHKERR fb->fieldCopy(1., "x_1", "x_2");
601 CHKERR fb->fieldAxpy(dt, "V", "x_2");
602 CHKERR fb->fieldCopy(1., "x_2", "x_1");
603
604 CHKERR fb->fieldCopy(-inv_num_step / dt, "F_0", "F_dot");
605 CHKERR fb->fieldAxpy(inv_num_step / dt, "F", "F_dot");
606 CHKERR fb->fieldCopy(1., "F", "F_0");
607 }
609}
610
613
614 if (auto ptr = tsPrePostProc.lock()) {
615 double dt;
616 CHKERR TSGetTimeStep(ts, &dt);
617 double time;
618 CHKERR TSGetTime(ts, &time);
619 }
621}
622
625
626 if (auto ptr = tsPrePostProc.lock()) {
627 double dt;
628 CHKERR TSGetTimeStep(ts, &dt);
629 double time;
630 CHKERR TSGetTime(ts, &time);
631 int step_num;
632 CHKERR TSGetStepNumber(ts, &step_num);
633 }
635}
636
637//! [Push operators to pipeline]
640 auto get_body_force = [this](const double, const double, const double) {
643 t_source(i) = 0.;
644 t_source(0) = 0.1;
645 t_source(1) = 1.;
646 return t_source;
647 };
648
649 // specific time scaling
650 auto get_time_scale = [this](const double time) {
651 return sin(time * omega * M_PI);
652 };
653
654 auto apply_rhs = [&](auto &pip) {
656
658 "GEOMETRY");
659
660 // Calculate Gradient of velocity
661 auto mat_v_grad_ptr = boost::make_shared<MatrixDouble>();
663 "V", mat_v_grad_ptr));
664
665 auto gravity_vector_ptr = boost::make_shared<MatrixDouble>();
666 gravity_vector_ptr->resize(SPACE_DIM, 1);
667 auto set_body_force = [&]() {
670 auto t_force = getFTensor1FromMat<SPACE_DIM, 0>(*gravity_vector_ptr);
671 double unit_weight = 0.;
672 CHKERR PetscOptionsGetReal(PETSC_NULLPTR, "", "-unit_weight", &unit_weight,
673 PETSC_NULLPTR);
674 t_force(i) = 0;
675 if (SPACE_DIM == 2) {
676 t_force(1) = -unit_weight;
677 } else if (SPACE_DIM == 3) {
678 t_force(2) = unit_weight;
679 }
681 };
682
683 CHKERR set_body_force();
684 pip.push_back(new OpBodyForce("V", gravity_vector_ptr,
685 [](double, double, double) { return 1.; }));
686
687 // Calculate unknown F
688 auto mat_H_tensor_ptr = boost::make_shared<MatrixDouble>();
690 "F", mat_H_tensor_ptr));
691
692 // // Calculate F
693 double tau = 0.2;
694 CHKERR PetscOptionsGetReal(PETSC_NULLPTR, "", "-tau", &tau, PETSC_NULLPTR);
695
696 double xi = 0.;
697 CHKERR PetscOptionsGetReal(PETSC_NULLPTR, "", "-xi", &xi, PETSC_NULLPTR);
698
699 // Calculate P stab
700 auto one = [&](const double, const double, const double) {
701 return 3. * bulk_modulus_K;
702 };
703 auto minus_one = [](const double, const double, const double) {
704 return -1.;
705 };
706
707 auto mat_dot_F_tensor_ptr = boost::make_shared<MatrixDouble>();
709 "F_dot", mat_dot_F_tensor_ptr));
710
711 // Calculate Gradient of Spatial Positions
712 auto mat_x_grad_ptr = boost::make_shared<MatrixDouble>();
714 "x_2", mat_x_grad_ptr));
715
716 auto mat_F_tensor_ptr = boost::make_shared<MatrixDouble>();
718 mat_F_tensor_ptr, mat_H_tensor_ptr));
719
720 auto mat_F_stab_ptr = boost::make_shared<MatrixDouble>();
722 mat_F_tensor_ptr, mat_F_stab_ptr, mat_dot_F_tensor_ptr, tau, xi,
723 mat_x_grad_ptr, mat_v_grad_ptr));
724
725 PetscBool is_linear_elasticity = PETSC_TRUE;
726 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-is_linear_elasticity",
727 &is_linear_elasticity, PETSC_NULLPTR);
728
729 auto mat_P_stab_ptr = boost::make_shared<MatrixDouble>();
730 if (is_linear_elasticity) {
733 mat_F_stab_ptr));
734 } else {
735 auto inv_F = boost::make_shared<MatrixDouble>();
736 auto det_ptr = boost::make_shared<VectorDouble>();
737
738 pip.push_back(
739 new OpInvertMatrix<SPACE_DIM>(mat_F_stab_ptr, det_ptr, inv_F));
740
741 // OpCalculatePiolaIncompressibleNH
744 mat_F_stab_ptr, inv_F, det_ptr));
745 }
746
747 pip.push_back(new OpGradTimesTensor2("V", mat_P_stab_ptr, minus_one));
748 pip.push_back(new OpRhsTestPiola("F", mat_v_grad_ptr, one));
749
751 };
752
753 auto *pipeline_mng = mField.getInterface<PipelineManager>();
754 CHKERR apply_rhs(pipeline_mng->getOpDomainExplicitRhsPipeline());
755
756 auto integration_rule = [](int, int, int approx_order) {
757 return 2 * approx_order;
758 };
759 CHKERR pipeline_mng->setDomainExplicitRhsIntegrationRule(integration_rule);
760
762}
763//! [Push operators to pipeline]
764
765/**
766 * @brief Monitor solution
767 *
768 * This functions is called by TS solver at the end of each step. It is used
769 * to output results to the hard drive.
770 */
771
772struct Monitor : public FEMethod {
775 boost::shared_ptr<PostProcEle> post_proc,
776 boost::shared_ptr<PostProcFaceEle> post_proc_bdry,
777 boost::shared_ptr<MatrixDouble> velocity_field_ptr,
778 boost::shared_ptr<MatrixDouble> x2_field_ptr,
779 boost::shared_ptr<MatrixDouble> geometry_field_ptr,
780 std::array<double, 3> pass_field_eval_coords,
781 boost::shared_ptr<SetPtsData> pass_field_eval_data)
782 : dM(dm), mField(m_field), postProc(post_proc),
783 postProcBdy(post_proc_bdry), velocityFieldPtr(velocity_field_ptr),
784 x2FieldPtr(x2_field_ptr), geometryFieldPtr(geometry_field_ptr),
785 fieldEvalCoords(pass_field_eval_coords),
786 fieldEvalData(pass_field_eval_data){};
789
790 auto *simple = mField.getInterface<Simple>();
791
793 ->evalFEAtThePoint<SPACE_DIM>(
794 fieldEvalCoords.data(), 1e-12, simple->getProblemName(),
795 simple->getDomainFEName(), fieldEvalData, mField.get_comm_rank(),
797
798 if (velocityFieldPtr->size1()) {
799 auto t_vel = getFTensor1FromMat<SPACE_DIM>(*velocityFieldPtr);
800 auto t_x2_field = getFTensor1FromMat<SPACE_DIM>(*x2FieldPtr);
801 auto t_geom = getFTensor1FromMat<SPACE_DIM>(*geometryFieldPtr);
802
803 double u_x = t_x2_field(0) - t_geom(0);
804 double u_y = t_x2_field(1) - t_geom(1);
805 double u_z = t_x2_field(2) - t_geom(2);
806
807 MOFEM_LOG("SYNC", Sev::inform)
808 << "Velocities x: " << t_vel(0) << " y: " << t_vel(1)
809 << " z: " << t_vel(2) << "\n";
810 MOFEM_LOG("SYNC", Sev::inform) << "Displacement x: " << u_x
811 << " y: " << u_y << " z: " << u_z << "\n";
812 }
813
814 for (auto m : mField.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(
815 std::regex((boost::format("%s(.*)") % "Data_Vertex").str()))) {
816 Range ents;
817 mField.get_moab().get_entities_by_dimension(m->getMeshset(), 0, ents,
818 true);
819 auto print_vets = [](boost::shared_ptr<FieldEntity> ent_ptr) {
821 if (!(ent_ptr->getPStatus() & PSTATUS_NOT_OWNED)) {
822 MOFEM_LOG("SYNC", Sev::inform)
823 << "Velocities: " << ent_ptr->getEntFieldData()[0] << " "
824 << ent_ptr->getEntFieldData()[1] << " "
825 << ent_ptr->getEntFieldData()[2] << "\n";
826 }
828 };
829 CHKERR mField.getInterface<FieldBlas>()->fieldLambdaOnEntities(
830 print_vets, "V", &ents);
831 }
833
834 PetscBool print_volume = PETSC_FALSE;
835 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-print_volume", &print_volume,
836 PETSC_NULLPTR);
837
838 PetscBool print_skin = PETSC_FALSE;
839 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-print_skin", &print_skin,
840 PETSC_NULLPTR);
841
842 int save_every_nth_step = 1;
843 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-save_step",
844 &save_every_nth_step, PETSC_NULLPTR);
845 if (ts_step % save_every_nth_step == 0) {
846
847 if (print_volume) {
849 CHKERR postProc->writeFile(
850 "out_step_" + boost::lexical_cast<std::string>(ts_step) + ".h5m");
851 }
852
853 if (print_skin) {
855 CHKERR postProcBdy->writeFile(
856 "out_boundary_" + boost::lexical_cast<std::string>(ts_step) +
857 ".h5m");
858 }
859 }
861 }
862
863private:
865 boost::shared_ptr<PostProcEle> postProc;
866 boost::shared_ptr<PostProcFaceEle> postProcBdy;
867 boost::shared_ptr<MatrixDouble> velocityFieldPtr;
868 boost::shared_ptr<MatrixDouble> x2FieldPtr;
869 boost::shared_ptr<MatrixDouble> geometryFieldPtr;
870 std::array<double, 3> fieldEvalCoords;
871 boost::shared_ptr<SetPtsData> fieldEvalData;
872};
873
874//! [Solve]
877 auto *simple = mField.getInterface<Simple>();
878 auto *pipeline_mng = mField.getInterface<PipelineManager>();
879
880 auto dm = simple->getDM();
881
882 auto calculate_stress_ops = [&](auto &pip) {
884
885 auto v_ptr = boost::make_shared<MatrixDouble>();
886 pip.push_back(new OpCalculateVectorFieldValues<SPACE_DIM>("V", v_ptr));
887 auto X_ptr = boost::make_shared<MatrixDouble>();
888 pip.push_back(
889 new OpCalculateVectorFieldValues<SPACE_DIM>("GEOMETRY", X_ptr));
890
891 auto x_ptr = boost::make_shared<MatrixDouble>();
892 pip.push_back(new OpCalculateVectorFieldValues<SPACE_DIM>("x_1", x_ptr));
893
894 // Calculate unknown F
895 auto mat_H_tensor_ptr = boost::make_shared<MatrixDouble>();
897 "F", mat_H_tensor_ptr));
898
899 auto u_ptr = boost::make_shared<MatrixDouble>();
900 pip.push_back(new OpCalculateDisplacement<SPACE_DIM>(x_ptr, X_ptr, u_ptr));
901 // Calculate P
902
903 auto mat_F_ptr = boost::make_shared<MatrixDouble>();
905 mat_F_ptr, mat_H_tensor_ptr));
906
907 PetscBool is_linear_elasticity = PETSC_TRUE;
908 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-is_linear_elasticity",
909 &is_linear_elasticity, PETSC_NULLPTR);
910
911 auto mat_P_ptr = boost::make_shared<MatrixDouble>();
912 if (is_linear_elasticity) {
915 mat_F_ptr));
916 } else {
917 auto inv_F = boost::make_shared<MatrixDouble>();
918 auto det_ptr = boost::make_shared<VectorDouble>();
919
920 pip.push_back(new OpInvertMatrix<SPACE_DIM>(mat_F_ptr, det_ptr, inv_F));
921
924 mat_F_ptr, inv_F, det_ptr));
925 }
926
927 auto mat_v_grad_ptr = boost::make_shared<MatrixDouble>();
929 "V", mat_v_grad_ptr));
930
931 return boost::make_tuple(v_ptr, X_ptr, x_ptr, mat_P_ptr, mat_F_ptr, u_ptr);
932 };
933
934 auto post_proc_boundary = [&]() {
935 auto boundary_post_proc_fe = boost::make_shared<PostProcFaceEle>(mField);
936
938 boundary_post_proc_fe->getOpPtrVector(), {}, "GEOMETRY");
939 auto op_loop_side =
940 new OpLoopSide<SideEle>(mField, simple->getDomainFEName(), SPACE_DIM);
941 // push ops to side element, through op_loop_side operator
942 auto [boundary_v_ptr, boundary_X_ptr, boundary_x_ptr, boundary_mat_P_ptr,
943 boundary_mat_F_ptr, boundary_u_ptr] =
944 calculate_stress_ops(op_loop_side->getOpPtrVector());
945 boundary_post_proc_fe->getOpPtrVector().push_back(op_loop_side);
946
948
949 boundary_post_proc_fe->getOpPtrVector().push_back(
950
951 new OpPPMap(
952
953 boundary_post_proc_fe->getPostProcMesh(),
954 boundary_post_proc_fe->getMapGaussPts(),
955
957
958 OpPPMap::DataMapMat{{"V", boundary_v_ptr},
959 {"GEOMETRY", boundary_X_ptr},
960 {"x", boundary_x_ptr},
961 {"U", boundary_u_ptr}},
962
963 OpPPMap::DataMapMat{{"FIRST_PIOLA", boundary_mat_P_ptr},
964 {"F", boundary_mat_F_ptr}},
965
967
968 )
969
970 );
971 return boundary_post_proc_fe;
972 };
973
974 // Add monitor to time solver
975
976 double rho = 1.;
977 CHKERR PetscOptionsGetReal(PETSC_NULLPTR, "", "-density", &rho, PETSC_NULLPTR);
978 auto get_rho = [rho](const double, const double, const double) {
979 return rho;
980 };
981
982 SmartPetscObj<Mat> M; ///< Mass matrix
983 SmartPetscObj<KSP> ksp; ///< Linear solver
984
985 auto ts_pre_post_proc = boost::make_shared<TSPrePostProc>();
986 tsPrePostProc = ts_pre_post_proc;
987
989 CHKERR MatZeroEntries(M);
990
991 boost::shared_ptr<DomainEle> vol_mass_ele(new DomainEle(mField));
992
993 vol_mass_ele->B = M;
994
995 auto integration_rule = [](int, int, int approx_order) {
996 return 2 * approx_order;
997 };
998
999 vol_mass_ele->getRuleHook = integration_rule;
1000
1001 vol_mass_ele->getOpPtrVector().push_back(new OpMassV("V", "V", get_rho));
1002 vol_mass_ele->getOpPtrVector().push_back(new OpMassF("F", "F"));
1003
1004 CHKERR DMoFEMLoopFiniteElements(dm, simple->getDomainFEName(), vol_mass_ele);
1005 CHKERR MatAssemblyBegin(M, MAT_FINAL_ASSEMBLY);
1006 CHKERR MatAssemblyEnd(M, MAT_FINAL_ASSEMBLY);
1007
1008 auto lumpVec = createDMVector(simple->getDM());
1009 CHKERR MatGetRowSum(M, lumpVec);
1010
1011 CHKERR MatZeroEntries(M);
1012 CHKERR MatDiagonalSet(M, lumpVec, INSERT_VALUES);
1013
1014 // Create and septup KSP (linear solver), we need this to calculate g(t,u) =
1015 // M^-1G(t,u)
1016 ksp = createKSP(mField.get_comm());
1017 CHKERR KSPSetOperators(ksp, M, M);
1018 CHKERR KSPSetFromOptions(ksp);
1019 CHKERR KSPSetUp(ksp);
1020
1021 auto solve_boundary_for_g = [&]() {
1023 if (*(pipeline_mng->getBoundaryExplicitRhsFE()->vecAssembleSwitch)) {
1024
1025 CHKERR VecGhostUpdateBegin(pipeline_mng->getBoundaryExplicitRhsFE()->ts_F,
1026 ADD_VALUES, SCATTER_REVERSE);
1027 CHKERR VecGhostUpdateEnd(pipeline_mng->getBoundaryExplicitRhsFE()->ts_F,
1028 ADD_VALUES, SCATTER_REVERSE);
1029 CHKERR VecAssemblyBegin(pipeline_mng->getBoundaryExplicitRhsFE()->ts_F);
1030 CHKERR VecAssemblyEnd(pipeline_mng->getBoundaryExplicitRhsFE()->ts_F);
1031 *(pipeline_mng->getBoundaryExplicitRhsFE()->vecAssembleSwitch) = false;
1032
1033 auto D =
1034 vectorDuplicate(pipeline_mng->getBoundaryExplicitRhsFE()->ts_F);
1035 CHKERR KSPSolve(ksp, pipeline_mng->getBoundaryExplicitRhsFE()->ts_F, D);
1036 CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
1037 CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
1038 CHKERR VecCopy(D, pipeline_mng->getBoundaryExplicitRhsFE()->ts_F);
1039 }
1040
1042 };
1043
1044 pipeline_mng->getBoundaryExplicitRhsFE()->postProcessHook =
1045 solve_boundary_for_g;
1046
1048 ts = pipeline_mng->createTSEX(dm);
1049
1050 // Field eval
1051 PetscBool field_eval_flag = PETSC_TRUE;
1052 boost::shared_ptr<MatrixDouble> velocity_field_ptr;
1053 boost::shared_ptr<MatrixDouble> geometry_field_ptr;
1054 boost::shared_ptr<MatrixDouble> spatial_position_field_ptr;
1055 boost::shared_ptr<SetPtsData> field_eval_data;
1056
1057 std::array<double, 3> field_eval_coords = {0.5, 0.5, 5.};
1058 int dim = 3;
1059 CHKERR PetscOptionsGetRealArray(NULL, NULL, "-field_eval_coords",
1060 field_eval_coords.data(), &dim,
1061 &field_eval_flag);
1062
1063 if (field_eval_flag) {
1064 field_eval_data =
1065 mField.getInterface<FieldEvaluatorInterface>()->getData<DomainEle>();
1066 CHKERR mField.getInterface<FieldEvaluatorInterface>()->buildTree<SPACE_DIM>(
1067 field_eval_data, simple->getDomainFEName());
1068
1069 field_eval_data->setEvalPoints(field_eval_coords.data(), 1);
1070
1071 auto no_rule = [](int, int, int) { return -1; };
1072
1073 auto fe_ptr = field_eval_data->feMethodPtr.lock();
1074 fe_ptr->getRuleHook = no_rule;
1075 velocity_field_ptr = boost::make_shared<MatrixDouble>();
1076 geometry_field_ptr = boost::make_shared<MatrixDouble>();
1077 spatial_position_field_ptr = boost::make_shared<MatrixDouble>();
1078 fe_ptr->getOpPtrVector().push_back(
1079 new OpCalculateVectorFieldValues<SPACE_DIM>("V", velocity_field_ptr));
1080 fe_ptr->getOpPtrVector().push_back(
1082 geometry_field_ptr));
1083 fe_ptr->getOpPtrVector().push_back(
1085 "x_2", spatial_position_field_ptr));
1086 }
1087
1088 auto post_proc_domain = [&]() {
1089 auto post_proc_fe_vol = boost::make_shared<PostProcEle>(mField);
1090
1092
1093 auto [boundary_v_ptr, boundary_X_ptr, boundary_x_ptr, boundary_mat_P_ptr,
1094 boundary_mat_F_ptr, boundary_u_ptr] =
1095 calculate_stress_ops(post_proc_fe_vol->getOpPtrVector());
1096
1097 post_proc_fe_vol->getOpPtrVector().push_back(
1098
1099 new OpPPMap(
1100
1101 post_proc_fe_vol->getPostProcMesh(),
1102 post_proc_fe_vol->getMapGaussPts(),
1103
1104 {},
1105
1106 {{"V", boundary_v_ptr},
1107 {"GEOMETRY", boundary_X_ptr},
1108 {"x", boundary_x_ptr},
1109 {"U", boundary_u_ptr}},
1110
1111 {{"FIRST_PIOLA", boundary_mat_P_ptr}, {"F", boundary_mat_F_ptr}},
1112
1113 {}
1114
1115 )
1116
1117 );
1118 return post_proc_fe_vol;
1119 };
1120
1121 boost::shared_ptr<FEMethod> null_fe;
1122 auto monitor_ptr = boost::make_shared<Monitor>(
1123 SmartPetscObj<DM>(dm, true), mField, post_proc_domain(),
1124 post_proc_boundary(), velocity_field_ptr, spatial_position_field_ptr,
1125 geometry_field_ptr, field_eval_coords, field_eval_data);
1126
1127 CHKERR DMMoFEMTSSetMonitor(dm, ts, simple->getDomainFEName(), null_fe,
1128 null_fe, monitor_ptr);
1129
1130 double ftime = 1;
1131 // CHKERR TSSetMaxTime(ts, ftime);
1132 CHKERR TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);
1133
1134 auto T = createDMVector(simple->getDM());
1135 CHKERR DMoFEMMeshToLocalVector(simple->getDM(), T, INSERT_VALUES,
1136 SCATTER_FORWARD);
1137 CHKERR TSSetSolution(ts, T);
1138 CHKERR TSSetFromOptions(ts);
1139
1140 CHKERR TSSetPostStage(ts, TSPrePostProc::tsPostStage);
1141 CHKERR TSSetPostStep(ts, TSPrePostProc::tsPostStep);
1142 CHKERR TSSetPreStep(ts, TSPrePostProc::tsPreStep);
1143
1144 boost::shared_ptr<ForcesAndSourcesCore> null;
1145
1146 if (auto ptr = tsPrePostProc.lock()) {
1147 ptr->fsRawPtr = this;
1148 CHKERR TSSetUp(ts);
1149 CHKERR TSSolve(ts, NULL);
1150 CHKERR TSGetTime(ts, &ftime);
1151 }
1152
1154}
1155//! [Solve]
1156
1157//! [Postprocess results]
1160 PetscBool test_flg = PETSC_FALSE;
1161 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-test", &test_flg, PETSC_NULLPTR);
1162 if (test_flg) {
1163 auto *simple = mField.getInterface<Simple>();
1164 auto T = createDMVector(simple->getDM());
1165 CHKERR DMoFEMMeshToLocalVector(simple->getDM(), T, INSERT_VALUES,
1166 SCATTER_FORWARD);
1167 double nrm2;
1168 CHKERR VecNorm(T, NORM_2, &nrm2);
1169 MOFEM_LOG("EXAMPLE", Sev::inform) << "Regression norm " << nrm2;
1170 constexpr double regression_value = 0.0194561;
1171 if (fabs(nrm2 - regression_value) > 1e-2)
1172 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
1173 "Regression test failed; wrong norm value.");
1174 }
1176}
1177//! [Postprocess results]
1178
1179//! [Check]
1184//! [Check]
1185
1186static char help[] = "...\n\n";
1187
1188int main(int argc, char *argv[]) {
1189
1190 // Initialisation of MoFEM/PETSc and MOAB data structures
1191 const char param_file[] = "param_file.petsc";
1192 MoFEM::Core::Initialize(&argc, &argv, param_file, help);
1193
1194 // Add logging channel for example
1195 auto core_log = logging::core::get();
1196 core_log->add_sink(
1198 LogManager::setLog("EXAMPLE");
1199 MOFEM_LOG_TAG("EXAMPLE", "example");
1200
1201 try {
1202
1203 //! [Register MoFEM discrete manager in PETSc]
1204 DMType dm_name = "DMMOFEM";
1205 CHKERR DMRegister_MoFEM(dm_name);
1206 //! [Register MoFEM discrete manager in PETSc
1207
1208 //! [Create MoAB]
1209 moab::Core mb_instance; ///< mesh database
1210 moab::Interface &moab = mb_instance; ///< mesh database interface
1211 //! [Create MoAB]
1212
1213 //! [Create MoFEM]
1214 MoFEM::Core core(moab); ///< finite element database
1215 MoFEM::Interface &m_field = core; ///< finite element database interface
1216 //! [Create MoFEM]
1217
1218 //! [Example]
1219 Example ex(m_field);
1220 CHKERR ex.runProblem();
1221 //! [Example]
1222 }
1224
1226}
constexpr double third
#define MOFEM_LOG_SEVERITY_SYNC(comm, severity)
Synchronise "SYNC" on curtain severity level.
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
int main()
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
Kronecker Delta class.
@ QUIET
#define CATCH_ERRORS
Catch errors.
@ MF_EXIST
FieldApproximationBase
approximation base
Definition definitions.h:58
@ LASTBASE
Definition definitions.h:69
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
@ DEMKOWICZ_JACOBI_BASE
Definition definitions.h:66
@ 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 ...
@ MOFEM_ATOM_TEST_INVALID
Definition definitions.h:40
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
constexpr int order
PostProcEleByDim< SPACE_DIM >::PostProcEleDomain PostProcEleDomain
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpMassV
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM *SPACE_DIM > OpMassF
static boost::weak_ptr< TSPrePostProc > tsPrePostProc
static char help[]
[Check]
constexpr int SPACE_DIM
FormsIntegrators< DomainEleOp >::Assembly< AssemblyType::PETSC >::LinearForm< IntegrationType::GAUSS >::OpBaseTimesVector< 1, SPACE_DIM *SPACE_DIM, 1 > OpRhsTestPiola
constexpr double poisson_ratio
constexpr double omega
Save field DOFS on vertices/tags.
PostProcEleByDim< SPACE_DIM >::PostProcEleBdy PostProcEleBdy
PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::DomainEle DomainEle
double trace(FTensor::Tensor2< T, 2, 2 > &t_stress)
FormsIntegrators< DomainEleOp >::Assembly< AssemblyType::PETSC >::LinearForm< IntegrationType::GAUSS >::OpGradTimesTensor< 1, SPACE_DIM, SPACE_DIM > OpGradTimesPiola
double bulk_modulus_K
FormsIntegrators< DomainEleOp >::Assembly< AssemblyType::PETSC >::LinearForm< IntegrationType::GAUSS >::OpGradTimesTensor< 1, SPACE_DIM, SPACE_DIM > OpGradTimesTensor2
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesVector< 1, SPACE_DIM, 1 > OpInertiaForce
double shear_modulus_G
constexpr double young_modulus
auto integration_rule
constexpr auto t_kd
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition DMMoFEM.cpp:514
PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M)
Definition DMMoFEM.cpp:1187
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition DMMoFEM.cpp:43
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition DMMoFEM.cpp:576
auto createDMVector(DM dm)
Get smart vector from DM.
Definition DMMoFEM.hpp:1234
@ PETSC
Standard PETSc assembly.
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
virtual MoFEMErrorCode loop_dofs(const Problem *problem_ptr, const std::string &field_name, RowColData rc, DofMethod &method, int lower_rank, int upper_rank, int verb=DEFAULT_VERBOSITY)=0
Make a loop over dofs.
FTensor::Index< 'i', SPACE_DIM > i
double dt
double D
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
auto createKSP(MPI_Comm comm)
PetscErrorCode DMMoFEMTSSetMonitor(DM dm, TS ts, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
Set Monitor To TS solver.
Definition DMMoFEM.cpp:1046
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
PetscErrorCode PetscOptionsGetReal(PetscOptions *, const char pre[], const char name[], PetscReal *dval, PetscBool *set)
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
PetscErrorCode PetscOptionsGetRealArray(PetscOptions *, const char pre[], const char name[], PetscReal dval[], PetscInt *nmax, PetscBool *set)
MoFEM::OpBrokenTopoBase GAUSS
Gaussian quadrature integration.
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
FTensor::Index< 'M', 3 > M
FormsIntegrators< DomainEleOp >::Assembly< A >::LinearForm< I >::OpGradTimesTensor< 1, FIELD_DIM, SPACE_DIM > OpGradTimesTensor
int save_every_nth_step
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
static constexpr int approx_order
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
SmartPetscObj< Mat > M
Mass matrix.
SmartPetscObj< KSP > ksp
Linear solver.
double getScale(const double time)
Get scaling at given time.
double getScale(const double time)
Get scaling at given time.
[Example]
Definition plastic.cpp:216
MoFEMErrorCode boundaryCondition()
[Set up problem]
MoFEMErrorCode assembleSystem()
[Push operators to pipeline]
MoFEMErrorCode readMesh()
[Run problem]
FieldApproximationBase base
Choice of finite element basis functions
Definition plot_base.cpp:68
MoFEMErrorCode checkResults()
[Postprocess results]
MoFEMErrorCode solveSystem()
[Solve]
Example(MoFEM::Interface &m_field)
MoFEMErrorCode runProblem()
MoFEM::Interface & mField
Reference to MoFEM interface.
Definition plastic.cpp:226
MoFEMErrorCode setupProblem()
MoFEMErrorCode outputResults()
[Solve]
double getScale(const double time)
Get scaling at given time.
Add operators pushing bases from local to physical configuration.
boost::weak_ptr< CacheTuple > getCacheWeakPtr() const
Get the cache weak pointer object.
Boundary condition manager for finite element problem setup.
virtual moab::Interface & get_moab()=0
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0
Core (interface) class.
Definition Core.hpp:82
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition Core.cpp:72
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:118
Deprecated interface functions.
Definition of the displacement bc data structure.
Definition BCData.hpp:72
Data on single entity (This is passed as argument to DataOperator::doWork)
Class (Function) to enforce essential constrains.
Definition Essential.hpp:25
Structure for user loop methods on finite elements.
Basic algebra on fields.
Definition FieldBlas.hpp:21
Field evaluator interface.
SetIntegrationPtsMethodData SetPtsData
structure to get information from mofem into EntitiesFieldData
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
Interface for managing meshsets containing materials and boundary conditions.
Assembly methods.
Definition Natural.hpp:65
Get values at integration pts for tensor field rank 2, i.e. matrix field.
Get field gradients at integration pts for scalar field rank 0, i.e. vector field.
Specialization for double precision vector field values calculation.
Operator for inverting matrices at integration points.
Element used to execute operators on side of the element.
Post post-proc data at points from hash maps.
std::map< std::string, boost::shared_ptr< MatrixDouble > > DataMapMat
std::map< std::string, boost::shared_ptr< VectorDouble > > DataMapVec
Template struct for dimension-specific finite element types.
PipelineManager interface.
Projection of edge entities with one mid-node on hierarchical basis.
Simple interface for fast problem set-up.
Definition Simple.hpp:27
MoFEMErrorCode getOptions()
get options
Definition Simple.cpp:180
intrusive_ptr for managing petsc objects
PetscInt ts_step
Current time step number.
Force scale operator for reading two columns.
double getScale(const double time)
Get scaling at a given time.
TimeScale(std::string file_name="", bool error_if_file_not_given=false, ScalingFun def_scaling_fun=[](double time) { return time;})
TimeScale constructor.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
[Push operators to pipeline]
boost::shared_ptr< PostProcFaceEle > postProcBdy
std::array< double, 3 > fieldEvalCoords
MoFEM::Interface & mField
Monitor(SmartPetscObj< DM > dm, MoFEM::Interface &m_field, boost::shared_ptr< PostProcEle > post_proc, boost::shared_ptr< PostProcFaceEle > post_proc_bdry, boost::shared_ptr< MatrixDouble > velocity_field_ptr, boost::shared_ptr< MatrixDouble > x2_field_ptr, boost::shared_ptr< MatrixDouble > geometry_field_ptr, std::array< double, 3 > pass_field_eval_coords, boost::shared_ptr< SetPtsData > pass_field_eval_data)
boost::shared_ptr< MatrixDouble > geometryFieldPtr
SmartPetscObj< DM > dM
MoFEMErrorCode postProcess()
Post-processing function executed at loop completion.
boost::shared_ptr< MatrixDouble > velocityFieldPtr
boost::shared_ptr< SetPtsData > fieldEvalData
boost::shared_ptr< MatrixDouble > x2FieldPtr
boost::shared_ptr< PostProcEle > postProc
MoFEMErrorCode doWork(int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
boost::shared_ptr< MatrixDouble > defGradPtr
boost::shared_ptr< MatrixDouble > gradTensorPtr
OpCalculateDeformationGradient(boost::shared_ptr< MatrixDouble > def_grad_ptr, boost::shared_ptr< MatrixDouble > grad_tensor_ptr)
boost::shared_ptr< MatrixDouble > XPtr
boost::shared_ptr< MatrixDouble > uPtr
MoFEMErrorCode doWork(int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
OpCalculateDisplacement(boost::shared_ptr< MatrixDouble > spatial_pos_ptr, boost::shared_ptr< MatrixDouble > reference_pos_ptr, boost::shared_ptr< MatrixDouble > u_ptr)
boost::shared_ptr< MatrixDouble > xPtr
boost::shared_ptr< MatrixDouble > gradxPtr
MoFEMErrorCode doWork(int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
OpCalculateFStab(boost::shared_ptr< MatrixDouble > def_grad_ptr, boost::shared_ptr< MatrixDouble > def_grad_stab_ptr, boost::shared_ptr< MatrixDouble > def_grad_dot_ptr, double tau_F_ptr, double xi_F_ptr, boost::shared_ptr< MatrixDouble > grad_x_ptr, boost::shared_ptr< MatrixDouble > grad_vel_ptr)
boost::shared_ptr< MatrixDouble > defGradStabPtr
boost::shared_ptr< MatrixDouble > gradVelPtr
boost::shared_ptr< MatrixDouble > defGradPtr
boost::shared_ptr< MatrixDouble > defGradDotPtr
OpCalculatePiolaIncompressibleNH(double shear_modulus, double bulk_modulus, double m_u, double lambda_lamme, boost::shared_ptr< MatrixDouble > first_piola_ptr, boost::shared_ptr< MatrixDouble > def_grad_ptr, boost::shared_ptr< MatrixDouble > inv_def_grad_ptr, boost::shared_ptr< VectorDouble > det)
boost::shared_ptr< VectorDouble > dEt
MoFEMErrorCode doWork(int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
boost::shared_ptr< MatrixDouble > invDefGradPtr
boost::shared_ptr< MatrixDouble > defGradPtr
boost::shared_ptr< MatrixDouble > firstPiolaPtr
OpCalculatePiola(double shear_modulus, double bulk_modulus, double m_u, double lambda_lamme, boost::shared_ptr< MatrixDouble > first_piola_ptr, boost::shared_ptr< MatrixDouble > def_grad_ptr)
MoFEMErrorCode doWork(int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
boost::shared_ptr< MatrixDouble > defGradPtr
boost::shared_ptr< MatrixDouble > firstPiolaPtr
PipelineManager::ElementsAndOpsByDim< 2 >::FaceSideEle SideEle
PipelineManager::ElementsAndOpsByDim< 3 >::FaceSideEle SideEle
Set of functions called by PETSc solver used to refine and update mesh.
static MoFEMErrorCode tsPostStep(TS ts)
virtual ~TSPrePostProc()=default
static MoFEMErrorCode tsPreStep(TS ts)
TSPrePostProc()=default
static MoFEMErrorCode tsPostStage(TS ts, PetscReal stagetime, PetscInt stageindex, Vec *Y)
[Boundary condition]
MoFEMErrorCode tsSetUp(TS ts)
Used to setup TS solver.
static boost::weak_ptr< TSPrePostProc > tsPrePostProc
double rho
Definition plastic.cpp:144
#define EXECUTABLE_DIMENSION
Definition plastic.cpp:13
ElementsAndOps< SPACE_DIM >::SideEle SideEle
Definition plastic.cpp:61
constexpr int SPACE_DIM
DomainNaturalBC::OpFlux< NaturalMeshsetType< BLOCKSET >, 1, SPACE_DIM > OpBodyForce