v0.15.0
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 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>
103struct OpCalculateFStab : public ForcesAndSourcesCore::UserDataOperator {
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
132
133 // tau_F = alpha deltaT
134 auto tau_F = tauFPtr;
135 double xi_F = xiF;
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>
167struct OpCalculatePiola : public ForcesAndSourcesCore::UserDataOperator {
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
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>
225struct OpCalculateDisplacement : public ForcesAndSourcesCore::UserDataOperator {
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>
267 : public ForcesAndSourcesCore::UserDataOperator {
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
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>
339 : public ForcesAndSourcesCore::UserDataOperator {
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
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", 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 auto &m_field = ptr->fsRawPtr->mField;
616 double dt;
617 CHKERR TSGetTimeStep(ts, &dt);
618 double time;
619 CHKERR TSGetTime(ts, &time);
620 }
622}
623
626
627 if (auto ptr = tsPrePostProc.lock()) {
628 auto &m_field = ptr->fsRawPtr->mField;
629 double dt;
630 CHKERR TSGetTimeStep(ts, &dt);
631 double time;
632 CHKERR TSGetTime(ts, &time);
633 int step_num;
634 CHKERR TSGetStepNumber(ts, &step_num);
635 }
637}
638
639//! [Push operators to pipeline]
642 auto get_body_force = [this](const double, const double, const double) {
645 t_source(i) = 0.;
646 t_source(0) = 0.1;
647 t_source(1) = 1.;
648 return t_source;
649 };
650
651 // specific time scaling
652 auto get_time_scale = [this](const double time) {
653 return sin(time * omega * M_PI);
654 };
655
656 auto apply_rhs = [&](auto &pip) {
658
660 "GEOMETRY");
661
662 // Calculate Gradient of velocity
663 auto mat_v_grad_ptr = boost::make_shared<MatrixDouble>();
665 "V", mat_v_grad_ptr));
666
667 auto gravity_vector_ptr = boost::make_shared<MatrixDouble>();
668 gravity_vector_ptr->resize(SPACE_DIM, 1);
669 auto set_body_force = [&]() {
672 auto t_force = getFTensor1FromMat<SPACE_DIM, 0>(*gravity_vector_ptr);
673 double unit_weight = 0.;
674 CHKERR PetscOptionsGetReal(PETSC_NULLPTR, "", "-unit_weight", &unit_weight,
675 PETSC_NULLPTR);
676 t_force(i) = 0;
677 if (SPACE_DIM == 2) {
678 t_force(1) = -unit_weight;
679 } else if (SPACE_DIM == 3) {
680 t_force(2) = unit_weight;
681 }
683 };
684
685 CHKERR set_body_force();
686 pip.push_back(new OpBodyForce("V", gravity_vector_ptr,
687 [](double, double, double) { return 1.; }));
688
689 // Calculate unknown F
690 auto mat_H_tensor_ptr = boost::make_shared<MatrixDouble>();
692 "F", mat_H_tensor_ptr));
693
694 // // Calculate F
695 double tau = 0.2;
696 CHKERR PetscOptionsGetReal(PETSC_NULLPTR, "", "-tau", &tau, PETSC_NULLPTR);
697
698 double xi = 0.;
699 CHKERR PetscOptionsGetReal(PETSC_NULLPTR, "", "-xi", &xi, PETSC_NULLPTR);
700
701 // Calculate P stab
702 auto one = [&](const double, const double, const double) {
703 return 3. * bulk_modulus_K;
704 };
705 auto minus_one = [](const double, const double, const double) {
706 return -1.;
707 };
708
709 auto mat_dot_F_tensor_ptr = boost::make_shared<MatrixDouble>();
711 "F_dot", mat_dot_F_tensor_ptr));
712
713 // Calculate Gradient of Spatial Positions
714 auto mat_x_grad_ptr = boost::make_shared<MatrixDouble>();
716 "x_2", mat_x_grad_ptr));
717
718 auto mat_F_tensor_ptr = boost::make_shared<MatrixDouble>();
720 mat_F_tensor_ptr, mat_H_tensor_ptr));
721
722 auto mat_F_stab_ptr = boost::make_shared<MatrixDouble>();
724 mat_F_tensor_ptr, mat_F_stab_ptr, mat_dot_F_tensor_ptr, tau, xi,
725 mat_x_grad_ptr, mat_v_grad_ptr));
726
727 PetscBool is_linear_elasticity = PETSC_TRUE;
728 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-is_linear_elasticity",
729 &is_linear_elasticity, PETSC_NULLPTR);
730
731 auto mat_P_stab_ptr = boost::make_shared<MatrixDouble>();
732 if (is_linear_elasticity) {
735 mat_F_stab_ptr));
736 } else {
737 auto inv_F = boost::make_shared<MatrixDouble>();
738 auto det_ptr = boost::make_shared<VectorDouble>();
739
740 pip.push_back(
741 new OpInvertMatrix<SPACE_DIM>(mat_F_stab_ptr, det_ptr, inv_F));
742
743 // OpCalculatePiolaIncompressibleNH
746 mat_F_stab_ptr, inv_F, det_ptr));
747 }
748
749 pip.push_back(new OpGradTimesTensor2("V", mat_P_stab_ptr, minus_one));
750 pip.push_back(new OpRhsTestPiola("F", mat_v_grad_ptr, one));
751
753 };
754
755 auto *pipeline_mng = mField.getInterface<PipelineManager>();
756 CHKERR apply_rhs(pipeline_mng->getOpDomainExplicitRhsPipeline());
757
758 auto integration_rule = [](int, int, int approx_order) {
759 return 2 * approx_order;
760 };
761 CHKERR pipeline_mng->setDomainExplicitRhsIntegrationRule(integration_rule);
762
764}
765//! [Push operators to pipeline]
766
767/**
768 * @brief Monitor solution
769 *
770 * This functions is called by TS solver at the end of each step. It is used
771 * to output results to the hard drive.
772 */
773
774struct Monitor : public FEMethod {
777 boost::shared_ptr<PostProcEle> post_proc,
778 boost::shared_ptr<PostProcFaceEle> post_proc_bdry,
779 boost::shared_ptr<MatrixDouble> velocity_field_ptr,
780 boost::shared_ptr<MatrixDouble> x2_field_ptr,
781 boost::shared_ptr<MatrixDouble> geometry_field_ptr,
782 std::array<double, 3> pass_field_eval_coords,
783 boost::shared_ptr<SetPtsData> pass_field_eval_data)
784 : dM(dm), mField(m_field), postProc(post_proc),
785 postProcBdy(post_proc_bdry), velocityFieldPtr(velocity_field_ptr),
786 x2FieldPtr(x2_field_ptr), geometryFieldPtr(geometry_field_ptr),
787 fieldEvalCoords(pass_field_eval_coords),
788 fieldEvalData(pass_field_eval_data){};
791
792 auto *simple = mField.getInterface<Simple>();
793
795 ->evalFEAtThePoint<SPACE_DIM>(
796 fieldEvalCoords.data(), 1e-12, simple->getProblemName(),
797 simple->getDomainFEName(), fieldEvalData, mField.get_comm_rank(),
799
800 if (velocityFieldPtr->size1()) {
802 auto t_x2_field = getFTensor1FromMat<SPACE_DIM>(*x2FieldPtr);
804
805 double u_x = t_x2_field(0) - t_geom(0);
806 double u_y = t_x2_field(1) - t_geom(1);
807 double u_z = t_x2_field(2) - t_geom(2);
808
809 MOFEM_LOG("SYNC", Sev::inform)
810 << "Velocities x: " << t_vel(0) << " y: " << t_vel(1)
811 << " z: " << t_vel(2) << "\n";
812 MOFEM_LOG("SYNC", Sev::inform) << "Displacement x: " << u_x
813 << " y: " << u_y << " z: " << u_z << "\n";
814 }
815
817 std::regex((boost::format("%s(.*)") % "Data_Vertex").str()))) {
818 Range ents;
819 mField.get_moab().get_entities_by_dimension(m->getMeshset(), 0, ents,
820 true);
821 auto print_vets = [](boost::shared_ptr<FieldEntity> ent_ptr) {
823 if (!(ent_ptr->getPStatus() & PSTATUS_NOT_OWNED)) {
824 MOFEM_LOG("SYNC", Sev::inform)
825 << "Velocities: " << ent_ptr->getEntFieldData()[0] << " "
826 << ent_ptr->getEntFieldData()[1] << " "
827 << ent_ptr->getEntFieldData()[2] << "\n";
828 }
830 };
831 CHKERR mField.getInterface<FieldBlas>()->fieldLambdaOnEntities(
832 print_vets, "V", &ents);
833 }
835
836 PetscBool print_volume = PETSC_FALSE;
837 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-print_volume", &print_volume,
838 PETSC_NULLPTR);
839
840 PetscBool print_skin = PETSC_FALSE;
841 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-print_skin", &print_skin,
842 PETSC_NULLPTR);
843
844 int save_every_nth_step = 1;
845 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-save_step",
846 &save_every_nth_step, PETSC_NULLPTR);
847 if (ts_step % save_every_nth_step == 0) {
848
849 if (print_volume) {
851 CHKERR postProc->writeFile(
852 "out_step_" + boost::lexical_cast<std::string>(ts_step) + ".h5m");
853 }
854
855 if (print_skin) {
857 CHKERR postProcBdy->writeFile(
858 "out_boundary_" + boost::lexical_cast<std::string>(ts_step) +
859 ".h5m");
860 }
861 }
863 }
864
865private:
867 boost::shared_ptr<PostProcEle> postProc;
868 boost::shared_ptr<PostProcFaceEle> postProcBdy;
869 boost::shared_ptr<MatrixDouble> velocityFieldPtr;
870 boost::shared_ptr<MatrixDouble> x2FieldPtr;
871 boost::shared_ptr<MatrixDouble> geometryFieldPtr;
872 std::array<double, 3> fieldEvalCoords;
873 boost::shared_ptr<SetPtsData> fieldEvalData;
874};
875
876//! [Solve]
879 auto *simple = mField.getInterface<Simple>();
880 auto *pipeline_mng = mField.getInterface<PipelineManager>();
881
882 auto dm = simple->getDM();
883
884 auto calculate_stress_ops = [&](auto &pip) {
886
887 auto v_ptr = boost::make_shared<MatrixDouble>();
888 pip.push_back(new OpCalculateVectorFieldValues<SPACE_DIM>("V", v_ptr));
889 auto X_ptr = boost::make_shared<MatrixDouble>();
890 pip.push_back(
891 new OpCalculateVectorFieldValues<SPACE_DIM>("GEOMETRY", X_ptr));
892
893 auto x_ptr = boost::make_shared<MatrixDouble>();
894 pip.push_back(new OpCalculateVectorFieldValues<SPACE_DIM>("x_1", x_ptr));
895
896 // Calculate unknown F
897 auto mat_H_tensor_ptr = boost::make_shared<MatrixDouble>();
899 "F", mat_H_tensor_ptr));
900
901 auto u_ptr = boost::make_shared<MatrixDouble>();
902 pip.push_back(new OpCalculateDisplacement<SPACE_DIM>(x_ptr, X_ptr, u_ptr));
903 // Calculate P
904
905 auto mat_F_ptr = boost::make_shared<MatrixDouble>();
907 mat_F_ptr, mat_H_tensor_ptr));
908
909 PetscBool is_linear_elasticity = PETSC_TRUE;
910 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-is_linear_elasticity",
911 &is_linear_elasticity, PETSC_NULLPTR);
912
913 auto mat_P_ptr = boost::make_shared<MatrixDouble>();
914 if (is_linear_elasticity) {
917 mat_F_ptr));
918 } else {
919 auto inv_F = boost::make_shared<MatrixDouble>();
920 auto det_ptr = boost::make_shared<VectorDouble>();
921
922 pip.push_back(new OpInvertMatrix<SPACE_DIM>(mat_F_ptr, det_ptr, inv_F));
923
926 mat_F_ptr, inv_F, det_ptr));
927 }
928
929 auto mat_v_grad_ptr = boost::make_shared<MatrixDouble>();
931 "V", mat_v_grad_ptr));
932
933 return boost::make_tuple(v_ptr, X_ptr, x_ptr, mat_P_ptr, mat_F_ptr, u_ptr);
934 };
935
936 auto post_proc_boundary = [&]() {
937 auto boundary_post_proc_fe = boost::make_shared<PostProcFaceEle>(mField);
938
940 boundary_post_proc_fe->getOpPtrVector(), {}, "GEOMETRY");
941 auto op_loop_side =
942 new OpLoopSide<SideEle>(mField, simple->getDomainFEName(), SPACE_DIM);
943 // push ops to side element, through op_loop_side operator
944 auto [boundary_v_ptr, boundary_X_ptr, boundary_x_ptr, boundary_mat_P_ptr,
945 boundary_mat_F_ptr, boundary_u_ptr] =
946 calculate_stress_ops(op_loop_side->getOpPtrVector());
947 boundary_post_proc_fe->getOpPtrVector().push_back(op_loop_side);
948
950
951 boundary_post_proc_fe->getOpPtrVector().push_back(
952
953 new OpPPMap(
954
955 boundary_post_proc_fe->getPostProcMesh(),
956 boundary_post_proc_fe->getMapGaussPts(),
957
959
960 OpPPMap::DataMapMat{{"V", boundary_v_ptr},
961 {"GEOMETRY", boundary_X_ptr},
962 {"x", boundary_x_ptr},
963 {"U", boundary_u_ptr}},
964
965 OpPPMap::DataMapMat{{"FIRST_PIOLA", boundary_mat_P_ptr},
966 {"F", boundary_mat_F_ptr}},
967
969
970 )
971
972 );
973 return boundary_post_proc_fe;
974 };
975
976 // Add monitor to time solver
977
978 double rho = 1.;
979 CHKERR PetscOptionsGetReal(PETSC_NULLPTR, "", "-density", &rho, PETSC_NULLPTR);
980 auto get_rho = [rho](const double, const double, const double) {
981 return rho;
982 };
983
984 SmartPetscObj<Mat> M; ///< Mass matrix
985 SmartPetscObj<KSP> ksp; ///< Linear solver
986
987 auto ts_pre_post_proc = boost::make_shared<TSPrePostProc>();
988 tsPrePostProc = ts_pre_post_proc;
989
991 CHKERR MatZeroEntries(M);
992
993 boost::shared_ptr<DomainEle> vol_mass_ele(new DomainEle(mField));
994
995 vol_mass_ele->B = M;
996
997 auto integration_rule = [](int, int, int approx_order) {
998 return 2 * approx_order;
999 };
1000
1001 vol_mass_ele->getRuleHook = integration_rule;
1002
1003 vol_mass_ele->getOpPtrVector().push_back(new OpMassV("V", "V", get_rho));
1004 vol_mass_ele->getOpPtrVector().push_back(new OpMassF("F", "F"));
1005
1006 CHKERR DMoFEMLoopFiniteElements(dm, simple->getDomainFEName(), vol_mass_ele);
1007 CHKERR MatAssemblyBegin(M, MAT_FINAL_ASSEMBLY);
1008 CHKERR MatAssemblyEnd(M, MAT_FINAL_ASSEMBLY);
1009
1010 auto lumpVec = createDMVector(simple->getDM());
1011 CHKERR MatGetRowSum(M, lumpVec);
1012
1013 CHKERR MatZeroEntries(M);
1014 CHKERR MatDiagonalSet(M, lumpVec, INSERT_VALUES);
1015
1016 // Create and septup KSP (linear solver), we need this to calculate g(t,u) =
1017 // M^-1G(t,u)
1018 ksp = createKSP(mField.get_comm());
1019 CHKERR KSPSetOperators(ksp, M, M);
1020 CHKERR KSPSetFromOptions(ksp);
1021 CHKERR KSPSetUp(ksp);
1022
1023 auto solve_boundary_for_g = [&]() {
1025 if (*(pipeline_mng->getBoundaryExplicitRhsFE()->vecAssembleSwitch)) {
1026
1027 CHKERR VecGhostUpdateBegin(pipeline_mng->getBoundaryExplicitRhsFE()->ts_F,
1028 ADD_VALUES, SCATTER_REVERSE);
1029 CHKERR VecGhostUpdateEnd(pipeline_mng->getBoundaryExplicitRhsFE()->ts_F,
1030 ADD_VALUES, SCATTER_REVERSE);
1031 CHKERR VecAssemblyBegin(pipeline_mng->getBoundaryExplicitRhsFE()->ts_F);
1032 CHKERR VecAssemblyEnd(pipeline_mng->getBoundaryExplicitRhsFE()->ts_F);
1033 *(pipeline_mng->getBoundaryExplicitRhsFE()->vecAssembleSwitch) = false;
1034
1035 auto D =
1036 smartVectorDuplicate(pipeline_mng->getBoundaryExplicitRhsFE()->ts_F);
1037 CHKERR KSPSolve(ksp, pipeline_mng->getBoundaryExplicitRhsFE()->ts_F, D);
1038 CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
1039 CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
1040 CHKERR VecCopy(D, pipeline_mng->getBoundaryExplicitRhsFE()->ts_F);
1041 }
1042
1044 };
1045
1046 pipeline_mng->getBoundaryExplicitRhsFE()->postProcessHook =
1047 solve_boundary_for_g;
1048
1050 ts = pipeline_mng->createTSEX(dm);
1051
1052 // Field eval
1053 PetscBool field_eval_flag = PETSC_TRUE;
1054 boost::shared_ptr<MatrixDouble> velocity_field_ptr;
1055 boost::shared_ptr<MatrixDouble> geometry_field_ptr;
1056 boost::shared_ptr<MatrixDouble> spatial_position_field_ptr;
1057 boost::shared_ptr<SetPtsData> field_eval_data;
1058
1059 std::array<double, 3> field_eval_coords = {0.5, 0.5, 5.};
1060 int dim = 3;
1061 CHKERR PetscOptionsGetRealArray(NULL, NULL, "-field_eval_coords",
1062 field_eval_coords.data(), &dim,
1063 &field_eval_flag);
1064
1065 if (field_eval_flag) {
1066 field_eval_data =
1067 mField.getInterface<FieldEvaluatorInterface>()->getData<DomainEle>();
1068 CHKERR mField.getInterface<FieldEvaluatorInterface>()->buildTree<SPACE_DIM>(
1069 field_eval_data, simple->getDomainFEName());
1070
1071 field_eval_data->setEvalPoints(field_eval_coords.data(), 1);
1072
1073 auto no_rule = [](int, int, int) { return -1; };
1074
1075 auto fe_ptr = field_eval_data->feMethodPtr.lock();
1076 fe_ptr->getRuleHook = no_rule;
1077 velocity_field_ptr = boost::make_shared<MatrixDouble>();
1078 geometry_field_ptr = boost::make_shared<MatrixDouble>();
1079 spatial_position_field_ptr = boost::make_shared<MatrixDouble>();
1080 fe_ptr->getOpPtrVector().push_back(
1081 new OpCalculateVectorFieldValues<SPACE_DIM>("V", velocity_field_ptr));
1082 fe_ptr->getOpPtrVector().push_back(
1084 geometry_field_ptr));
1085 fe_ptr->getOpPtrVector().push_back(
1087 "x_2", spatial_position_field_ptr));
1088 }
1089
1090 auto post_proc_domain = [&]() {
1091 auto post_proc_fe_vol = boost::make_shared<PostProcEle>(mField);
1092
1094
1095 auto [boundary_v_ptr, boundary_X_ptr, boundary_x_ptr, boundary_mat_P_ptr,
1096 boundary_mat_F_ptr, boundary_u_ptr] =
1097 calculate_stress_ops(post_proc_fe_vol->getOpPtrVector());
1098
1099 post_proc_fe_vol->getOpPtrVector().push_back(
1100
1101 new OpPPMap(
1102
1103 post_proc_fe_vol->getPostProcMesh(),
1104 post_proc_fe_vol->getMapGaussPts(),
1105
1106 {},
1107
1108 {{"V", boundary_v_ptr},
1109 {"GEOMETRY", boundary_X_ptr},
1110 {"x", boundary_x_ptr},
1111 {"U", boundary_u_ptr}},
1112
1113 {{"FIRST_PIOLA", boundary_mat_P_ptr}, {"F", boundary_mat_F_ptr}},
1114
1115 {}
1116
1117 )
1118
1119 );
1120 return post_proc_fe_vol;
1121 };
1122
1123 boost::shared_ptr<FEMethod> null_fe;
1124 auto monitor_ptr = boost::make_shared<Monitor>(
1125 SmartPetscObj<DM>(dm, true), mField, post_proc_domain(),
1126 post_proc_boundary(), velocity_field_ptr, spatial_position_field_ptr,
1127 geometry_field_ptr, field_eval_coords, field_eval_data);
1128
1129 CHKERR DMMoFEMTSSetMonitor(dm, ts, simple->getDomainFEName(), null_fe,
1130 null_fe, monitor_ptr);
1131
1132 double ftime = 1;
1133 // CHKERR TSSetMaxTime(ts, ftime);
1134 CHKERR TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);
1135
1136 auto T = createDMVector(simple->getDM());
1137 CHKERR DMoFEMMeshToLocalVector(simple->getDM(), T, INSERT_VALUES,
1138 SCATTER_FORWARD);
1139 CHKERR TSSetSolution(ts, T);
1140 CHKERR TSSetFromOptions(ts);
1141
1142 CHKERR TSSetPostStage(ts, TSPrePostProc::tsPostStage);
1143 CHKERR TSSetPostStep(ts, TSPrePostProc::tsPostStep);
1144 CHKERR TSSetPreStep(ts, TSPrePostProc::tsPreStep);
1145
1146 boost::shared_ptr<ForcesAndSourcesCore> null;
1147
1148 if (auto ptr = tsPrePostProc.lock()) {
1149 ptr->fsRawPtr = this;
1150 CHKERR TSSetUp(ts);
1151 CHKERR TSSolve(ts, NULL);
1152 CHKERR TSGetTime(ts, &ftime);
1153 }
1154
1156}
1157//! [Solve]
1158
1159//! [Postprocess results]
1162 PetscBool test_flg = PETSC_FALSE;
1163 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-test", &test_flg, PETSC_NULLPTR);
1164 if (test_flg) {
1165 auto *simple = mField.getInterface<Simple>();
1166 auto T = createDMVector(simple->getDM());
1167 CHKERR DMoFEMMeshToLocalVector(simple->getDM(), T, INSERT_VALUES,
1168 SCATTER_FORWARD);
1169 double nrm2;
1170 CHKERR VecNorm(T, NORM_2, &nrm2);
1171 MOFEM_LOG("EXAMPLE", Sev::inform) << "Regression norm " << nrm2;
1172 constexpr double regression_value = 0.0194561;
1173 if (fabs(nrm2 - regression_value) > 1e-2)
1174 SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
1175 "Regression test failed; wrong norm value.");
1176 }
1178}
1179//! [Postprocess results]
1180
1181//! [Check]
1186//! [Check]
1187
1188static char help[] = "...\n\n";
1189
1190int main(int argc, char *argv[]) {
1191
1192 // Initialisation of MoFEM/PETSc and MOAB data structures
1193 const char param_file[] = "param_file.petsc";
1194 MoFEM::Core::Initialize(&argc, &argv, param_file, help);
1195
1196 // Add logging channel for example
1197 auto core_log = logging::core::get();
1198 core_log->add_sink(
1200 LogManager::setLog("EXAMPLE");
1201 MOFEM_LOG_TAG("EXAMPLE", "example");
1202
1203 try {
1204
1205 //! [Register MoFEM discrete manager in PETSc]
1206 DMType dm_name = "DMMOFEM";
1207 CHKERR DMRegister_MoFEM(dm_name);
1208 //! [Register MoFEM discrete manager in PETSc
1209
1210 //! [Create MoAB]
1211 moab::Core mb_instance; ///< mesh database
1212 moab::Interface &moab = mb_instance; ///< mesh database interface
1213 //! [Create MoAB]
1214
1215 //! [Create MoFEM]
1216 MoFEM::Core core(moab); ///< finite element database
1217 MoFEM::Interface &m_field = core; ///< finite element database interface
1218 //! [Create MoFEM]
1219
1220 //! [Example]
1221 Example ex(m_field);
1222 CHKERR ex.runProblem();
1223 //! [Example]
1224 }
1226
1228}
constexpr double third
ForcesAndSourcesCore::UserDataOperator UserDataOperator
#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
constexpr int SPACE_DIM
[Define dimension]
Definition adjoint.cpp:22
constexpr double shear_modulus_G
Definition adjoint.cpp:76
constexpr double bulk_modulus_K
Definition adjoint.cpp:75
int main()
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
[Define dimension]
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 nme:nme847.
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
static boost::weak_ptr< TSPrePostProc > tsPrePostProc
static char help[]
[Check]
FormsIntegrators< DomainEleOp >::Assembly< AssemblyType::PETSC >::LinearForm< IntegrationType::GAUSS >::OpBaseTimesVector< 1, SPACE_DIM *SPACE_DIM, 1 > OpRhsTestPiola
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpMassV
constexpr int SPACE_DIM
FormsIntegrators< DomainEleOp >::Assembly< AssemblyType::PETSC >::LinearForm< IntegrationType::GAUSS >::OpGradTimesTensor< 1, SPACE_DIM, SPACE_DIM > OpGradTimesTensor2
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM *SPACE_DIM > OpMassF
constexpr double poisson_ratio
constexpr double omega
Save field DOFS on vertices/tags.
PostProcEleByDim< SPACE_DIM >::PostProcEleBdy PostProcEleBdy
PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::DomainEle DomainEle
FormsIntegrators< DomainEleOp >::Assembly< AssemblyType::PETSC >::LinearForm< IntegrationType::GAUSS >::OpGradTimesTensor< 1, SPACE_DIM, SPACE_DIM > OpGradTimesPiola
double trace(FTensor::Tensor2< T, 2, 2 > &t_stress)
double bulk_modulus_K
double shear_modulus_G
constexpr double young_modulus
static boost::weak_ptr< TSPrePostProc > tsPrePostProc
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:1102
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.
MoFEMErrorCode getCubitMeshsetPtr(const int ms_id, const CubitBCType cubit_bc_type, const CubitMeshSets **cubit_meshset_ptr) const
get cubit meshset
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)
DEPRECATED SmartPetscObj< Vec > smartVectorDuplicate(Vec vec)
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)
PetscErrorCode PetscOptionsGetRealArray(PetscOptions *, const char pre[], const char name[], PetscReal dval[], PetscInt *nmax, PetscBool *set)
FTensor::Tensor1< FTensor::PackPtr< T *, S >, Tensor_Dim > getFTensor1FromMat(ublas::matrix< T, L, A > &data)
Get tensor rank 1 (vector) form data matrix.
FTensor::Tensor2< FTensor::PackPtr< double *, 1 >, Tensor_Dim1, Tensor_Dim2 > getFTensor2FromMat(MatrixDouble &data)
Get tensor rank 2 (matrix) form data matrix.
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
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
double rho
Definition plastic.cpp:144
#define EXECUTABLE_DIMENSION
Definition plastic.cpp:13
ElementsAndOps< SPACE_DIM >::SideEle SideEle
Definition plastic.cpp:61
static constexpr int approx_order
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
Definition plot_base.cpp:68
MoFEMErrorCode checkResults()
[Postprocess results]
MoFEMErrorCode solveSystem()
[Solve]
Example(MoFEM::Interface &m_field)
MoFEMErrorCode runProblem()
MoFEM::Interface & mField
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 ptr object.
Simple interface for fast problem set-up.
Definition BcManager.hpp:29
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.
structure to get information form 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.
Get values at integration pts for tensor field rank 1, i.e. vector field.
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
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
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()
function is run at the end of loop
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.
DomainNaturalBCRhs::OpFlux< NaturalMeshsetType< BLOCKSET >, 1, SPACE_DIM > OpBodyForce