v0.14.0
nonlinear_elastic.cpp
Go to the documentation of this file.
1 /**
2  * \file nonlinear_elastic.cpp
3  * \example nonlinear_elastic.cpp
4  *
5  * Plane stress elastic dynamic problem
6  *
7  */
8 
9 #include <MoFEM.hpp>
10 #include <MatrixFunction.hpp>
11 
12 using namespace MoFEM;
13 
14 constexpr int SPACE_DIM =
15  EXECUTABLE_DIMENSION; //< Space dimension of problem, mesh
16 
19 using BoundaryEle =
23 
24 using DomainNaturalBC =
26 using OpBodyForce =
28 
29 using BoundaryNaturalBC =
32 
33 template <int DIM> struct PostProcEleByDim;
34 
35 template <> struct PostProcEleByDim<2> {
39 };
40 
41 template <> struct PostProcEleByDim<3> {
45 };
46 
49 using PostProcEdges =
51 
52 #include <HenckyOps.hpp>
53 using namespace HenckyOps;
54 
55 struct Example {
56 
57  Example(MoFEM::Interface &m_field) : mField(m_field) {}
58 
59  MoFEMErrorCode runProblem();
60 
61 private:
62  MoFEM::Interface &mField;
63 
64  MoFEMErrorCode readMesh();
65  MoFEMErrorCode setupProblem();
66  MoFEMErrorCode boundaryCondition();
67  MoFEMErrorCode assembleSystem();
68  MoFEMErrorCode solveSystem();
69  MoFEMErrorCode gettingNorms();
70  MoFEMErrorCode outputResults();
71  MoFEMErrorCode checkResults();
72 };
73 
74 //! [Run problem]
77  CHKERR readMesh();
78  CHKERR setupProblem();
79  CHKERR boundaryCondition();
80  CHKERR assembleSystem();
81  CHKERR solveSystem();
82  CHKERR gettingNorms();
83  CHKERR outputResults();
84  CHKERR checkResults();
86 }
87 //! [Run problem]
88 
89 //! [Read mesh]
92  auto simple = mField.getInterface<Simple>();
93  CHKERR simple->getOptions();
94  char meshFileName[255];
95  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-file_name",
96  meshFileName, 255, PETSC_NULL);
97  CHKERR simple->loadFile("", meshFileName,
100 }
101 //! [Read mesh]
102 
103 //! [Set up problem]
106  Simple *simple = mField.getInterface<Simple>();
107 
108  enum bases { AINSWORTH, DEMKOWICZ, LASBASETOPT };
109  const char *list_bases[LASBASETOPT] = {"ainsworth", "demkowicz"};
110  PetscInt choice_base_value = AINSWORTH;
111  CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases,
112  LASBASETOPT, &choice_base_value, PETSC_NULL);
113 
115  switch (choice_base_value) {
116  case AINSWORTH:
118  MOFEM_LOG("WORLD", Sev::inform)
119  << "Set AINSWORTH_LEGENDRE_BASE for displacements";
120  break;
121  case DEMKOWICZ:
122  base = DEMKOWICZ_JACOBI_BASE;
123  MOFEM_LOG("WORLD", Sev::inform)
124  << "Set DEMKOWICZ_JACOBI_BASE for displacements";
125  break;
126  default:
127  base = LASTBASE;
128  break;
129  }
130 
131  // Add field
132  CHKERR simple->addDomainField("U", H1, base, SPACE_DIM);
133  CHKERR simple->addBoundaryField("U", H1, base, SPACE_DIM);
134  int order = 2;
135  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
136  CHKERR simple->setFieldOrder("U", order);
137  CHKERR simple->setUp();
139 }
140 //! [Set up problem]
141 
142 //! [Boundary condition]
145  auto *pipeline_mng = mField.getInterface<PipelineManager>();
146  auto simple = mField.getInterface<Simple>();
147  auto bc_mng = mField.getInterface<BcManager>();
148  auto time_scale = boost::make_shared<TimeScale>();
149 
150  auto integration_rule = [](int, int, int approx_order) {
151  return 2 * (approx_order - 1);
152  };
153 
154  CHKERR pipeline_mng->setDomainRhsIntegrationRule(integration_rule);
155  CHKERR pipeline_mng->setDomainLhsIntegrationRule(integration_rule);
156  CHKERR pipeline_mng->setBoundaryRhsIntegrationRule(integration_rule);
157 
159  pipeline_mng->getOpBoundaryRhsPipeline(), mField, "U", {time_scale},
160  "FORCE", Sev::inform);
161 
162  //! [Define gravity vector]
164  pipeline_mng->getOpDomainRhsPipeline(), mField, "U", {time_scale},
165  "BODY_FORCE", Sev::inform);
166 
167  // Essential BC
168  CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_X",
169  "U", 0, 0);
170  CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_Y",
171  "U", 1, 1);
172  CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_Z",
173  "U", 2, 2);
174  CHKERR bc_mng->pushMarkDOFsOnEntities<DisplacementCubitBcData>(
175  simple->getProblemName(), "U");
176 
177  // Essential MPCs BC
178  CHKERR bc_mng->addBlockDOFsToMPCs(simple->getProblemName(), "U");
179 
181 }
182 //! [Boundary condition]
183 
184 //! [Push operators to pipeline]
187  auto *simple = mField.getInterface<Simple>();
188  auto *pipeline_mng = mField.getInterface<PipelineManager>();
189 
190  auto add_domain_ops_lhs = [&](auto &pip) {
193  CHKERR HenckyOps::opFactoryDomainLhs<SPACE_DIM, PETSC, GAUSS, DomainEleOp>(
194  mField, pip, "U", "MAT_ELASTIC", Sev::inform);
196  };
197 
198  auto add_domain_ops_rhs = [&](auto &pip) {
201  CHKERR HenckyOps::opFactoryDomainRhs<SPACE_DIM, PETSC, GAUSS, DomainEleOp>(
202  mField, pip, "U", "MAT_ELASTIC", Sev::inform);
204  };
205 
206  CHKERR add_domain_ops_lhs(pipeline_mng->getOpDomainLhsPipeline());
207  CHKERR add_domain_ops_rhs(pipeline_mng->getOpDomainRhsPipeline());
208 
210 }
211 //! [Push operators to pipeline]
212 
213 /**
214  * @brief Monitor solution
215  *
216  * This functions is called by TS solver at the end of each step. It is used
217  * to output results to the hard drive.
218  */
219 struct Monitor : public FEMethod {
220  Monitor(SmartPetscObj<DM> dm, boost::shared_ptr<PostProcEleBdy> post_proc)
221  : dM(dm), postProc(post_proc){};
224  constexpr int save_every_nth_step = 1;
225  if (ts_step % save_every_nth_step == 0) {
226  CHKERR DMoFEMLoopFiniteElements(dM, "bFE", postProc);
228  postProc->mField, postProc->getPostProcMesh(), {"U"});
229 
230  CHKERR postProc->writeFile(
231  "out_step_" + boost::lexical_cast<std::string>(ts_step) + ".h5m");
232  }
234  }
235 
236 private:
238  boost::shared_ptr<PostProcEleBdy> postProc;
239 };
240 
241 //! [Solve]
244  auto *simple = mField.getInterface<Simple>();
245  auto *pipeline_mng = mField.getInterface<PipelineManager>();
246 
247  auto dm = simple->getDM();
248  auto ts = pipeline_mng->createTSIM();
249 
250  // Setup postprocessing
251  auto create_post_proc_fe = [dm, this, simple]() {
252  auto post_proc_ele_domain = [dm, this](auto &pip_domain) {
254  auto common_ptr = commonDataFactory<SPACE_DIM, GAUSS, DomainEleOp>(
255  mField, pip_domain, "U", "MAT_ELASTIC", Sev::inform);
256  return common_ptr;
257  };
258 
259  auto post_proc_fe_bdy = boost::make_shared<PostProcEleBdy>(mField);
260  auto u_ptr = boost::make_shared<MatrixDouble>();
261  post_proc_fe_bdy->getOpPtrVector().push_back(
263  auto op_loop_side =
264  new OpLoopSide<SideEle>(mField, simple->getDomainFEName(), SPACE_DIM);
265  auto common_ptr = post_proc_ele_domain(op_loop_side->getOpPtrVector());
266  post_proc_fe_bdy->getOpPtrVector().push_back(op_loop_side);
267 
269 
270  post_proc_fe_bdy->getOpPtrVector().push_back(
271 
272  new OpPPMap(
273 
274  post_proc_fe_bdy->getPostProcMesh(),
275  post_proc_fe_bdy->getMapGaussPts(),
276 
277  {},
278 
279  {{"U", u_ptr}},
280 
281  {{"GRAD", common_ptr->matGradPtr},
282  {"FIRST_PIOLA", common_ptr->getMatFirstPiolaStress()}},
283 
284  {}
285 
286  )
287 
288  );
289  return post_proc_fe_bdy;
290  };
291 
292  auto add_extra_finite_elements_to_ksp_solver_pipelines = [&]() {
294 
295  auto pre_proc_ptr = boost::make_shared<FEMethod>();
296  auto post_proc_rhs_ptr = boost::make_shared<FEMethod>();
297  auto post_proc_lhs_ptr = boost::make_shared<FEMethod>();
298 
299  auto time_scale = boost::make_shared<TimeScale>();
300 
301  auto get_bc_hook_rhs = [this, pre_proc_ptr, time_scale]() {
304  {time_scale}, false)();
306  };
307 
308  pre_proc_ptr->preProcessHook = get_bc_hook_rhs;
309 
310  auto get_post_proc_hook_rhs = [this, post_proc_rhs_ptr]() {
313  mField, post_proc_rhs_ptr, nullptr, Sev::verbose)();
315  mField, post_proc_rhs_ptr, 1.)();
316  CHKERR EssentialPostProcRhs<MPCsType>(mField, post_proc_rhs_ptr)();
318  };
319  auto get_post_proc_hook_lhs = [this, post_proc_lhs_ptr]() {
322  mField, post_proc_lhs_ptr, 1.)();
323  CHKERR EssentialPostProcLhs<MPCsType>(mField, post_proc_lhs_ptr)();
325  };
326  post_proc_rhs_ptr->postProcessHook = get_post_proc_hook_rhs;
327  post_proc_lhs_ptr->postProcessHook = get_post_proc_hook_lhs;
328 
329  // This is low level pushing finite elements (pipelines) to solver
330  auto ts_ctx_ptr = getDMTsCtx(simple->getDM());
331  ts_ctx_ptr->getPreProcessIFunction().push_front(pre_proc_ptr);
332  ts_ctx_ptr->getPreProcessIJacobian().push_front(pre_proc_ptr);
333  ts_ctx_ptr->getPostProcessIFunction().push_back(post_proc_rhs_ptr);
334  ts_ctx_ptr->getPostProcessIJacobian().push_back(post_proc_lhs_ptr);
336  };
337 
338  // Add extra finite elements to SNES solver pipelines to resolve essential
339  // boundary conditions
340  CHKERR add_extra_finite_elements_to_ksp_solver_pipelines();
341 
342  auto create_monitor_fe = [dm](auto &&post_proc_fe) {
343  return boost::make_shared<Monitor>(dm, post_proc_fe);
344  };
345 
346  // Set monitor which postprocessing results and saves them to the hard drive
347  boost::shared_ptr<FEMethod> null_fe;
348  auto monitor_ptr = create_monitor_fe(create_post_proc_fe());
349  CHKERR DMMoFEMTSSetMonitor(dm, ts, simple->getDomainFEName(), null_fe,
350  null_fe, monitor_ptr);
351 
352  // Set time solver
353  double ftime = 1;
354  CHKERR TSSetDuration(ts, PETSC_DEFAULT, ftime);
355  CHKERR TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);
356 
357  auto D = createDMVector(simple->getDM());
358 
359  CHKERR TSSetSolution(ts, D);
360  CHKERR TSSetFromOptions(ts);
361 
362  CHKERR TSSolve(ts, NULL);
363  CHKERR TSGetTime(ts, &ftime);
364 
365  PetscInt steps, snesfails, rejects, nonlinits, linits;
366  CHKERR TSGetStepNumber(ts, &steps);
367  CHKERR TSGetSNESFailures(ts, &snesfails);
368  CHKERR TSGetStepRejections(ts, &rejects);
369  CHKERR TSGetSNESIterations(ts, &nonlinits);
370  CHKERR TSGetKSPIterations(ts, &linits);
371  MOFEM_LOG_C("EXAMPLE", Sev::inform,
372  "steps %d (%d rejected, %d SNES fails), ftime %g, nonlinits "
373  "%d, linits %d",
374  steps, rejects, snesfails, ftime, nonlinits, linits);
375 
377 }
378 //! [Solve]
379 
380 //! [Getting norms]
383 
384  auto simple = mField.getInterface<Simple>();
385  auto dm = simple->getDM();
386 
387  auto T = createDMVector(simple->getDM());
388  CHKERR DMoFEMMeshToLocalVector(simple->getDM(), T, INSERT_VALUES,
389  SCATTER_FORWARD);
390  double nrm2;
391  CHKERR VecNorm(T, NORM_2, &nrm2);
392  MOFEM_LOG("EXAMPLE", Sev::inform) << "Solution norm " << nrm2;
393 
394  auto post_proc_norm_fe = boost::make_shared<DomainEle>(mField);
395 
396  auto post_proc_norm_rule_hook = [](int, int, int p) -> int { return 2 * p; };
397  post_proc_norm_fe->getRuleHook = post_proc_norm_rule_hook;
398 
400  post_proc_norm_fe->getOpPtrVector(), {H1});
401 
402  enum NORMS { U_NORM_L2 = 0, PIOLA_NORM, LAST_NORM };
403  auto norms_vec =
404  createVectorMPI(mField.get_comm(),
405  (mField.get_comm_rank() == 0) ? LAST_NORM : 0, LAST_NORM);
406  CHKERR VecZeroEntries(norms_vec);
407 
408  auto u_ptr = boost::make_shared<MatrixDouble>();
409  post_proc_norm_fe->getOpPtrVector().push_back(
411 
412  post_proc_norm_fe->getOpPtrVector().push_back(
413  new OpCalcNormL2Tensor1<SPACE_DIM>(u_ptr, norms_vec, U_NORM_L2));
414 
415  auto common_ptr = commonDataFactory<SPACE_DIM, GAUSS, DomainEleOp>(
416  mField, post_proc_norm_fe->getOpPtrVector(), "U", "MAT_ELASTIC",
417  Sev::inform);
418 
419  post_proc_norm_fe->getOpPtrVector().push_back(
421  common_ptr->getMatFirstPiolaStress(), norms_vec, PIOLA_NORM));
422 
423  CHKERR DMoFEMLoopFiniteElements(dm, simple->getDomainFEName(),
424  post_proc_norm_fe);
425 
426  CHKERR VecAssemblyBegin(norms_vec);
427  CHKERR VecAssemblyEnd(norms_vec);
428 
429  MOFEM_LOG_CHANNEL("SELF"); // Clear channel from old tags
430  if (mField.get_comm_rank() == 0) {
431  const double *norms;
432  CHKERR VecGetArrayRead(norms_vec, &norms);
433  MOFEM_TAG_AND_LOG("SELF", Sev::inform, "example")
434  << "norm_u: " << std::scientific << std::sqrt(norms[U_NORM_L2]);
435  MOFEM_TAG_AND_LOG("SELF", Sev::inform, "example")
436  << "norm_piola: " << std::scientific << std::sqrt(norms[PIOLA_NORM]);
437  CHKERR VecRestoreArrayRead(norms_vec, &norms);
438  }
439 
441 }
442 //! [Getting norms]
443 
444 //! [Postprocessing results]
447  PetscInt test_nb = 0;
448  PetscBool test_flg = PETSC_FALSE;
449  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-test", &test_nb, &test_flg);
450 
451  if (test_flg) {
452  auto simple = mField.getInterface<Simple>();
453  auto T = createDMVector(simple->getDM());
454  CHKERR DMoFEMMeshToLocalVector(simple->getDM(), T, INSERT_VALUES,
455  SCATTER_FORWARD);
456  double nrm2;
457  CHKERR VecNorm(T, NORM_2, &nrm2);
458  MOFEM_LOG("EXAMPLE", Sev::verbose) << "Regression norm " << nrm2;
459  double regression_value = 0;
460  switch (test_nb) {
461  case 1:
462  regression_value = 1.02789;
463  break;
464  case 2:
465  regression_value = 1.8841e+00;
466  break;
467  case 3:
468  regression_value = 1.8841e+00;
469  break;
470  case 4: // just links
471  regression_value = 4.9625e+00;
472  break;
473  case 5: // link master
474  regression_value = 6.6394e+00;
475  break;
476  case 6: // link master swap
477  regression_value = 4.98764e+00;
478  break;
479  case 7: // link Y
480  regression_value = 4.9473e+00;
481  break;
482  case 8: // link_3D_repr
483  regression_value = 2.5749e-01;
484  break;
485 
486  default:
487  SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Wrong test number.");
488  break;
489  }
490  if (fabs(nrm2 - regression_value) > 1e-2)
491  SETERRQ2(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
492  "Regression test field; wrong norm value. %6.4e != %6.4e", nrm2,
493  regression_value);
494  }
496 }
497 //! [Postprocessing results]
498 
499 //! [Check]
503 }
504 //! [Check]
505 
506 static char help[] = "...\n\n";
507 
508 int main(int argc, char *argv[]) {
509 
510  // Initialisation of MoFEM/PETSc and MOAB data structures
511  const char param_file[] = "param_file.petsc";
512  MoFEM::Core::Initialize(&argc, &argv, param_file, help);
513 
514  // Add logging channel for example
515  auto core_log = logging::core::get();
516  core_log->add_sink(
518  LogManager::setLog("EXAMPLE");
519  MOFEM_LOG_TAG("EXAMPLE", "example");
520 
521  try {
522 
523  //! [Register MoFEM discrete manager in PETSc]
524  DMType dm_name = "DMMOFEM";
525  CHKERR DMRegister_MoFEM(dm_name);
526  //! [Register MoFEM discrete manager in PETSc
527 
528  //! [Create MoAB]
529  moab::Core mb_instance; ///< mesh database
530  moab::Interface &moab = mb_instance; ///< mesh database interface
531  //! [Create MoAB]
532 
533  //! [Create MoFEM]
534  MoFEM::Core core(moab); ///< finite element database
535  MoFEM::Interface &m_field = core; ///< finite element database interface
536  //! [Create MoFEM]
537 
538  //! [Example]
539  Example ex(m_field);
540  CHKERR ex.runProblem();
541  //! [Example]
542  }
543  CATCH_ERRORS;
544 
546 }
MoFEM::NaturalBC::Assembly::LinearForm
Definition: Natural.hpp:67
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:128
Example::checkResults
MoFEMErrorCode checkResults()
[Postprocess results]
Definition: dynamic_first_order_con_law.cpp:1205
EXECUTABLE_DIMENSION
#define EXECUTABLE_DIMENSION
Definition: plastic.cpp:13
BoundaryEle
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
Definition: child_and_parent.cpp:39
MoFEM::EssentialPreProcReaction< DisplacementCubitBcData >
Specialization for DisplacementCubitBcData.
Definition: EssentialDisplacementCubitBcData.hpp:157
Example::assembleSystem
MoFEMErrorCode assembleSystem()
[Push operators to pipeline]
Definition: dynamic_first_order_con_law.cpp:647
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
LASTBASE
@ LASTBASE
Definition: definitions.h:69
Example::Example
Example(MoFEM::Interface &m_field)
Definition: nonlinear_elastic.cpp:57
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
MoFEM::FEMethod
structure for User Loop Methods on finite elements
Definition: LoopMethods.hpp:369
MoFEM::OpFluxRhsImpl
Definition: Natural.hpp:39
MoFEM::PipelineManager::ElementsAndOpsByDim
Definition: PipelineManager.hpp:38
MoFEM::EssentialPostProcLhs< DisplacementCubitBcData >
Specialization for DisplacementCubitBcData.
Definition: EssentialDisplacementCubitBcData.hpp:134
Monitor::Monitor
Monitor(SmartPetscObj< DM > dm, boost::shared_ptr< PostProcEleBdy > post_proc)
Definition: nonlinear_elastic.cpp:220
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
PostProcEleByDim< 2 >::SideEle
PipelineManager::ElementsAndOpsByDim< 2 >::FaceSideEle SideEle
Definition: dynamic_first_order_con_law.cpp:44
MoFEM::EssentialPostProcRhs< DisplacementCubitBcData >
Specialization for DisplacementCubitBcData.
Definition: EssentialDisplacementCubitBcData.hpp:113
MoFEM::PETSC
@ PETSC
Definition: FormsIntegrators.hpp:105
MoFEM::PipelineManager
PipelineManager interface.
Definition: PipelineManager.hpp:24
MoFEM.hpp
MoFEM::DMoFEMMeshToLocalVector
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:523
MoFEM::DisplacementCubitBcData
Definition of the displacement bc data structure.
Definition: BCData.hpp:76
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::EssentialPostProcRhs< MPCsType >
Specialization for DisplacementCubitBcData.
Definition: EssentialMPCsData.hpp:81
HenckyOps
Definition: HenckyOps.hpp:11
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
Monitor::postProcess
MoFEMErrorCode postProcess()
function is run at the end of loop
Definition: nonlinear_elastic.cpp:222
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
EntData
EntitiesFieldData::EntData EntData
Definition: nonlinear_elastic.cpp:17
Example::readMesh
MoFEMErrorCode readMesh()
[Run problem]
Definition: dynamic_first_order_con_law.cpp:463
PostProcEleByDim
Definition: dynamic_first_order_con_law.cpp:39
MoFEM::OpCalcNormL2Tensor2
Get norm of input MatrixDouble for Tensor2.
Definition: NormsOperators.hpp:72
save_every_nth_step
int save_every_nth_step
Definition: photon_diffusion.cpp:67
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
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::BcManager
Simple interface for fast problem set-up.
Definition: BcManager.hpp:25
BoundaryEleOp
MOFEM_LOG_C
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
simple
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69
Example
[Example]
Definition: plastic.cpp:177
main
int main(int argc, char *argv[])
Definition: nonlinear_elastic.cpp:508
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
MoFEM::PostProcBrokenMeshInMoabBaseContImpl
Definition: PostProcBrokenMeshInMoabBase.hpp:880
Example::solveSystem
MoFEMErrorCode solveSystem()
[Solve]
Definition: dynamic_first_order_con_law.cpp:893
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::EssentialPostProcLhs< MPCsType >
Specialization for MPCsType.
Definition: EssentialMPCsData.hpp:102
SPACE_DIM
constexpr int SPACE_DIM
Definition: nonlinear_elastic.cpp:14
MatrixFunction.hpp
MoFEM::OpCalcNormL2Tensor1
Get norm of input MatrixDouble for Tensor1.
Definition: NormsOperators.hpp:44
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
MoFEM::AddFluxToRhsPipelineImpl
Definition: Natural.hpp:46
PostProcEleBdy
PostProcEleByDim< SPACE_DIM >::PostProcEleBdy PostProcEleBdy
Definition: dynamic_first_order_con_law.cpp:55
MoFEM::PostProcBrokenMeshInMoabBase
Definition: PostProcBrokenMeshInMoabBase.hpp:94
MoFEM::AddHOOps
Add operators pushing bases from local to physical configuration.
Definition: HODataOperators.hpp:413
integration_rule
auto integration_rule
Definition: free_surface.cpp:185
DomainEleOp
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
MOFEM_TAG_AND_LOG
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
Definition: LogManager.hpp:362
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::NaturalBC::Assembly
Assembly methods.
Definition: Natural.hpp:65
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
Example::setupProblem
MoFEMErrorCode setupProblem()
[Run problem]
Definition: plastic.cpp:233
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
Example::boundaryCondition
MoFEMErrorCode boundaryCondition()
[Set up problem]
Definition: dynamic_first_order_con_law.cpp:536
Monitor
[Push operators to pipeline]
Definition: dynamic_first_order_con_law.cpp:783
approx_order
int approx_order
Definition: test_broken_space.cpp:50
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
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
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
help
static char help[]
[Check]
Definition: nonlinear_elastic.cpp:506
MoFEM::createVectorMPI
auto createVectorMPI(MPI_Comm comm, PetscInt n, PetscInt N)
Create MPI Vector.
Definition: PetscSmartObj.hpp:202
Example::gettingNorms
MoFEMErrorCode gettingNorms()
[Solve]
Definition: nonlinear_elastic.cpp:381
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
Example::runProblem
MoFEMErrorCode runProblem()
[Run problem]
Definition: plastic.cpp:214
Monitor::postProc
boost::shared_ptr< PostProcEleBdy > postProc
Definition: nonlinear_elastic.cpp:238
ReactionDiffusionEquation::D
const double D
diffusivity
Definition: reaction_diffusion.cpp:20
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
DomainEle
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
Definition: child_and_parent.cpp:34
MoFEM::LogManager::setLog
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
Definition: LogManager.cpp:389
MoFEM::EssentialPreProc
Class (Function) to enforce essential constrains.
Definition: Essential.hpp:25
MoFEM::SmartPetscObj< DM >
Example::outputResults
MoFEMErrorCode outputResults()
[Solve]
Definition: dynamic_first_order_con_law.cpp:1183
MoFEM::DMoFEMLoopFiniteElements
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:586
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
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEM::getDMTsCtx
auto getDMTsCtx(DM dm)
Get TS context data structure used by DM.
Definition: DMMoFEM.hpp:1141
PostProcEleByDim< 3 >::SideEle
PipelineManager::ElementsAndOpsByDim< 3 >::FaceSideEle SideEle
Definition: dynamic_first_order_con_law.cpp:50
MoFEM::OpLoopSide
Element used to execute operators on side of the element.
Definition: ForcesAndSourcesCore.hpp:1290
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::OpPostProcMapInMoab
Post post-proc data at points from hash maps.
Definition: PostProcBrokenMeshInMoabBase.hpp:698