12 using namespace MoFEM;
15 constexpr
double third = boost::math::constants::third<double>();
16 return (t_stress(0, 0) + t_stress(1, 1));
20 constexpr
double third = boost::math::constants::third<double>();
21 return (t_stress(0, 0) + t_stress(1, 1) + t_stress(2, 2));
66 GAUSS>::OpBaseTimesVector<1, SPACE_DIM, 0>;
103 template <
int DIM_0,
int DIM_1>
106 boost::shared_ptr<MatrixDouble> def_grad_stab_ptr,
107 boost::shared_ptr<MatrixDouble> def_grad_dot_ptr,
108 double tau_F_ptr,
double xi_F_ptr,
109 boost::shared_ptr<MatrixDouble> grad_x_ptr,
110 boost::shared_ptr<MatrixDouble> grad_vel_ptr)
112 defGradPtr(def_grad_ptr), defGradStabPtr(def_grad_stab_ptr),
113 defGradDotPtr(def_grad_dot_ptr), tauFPtr(tau_F_ptr), xiF(xi_F_ptr),
114 gradxPtr(grad_x_ptr), gradVelPtr(grad_vel_ptr) {}
124 const size_t nb_gauss_pts = getGaussPts().size2();
126 defGradStabPtr->resize(DIM_0 * DIM_1, nb_gauss_pts,
false);
127 defGradStabPtr->clear();
130 auto t_F = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*defGradPtr);
131 auto t_Fstab = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*defGradStabPtr);
132 auto t_F_dot = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*defGradDotPtr);
135 auto tau_F = tauFPtr;
137 auto t_gradx = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*gradxPtr);
138 auto t_gradVel = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*gradVelPtr);
140 for (
auto gg = 0; gg != nb_gauss_pts; ++gg) {
142 t_Fstab(
i,
j) = t_F(
i,
j) + tau_F * (t_gradVel(
i,
j) - t_F_dot(
i,
j)) +
143 xi_F * (t_gradx(
i,
j) - t_F(
i,
j));
167 template <
int DIM_0,
int DIM_1>
171 boost::shared_ptr<MatrixDouble> first_piola_ptr,
172 boost::shared_ptr<MatrixDouble> def_grad_ptr)
174 shearModulus(shear_modulus), bulkModulus(bulk_modulus), mU(m_u),
175 lammeLambda(lambda_lamme), firstPiolaPtr(first_piola_ptr),
176 defGradPtr(def_grad_ptr) {}
190 const size_t nb_gauss_pts = getGaussPts().size2();
193 firstPiolaPtr->resize(DIM_0 * DIM_1, nb_gauss_pts,
false);
194 firstPiolaPtr->clear();
197 auto t_P = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*firstPiolaPtr);
198 auto t_F = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*defGradPtr);
199 const double two_o_three = 2. / 3.;
200 const double trace_t_dk = DIM_0;
201 for (
auto gg = 0; gg != nb_gauss_pts; ++gg) {
203 t_P(
i,
j) = shearModulus * (t_F(
i,
j) + t_F(
j,
i) - 2. *
t_kd(
i,
j) -
205 two_o_three * trace_t_dk *
t_kd(
i,
j)) +
207 bulkModulus * trace_t_dk *
t_kd(
i,
j);
228 boost::shared_ptr<MatrixDouble> reference_pos_ptr,
229 boost::shared_ptr<MatrixDouble> u_ptr)
231 xPtr(spatial_pos_ptr), XPtr(reference_pos_ptr), uPtr(u_ptr) {}
240 const size_t nb_gauss_pts = getGaussPts().size2();
242 uPtr->resize(DIM, nb_gauss_pts,
false);
246 auto t_x = getFTensor1FromMat<DIM>(*xPtr);
247 auto t_X = getFTensor1FromMat<DIM>(*XPtr);
248 auto t_u = getFTensor1FromMat<DIM>(*uPtr);
249 for (
auto gg = 0; gg != nb_gauss_pts; ++gg) {
251 t_u(
i) = t_x(
i) - t_X(
i);
261 boost::shared_ptr<MatrixDouble>
xPtr;
262 boost::shared_ptr<MatrixDouble>
XPtr;
263 boost::shared_ptr<MatrixDouble>
uPtr;
266 template <
int DIM_0,
int DIM_1>
270 double shear_modulus,
double bulk_modulus,
double m_u,
271 double lambda_lamme, boost::shared_ptr<MatrixDouble> first_piola_ptr,
272 boost::shared_ptr<MatrixDouble> def_grad_ptr,
273 boost::shared_ptr<MatrixDouble> inv_def_grad_ptr,
274 boost::shared_ptr<VectorDouble> det)
276 shearModulus(shear_modulus), bulkModulus(bulk_modulus), mU(m_u),
277 lammeLambda(lambda_lamme), firstPiolaPtr(first_piola_ptr),
278 defGradPtr(def_grad_ptr), invDefGradPtr(inv_def_grad_ptr), dEt(det) {}
293 const size_t nb_gauss_pts = getGaussPts().size2();
296 firstPiolaPtr->resize(DIM_0 * DIM_1, nb_gauss_pts,
false);
297 firstPiolaPtr->clear();
300 auto t_P = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*firstPiolaPtr);
301 auto t_F = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*defGradPtr);
302 auto t_inv_F = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*invDefGradPtr);
303 auto t_det = getFTensor0FromVec<1>(*dEt);
304 const double two_o_three = 2. / 3.;
305 const double one_o_three = 1. / 3.;
306 const double bulk_mod = bulkModulus;
307 const double shear_mod = shearModulus;
308 for (
auto gg = 0; gg != nb_gauss_pts; ++gg) {
312 t_P(
i,
j) = bulk_mod * (t_det - 1.) * t_det * t_inv_F(
j,
i);
315 shear_mod * pow(t_det, two_o_three) *
316 (t_F(
i,
j) - one_o_three * (t_F(
l,
k) * t_F(
l,
k)) * t_inv_F(
j,
i));
335 boost::shared_ptr<VectorDouble>
dEt;
338 template <
int DIM_0,
int DIM_1>
342 boost::shared_ptr<MatrixDouble> def_grad_ptr,
343 boost::shared_ptr<MatrixDouble> grad_tensor_ptr)
345 defGradPtr(def_grad_ptr), gradTensorPtr(grad_tensor_ptr) {}
358 const size_t nb_gauss_pts = getGaussPts().size2();
361 defGradPtr->resize(DIM_0 * DIM_1, nb_gauss_pts,
false);
365 auto t_F = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*defGradPtr);
366 auto t_H = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*gradTensorPtr);
367 for (
auto gg = 0; gg != nb_gauss_pts; ++gg) {
400 PetscInt stageindex,
Vec *Y);
439 double getScale(
const double time) {
return 0.001 * sin(0.1 * time); };
444 double getScale(
const double time) {
return 0.001; };
453 CHKERR boundaryCondition();
477 enum bases { AINSWORTH, DEMKOWICZ, LASBASETOPT };
478 const char *list_bases[LASBASETOPT] = {
"ainsworth",
"demkowicz"};
479 PetscInt choice_base_value = AINSWORTH;
481 LASBASETOPT, &choice_base_value, PETSC_NULL);
484 switch (choice_base_value) {
488 <<
"Set AINSWORTH_LEGENDRE_BASE for displacements";
493 <<
"Set DEMKOWICZ_JACOBI_BASE for displacements";
520 auto project_ho_geometry = [&]() {
522 CHKERR mField.loop_dofs(
"x_1", ent_method_x);
524 CHKERR mField.loop_dofs(
"x_2", ent_method_x_2);
527 return mField.loop_dofs(
"GEOMETRY", ent_method);
529 CHKERR project_ho_geometry();
542 auto time_scale = boost::make_shared<TimeScale>();
544 PetscBool sin_time_function = PETSC_FALSE;
546 &sin_time_function, PETSC_NULL);
548 if (sin_time_function)
549 time_scale = boost::make_shared<DynamicFirstOrderConsSinusTimeScale>();
551 time_scale = boost::make_shared<DynamicFirstOrderConsConstantTimeScale>();
553 pipeline_mng->getBoundaryExplicitRhsFE().reset();
555 pipeline_mng->getOpBoundaryExplicitRhsPipeline(), {NOSPACE},
"GEOMETRY");
558 pipeline_mng->getOpBoundaryExplicitRhsPipeline(), mField,
"V",
559 {time_scale},
"FORCE", Sev::inform);
569 simple->getProblemName(),
"V");
571 auto get_pre_proc_hook = [&]() {
573 mField, pipeline_mng->getDomainExplicitRhsFE(), {time_scale});
575 pipeline_mng->getDomainExplicitRhsFE()->preProcessHook = get_pre_proc_hook();
582 PetscInt stageindex,
Vec *Y) {
586 auto &m_field = ptr->fsRawPtr->mField;
594 CHKERR TSGetTime(ts, &time);
596 Vec *stage_solutions;
598 CHKERR TSGetStages(ts, &num_stages, &stage_solutions);
599 PetscPrintf(PETSC_COMM_WORLD,
"Check timestep %d time %e dt %e\n",
600 num_stages, time,
dt);
602 const double inv_num_step = (
double)num_stages;
603 CHKERR fb->fieldCopy(1.,
"x_1",
"x_2");
604 CHKERR fb->fieldAxpy(
dt,
"V",
"x_2");
605 CHKERR fb->fieldCopy(1.,
"x_2",
"x_1");
607 CHKERR fb->fieldCopy(-inv_num_step /
dt,
"F_0",
"F_dot");
608 CHKERR fb->fieldAxpy(inv_num_step /
dt,
"F",
"F_dot");
609 CHKERR fb->fieldCopy(1.,
"F",
"F_0");
618 auto &m_field = ptr->fsRawPtr->mField;
619 auto fb = m_field.getInterface<
FieldBlas>();
623 CHKERR TSGetTime(ts, &time);
632 auto &m_field = ptr->fsRawPtr->mField;
639 CHKERR TSGetTime(ts, &time);
641 CHKERR TSGetStepNumber(ts, &step_num);
662 auto get_time_scale = [
this](
const double time) {
663 return sin(time *
omega * M_PI);
666 auto apply_rhs = [&](
auto &pip) {
673 auto mat_v_grad_ptr = boost::make_shared<MatrixDouble>();
675 "V", mat_v_grad_ptr));
677 auto gravity_vector_ptr = boost::make_shared<MatrixDouble>();
678 gravity_vector_ptr->resize(
SPACE_DIM, 1);
679 auto set_body_force = [&]() {
682 auto t_force = getFTensor1FromMat<SPACE_DIM, 0>(*gravity_vector_ptr);
683 double unit_weight = 0.;
688 t_force(1) = -unit_weight;
690 t_force(2) = unit_weight;
696 pip.push_back(
new OpBodyForce(
"V", gravity_vector_ptr,
697 [](
double,
double,
double) {
return 1.; }));
700 auto mat_H_tensor_ptr = boost::make_shared<MatrixDouble>();
702 "F", mat_H_tensor_ptr));
719 auto mat_dot_F_tensor_ptr = boost::make_shared<MatrixDouble>();
721 "F_dot", mat_dot_F_tensor_ptr));
724 auto mat_x_grad_ptr = boost::make_shared<MatrixDouble>();
726 "x_2", mat_x_grad_ptr));
728 auto mat_F_tensor_ptr = boost::make_shared<MatrixDouble>();
730 mat_F_tensor_ptr, mat_H_tensor_ptr));
732 auto mat_F_stab_ptr = boost::make_shared<MatrixDouble>();
734 mat_F_tensor_ptr, mat_F_stab_ptr, mat_dot_F_tensor_ptr, tau, xi,
735 mat_x_grad_ptr, mat_v_grad_ptr));
737 PetscBool is_linear_elasticity = PETSC_TRUE;
739 &is_linear_elasticity, PETSC_NULL);
741 auto mat_P_stab_ptr = boost::make_shared<MatrixDouble>();
742 if (is_linear_elasticity) {
747 auto inv_F = boost::make_shared<MatrixDouble>();
748 auto det_ptr = boost::make_shared<VectorDouble>();
756 mat_F_stab_ptr, inv_F, det_ptr));
765 CHKERR apply_rhs(pipeline_mng->getOpDomainExplicitRhsPipeline());
786 boost::shared_ptr<PostProcEle> post_proc,
787 boost::shared_ptr<PostProcFaceEle> post_proc_bdry,
788 boost::shared_ptr<MatrixDouble> velocity_field_ptr,
789 boost::shared_ptr<MatrixDouble> x2_field_ptr,
790 boost::shared_ptr<MatrixDouble> geometry_field_ptr,
791 std::array<double, 3> pass_field_eval_coords,
792 boost::shared_ptr<SetPtsData> pass_field_eval_data)
793 : dM(dm), mField(m_field), postProc(post_proc),
794 postProcBdy(post_proc_bdry), velocityFieldPtr(velocity_field_ptr),
795 x2FieldPtr(x2_field_ptr), geometryFieldPtr(geometry_field_ptr),
796 fieldEvalCoords(pass_field_eval_coords),
797 fieldEvalData(pass_field_eval_data){};
806 fieldEvalCoords.data(), 1e-12,
simple->getProblemName(),
811 fieldEvalCoords.data(), 1e-12,
simple->getProblemName(),
816 if (velocityFieldPtr->size1()) {
817 auto t_vel = getFTensor1FromMat<SPACE_DIM>(*velocityFieldPtr);
818 auto t_x2_field = getFTensor1FromMat<SPACE_DIM>(*x2FieldPtr);
819 auto t_geom = getFTensor1FromMat<SPACE_DIM>(*geometryFieldPtr);
821 double u_x = t_x2_field(0) - t_geom(0);
822 double u_y = t_x2_field(1) - t_geom(1);
823 double u_z = t_x2_field(2) - t_geom(2);
826 <<
"Velocities x: " << t_vel(0) <<
" y: " << t_vel(1)
827 <<
" z: " << t_vel(2) <<
"\n";
828 MOFEM_LOG(
"SYNC", Sev::inform) <<
"Displacement x: " << u_x
829 <<
" y: " << u_y <<
" z: " << u_z <<
"\n";
833 std::regex((boost::format(
"%s(.*)") %
"Data_Vertex").str()))) {
835 mField.
get_moab().get_entities_by_dimension(
m->getMeshset(), 0, ents,
837 auto print_vets = [](boost::shared_ptr<FieldEntity> ent_ptr) {
839 if (!(ent_ptr->getPStatus() & PSTATUS_NOT_OWNED)) {
841 <<
"Velocities: " << ent_ptr->getEntFieldData()[0] <<
" "
842 << ent_ptr->getEntFieldData()[1] <<
" "
843 << ent_ptr->getEntFieldData()[2] <<
"\n";
848 print_vets,
"V", &ents);
852 PetscBool print_volume = PETSC_FALSE;
856 PetscBool print_skin = PETSC_FALSE;
867 CHKERR postProc->writeFile(
868 "out_step_" + boost::lexical_cast<std::string>(ts_step) +
".h5m");
873 CHKERR postProcBdy->writeFile(
874 "out_boundary_" + boost::lexical_cast<std::string>(ts_step) +
898 auto dm =
simple->getDM();
900 auto calculate_stress_ops = [&](
auto &pip) {
903 auto v_ptr = boost::make_shared<MatrixDouble>();
905 auto X_ptr = boost::make_shared<MatrixDouble>();
909 auto x_ptr = boost::make_shared<MatrixDouble>();
913 auto mat_H_tensor_ptr = boost::make_shared<MatrixDouble>();
915 "F", mat_H_tensor_ptr));
917 auto u_ptr = boost::make_shared<MatrixDouble>();
921 auto mat_F_ptr = boost::make_shared<MatrixDouble>();
923 mat_F_ptr, mat_H_tensor_ptr));
925 PetscBool is_linear_elasticity = PETSC_TRUE;
927 &is_linear_elasticity, PETSC_NULL);
929 auto mat_P_ptr = boost::make_shared<MatrixDouble>();
930 if (is_linear_elasticity) {
935 auto inv_F = boost::make_shared<MatrixDouble>();
936 auto det_ptr = boost::make_shared<VectorDouble>();
942 mat_F_ptr, inv_F, det_ptr));
945 auto mat_v_grad_ptr = boost::make_shared<MatrixDouble>();
947 "V", mat_v_grad_ptr));
949 return boost::make_tuple(v_ptr, X_ptr, x_ptr, mat_P_ptr, mat_F_ptr, u_ptr);
952 auto post_proc_boundary = [&]() {
953 auto boundary_post_proc_fe = boost::make_shared<PostProcFaceEle>(mField);
956 boundary_post_proc_fe->getOpPtrVector(), {},
"GEOMETRY");
960 auto [boundary_v_ptr, boundary_X_ptr, boundary_x_ptr, boundary_mat_P_ptr,
961 boundary_mat_F_ptr, boundary_u_ptr] =
962 calculate_stress_ops(op_loop_side->getOpPtrVector());
963 boundary_post_proc_fe->getOpPtrVector().push_back(op_loop_side);
967 boundary_post_proc_fe->getOpPtrVector().push_back(
971 boundary_post_proc_fe->getPostProcMesh(),
972 boundary_post_proc_fe->getMapGaussPts(),
977 {
"GEOMETRY", boundary_X_ptr},
978 {
"x", boundary_x_ptr},
979 {
"U", boundary_u_ptr}},
982 {
"F", boundary_mat_F_ptr}},
989 return boundary_post_proc_fe;
1003 auto ts_pre_post_proc = boost::make_shared<TSPrePostProc>();
1009 boost::shared_ptr<DomainEle> vol_mass_ele(
new DomainEle(mField));
1011 vol_mass_ele->B =
M;
1019 vol_mass_ele->getOpPtrVector().push_back(
new OpMassV(
"V",
"V", get_rho));
1020 vol_mass_ele->getOpPtrVector().push_back(
new OpMassF(
"F",
"F"));
1023 CHKERR MatAssemblyBegin(
M, MAT_FINAL_ASSEMBLY);
1024 CHKERR MatAssemblyEnd(
M, MAT_FINAL_ASSEMBLY);
1027 CHKERR MatGetRowSum(
M, lumpVec);
1030 CHKERR MatDiagonalSet(
M, lumpVec, INSERT_VALUES);
1036 CHKERR KSPSetFromOptions(ksp);
1039 auto solve_boundary_for_g = [&]() {
1041 if (*(pipeline_mng->getBoundaryExplicitRhsFE()->vecAssembleSwitch)) {
1043 CHKERR VecGhostUpdateBegin(pipeline_mng->getBoundaryExplicitRhsFE()->ts_F,
1044 ADD_VALUES, SCATTER_REVERSE);
1045 CHKERR VecGhostUpdateEnd(pipeline_mng->getBoundaryExplicitRhsFE()->ts_F,
1046 ADD_VALUES, SCATTER_REVERSE);
1047 CHKERR VecAssemblyBegin(pipeline_mng->getBoundaryExplicitRhsFE()->ts_F);
1048 CHKERR VecAssemblyEnd(pipeline_mng->getBoundaryExplicitRhsFE()->ts_F);
1049 *(pipeline_mng->getBoundaryExplicitRhsFE()->vecAssembleSwitch) =
false;
1053 CHKERR KSPSolve(ksp, pipeline_mng->getBoundaryExplicitRhsFE()->ts_F,
D);
1054 CHKERR VecGhostUpdateBegin(
D, INSERT_VALUES, SCATTER_FORWARD);
1055 CHKERR VecGhostUpdateEnd(
D, INSERT_VALUES, SCATTER_FORWARD);
1056 CHKERR VecCopy(
D, pipeline_mng->getBoundaryExplicitRhsFE()->ts_F);
1062 pipeline_mng->getBoundaryExplicitRhsFE()->postProcessHook =
1063 solve_boundary_for_g;
1066 ts = pipeline_mng->createTSEX(dm);
1069 PetscBool field_eval_flag = PETSC_TRUE;
1070 boost::shared_ptr<MatrixDouble> velocity_field_ptr;
1071 boost::shared_ptr<MatrixDouble> geometry_field_ptr;
1072 boost::shared_ptr<MatrixDouble> spatial_position_field_ptr;
1073 boost::shared_ptr<SetPtsData> field_eval_data;
1075 std::array<double, 3> field_eval_coords = {0.5, 0.5, 5.};
1078 field_eval_coords.data(), &dim,
1081 if (field_eval_flag) {
1086 field_eval_data,
simple->getDomainFEName());
1089 field_eval_data,
simple->getDomainFEName());
1092 field_eval_data->setEvalPoints(field_eval_coords.data(), 1);
1094 auto no_rule = [](
int,
int,
int) {
return -1; };
1096 auto fe_ptr = field_eval_data->feMethodPtr.lock();
1097 fe_ptr->getRuleHook = no_rule;
1098 velocity_field_ptr = boost::make_shared<MatrixDouble>();
1099 geometry_field_ptr = boost::make_shared<MatrixDouble>();
1100 spatial_position_field_ptr = boost::make_shared<MatrixDouble>();
1101 fe_ptr->getOpPtrVector().push_back(
1103 fe_ptr->getOpPtrVector().push_back(
1105 geometry_field_ptr));
1106 fe_ptr->getOpPtrVector().push_back(
1108 "x_2", spatial_position_field_ptr));
1111 auto post_proc_domain = [&]() {
1112 auto post_proc_fe_vol = boost::make_shared<PostProcEle>(mField);
1116 auto [boundary_v_ptr, boundary_X_ptr, boundary_x_ptr, boundary_mat_P_ptr,
1117 boundary_mat_F_ptr, boundary_u_ptr] =
1118 calculate_stress_ops(post_proc_fe_vol->getOpPtrVector());
1120 post_proc_fe_vol->getOpPtrVector().push_back(
1124 post_proc_fe_vol->getPostProcMesh(),
1125 post_proc_fe_vol->getMapGaussPts(),
1129 {{
"V", boundary_v_ptr},
1130 {
"GEOMETRY", boundary_X_ptr},
1131 {
"x", boundary_x_ptr},
1132 {
"U", boundary_u_ptr}},
1134 {{
"FIRST_PIOLA", boundary_mat_P_ptr}, {
"F", boundary_mat_F_ptr}},
1141 return post_proc_fe_vol;
1144 boost::shared_ptr<FEMethod> null_fe;
1145 auto monitor_ptr = boost::make_shared<Monitor>(
1147 post_proc_boundary(), velocity_field_ptr, spatial_position_field_ptr,
1148 geometry_field_ptr, field_eval_coords, field_eval_data);
1151 null_fe, monitor_ptr);
1155 CHKERR TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);
1160 CHKERR TSSetSolution(ts, T);
1161 CHKERR TSSetFromOptions(ts);
1163 auto fb = mField.getInterface<
FieldBlas>();
1169 boost::shared_ptr<ForcesAndSourcesCore>
null;
1172 ptr->fsRawPtr =
this;
1174 CHKERR TSSolve(ts, NULL);
1175 CHKERR TSGetTime(ts, &ftime);
1185 PetscBool test_flg = PETSC_FALSE;
1193 CHKERR VecNorm(T, NORM_2, &nrm2);
1194 MOFEM_LOG(
"EXAMPLE", Sev::inform) <<
"Regression norm " << nrm2;
1195 constexpr
double regression_value = 0.0194561;
1196 if (fabs(nrm2 - regression_value) > 1e-2)
1198 "Regression test failed; wrong norm value.");
1216 const char param_file[] =
"param_file.petsc";
1220 auto core_log = logging::core::get();
1222 LogManager::createSink(LogManager::getStrmWorld(),
"EXAMPLE"));
1223 LogManager::setLog(
"EXAMPLE");
1229 DMType dm_name =
"DMMOFEM";