13 #ifndef EXECUTABLE_DIMENSION
14 #define EXECUTABLE_DIMENSION 3
17 #ifndef SCHUR_ASSEMBLE
18 #define SCHUR_ASSEMBLE 0
24 using namespace MoFEM;
26 #include <GenericElementInterface.hpp>
29 #include <boost/python.hpp>
30 #include <boost/python/def.hpp>
31 #include <boost/python/numpy.hpp>
32 namespace bp = boost::python;
33 namespace np = boost::python::numpy;
75 IT>::OpBaseTimesVector<1, SPACE_DIM, 1>;
109 #ifdef WITH_MODULE_MFRONT_INTERFACE
110 #include <MFrontMoFEMInterface.hpp>
149 boost::shared_ptr<SDFPython> sdfPythonPtr;
164 CHKERR createCommonData();
193 enum bases { AINSWORTH, DEMKOWICZ, LASBASETOPT };
194 const char *list_bases[LASBASETOPT] = {
"ainsworth",
"demkowicz"};
195 PetscInt choice_base_value = AINSWORTH;
197 LASBASETOPT, &choice_base_value, PETSC_NULL);
200 switch (choice_base_value) {
204 <<
"Set AINSWORTH_LEGENDRE_BASE for displacements";
209 <<
"Set DEMKOWICZ_JACOBI_BASE for displacements";
229 auto get_skin = [&]() {
231 CHKERR mField.get_moab().get_entities_by_dimension(0,
SPACE_DIM, body_ents);
232 Skinner skin(&mField.get_moab());
234 CHKERR skin.find_skin(0, body_ents,
false, skin_ents);
238 auto filter_blocks = [&](
auto skin) {
239 bool is_contact_block =
false;
244 (boost::format(
"%s(.*)") %
"CONTACT").str()
253 <<
"Find contact block set: " <<
m->getName();
254 auto meshset =
m->getMeshset();
255 Range contact_meshset_range;
256 CHKERR mField.get_moab().get_entities_by_dimension(
257 meshset,
SPACE_DIM - 1, contact_meshset_range,
true);
260 contact_meshset_range);
261 contact_range.merge(contact_meshset_range);
263 if (is_contact_block) {
265 <<
"Nb entities in contact surface: " << contact_range.size();
267 skin = intersect(skin, contact_range);
272 auto filter_true_skin = [&](
auto skin) {
274 ParallelComm *pcomm =
276 CHKERR pcomm->filter_pstatus(skin, PSTATUS_SHARED | PSTATUS_MULTISHARED,
277 PSTATUS_NOT, -1, &boundary_ents);
278 return boundary_ents;
281 auto boundary_ents = filter_true_skin(filter_blocks(get_skin()));
288 CHKERR mField.get_moab().get_adjacencies(boundary_ents, 1,
false, ho_ents,
289 moab::Interface::UNION);
298 auto project_ho_geometry = [&]() {
300 return mField.loop_dofs(
"GEOMETRY", ent_method);
303 PetscBool project_geometry = PETSC_TRUE;
305 &project_geometry, PETSC_NULL);
306 if (project_geometry) {
307 CHKERR project_ho_geometry();
317 PetscBool use_mfront = PETSC_FALSE;
339 if (!mfrontInterface) {
343 MOFEM_LOG(
"CONTACT", Sev::inform) <<
"Using MFront for material model";
354 PetscBool use_scale = PETSC_FALSE;
369 char sdf_file_name[255];
371 sdf_file_name, 255, PETSC_NULL);
373 sdfPythonPtr = boost::make_shared<SDFPython>();
374 CHKERR sdfPythonPtr->sdfInit(sdf_file_name);
375 sdfPythonWeakPtr = sdfPythonPtr;
381 "Use executable contact_2d with axisymmetric model");
385 "Axisymmetric model is only available with MFront (set "
388 MOFEM_LOG(
"CONTACT", Sev::inform) <<
"Using axisymmetric model";
393 MOFEM_LOG(
"CONTACT", Sev::inform) <<
"Using plane strain model";
398 #ifndef WITH_MODULE_MFRONT_INTERFACE
401 "MFrontInterface module was not found while use_mfront was set to 1");
405 boost::make_shared<MFrontMoFEMInterface<TRIDIMENSIONAL>>(
410 boost::make_shared<MFrontMoFEMInterface<AXISYMMETRICAL>>(
413 mfrontInterface = boost::make_shared<MFrontMoFEMInterface<PLANESTRAIN>>(
418 CHKERR mfrontInterface->getCommandLineParameters();
422 auto dm =
simple->getDM();
427 mfrontInterface->setMonitorPtr(monitorPtr);
437 auto bc_mng = mField.getInterface<
BcManager>();
440 for (
auto f : {
"U",
"SIGMA"}) {
441 CHKERR bc_mng->removeBlockDOFsOnEntities(
simple->getProblemName(),
442 "REMOVE_X",
f, 0, 0);
443 CHKERR bc_mng->removeBlockDOFsOnEntities(
simple->getProblemName(),
444 "REMOVE_Y",
f, 1, 1);
445 CHKERR bc_mng->removeBlockDOFsOnEntities(
simple->getProblemName(),
446 "REMOVE_Z",
f, 2, 2);
447 CHKERR bc_mng->removeBlockDOFsOnEntities(
simple->getProblemName(),
448 "REMOVE_ALL",
f, 0, 3);
451 CHKERR bc_mng->removeBlockDOFsOnEntities(
simple->getProblemName(),
"FIX_X",
452 "SIGMA", 0, 0,
false,
true);
453 CHKERR bc_mng->removeBlockDOFsOnEntities(
simple->getProblemName(),
"FIX_Y",
454 "SIGMA", 1, 1,
false,
true);
455 CHKERR bc_mng->removeBlockDOFsOnEntities(
simple->getProblemName(),
"FIX_Z",
456 "SIGMA", 2, 2,
false,
true);
457 CHKERR bc_mng->removeBlockDOFsOnEntities(
simple->getProblemName(),
"FIX_ALL",
458 "SIGMA", 0, 3,
false,
true);
459 CHKERR bc_mng->removeBlockDOFsOnEntities(
460 simple->getProblemName(),
"NO_CONTACT",
"SIGMA", 0, 3,
false,
true);
465 simple->getProblemName(),
"U");
477 auto time_scale = boost::make_shared<ScaledTimeScale>();
478 auto body_force_time_scale =
479 boost::make_shared<ScaledTimeScale>(
"body_force_hist.txt");
488 auto add_domain_base_ops = [&](
auto &pip) {
495 auto hencky_common_data_ptr = boost::make_shared<HenckyOps::CommonData>();
496 hencky_common_data_ptr->matDPtr = boost::make_shared<MatrixDouble>();
497 hencky_common_data_ptr->matGradPtr = boost::make_shared<MatrixDouble>();
499 auto add_domain_ops_lhs = [&](
auto &pip) {
509 auto fe_domain_lhs = pip_mng->getDomainLhsFE();
511 auto get_inertia_and_mass_damping =
513 return (
rho *
scale) * fe_domain_lhs->ts_aa +
516 pip.push_back(
new OpMass(
"U",
"U", get_inertia_and_mass_damping));
520 auto fe_domain_lhs = pip_mng->getDomainLhsFE();
522 auto get_inertia_and_mass_damping =
526 pip.push_back(
new OpMass(
"U",
"U", get_inertia_and_mass_damping));
529 if (!mfrontInterface) {
530 CHKERR HenckyOps::opFactoryDomainLhs<SPACE_DIM, AT, IT, DomainEleOp>(
531 mField, pip,
"U",
"MAT_ELASTIC", Sev::verbose,
scale);
533 CHKERR mfrontInterface->opFactoryDomainLhs(pip);
539 auto add_domain_ops_rhs = [&](
auto &pip) {
543 pip, mField,
"U", {body_force_time_scale}, Sev::inform);
547 AT>::LinearForm<IT>::OpBaseTimesVector<1,
SPACE_DIM, 1>;
552 auto mat_acceleration = boost::make_shared<MatrixDouble>();
554 "U", mat_acceleration));
556 new OpInertiaForce(
"U", mat_acceleration, [](
double,
double,
double) {
563 auto mat_velocity = boost::make_shared<MatrixDouble>();
567 new OpInertiaForce(
"U", mat_velocity, [](
double,
double,
double) {
572 if (!mfrontInterface) {
573 CHKERR HenckyOps::opFactoryDomainRhs<SPACE_DIM, AT, IT, DomainEleOp>(
574 mField, pip,
"U",
"MAT_ELASTIC", Sev::inform,
scale);
576 CHKERR mfrontInterface->opFactoryDomainRhs(pip);
579 CHKERR ContactOps::opFactoryDomainRhs<SPACE_DIM, AT, IT, DomainEleOp>(
585 auto add_boundary_base_ops = [&](
auto &pip) {
595 auto add_boundary_ops_lhs = [&](
auto &pip) {
605 pip, mField,
"U", Sev::inform);
610 auto fe_boundary_lhs = pip_mng->getBoundaryLhsFE();
615 [
this, fe_boundary_lhs](
double,
double,
double) {
624 ContactOps::opFactoryBoundaryLhs<SPACE_DIM, AT, GAUSS, BoundaryEleOp>(
628 mField, pip,
simple->getDomainFEName(),
"SIGMA",
"U",
"GEOMETRY",
634 auto add_boundary_ops_rhs = [&](
auto &pip) {
639 AT>::LinearForm<IT>::OpBaseTimesVector<1,
SPACE_DIM, 1>;
644 pip, mField,
"U", {time_scale}, Sev::inform);
647 auto u_disp = boost::make_shared<MatrixDouble>();
648 auto dot_u_disp = boost::make_shared<MatrixDouble>();
653 new OpSpringRhs(
"U", u_disp, [
this](
double,
double,
double) {
657 new OpSpringRhs(
"U", dot_u_disp, [
this](
double,
double,
double) {
663 ContactOps::opFactoryBoundaryRhs<SPACE_DIM, AT, GAUSS, BoundaryEleOp>(
669 CHKERR add_domain_base_ops(pip_mng->getOpDomainLhsPipeline());
670 CHKERR add_domain_base_ops(pip_mng->getOpDomainRhsPipeline());
671 CHKERR add_domain_ops_lhs(pip_mng->getOpDomainLhsPipeline());
672 CHKERR add_domain_ops_rhs(pip_mng->getOpDomainRhsPipeline());
674 CHKERR add_boundary_base_ops(pip_mng->getOpBoundaryLhsPipeline());
675 CHKERR add_boundary_base_ops(pip_mng->getOpBoundaryRhsPipeline());
676 CHKERR add_boundary_ops_lhs(pip_mng->getOpBoundaryLhsPipeline());
677 CHKERR add_boundary_ops_rhs(pip_mng->getOpBoundaryRhsPipeline());
679 if (mfrontInterface) {
680 CHKERR mfrontInterface->setUpdateElementVariablesOperators();
683 CHKERR pip_mng->setDomainRhsIntegrationRule(integration_rule_vol);
684 CHKERR pip_mng->setDomainLhsIntegrationRule(integration_rule_vol);
685 CHKERR pip_mng->setBoundaryRhsIntegrationRule(integration_rule_boundary);
686 CHKERR pip_mng->setBoundaryLhsIntegrationRule(integration_rule_boundary);
694 static boost::shared_ptr<SetUpSchur>
709 auto set_section_monitor = [&](
auto solver) {
712 CHKERR TSGetSNES(solver, &snes);
713 PetscViewerAndFormat *vf;
714 CHKERR PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,
715 PETSC_VIEWER_DEFAULT, &vf);
718 (
MoFEMErrorCode(*)(SNES, PetscInt, PetscReal,
void *))SNESMonitorFields,
723 auto scatter_create = [&](
auto D,
auto coeff) {
725 CHKERR is_manager->isCreateProblemFieldAndRank(
simple->getProblemName(),
726 ROW,
"U", coeff, coeff, is);
728 CHKERR ISGetLocalSize(is, &loc_size);
730 CHKERR VecCreateMPI(mField.get_comm(), loc_size, PETSC_DETERMINE, &
v);
732 CHKERR VecScatterCreate(
D, is,
v, PETSC_NULL, &scatter);
737 auto set_time_monitor = [&](
auto dm,
auto solver) {
739 monitorPtr->setScatterVectors(uXScatter, uYScatter, uZScatter);
740 boost::shared_ptr<ForcesAndSourcesCore>
null;
742 monitorPtr,
null,
null);
746 auto set_essential_bc = [&]() {
750 auto pre_proc_ptr = boost::make_shared<FEMethod>();
751 auto post_proc_rhs_ptr = boost::make_shared<FEMethod>();
752 auto post_proc_lhs_ptr = boost::make_shared<FEMethod>();
755 auto time_scale = boost::make_shared<TimeScale>();
757 auto get_bc_hook_rhs = [&]() {
759 {time_scale},
false);
762 pre_proc_ptr->preProcessHook = get_bc_hook_rhs();
764 auto get_post_proc_hook_rhs = [&]() {
766 mField, post_proc_rhs_ptr, 1.);
768 auto get_post_proc_hook_lhs = [&]() {
770 mField, post_proc_lhs_ptr, 1.);
772 post_proc_rhs_ptr->postProcessHook = get_post_proc_hook_rhs();
774 ts_ctx_ptr->getPreProcessIFunction().push_front(pre_proc_ptr);
775 ts_ctx_ptr->getPreProcessIJacobian().push_front(pre_proc_ptr);
776 ts_ctx_ptr->getPostProcessIFunction().push_back(post_proc_rhs_ptr);
777 post_proc_lhs_ptr->postProcessHook = get_post_proc_hook_lhs();
778 ts_ctx_ptr->getPostProcessIJacobian().push_back(post_proc_lhs_ptr);
783 auto set_schur_pc = [&](
auto solver) {
784 boost::shared_ptr<SetUpSchur> schur_ptr;
793 auto dm =
simple->getDM();
798 uXScatter = scatter_create(
D, 0);
799 uYScatter = scatter_create(
D, 1);
801 uZScatter = scatter_create(
D, 2);
805 CHKERR set_essential_bc();
808 auto solver = pip_mng->createTSIM();
809 CHKERR TSSetFromOptions(solver);
810 auto schur_pc_ptr = set_schur_pc(solver);
813 CHKERR set_section_monitor(solver);
814 CHKERR set_time_monitor(dm, solver);
815 CHKERR TSSetSolution(solver,
D);
817 CHKERR TSSolve(solver, NULL);
819 auto solver = pip_mng->createTSIM2();
820 CHKERR TSSetFromOptions(solver);
821 auto schur_pc_ptr = set_schur_pc(solver);
823 auto dm =
simple->getDM();
827 CHKERR set_section_monitor(solver);
828 CHKERR set_time_monitor(dm, solver);
829 CHKERR TS2SetSolution(solver,
D, DD);
831 CHKERR TSSolve(solver, NULL);
841 if (
atom_test && !mField.get_comm_rank()) {
846 double analytical_active_area = 1.0;
848 double tol_force = 1e-3;
849 double tol_norm = 7.5;
850 double tol_area = 3e-2;
851 double fem_active_area = t_ptr[3];
856 fem_force = t_ptr[1];
861 fem_force = t_ptr[1];
862 norm = monitorPtr->getErrorNorm(1);
868 fem_force = t_ptr[2];
869 analytical_active_area = M_PI / 4;
879 hertz_force = 15.873;
881 fem_force = t_ptr[1];
882 norm = monitorPtr->getErrorNorm(1);
883 analytical_active_area = M_PI;
888 fem_force = t_ptr[1];
892 hertz_force = 0.5289;
893 fem_force = t_ptr[2];
898 "atom test %d does not exist",
atom_test);
900 if (fabs(fem_force - hertz_force) / hertz_force > tol_force) {
902 "atom test %d failed: Wrong FORCE output: %3.4e != %3.4e",
905 if (norm > tol_norm) {
907 "atom test %d failed: Wrong NORM output: %3.4e > %3.4e",
910 if (fabs(fem_active_area - analytical_active_area) > tol_area) {
912 "atom test %d failed: AREA computed %3.4e but should be %3.4e",
913 atom_test, fem_active_area, analytical_active_area);
924 static char help[] =
"...\n\n";
926 int main(
int argc,
char *argv[]) {
934 const char param_file[] =
"param_file.petsc";
938 auto core_log = logging::core::get();
947 DMType dm_name =
"DMMOFEM";
949 DMType dm_name_mg =
"DMMOFEM_MG";
979 if (Py_FinalizeEx() < 0) {
1014 CHKERR TSGetSNES(solver, &snes);
1016 CHKERR SNESGetKSP(snes, &ksp);
1017 CHKERR KSPSetFromOptions(ksp);
1020 CHKERR KSPGetPC(ksp, &pc);
1022 PetscBool is_pcfs = PETSC_FALSE;
1023 PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &is_pcfs);
1026 MOFEM_LOG(
"CONTACT", Sev::inform) <<
"Setup Schur pc";
1031 "It is expected that Schur matrix is not allocated. This is "
1032 "possible only if PC is set up twice");
1044 CHKERR TSGetDM(solver, &solver_dm);
1045 CHKERR DMSetMatType(solver_dm, MATSHELL);
1052 auto swap_assemble = [](TS ts, PetscReal
t,
Vec u,
Vec u_t, PetscReal
a,
1053 Mat
A, Mat B,
void *ctx) {
1056 CHKERR TSSetIJacobian(solver,
A,
P, swap_assemble, ts_ctx_ptr.get());
1058 auto swap_assemble = [](TS ts, PetscReal
t,
Vec u,
Vec u_t,
Vec utt,
1059 PetscReal
a, PetscReal aa, Mat
A, Mat B,
1063 CHKERR TSSetI2Jacobian(solver,
A,
P, swap_assemble, ts_ctx_ptr.get());
1071 CHKERR setDiagonalPC(pc);
1074 MOFEM_LOG(
"CONTACT", Sev::inform) <<
"No Schur PC";
1087 auto create_dm = [&](
const char *name,
const char *
field_name,
auto dm_type) {
1088 auto dm =
createDM(mField.get_comm(), dm_type);
1089 auto create_dm_imp = [&]() {
1101 "Error in creating schurDM. It is possible that schurDM is "
1109 schurDM = create_dm(
"SCHUR",
"U",
"DMMOFEM_MG");
1110 blockDM = create_dm(
"BLOCK",
"SIGMA",
"DMMOFEM");
1114 auto get_nested_mat_data = [&](
auto schur_dm,
auto block_dm) {
1120 simple->getDomainFEName(),
1124 {
"U",
"U"}, {
"SIGMA",
"U"}, {
"U",
"SIGMA"}, {
"SIGMA",
"SIGMA"}
1132 {schur_dm, block_dm}, block_mat_data,
1134 {
"SIGMA"}, {
nullptr}, true
1139 auto nested_mat_data = get_nested_mat_data(schurDM, blockDM);
1144 "Only BLOCK_SCHUR is implemented");
1153 double eps_stab = 1e-4;
1159 using OpMassStab = B::OpMass<3, SPACE_DIM * SPACE_DIM>;
1169 pip->getOpBoundaryLhsPipeline().push_back(
1170 new OpMassStab(
"SIGMA",
"SIGMA",
1171 [eps_stab](
double,
double,
double) {
return eps_stab; }));
1172 pip->getOpBoundaryLhsPipeline().push_back(
1180 pip->getOpDomainLhsPipeline().push_back(
1187 auto pre_proc_schur_lhs_ptr = boost::make_shared<FEMethod>();
1188 auto post_proc_schur_lhs_ptr = boost::make_shared<FEMethod>();
1190 pre_proc_schur_lhs_ptr->preProcessHook = [
this]() {
1192 CHKERR MatZeroEntries(S);
1193 MOFEM_LOG(
"CONTACT", Sev::verbose) <<
"Lhs Assemble Begin";
1197 post_proc_schur_lhs_ptr->postProcessHook = [
this, ao_up,
1198 post_proc_schur_lhs_ptr]() {
1200 MOFEM_LOG(
"CONTACT", Sev::verbose) <<
"Lhs Assemble End";
1201 auto print_mat_norm = [
this](
auto a, std::string prefix) {
1204 CHKERR MatNorm(
a, NORM_FROBENIUS, &nrm);
1205 MOFEM_LOG(
"CONTACT", Sev::noisy) << prefix <<
" norm = " << nrm;
1208 CHKERR MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY);
1209 CHKERR MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY);
1211 mField, post_proc_schur_lhs_ptr, 1, S, ao_up)();
1213 CHKERR print_mat_norm(S,
"S");
1215 MOFEM_LOG(
"CONTACT", Sev::verbose) <<
"Lhs Assemble Finish";
1220 ts_ctx_ptr->getPreProcessIJacobian().push_front(pre_proc_schur_lhs_ptr);
1221 ts_ctx_ptr->getPostProcessIJacobian().push_back(post_proc_schur_lhs_ptr);
1229 auto block_is =
getDMSubData(blockDM)->getSmartRowIs();
1230 CHKERR PCFieldSplitSetIS(pc, NULL, block_is);
1231 CHKERR PCFieldSplitSetSchurPre(pc, PC_FIELDSPLIT_SCHUR_PRE_USER, S);
1238 CHKERR PCFieldSplitSchurGetSubKSP(pc, PETSC_NULL, &subksp);
1239 auto get_pc = [](
auto ksp) {
1241 CHKERR KSPGetPC(ksp, &pc_raw);
1246 auto set_pc_p_mg = [](
auto dm,
auto pc,
auto S) {
1249 PetscBool same = PETSC_FALSE;
1250 PetscObjectTypeCompare((PetscObject)pc, PCMG, &same);
1254 CHKERR PCSetFromOptions(pc);
1259 CHKERR set_pc_p_mg(schurDM, get_pc(subksp[1]), S);
1261 CHKERR PetscFree(subksp);
1265 boost::shared_ptr<SetUpSchur>
1267 return boost::shared_ptr<SetUpSchur>(
new SetUpSchurImpl(m_field));