v0.14.0
contact.cpp
Go to the documentation of this file.
1 /**
2  * \file contact.cpp
3  * \CONTACT contact.cpp
4  *
5  * CONTACT of contact problem
6  *
7  * @copyright Copyright (c) 2023
8  */
9 
10 /* The above code is a preprocessor directive in C++ that checks if the macro
11 "EXECUTABLE_DIMENSION" has been defined. If it has not been defined, it is set
12 to 3" */
13 #ifndef EXECUTABLE_DIMENSION
14 #define EXECUTABLE_DIMENSION 3
15 #endif
16 
17 #ifndef SCHUR_ASSEMBLE
18 #define SCHUR_ASSEMBLE 0
19 #endif
20 
21 #include <MoFEM.hpp>
22 #include <MatrixFunction.hpp>
23 
24 using namespace MoFEM;
25 
26 #include <GenericElementInterface.hpp>
27 
28 #ifdef PYTHON_SDF
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;
34 #endif
35 
36 constexpr AssemblyType AT =
38  : AssemblyType::PETSC; //< selected assembly type
39 constexpr IntegrationType IT =
40  IntegrationType::GAUSS; //< selected integration type
41 
42 template <int DIM> struct ElementsAndOps;
43 
44 template <> struct ElementsAndOps<2> : PipelineManager::ElementsAndOpsByDim<2> {
45  static constexpr FieldSpace CONTACT_SPACE = HCURL;
46 };
47 
48 template <> struct ElementsAndOps<3> : PipelineManager::ElementsAndOpsByDim<3> {
49  static constexpr FieldSpace CONTACT_SPACE = HDIV;
50 };
51 
54 
55 constexpr int SPACE_DIM =
56  EXECUTABLE_DIMENSION; //< Space dimension of problem, mesh
57 
58 /* The above code is defining an alias `EntData` for the type
59 `EntitiesFieldData::EntData`. This is a C++ syntax for creating a new name for
60 an existing type. */
67 //! [Specialisation for assembly]
68 
70 
71 //! [Operators used for contact]
75  IT>::OpBaseTimesVector<1, SPACE_DIM, 1>;
76 //! [Operators used for contact]
77 
78 PetscBool is_quasi_static = PETSC_TRUE;
79 
80 int order = 2; //< Order of displacements in the domain
81 int contact_order = 2; //< Order of displacements in boundary and side elements
82 int sigma_order = 1; //< Order of Lagrange multiplier in side elements
83 int geom_order = 1;
84 double young_modulus = 100;
85 double poisson_ratio = 0.25;
86 double rho = 0.0;
87 double spring_stiffness = 0.0;
88 double vis_spring_stiffness = 0.0;
89 double alpha_damping = 0;
90 
91 double scale = 1.;
92 
93 PetscBool is_axisymmetric = PETSC_FALSE; //< Axisymmetric model
94 
95 // #define HENCKY_SMALL_STRAIN
96 
97 int atom_test = 0;
98 
99 namespace ContactOps {
100 double cn_contact = 0.1;
101 }; // namespace ContactOps
102 
103 #include <HenckyOps.hpp>
104 using namespace HenckyOps;
105 #include <ContactOps.hpp>
106 #include <PostProcContact.hpp>
107 #include <ContactNaturalBC.hpp>
108 
109 #ifdef WITH_MODULE_MFRONT_INTERFACE
110 #include <MFrontMoFEMInterface.hpp>
111 #endif
112 
114 using OpDomainRhsBCs =
117 using OpBoundaryRhsBCs =
120 using OpBoundaryLhsBCs =
122 
123 using namespace ContactOps;
124 
125 struct Contact {
126 
127  Contact(MoFEM::Interface &m_field) : mField(m_field) {}
128 
129  MoFEMErrorCode runProblem();
130 
131 private:
133 
134  MoFEMErrorCode setupProblem();
135  MoFEMErrorCode createCommonData();
136  MoFEMErrorCode bC();
137  MoFEMErrorCode OPs();
138  MoFEMErrorCode tsSolve();
139  MoFEMErrorCode checkResults();
140 
141  std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uXScatter;
142  std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uYScatter;
143  std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uZScatter;
144 
145  boost::shared_ptr<GenericElementInterface> mfrontInterface;
146  boost::shared_ptr<Monitor> monitorPtr;
147 
148 #ifdef PYTHON_SDF
149  boost::shared_ptr<SDFPython> sdfPythonPtr;
150 #endif
151 
154  double getScale(const double time) {
155  return scale * MoFEM::TimeScale::getScale(time);
156  };
157  };
158 };
159 
160 //! [Run problem]
163  CHKERR setupProblem();
164  CHKERR createCommonData();
165  CHKERR bC();
166  CHKERR OPs();
167  CHKERR tsSolve();
168  CHKERR checkResults();
170 }
171 //! [Run problem]
172 
173 //! [Set up problem]
176  Simple *simple = mField.getInterface<Simple>();
177 
178  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
179  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-contact_order", &contact_order,
180  PETSC_NULL);
181  sigma_order = std::max(order, contact_order) - 1;
182  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-sigma_order", &sigma_order,
183  PETSC_NULL);
184  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-geom_order", &geom_order,
185  PETSC_NULL);
186 
187  MOFEM_LOG("CONTACT", Sev::inform) << "Order " << order;
188  MOFEM_LOG("CONTACT", Sev::inform) << "Contact order " << contact_order;
189  MOFEM_LOG("CONTACT", Sev::inform) << "Sigma order " << sigma_order;
190  MOFEM_LOG("CONTACT", Sev::inform) << "Geom order " << geom_order;
191 
192  // Select base
193  enum bases { AINSWORTH, DEMKOWICZ, LASBASETOPT };
194  const char *list_bases[LASBASETOPT] = {"ainsworth", "demkowicz"};
195  PetscInt choice_base_value = AINSWORTH;
196  CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases,
197  LASBASETOPT, &choice_base_value, PETSC_NULL);
198 
200  switch (choice_base_value) {
201  case AINSWORTH:
203  MOFEM_LOG("CONTACT", Sev::inform)
204  << "Set AINSWORTH_LEGENDRE_BASE for displacements";
205  break;
206  case DEMKOWICZ:
207  base = DEMKOWICZ_JACOBI_BASE;
208  MOFEM_LOG("CONTACT", Sev::inform)
209  << "Set DEMKOWICZ_JACOBI_BASE for displacements";
210  break;
211  default:
212  base = LASTBASE;
213  break;
214  }
215 
216  // Note: For tets we have only H1 Ainsworth base, for Hex we have only H1
217  // Demkowicz base. We need to implement Demkowicz H1 base on tet.
218  CHKERR simple->addDomainField("U", H1, base, SPACE_DIM);
219  CHKERR simple->addBoundaryField("U", H1, base, SPACE_DIM);
220  CHKERR simple->addDomainField("SIGMA", CONTACT_SPACE, DEMKOWICZ_JACOBI_BASE,
221  SPACE_DIM);
222  CHKERR simple->addBoundaryField("SIGMA", CONTACT_SPACE, DEMKOWICZ_JACOBI_BASE,
223  SPACE_DIM);
224  CHKERR simple->addDataField("GEOMETRY", H1, base, SPACE_DIM);
225 
226  CHKERR simple->setFieldOrder("U", order);
227  CHKERR simple->setFieldOrder("GEOMETRY", geom_order);
228 
229  auto get_skin = [&]() {
230  Range body_ents;
231  CHKERR mField.get_moab().get_entities_by_dimension(0, SPACE_DIM, body_ents);
232  Skinner skin(&mField.get_moab());
233  Range skin_ents;
234  CHKERR skin.find_skin(0, body_ents, false, skin_ents);
235  return skin_ents;
236  };
237 
238  auto filter_blocks = [&](auto skin) {
239  bool is_contact_block = false;
240  Range contact_range;
241  for (auto m :
242  mField.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
243 
244  (boost::format("%s(.*)") % "CONTACT").str()
245 
246  ))
247 
248  ) {
249  is_contact_block =
250  true; ///< blocs interation is collective, so that is set irrespective
251  ///< if there are entities in given rank or not in the block
252  MOFEM_LOG("CONTACT", Sev::inform)
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);
258 
259  CHKERR mField.getInterface<CommInterface>()->synchroniseEntities(
260  contact_meshset_range);
261  contact_range.merge(contact_meshset_range);
262  }
263  if (is_contact_block) {
264  MOFEM_LOG("SYNC", Sev::inform)
265  << "Nb entities in contact surface: " << contact_range.size();
266  MOFEM_LOG_SYNCHRONISE(mField.get_comm());
267  skin = intersect(skin, contact_range);
268  }
269  return skin;
270  };
271 
272  auto filter_true_skin = [&](auto skin) {
273  Range boundary_ents;
274  ParallelComm *pcomm =
275  ParallelComm::get_pcomm(&mField.get_moab(), MYPCOMM_INDEX);
276  CHKERR pcomm->filter_pstatus(skin, PSTATUS_SHARED | PSTATUS_MULTISHARED,
277  PSTATUS_NOT, -1, &boundary_ents);
278  return boundary_ents;
279  };
280 
281  auto boundary_ents = filter_true_skin(filter_blocks(get_skin()));
282  CHKERR simple->setFieldOrder("SIGMA", 0);
283  CHKERR simple->setFieldOrder("SIGMA", sigma_order, &boundary_ents);
284 
285  if (contact_order > order) {
286  Range ho_ents;
287 
288  CHKERR mField.get_moab().get_adjacencies(boundary_ents, 1, false, ho_ents,
289  moab::Interface::UNION);
290 
291  CHKERR mField.getInterface<CommInterface>()->synchroniseEntities(ho_ents);
292  CHKERR simple->setFieldOrder("U", contact_order, &ho_ents);
293  CHKERR mField.getInterface<CommInterface>()->synchroniseFieldEntities("U");
294  }
295 
296  CHKERR simple->setUp();
297 
298  auto project_ho_geometry = [&]() {
299  Projection10NodeCoordsOnField ent_method(mField, "GEOMETRY");
300  return mField.loop_dofs("GEOMETRY", ent_method);
301  };
302 
303  PetscBool project_geometry = PETSC_TRUE;
304  CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-project_geometry",
305  &project_geometry, PETSC_NULL);
306  if (project_geometry) {
307  CHKERR project_ho_geometry();
308  }
309 
311 } //! [Set up problem]
312 
313 //! [Create common data]
316 
317  PetscBool use_mfront = PETSC_FALSE;
318  CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-use_mfront", &use_mfront,
319  PETSC_NULL);
320  CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-is_axisymmetric",
321  &is_axisymmetric, PETSC_NULL);
322  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-atom_test", &atom_test,
323  PETSC_NULL);
324 
325  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-scale", &scale, PETSC_NULL);
326  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-young_modulus", &young_modulus,
327  PETSC_NULL);
328  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-poisson_ratio", &poisson_ratio,
329  PETSC_NULL);
330  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-rho", &rho, PETSC_NULL);
331  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-cn", &cn_contact, PETSC_NULL);
332  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-spring_stiffness",
333  &spring_stiffness, PETSC_NULL);
334  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-vis_spring_stiffness",
335  &vis_spring_stiffness, PETSC_NULL);
336  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-alpha_damping", &alpha_damping,
337  PETSC_NULL);
338 
339  if (!mfrontInterface) {
340  MOFEM_LOG("CONTACT", Sev::inform) << "Young modulus " << young_modulus;
341  MOFEM_LOG("CONTACT", Sev::inform) << "Poisson_ratio " << poisson_ratio;
342  } else {
343  MOFEM_LOG("CONTACT", Sev::inform) << "Using MFront for material model";
344  }
345 
346  MOFEM_LOG("CONTACT", Sev::inform) << "Density " << rho;
347  MOFEM_LOG("CONTACT", Sev::inform) << "cn_contact " << cn_contact;
348  MOFEM_LOG("CONTACT", Sev::inform) << "Spring stiffness " << spring_stiffness;
349  MOFEM_LOG("CONTACT", Sev::inform)
350  << "Vis spring_stiffness " << vis_spring_stiffness;
351 
352  MOFEM_LOG("CONTACT", Sev::inform) << "alpha_damping " << alpha_damping;
353 
354  PetscBool use_scale = PETSC_FALSE;
355  CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-use_scale", &use_scale,
356  PETSC_NULL);
357  if (use_scale) {
358  scale /= young_modulus;
359  }
360 
361  MOFEM_LOG("CONTACT", Sev::inform) << "Scale " << scale;
362 
363  CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-is_quasi_static",
364  &is_quasi_static, PETSC_NULL);
365  MOFEM_LOG("CONTACT", Sev::inform)
366  << "Is quasi-static: " << (is_quasi_static ? "true" : "false");
367 
368 #ifdef PYTHON_SDF
369  char sdf_file_name[255];
370  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-sdf_file",
371  sdf_file_name, 255, PETSC_NULL);
372 
373  sdfPythonPtr = boost::make_shared<SDFPython>();
374  CHKERR sdfPythonPtr->sdfInit(sdf_file_name);
375  sdfPythonWeakPtr = sdfPythonPtr;
376 #endif
377 
378  if (is_axisymmetric) {
379  if (SPACE_DIM == 3) {
380  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
381  "Use executable contact_2d with axisymmetric model");
382  } else {
383  if (!use_mfront) {
384  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
385  "Axisymmetric model is only available with MFront (set "
386  "use_mfront to 1)");
387  } else {
388  MOFEM_LOG("CONTACT", Sev::inform) << "Using axisymmetric model";
389  }
390  }
391  } else {
392  if (SPACE_DIM == 2) {
393  MOFEM_LOG("CONTACT", Sev::inform) << "Using plane strain model";
394  }
395  }
396 
397  if (use_mfront) {
398 #ifndef WITH_MODULE_MFRONT_INTERFACE
399  SETERRQ(
400  PETSC_COMM_SELF, MOFEM_NOT_FOUND,
401  "MFrontInterface module was not found while use_mfront was set to 1");
402 #else
403  if (SPACE_DIM == 3) {
404  mfrontInterface =
405  boost::make_shared<MFrontMoFEMInterface<TRIDIMENSIONAL>>(
406  mField, "U", "GEOMETRY", true, is_quasi_static);
407  } else if (SPACE_DIM == 2) {
408  if (is_axisymmetric) {
409  mfrontInterface =
410  boost::make_shared<MFrontMoFEMInterface<AXISYMMETRICAL>>(
411  mField, "U", "GEOMETRY", true, is_quasi_static);
412  } else {
413  mfrontInterface = boost::make_shared<MFrontMoFEMInterface<PLANESTRAIN>>(
414  mField, "U", "GEOMETRY", true, is_quasi_static);
415  }
416  }
417 #endif
418  CHKERR mfrontInterface->getCommandLineParameters();
419  }
420 
421  Simple *simple = mField.getInterface<Simple>();
422  auto dm = simple->getDM();
423  monitorPtr =
424  boost::make_shared<Monitor>(dm, scale, mfrontInterface, is_axisymmetric);
425 
426  if (use_mfront) {
427  mfrontInterface->setMonitorPtr(monitorPtr);
428  }
429 
431 }
432 //! [Create common data]
433 
434 //! [Boundary condition]
437  auto bc_mng = mField.getInterface<BcManager>();
438  auto simple = mField.getInterface<Simple>();
439 
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);
449  }
450 
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);
461 
462  // Note remove has to be always before push. Then node marking will be
463  // corrupted.
464  CHKERR bc_mng->pushMarkDOFsOnEntities<DisplacementCubitBcData>(
465  simple->getProblemName(), "U");
466 
468 }
469 //! [Boundary condition]
470 
471 //! [Push operators to pip]
474  auto simple = mField.getInterface<Simple>();
475  auto *pip_mng = mField.getInterface<PipelineManager>();
476  auto bc_mng = mField.getInterface<BcManager>();
477  auto time_scale = boost::make_shared<ScaledTimeScale>();
478  auto body_force_time_scale =
479  boost::make_shared<ScaledTimeScale>("body_force_hist.txt");
480 
481  auto integration_rule_vol = [](int, int, int approx_order) {
482  return 2 * approx_order + geom_order - 1;
483  };
484  auto integration_rule_boundary = [](int, int, int approx_order) {
485  return 2 * approx_order + geom_order - 1;
486  };
487 
488  auto add_domain_base_ops = [&](auto &pip) {
491  "GEOMETRY");
493  };
494 
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>();
498 
499  auto add_domain_ops_lhs = [&](auto &pip) {
501 
502  //! [Only used for dynamics]
505  //! [Only used for dynamics]
506  if (is_quasi_static == PETSC_FALSE) {
507 
508  auto *pip_mng = mField.getInterface<PipelineManager>();
509  auto fe_domain_lhs = pip_mng->getDomainLhsFE();
510 
511  auto get_inertia_and_mass_damping =
512  [this, fe_domain_lhs](const double, const double, const double) {
513  return (rho * scale) * fe_domain_lhs->ts_aa +
514  (alpha_damping * scale) * fe_domain_lhs->ts_a;
515  };
516  pip.push_back(new OpMass("U", "U", get_inertia_and_mass_damping));
517  } else {
518 
519  auto *pip_mng = mField.getInterface<PipelineManager>();
520  auto fe_domain_lhs = pip_mng->getDomainLhsFE();
521 
522  auto get_inertia_and_mass_damping =
523  [this, fe_domain_lhs](const double, const double, const double) {
524  return (alpha_damping * scale) * fe_domain_lhs->ts_a;
525  };
526  pip.push_back(new OpMass("U", "U", get_inertia_and_mass_damping));
527  }
528 
529  if (!mfrontInterface) {
530  CHKERR HenckyOps::opFactoryDomainLhs<SPACE_DIM, AT, IT, DomainEleOp>(
531  mField, pip, "U", "MAT_ELASTIC", Sev::verbose, scale);
532  } else {
533  CHKERR mfrontInterface->opFactoryDomainLhs(pip);
534  }
535 
537  };
538 
539  auto add_domain_ops_rhs = [&](auto &pip) {
541 
543  pip, mField, "U", {body_force_time_scale}, Sev::inform);
544 
545  //! [Only used for dynamics]
547  AT>::LinearForm<IT>::OpBaseTimesVector<1, SPACE_DIM, 1>;
548  //! [Only used for dynamics]
549 
550  // only in case of dynamics
551  if (is_quasi_static == PETSC_FALSE) {
552  auto mat_acceleration = boost::make_shared<MatrixDouble>();
554  "U", mat_acceleration));
555  pip.push_back(
556  new OpInertiaForce("U", mat_acceleration, [](double, double, double) {
557  return rho * scale;
558  }));
559  }
560 
561  // only in case of viscosity
562  if (alpha_damping > 0) {
563  auto mat_velocity = boost::make_shared<MatrixDouble>();
564  pip.push_back(
565  new OpCalculateVectorFieldValuesDot<SPACE_DIM>("U", mat_velocity));
566  pip.push_back(
567  new OpInertiaForce("U", mat_velocity, [](double, double, double) {
568  return alpha_damping * scale;
569  }));
570  }
571 
572  if (!mfrontInterface) {
573  CHKERR HenckyOps::opFactoryDomainRhs<SPACE_DIM, AT, IT, DomainEleOp>(
574  mField, pip, "U", "MAT_ELASTIC", Sev::inform, scale);
575  } else {
576  CHKERR mfrontInterface->opFactoryDomainRhs(pip);
577  }
578 
579  CHKERR ContactOps::opFactoryDomainRhs<SPACE_DIM, AT, IT, DomainEleOp>(
580  pip, "SIGMA", "U", is_axisymmetric);
581 
583  };
584 
585  auto add_boundary_base_ops = [&](auto &pip) {
588  "GEOMETRY");
589  // We have to integrate on curved face geometry, thus integration weight
590  // have to adjusted.
591  pip.push_back(new OpSetHOWeightsOnSubDim<SPACE_DIM>());
593  };
594 
595  auto add_boundary_ops_lhs = [&](auto &pip) {
597 
598  //! [Operators used for contact]
601  //! [Operators used for contact]
602 
603  // Add Natural BCs to LHS
605  pip, mField, "U", Sev::inform);
606 
607  if (spring_stiffness > 0 || vis_spring_stiffness > 0) {
608 
609  auto *pip_mng = mField.getInterface<PipelineManager>();
610  auto fe_boundary_lhs = pip_mng->getBoundaryLhsFE();
611 
612  pip.push_back(new OpSpringLhs(
613  "U", "U",
614 
615  [this, fe_boundary_lhs](double, double, double) {
616  return spring_stiffness * scale +
617  (vis_spring_stiffness * scale) * fe_boundary_lhs->ts_a;
618  }
619 
620  ));
621  }
622 
623  CHKERR
624  ContactOps::opFactoryBoundaryLhs<SPACE_DIM, AT, GAUSS, BoundaryEleOp>(
625  pip, "SIGMA", "U", is_axisymmetric);
627  DomainEle>(
628  mField, pip, simple->getDomainFEName(), "SIGMA", "U", "GEOMETRY",
629  integration_rule_vol, is_axisymmetric);
630 
632  };
633 
634  auto add_boundary_ops_rhs = [&](auto &pip) {
636 
637  //! [Operators used for contact]
639  AT>::LinearForm<IT>::OpBaseTimesVector<1, SPACE_DIM, 1>;
640  //! [Operators used for contact]
641 
642  // Add Natural BCs to RHS
644  pip, mField, "U", {time_scale}, Sev::inform);
645 
646  if (spring_stiffness > 0 || vis_spring_stiffness > 0) {
647  auto u_disp = boost::make_shared<MatrixDouble>();
648  auto dot_u_disp = boost::make_shared<MatrixDouble>();
649  pip.push_back(new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_disp));
650  pip.push_back(
651  new OpCalculateVectorFieldValuesDot<SPACE_DIM>("U", dot_u_disp));
652  pip.push_back(
653  new OpSpringRhs("U", u_disp, [this](double, double, double) {
654  return spring_stiffness * scale;
655  }));
656  pip.push_back(
657  new OpSpringRhs("U", dot_u_disp, [this](double, double, double) {
658  return vis_spring_stiffness * scale;
659  }));
660  }
661 
662  CHKERR
663  ContactOps::opFactoryBoundaryRhs<SPACE_DIM, AT, GAUSS, BoundaryEleOp>(
664  pip, "SIGMA", "U", is_axisymmetric);
665 
667  };
668 
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());
673 
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());
678 
679  if (mfrontInterface) {
680  CHKERR mfrontInterface->setUpdateElementVariablesOperators();
681  }
682 
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);
687 
689 }
690 //! [Push operators to pip]
691 
692 //! [Solve]
693 struct SetUpSchur {
694  static boost::shared_ptr<SetUpSchur>
695  createSetUpSchur(MoFEM::Interface &m_field);
696  virtual MoFEMErrorCode setUp(SmartPetscObj<TS> solver) = 0;
697 
698 protected:
699  SetUpSchur() = default;
700 };
701 
704 
705  Simple *simple = mField.getInterface<Simple>();
706  PipelineManager *pip_mng = mField.getInterface<PipelineManager>();
707  ISManager *is_manager = mField.getInterface<ISManager>();
708 
709  auto set_section_monitor = [&](auto solver) {
711  SNES snes;
712  CHKERR TSGetSNES(solver, &snes);
713  PetscViewerAndFormat *vf;
714  CHKERR PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,
715  PETSC_VIEWER_DEFAULT, &vf);
716  CHKERR SNESMonitorSet(
717  snes,
718  (MoFEMErrorCode(*)(SNES, PetscInt, PetscReal, void *))SNESMonitorFields,
719  vf, (MoFEMErrorCode(*)(void **))PetscViewerAndFormatDestroy);
721  };
722 
723  auto scatter_create = [&](auto D, auto coeff) {
725  CHKERR is_manager->isCreateProblemFieldAndRank(simple->getProblemName(),
726  ROW, "U", coeff, coeff, is);
727  int loc_size;
728  CHKERR ISGetLocalSize(is, &loc_size);
729  Vec v;
730  CHKERR VecCreateMPI(mField.get_comm(), loc_size, PETSC_DETERMINE, &v);
731  VecScatter scatter;
732  CHKERR VecScatterCreate(D, is, v, PETSC_NULL, &scatter);
733  return std::make_tuple(SmartPetscObj<Vec>(v),
734  SmartPetscObj<VecScatter>(scatter));
735  };
736 
737  auto set_time_monitor = [&](auto dm, auto solver) {
739  monitorPtr->setScatterVectors(uXScatter, uYScatter, uZScatter);
740  boost::shared_ptr<ForcesAndSourcesCore> null;
741  CHKERR DMMoFEMTSSetMonitor(dm, solver, simple->getDomainFEName(),
742  monitorPtr, null, null);
744  };
745 
746  auto set_essential_bc = [&]() {
748  // This is low level pushing finite elements (pipelines) to solver
749  auto ts_ctx_ptr = getDMTsCtx(simple->getDM());
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>();
753 
754  // Add boundary condition scaling
755  auto time_scale = boost::make_shared<TimeScale>();
756 
757  auto get_bc_hook_rhs = [&]() {
758  EssentialPreProc<DisplacementCubitBcData> hook(mField, pre_proc_ptr,
759  {time_scale}, false);
760  return hook;
761  };
762  pre_proc_ptr->preProcessHook = get_bc_hook_rhs();
763 
764  auto get_post_proc_hook_rhs = [&]() {
766  mField, post_proc_rhs_ptr, 1.);
767  };
768  auto get_post_proc_hook_lhs = [&]() {
770  mField, post_proc_lhs_ptr, 1.);
771  };
772  post_proc_rhs_ptr->postProcessHook = get_post_proc_hook_rhs();
773 
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);
780  };
781 
782  // Set up Schur preconditioner
783  auto set_schur_pc = [&](auto solver) {
784  boost::shared_ptr<SetUpSchur> schur_ptr;
785  if (AT == AssemblyType::BLOCK_SCHUR) {
786  // Set up Schur preconditioner
787  schur_ptr = SetUpSchur::createSetUpSchur(mField);
788  CHK_MOAB_THROW(schur_ptr->setUp(solver), "SetUpSchur::setUp");
789  }
790  return schur_ptr;
791  };
792 
793  auto dm = simple->getDM();
794  auto D = createDMVector(dm);
795 
797 
798  uXScatter = scatter_create(D, 0);
799  uYScatter = scatter_create(D, 1);
800  if (SPACE_DIM == 3)
801  uZScatter = scatter_create(D, 2);
802 
803  // Add extra finite elements to SNES solver pipelines to resolve essential
804  // boundary conditions
805  CHKERR set_essential_bc();
806 
807  if (is_quasi_static == PETSC_TRUE) {
808  auto solver = pip_mng->createTSIM();
809  CHKERR TSSetFromOptions(solver);
810  auto schur_pc_ptr = set_schur_pc(solver);
811 
812  auto D = createDMVector(dm);
813  CHKERR set_section_monitor(solver);
814  CHKERR set_time_monitor(dm, solver);
815  CHKERR TSSetSolution(solver, D);
816  CHKERR TSSetUp(solver);
817  CHKERR TSSolve(solver, NULL);
818  } else {
819  auto solver = pip_mng->createTSIM2();
820  CHKERR TSSetFromOptions(solver);
821  auto schur_pc_ptr = set_schur_pc(solver);
822 
823  auto dm = simple->getDM();
824  auto D = createDMVector(dm);
825  auto DD = vectorDuplicate(D);
826 
827  CHKERR set_section_monitor(solver);
828  CHKERR set_time_monitor(dm, solver);
829  CHKERR TS2SetSolution(solver, D, DD);
830  CHKERR TSSetUp(solver);
831  CHKERR TSSolve(solver, NULL);
832  }
833 
835 }
836 //! [Solve]
837 
838 //! [Check]
841  if (atom_test && !mField.get_comm_rank()) {
842  const double *t_ptr;
843  CHKERR VecGetArrayRead(ContactOps::CommonData::totalTraction, &t_ptr);
844  double hertz_force;
845  double fem_force;
846  double analytical_active_area = 1.0;
847  double norm = 1e-5;
848  double tol_force = 1e-3;
849  double tol_norm = 7.5; // change when analytical functions are updated
850  double tol_area = 3e-2;
851  double fem_active_area = t_ptr[3];
852 
853  switch (atom_test) {
854  case 1: // plane stress
855  hertz_force = 3.927;
856  fem_force = t_ptr[1];
857  break;
858 
859  case 2: // plane strain
860  hertz_force = 4.675;
861  fem_force = t_ptr[1];
862  norm = monitorPtr->getErrorNorm(1);
863  break;
864 
865  case 3: // Hertz 3D
866  hertz_force = 3.968;
867  tol_force = 2e-3;
868  fem_force = t_ptr[2];
869  analytical_active_area = M_PI / 4;
870  tol_area = 0.2;
871  break;
872 
873  case 4: // axisymmetric
874  tol_force = 5e-3;
875  tol_area = 0.2;
876  // analytical_active_area = M_PI;
877 
878  case 5: // axisymmetric
879  hertz_force = 15.873;
880  tol_force = 5e-3;
881  fem_force = t_ptr[1];
882  norm = monitorPtr->getErrorNorm(1);
883  analytical_active_area = M_PI;
884  break;
885 
886  case 6: // wavy 2d
887  hertz_force = 0.374;
888  fem_force = t_ptr[1];
889  break;
890 
891  case 7: // wavy 3d
892  hertz_force = 0.5289;
893  fem_force = t_ptr[2];
894  break;
895 
896  default:
897  SETERRQ1(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
898  "atom test %d does not exist", atom_test);
899  }
900  if (fabs(fem_force - hertz_force) / hertz_force > tol_force) {
901  SETERRQ3(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
902  "atom test %d failed: Wrong FORCE output: %3.4e != %3.4e",
903  atom_test, fem_force, hertz_force);
904  }
905  if (norm > tol_norm) {
906  SETERRQ3(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
907  "atom test %d failed: Wrong NORM output: %3.4e > %3.4e",
908  atom_test, norm, tol_norm);
909  }
910  if (fabs(fem_active_area - analytical_active_area) > tol_area) {
911  SETERRQ3(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
912  "atom test %d failed: AREA computed %3.4e but should be %3.4e",
913  atom_test, fem_active_area, analytical_active_area);
914  }
915  CHKERR VecRestoreArrayRead(ContactOps::CommonData::totalTraction, &t_ptr);
916  }
917 
919 
921 }
922 //! [Check]
923 
924 static char help[] = "...\n\n";
925 
926 int main(int argc, char *argv[]) {
927 
928 #ifdef PYTHON_SDF
929  Py_Initialize();
930  np::initialize();
931 #endif
932 
933  // Initialisation of MoFEM/PETSc and MOAB data structures
934  const char param_file[] = "param_file.petsc";
935  MoFEM::Core::Initialize(&argc, &argv, param_file, help);
936 
937  // Add logging channel for CONTACT
938  auto core_log = logging::core::get();
939  core_log->add_sink(
941  LogManager::setLog("CONTACT");
942  MOFEM_LOG_TAG("CONTACT", "Indent");
943 
944  try {
945 
946  //! [Register MoFEM discrete manager in PETSc]
947  DMType dm_name = "DMMOFEM";
948  CHKERR DMRegister_MoFEM(dm_name);
949  DMType dm_name_mg = "DMMOFEM_MG";
951  //! [Register MoFEM discrete manager in PETSc
952 
953  //! [Create MoAB]
954  moab::Core mb_instance; ///< mesh database
955  moab::Interface &moab = mb_instance; ///< mesh database interface
956  //! [Create MoAB]
957 
958  //! [Create MoFEM]
959  MoFEM::Core core(moab); ///< finite element database
960  MoFEM::Interface &m_field = core; ///< finite element database interface
961  //! [Create MoFEM]
962 
963  //! [Load mesh]
964  Simple *simple = m_field.getInterface<Simple>();
965  CHKERR simple->getOptions();
966  CHKERR simple->loadFile("");
967  //! [Load mesh]
968 
969  //! [CONTACT]
970  Contact ex(m_field);
971  CHKERR ex.runProblem();
972  //! [CONTACT]
973  }
974  CATCH_ERRORS;
975 
977 
978 #ifdef PYTHON_SDF
979  if (Py_FinalizeEx() < 0) {
980  exit(120);
981  }
982 #endif
983 
984  return 0;
985 }
986 
987 struct SetUpSchurImpl : public SetUpSchur {
988 
989  SetUpSchurImpl(MoFEM::Interface &m_field) : SetUpSchur(), mField(m_field) {}
990 
991  virtual ~SetUpSchurImpl() {}
992 
993  MoFEMErrorCode setUp(SmartPetscObj<TS> solver);
994 
995 private:
996  MoFEMErrorCode createSubDM();
997  MoFEMErrorCode setOperator();
998  MoFEMErrorCode setPC(PC pc);
999  MoFEMErrorCode setDiagonalPC(PC pc);
1000 
1001  MoFEM::Interface &mField;
1002 
1006 };
1007 
1010  auto simple = mField.getInterface<Simple>();
1011  auto pip = mField.getInterface<PipelineManager>();
1012 
1013  SNES snes;
1014  CHKERR TSGetSNES(solver, &snes);
1015  KSP ksp;
1016  CHKERR SNESGetKSP(snes, &ksp);
1017  CHKERR KSPSetFromOptions(ksp);
1018 
1019  PC pc;
1020  CHKERR KSPGetPC(ksp, &pc);
1021 
1022  PetscBool is_pcfs = PETSC_FALSE;
1023  PetscObjectTypeCompare((PetscObject)pc, PCFIELDSPLIT, &is_pcfs);
1024  if (is_pcfs) {
1025 
1026  MOFEM_LOG("CONTACT", Sev::inform) << "Setup Schur pc";
1027 
1028  if (S) {
1031  "It is expected that Schur matrix is not allocated. This is "
1032  "possible only if PC is set up twice");
1033  }
1034 
1035  CHKERR createSubDM();
1036 
1037  // Add data to DM storage
1038  S = createDMMatrix(schurDM);
1039  CHKERR MatSetBlockSize(S, SPACE_DIM);
1040  // CHKERR MatSetOption(S, MAT_SYMMETRIC, PETSC_TRUE);
1041 
1042  // Set DM to use shell block matrix
1043  DM solver_dm;
1044  CHKERR TSGetDM(solver, &solver_dm);
1045  CHKERR DMSetMatType(solver_dm, MATSHELL);
1046 
1047  auto ts_ctx_ptr = getDMTsCtx(solver_dm);
1048  auto A = createDMBlockMat(simple->getDM());
1049  auto P = createDMNestSchurMat(simple->getDM());
1050 
1051  if (is_quasi_static == PETSC_TRUE) {
1052  auto swap_assemble = [](TS ts, PetscReal t, Vec u, Vec u_t, PetscReal a,
1053  Mat A, Mat B, void *ctx) {
1054  return TsSetIJacobian(ts, t, u, u_t, a, B, A, ctx);
1055  };
1056  CHKERR TSSetIJacobian(solver, A, P, swap_assemble, ts_ctx_ptr.get());
1057  } else {
1058  auto swap_assemble = [](TS ts, PetscReal t, Vec u, Vec u_t, Vec utt,
1059  PetscReal a, PetscReal aa, Mat A, Mat B,
1060  void *ctx) {
1061  return TsSetI2Jacobian(ts, t, u, u_t, utt, a, aa, B, A, ctx);
1062  };
1063  CHKERR TSSetI2Jacobian(solver, A, P, swap_assemble, ts_ctx_ptr.get());
1064  }
1065  CHKERR KSPSetOperators(ksp, A, P);
1066 
1067  CHKERR setOperator();
1068  CHKERR setPC(pc);
1069  CHKERR TSSetUp(solver);
1070  CHKERR KSPSetUp(ksp);
1071  CHKERR setDiagonalPC(pc);
1072 
1073  } else {
1074  MOFEM_LOG("CONTACT", Sev::inform) << "No Schur PC";
1075  pip->getOpBoundaryLhsPipeline().push_front(createOpSchurAssembleBegin());
1076  pip->getOpBoundaryLhsPipeline().push_back(createOpSchurAssembleEnd({}, {}));
1077  pip->getOpDomainLhsPipeline().push_front(createOpSchurAssembleBegin());
1078  pip->getOpDomainLhsPipeline().push_back(createOpSchurAssembleEnd({}, {}));
1079  }
1081 }
1082 
1085  auto simple = mField.getInterface<Simple>();
1086 
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 = [&]() {
1091  CHKERR DMMoFEMCreateSubDM(dm, simple->getDM(), name);
1092  CHKERR DMMoFEMSetSquareProblem(dm, PETSC_TRUE);
1093  CHKERR DMMoFEMAddElement(dm, simple->getDomainFEName());
1096  CHKERR DMSetUp(dm);
1098  };
1100  create_dm_imp(),
1101  "Error in creating schurDM. It is possible that schurDM is "
1102  "already created");
1103  return dm;
1104  };
1105 
1106  // Note: here we can make block with bubbles of "U" and "SIGMA" fields. See
1107  // vec-0 where bubbles are added.
1108 
1109  schurDM = create_dm("SCHUR", "U", "DMMOFEM_MG");
1110  blockDM = create_dm("BLOCK", "SIGMA", "DMMOFEM");
1111 
1112  if constexpr (AT == AssemblyType::BLOCK_SCHUR) {
1113 
1114  auto get_nested_mat_data = [&](auto schur_dm, auto block_dm) {
1115  auto block_mat_data = createBlockMatStructure(
1116  simple->getDM(),
1117 
1118  {{
1119 
1120  simple->getDomainFEName(),
1121 
1122  {
1123 
1124  {"U", "U"}, {"SIGMA", "U"}, {"U", "SIGMA"}, {"SIGMA", "SIGMA"}
1125 
1126  }}}
1127 
1128  );
1129 
1131 
1132  {schur_dm, block_dm}, block_mat_data,
1133 
1134  {"SIGMA"}, {nullptr}, true
1135 
1136  );
1137  };
1138 
1139  auto nested_mat_data = get_nested_mat_data(schurDM, blockDM);
1140  CHKERR DMMoFEMSetNestSchurData(simple->getDM(), nested_mat_data);
1141 
1142  } else {
1143  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
1144  "Only BLOCK_SCHUR is implemented");
1145  }
1146 
1148 }
1149 
1152 
1153  double eps_stab = 1e-4;
1154  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-eps_stab", &eps_stab,
1155  PETSC_NULL);
1156 
1159  using OpMassStab = B::OpMass<3, SPACE_DIM * SPACE_DIM>;
1160 
1161  auto simple = mField.getInterface<Simple>();
1162  auto pip = mField.getInterface<PipelineManager>();
1163 
1164  auto dm_is = getDMSubData(schurDM)->getSmartRowIs();
1165  auto ao_up = createAOMappingIS(dm_is, PETSC_NULL);
1166 
1167  // Boundary
1168  pip->getOpBoundaryLhsPipeline().push_front(createOpSchurAssembleBegin());
1169  pip->getOpBoundaryLhsPipeline().push_back(
1170  new OpMassStab("SIGMA", "SIGMA",
1171  [eps_stab](double, double, double) { return eps_stab; }));
1172  pip->getOpBoundaryLhsPipeline().push_back(
1173 
1174  createOpSchurAssembleEnd({"SIGMA"}, {nullptr}, ao_up, S, false, false)
1175 
1176  );
1177 
1178  // Domain
1179  pip->getOpDomainLhsPipeline().push_front(createOpSchurAssembleBegin());
1180  pip->getOpDomainLhsPipeline().push_back(
1181 
1182  createOpSchurAssembleEnd({"SIGMA"}, {nullptr}, ao_up, S, false,
1183  false)
1184 
1185  );
1186 
1187  auto pre_proc_schur_lhs_ptr = boost::make_shared<FEMethod>();
1188  auto post_proc_schur_lhs_ptr = boost::make_shared<FEMethod>();
1189 
1190  pre_proc_schur_lhs_ptr->preProcessHook = [this]() {
1192  CHKERR MatZeroEntries(S);
1193  MOFEM_LOG("CONTACT", Sev::verbose) << "Lhs Assemble Begin";
1195  };
1196 
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) {
1203  double nrm;
1204  CHKERR MatNorm(a, NORM_FROBENIUS, &nrm);
1205  MOFEM_LOG("CONTACT", Sev::noisy) << prefix << " norm = " << nrm;
1207  };
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)();
1212 #ifndef NDEBUG
1213  CHKERR print_mat_norm(S, "S");
1214 #endif // NDEBUG
1215  MOFEM_LOG("CONTACT", Sev::verbose) << "Lhs Assemble Finish";
1217  };
1218 
1219  auto ts_ctx_ptr = getDMTsCtx(simple->getDM());
1220  ts_ctx_ptr->getPreProcessIJacobian().push_front(pre_proc_schur_lhs_ptr);
1221  ts_ctx_ptr->getPostProcessIJacobian().push_back(post_proc_schur_lhs_ptr);
1222 
1224 }
1225 
1228  auto simple = mField.getInterface<Simple>();
1229  auto block_is = getDMSubData(blockDM)->getSmartRowIs();
1230  CHKERR PCFieldSplitSetIS(pc, NULL, block_is);
1231  CHKERR PCFieldSplitSetSchurPre(pc, PC_FIELDSPLIT_SCHUR_PRE_USER, S);
1233 }
1234 
1237  KSP *subksp;
1238  CHKERR PCFieldSplitSchurGetSubKSP(pc, PETSC_NULL, &subksp);
1239  auto get_pc = [](auto ksp) {
1240  PC pc_raw;
1241  CHKERR KSPGetPC(ksp, &pc_raw);
1242  return SmartPetscObj<PC>(pc_raw, true); // bump reference
1243  };
1244  CHKERR setSchurA00MatSolvePC(get_pc(subksp[0]));
1245 
1246  auto set_pc_p_mg = [](auto dm, auto pc, auto S) {
1248  CHKERR PCSetDM(pc, dm);
1249  PetscBool same = PETSC_FALSE;
1250  PetscObjectTypeCompare((PetscObject)pc, PCMG, &same);
1251  if (same) {
1253  pc, createPCMGSetUpViaApproxOrdersCtx(dm, S, true));
1254  CHKERR PCSetFromOptions(pc);
1255  }
1257  };
1258 
1259  CHKERR set_pc_p_mg(schurDM, get_pc(subksp[1]), S);
1260 
1261  CHKERR PetscFree(subksp);
1263 }
1264 
1265 boost::shared_ptr<SetUpSchur>
1267  return boost::shared_ptr<SetUpSchur>(new SetUpSchurImpl(m_field));
1268 }
MoFEM::NaturalBC::Assembly::LinearForm
Definition: Natural.hpp:67
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
CHK_MOAB_THROW
#define CHK_MOAB_THROW(err, msg)
Check error code of MoAB function and throw MoFEM exception.
Definition: definitions.h:589
MoFEM::createAOMappingIS
auto createAOMappingIS(IS isapp, IS ispetsc)
Creates an application mapping using two index sets.
Definition: PetscSmartObj.hpp:318
SetUpSchurImpl
Definition: test_broken_space.cpp:515
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:128
BoundaryEle
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
Definition: child_and_parent.cpp:39
MYPCOMM_INDEX
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:228
spring_stiffness
double spring_stiffness
Definition: contact.cpp:87
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
Contact::monitorPtr
boost::shared_ptr< Monitor > monitorPtr
Definition: contact.cpp:146
MoFEM::DMMoFEMAddSubFieldCol
PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[])
Definition: DMMoFEM.cpp:280
Contact::bC
MoFEMErrorCode bC()
[Create common data]
Definition: contact.cpp:435
LASTBASE
@ LASTBASE
Definition: definitions.h:69
ContactOps
Definition: contact.cpp:99
SetUpSchurImpl::setUp
MoFEMErrorCode setUp(SmartPetscObj< KSP >)
Definition: test_broken_space.cpp:528
IT
constexpr IntegrationType IT
Definition: contact.cpp:39
MoFEM::OpFluxRhsImpl
Definition: Natural.hpp:39
main
int main(int argc, char *argv[])
Definition: contact.cpp:926
MoFEM::PipelineManager::ElementsAndOpsByDim
Definition: PipelineManager.hpp:38
order
int order
Definition: contact.cpp:80
is_axisymmetric
PetscBool is_axisymmetric
Definition: contact.cpp:93
MoFEM::BLOCK_PRECONDITIONER_SCHUR
@ BLOCK_PRECONDITIONER_SCHUR
Definition: FormsIntegrators.hpp:109
MoFEM::EssentialPostProcLhs< DisplacementCubitBcData >
Specialization for DisplacementCubitBcData.
Definition: EssentialDisplacementCubitBcData.hpp:134
MoFEM::NaturalBC::Assembly::BiLinearForm
Definition: Natural.hpp:74
MoFEM::createPCMGSetUpViaApproxOrdersCtx
boost::shared_ptr< PCMGSetUpViaApproxOrdersCtx > createPCMGSetUpViaApproxOrdersCtx(DM dm, Mat A, bool use_shell_mat)
createPCMGSetUpViaApproxOrdersCtx
Definition: PCMGSetUpViaApproxOrders.cpp:630
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:609
MoFEM::DMMoFEMSetSquareProblem
PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem)
set squared problem
Definition: DMMoFEM.cpp:456
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::OpCalculateVectorFieldValues
Get values at integration pts for tensor filed rank 1, i.e. vector field.
Definition: UserDataOperators.hpp:466
SetUpSchurImpl::~SetUpSchurImpl
virtual ~SetUpSchurImpl()
Definition: contact.cpp:991
SetUpSchurImpl::setPC
MoFEMErrorCode setPC(PC pc)
Definition: contact.cpp:1226
MoFEM::EssentialPostProcRhs< DisplacementCubitBcData >
Specialization for DisplacementCubitBcData.
Definition: EssentialDisplacementCubitBcData.hpp:113
MoFEM::PETSC
@ PETSC
Definition: FormsIntegrators.hpp:105
Contact::uXScatter
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uXScatter
Definition: contact.cpp:141
MoFEM::PipelineManager
PipelineManager interface.
Definition: PipelineManager.hpp:24
MoFEM.hpp
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
MoFEM::DisplacementCubitBcData
Definition of the displacement bc data structure.
Definition: BCData.hpp:76
MoFEM::TsSetI2Jacobian
PetscErrorCode TsSetI2Jacobian(TS ts, PetscReal t, Vec u, Vec u_t, Vec u_tt, PetscReal a, PetscReal aa, Mat A, Mat B, void *ctx)
Calculation Jacobian for second order PDE in time.
Definition: TsCtx.cpp:511
MoFEM::Projection10NodeCoordsOnField
Projection of edge entities with one mid-node on hierarchical basis.
Definition: Projection10NodeCoordsOnField.hpp:24
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
rho
double rho
Definition: contact.cpp:86
HenckyOps
Definition: HenckyOps.hpp:11
Contact::mField
MoFEM::Interface & mField
Definition: contact.cpp:132
MoFEM::LogManager::createSink
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
Definition: LogManager.cpp:298
HenckyOps.hpp
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
ContactOps::opFactoryBoundaryToDomainLhs
MoFEMErrorCode opFactoryBoundaryToDomainLhs(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string fe_domain_name, std::string sigma, std::string u, std::string geom, ForcesAndSourcesCore::RuleHookFun rule, bool is_axisymmetric=false)
Definition: ContactOps.hpp:1265
MoFEM::getDMSubData
auto getDMSubData(DM dm)
Get sub problem data structure.
Definition: DMMoFEM.hpp:1157
SetUpSchurImpl::schurDM
SmartPetscObj< DM > schurDM
Definition: contact.cpp:1004
SetUpSchurImpl::createSubDM
MoFEMErrorCode createSubDM()
Definition: contact.cpp:1083
MoFEM::createDMMatrix
auto createDMMatrix(DM dm)
Get smart matrix from DM.
Definition: DMMoFEM.hpp:1056
MoFEM::DMMoFEMAddElement
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition: DMMoFEM.cpp:497
MoFEM::createDMBlockMat
auto createDMBlockMat(DM dm)
Definition: DMMoFEM.hpp:1076
SCHUR_ASSEMBLE
#define SCHUR_ASSEMBLE
Definition: contact.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
SPACE_DIM
constexpr int SPACE_DIM
Definition: contact.cpp:55
scale
double scale
Definition: contact.cpp:91
MoFEM::OpFluxLhsImpl
Definition: Natural.hpp:43
ROW
@ ROW
Definition: definitions.h:136
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
SetUpSchurImpl::blockDM
SmartPetscObj< DM > blockDM
Definition: contact.cpp:1005
MoFEM::OpCalculateVectorFieldValuesFromPetscVecImpl
Approximate field values for given petsc vector.
Definition: UserDataOperators.hpp:595
EntData
EntitiesFieldData::EntData EntData
Definition: contact.cpp:61
OpInertiaForce
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpBaseTimesVector< 1, SPACE_DIM, 1 > OpInertiaForce
Definition: dynamic_first_order_con_law.cpp:63
atom_test
int atom_test
Definition: contact.cpp:97
FieldSpace
FieldSpace
approximation spaces
Definition: definitions.h:82
MoFEM::createOpSchurAssembleEnd
OpSchurAssembleBase * createOpSchurAssembleEnd(std::vector< std::string > fields_name, std::vector< boost::shared_ptr< Range >> field_ents, SmartPetscObj< AO > ao, SmartPetscObj< Mat > schur, bool sym_schur, bool symm_op)
Construct a new Op Schur Assemble End object.
Definition: Schur.cpp:2186
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::createDMVector
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1099
MoFEM::EssentialPreProc< DisplacementCubitBcData >
Specialization for DisplacementCubitBcData.
Definition: EssentialDisplacementCubitBcData.hpp:91
vis_spring_stiffness
double vis_spring_stiffness
Definition: contact.cpp:88
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::ISManager
Section manager is used to create indexes and sections.
Definition: ISManager.hpp:23
MoFEM::createDM
auto createDM(MPI_Comm comm, const std::string dm_type_name)
Creates smart DM object.
Definition: PetscSmartObj.hpp:141
a
constexpr double a
Definition: approx_sphere.cpp:30
MoFEM::TimeScale::TimeScale
TimeScale(std::string file_name="", bool error_if_file_not_given=false)
TimeScale constructor.
Definition: ScalingMethod.cpp:22
MoFEM::BcManager
Simple interface for fast problem set-up.
Definition: BcManager.hpp:25
OpMass
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpMass
[Only used with Hooke equation (linear material model)]
Definition: seepage.cpp:57
BoundaryEleOp
ContactOps::CommonData::totalTraction
static SmartPetscObj< Vec > totalTraction
Definition: ContactOps.hpp:29
OpSpringRhs
FormsIntegrators< BoundaryEleOp >::Assembly< AT >::LinearForm< IT >::OpBaseTimesVector< 1, SPACE_DIM, 1 > OpSpringRhs
Definition: contact.cpp:75
alpha_damping
double alpha_damping
Definition: contact.cpp:89
simple
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69
MoFEM::DMMoFEMCreateSubDM
PetscErrorCode DMMoFEMCreateSubDM(DM subdm, DM dm, const char problem_name[])
Must be called by user to set Sub DM MoFEM data structures.
Definition: DMMoFEM.cpp:215
double
MoFEM::FormsIntegrators::Assembly
Assembly methods.
Definition: FormsIntegrators.hpp:317
MoFEM::DMRegister_MGViaApproxOrders
MoFEMErrorCode DMRegister_MGViaApproxOrders(const char sname[])
Register DM for Multi-Grid via approximation orders.
Definition: PCMGSetUpViaApproxOrders.cpp:302
poisson_ratio
double poisson_ratio
Definition: contact.cpp:85
MoFEM::TimeScale
Force scale operator for reading two columns.
Definition: ScalingMethod.hpp:32
MOFEM_LOG_SYNCHRONISE
#define MOFEM_LOG_SYNCHRONISE(comm)
Synchronise "SYNC" channel.
Definition: LogManager.hpp:345
Contact
Definition: contact.cpp:125
EshelbianPlasticity::P
@ P
Definition: EshelbianContact.cpp:197
Contact::uZScatter
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uZScatter
Definition: contact.cpp:143
MoFEM::LogManager::getStrmWorld
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
Definition: LogManager.cpp:344
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
SideEle
ElementsAndOps< SPACE_DIM >::SideEle SideEle
Definition: plastic.cpp:61
MoFEM::createOpSchurAssembleBegin
OpSchurAssembleBase * createOpSchurAssembleBegin()
Definition: Schur.cpp:2181
MoFEM::GAUSS
@ GAUSS
Definition: FormsIntegrators.hpp:136
contact_order
int contact_order
Definition: contact.cpp:81
MatrixFunction.hpp
Contact::mfrontInterface
boost::shared_ptr< GenericElementInterface > mfrontInterface
Definition: contact.cpp:145
MoFEM::AddFluxToLhsPipelineImpl
Definition: Natural.hpp:49
MoFEM::createDMNestSchurMat
auto createDMNestSchurMat(DM dm)
Definition: DMMoFEM.hpp:1083
MoFEM::TsSetIJacobian
PetscErrorCode TsSetIJacobian(TS ts, PetscReal t, Vec u, Vec u_t, PetscReal a, Mat A, Mat B, void *ctx)
Set function evaluating jacobian in TS solver.
Definition: TsCtx.cpp:165
Contact::Contact
Contact(MoFEM::Interface &m_field)
Definition: contact.cpp:127
EXECUTABLE_DIMENSION
#define EXECUTABLE_DIMENSION
Definition: contact.cpp:14
SetUpSchurImpl::setOperator
MoFEMErrorCode setOperator()
Definition: contact.cpp:1150
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
Contact::createCommonData
MoFEMErrorCode createCommonData()
[Set up problem]
Definition: contact.cpp:314
MoFEM::AddFluxToRhsPipelineImpl
Definition: Natural.hpp:46
t
constexpr double t
plate stiffness
Definition: plate.cpp:58
BiLinearForm
Contact::ScaledTimeScale
Definition: contact.cpp:152
MOFEM_NOT_FOUND
@ MOFEM_NOT_FOUND
Definition: definitions.h:33
MoFEM::OpSetHOWeightsOnSubDim
Definition: HODataOperators.hpp:145
Contact::ScaledTimeScale::getScale
double getScale(const double time)
Get scaling at a given time.
Definition: contact.cpp:154
MoFEM::TimeScale::getScale
double getScale(const double time)
Get scaling at a given time.
Definition: ScalingMethod.cpp:133
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
ContactOps.hpp
MoFEM::IntegrationType
IntegrationType
Form integrator integration types.
Definition: FormsIntegrators.hpp:136
MoFEM::AddHOOps
Add operators pushing bases from local to physical configuration.
Definition: HODataOperators.hpp:413
ContactOps::CommonData::createTotalTraction
static auto createTotalTraction(MoFEM::Interface &m_field)
Definition: ContactOps.hpp:31
ContactOps::cn_contact
double cn_contact
Definition: contact.cpp:100
MoFEM::setSchurA00MatSolvePC
MoFEMErrorCode setSchurA00MatSolvePC(SmartPetscObj< PC > pc)
Set PC for A00 block.
Definition: Schur.cpp:2223
Contact::OPs
MoFEMErrorCode OPs()
[Boundary condition]
Definition: contact.cpp:472
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
ElementsAndOps
Definition: child_and_parent.cpp:18
Range
DomainEleOp
is_quasi_static
PetscBool is_quasi_static
[Operators used for contact]
Definition: contact.cpp:78
MoFEM::CoreTmp< 0 >::Initialize
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
Contact::setupProblem
MoFEMErrorCode setupProblem()
[Run problem]
Definition: contact.cpp:174
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::NaturalBC::Assembly
Assembly methods.
Definition: Natural.hpp:65
MoFEM::vectorDuplicate
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
Definition: PetscSmartObj.hpp:221
young_modulus
double young_modulus
Definition: contact.cpp:84
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
MoFEM::AssemblyType
AssemblyType
[Storage and set boundary conditions]
Definition: FormsIntegrators.hpp:104
MoFEM::CommInterface
Managing BitRefLevels.
Definition: CommInterface.hpp:21
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
approx_order
int approx_order
Definition: test_broken_space.cpp:50
MoFEM::createSchurNestedMatrixStruture
boost::shared_ptr< NestSchurData > createSchurNestedMatrixStruture(std::pair< SmartPetscObj< DM >, SmartPetscObj< DM >> dms, boost::shared_ptr< BlockStructure > block_mat_data_ptr, std::vector< std::string > fields_names, std::vector< boost::shared_ptr< Range >> field_ents, bool add_preconditioner_block)
Get the Schur Nest Mat Array object.
Definition: Schur.cpp:1944
MoFEM::DMMoFEMTSSetMonitor
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:1056
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::MeshsetsManager
Interface for managing meshsets containing materials and boundary conditions.
Definition: MeshsetsManager.hpp:104
MoFEM::PetscOptionsGetEList
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
Definition: DeprecatedPetsc.hpp:203
Contact::tsSolve
MoFEMErrorCode tsSolve()
Definition: contact.cpp:702
Contact::runProblem
MoFEMErrorCode runProblem()
[Run problem]
Definition: contact.cpp:161
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
SetUpSchurImpl::SetUpSchurImpl
SetUpSchurImpl(MoFEM::Interface &m_field)
Definition: contact.cpp:989
geom_order
int geom_order
Definition: contact.cpp:83
ReactionDiffusionEquation::D
const double D
diffusivity
Definition: reaction_diffusion.cpp:20
MoFEM::BLOCK_SCHUR
@ BLOCK_SCHUR
Definition: FormsIntegrators.hpp:108
ContactNaturalBC.hpp
SetUpSchur
[Push operators to pipeline]
Definition: test_broken_space.cpp:40
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
MoFEM::DMMoFEMSetNestSchurData
MoFEMErrorCode DMMoFEMSetNestSchurData(DM dm, boost::shared_ptr< NestSchurData >)
Definition: DMMoFEM.cpp:1562
MoFEM::createBlockMatStructure
boost::shared_ptr< BlockStructure > createBlockMatStructure(DM dm, SchurFEOpsFEandFields schur_fe_op_vec)
Create a Mat Diag Blocks object.
Definition: Schur.cpp:1009
MoFEM::PCMGSetUpViaApproxOrders
MoFEMErrorCode PCMGSetUpViaApproxOrders(PC pc, boost::shared_ptr< PCMGSetUpViaApproxOrdersCtx > ctx, int verb)
Function build MG structure.
Definition: PCMGSetUpViaApproxOrders.cpp:634
MoFEM::PetscOptionsGetScalar
PetscErrorCode PetscOptionsGetScalar(PetscOptions *, const char pre[], const char name[], PetscScalar *dval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:162
MoFEM::DMMoFEMAddSubFieldRow
PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[])
Definition: DMMoFEM.cpp:238
DomainEle
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
Definition: child_and_parent.cpp:34
MoFEM::MeshsetsManager::getCubitMeshsetPtr
MoFEMErrorCode getCubitMeshsetPtr(const int ms_id, const CubitBCType cubit_bc_type, const CubitMeshSets **cubit_meshset_ptr) const
get cubit meshset
Definition: MeshsetsManager.cpp:578
help
static char help[]
[Check]
Definition: contact.cpp:924
MoFEM::LogManager::setLog
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:389
MoFEM::SmartPetscObj< VecScatter >
SetUpSchur::createSetUpSchur
static boost::shared_ptr< SetUpSchur > createSetUpSchur(MoFEM::Interface &m_field)
Definition: test_broken_space.cpp:768
OpSpringLhs
FormsIntegrators< BoundaryEleOp >::Assembly< AT >::BiLinearForm< IT >::OpMass< 1, SPACE_DIM > OpSpringLhs
[Operators used for contact]
Definition: contact.cpp:73
sigma_order
int sigma_order
Definition: contact.cpp:82
convert.int
int
Definition: convert.py:64
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
SetUpSchurImpl::setDiagonalPC
MoFEMErrorCode setDiagonalPC(PC pc)
Definition: contact.cpp:1235
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
AT
constexpr AssemblyType AT
Definition: contact.cpp:36
HDIV
@ HDIV
field with continuous normal traction
Definition: definitions.h:87
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
PostProcContact.hpp
MoFEM::getDMTsCtx
auto getDMTsCtx(DM dm)
Get TS context data structure used by DM.
Definition: DMMoFEM.hpp:1141
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
Contact::checkResults
MoFEMErrorCode checkResults()
[Solve]
Definition: contact.cpp:839
MOFEM_INVALID_DATA
@ MOFEM_INVALID_DATA
Definition: definitions.h:36
CONTACT_SPACE
constexpr FieldSpace CONTACT_SPACE
[Specialisation for assembly]
Definition: contact.cpp:69
MoFEM::PetscOptionsGetBool
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:182
Contact::uYScatter
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uYScatter
Definition: contact.cpp:142