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 
50 
51 constexpr double young_modulus = 100;
52 constexpr double poisson_ratio = 0.3;
53 
54 #include <HenckyOps.hpp>
55 using namespace HenckyOps;
56 
57 struct Example {
58 
59  Example(MoFEM::Interface &m_field) : mField(m_field) {}
60 
61  MoFEMErrorCode runProblem();
62 
63 private:
64  MoFEM::Interface &mField;
65 
66  MoFEMErrorCode readMesh();
67  MoFEMErrorCode setupProblem();
68  MoFEMErrorCode boundaryCondition();
69  MoFEMErrorCode assembleSystem();
70  MoFEMErrorCode solveSystem();
71  MoFEMErrorCode gettingNorms();
72  MoFEMErrorCode outputResults();
73  MoFEMErrorCode checkResults();
74 };
75 
76 //! [Run problem]
79  CHKERR readMesh();
80  CHKERR setupProblem();
81  CHKERR boundaryCondition();
82  CHKERR assembleSystem();
83  CHKERR solveSystem();
84  CHKERR gettingNorms();
85  CHKERR outputResults();
86  CHKERR checkResults();
88 }
89 //! [Run problem]
90 
91 //! [Read mesh]
94  auto simple = mField.getInterface<Simple>();
95  CHKERR simple->getOptions();
96  char meshFileName[255];
97  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-file_name",
98  meshFileName, 255, PETSC_NULL);
99  CHKERR simple->loadFile("", meshFileName,
102 }
103 //! [Read mesh]
104 
105 //! [Set up problem]
108  Simple *simple = mField.getInterface<Simple>();
109 
110  enum bases { AINSWORTH, DEMKOWICZ, LASBASETOPT };
111  const char *list_bases[LASBASETOPT] = {"ainsworth", "demkowicz"};
112  PetscInt choice_base_value = AINSWORTH;
113  CHKERR PetscOptionsGetEList(PETSC_NULL, NULL, "-base", list_bases,
114  LASBASETOPT, &choice_base_value, PETSC_NULL);
115 
117  switch (choice_base_value) {
118  case AINSWORTH:
120  MOFEM_LOG("WORLD", Sev::inform)
121  << "Set AINSWORTH_LEGENDRE_BASE for displacements";
122  break;
123  case DEMKOWICZ:
124  base = DEMKOWICZ_JACOBI_BASE;
125  MOFEM_LOG("WORLD", Sev::inform)
126  << "Set DEMKOWICZ_JACOBI_BASE for displacements";
127  break;
128  default:
129  base = LASTBASE;
130  break;
131  }
132 
133  // Add field
134  CHKERR simple->addDomainField("U", H1, base, SPACE_DIM);
135  CHKERR simple->addBoundaryField("U", H1, base, SPACE_DIM);
136  int order = 2;
137  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
138  CHKERR simple->setFieldOrder("U", order);
139  CHKERR simple->setUp();
141 }
142 //! [Set up problem]
143 
144 //! [Boundary condition]
147  auto *pipeline_mng = mField.getInterface<PipelineManager>();
148  auto simple = mField.getInterface<Simple>();
149  auto bc_mng = mField.getInterface<BcManager>();
150  auto time_scale = boost::make_shared<TimeScale>();
151 
152  auto integration_rule = [](int, int, int approx_order) {
153  return 2 * (approx_order - 1);
154  };
155 
156  CHKERR pipeline_mng->setDomainRhsIntegrationRule(integration_rule);
157  CHKERR pipeline_mng->setDomainLhsIntegrationRule(integration_rule);
158  CHKERR pipeline_mng->setBoundaryRhsIntegrationRule(integration_rule);
159 
161  pipeline_mng->getOpBoundaryRhsPipeline(), mField, "U", {time_scale},
162  "FORCE", Sev::inform);
163 
164  //! [Define gravity vector]
166  pipeline_mng->getOpDomainRhsPipeline(), mField, "U", {time_scale},
167  "BODY_FORCE", Sev::inform);
168 
169  // Essential BC
170  CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_X",
171  "U", 0, 0);
172  CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_Y",
173  "U", 1, 1);
174  CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_Z",
175  "U", 2, 2);
176  CHKERR bc_mng->pushMarkDOFsOnEntities<DisplacementCubitBcData>(
177  simple->getProblemName(), "U");
178 
179  // Essential MPCs BC
180  CHKERR bc_mng->addBlockDOFsToMPCs(simple->getProblemName(), "U");
181 
183 }
184 //! [Boundary condition]
185 
186 //! [Push operators to pipeline]
189  auto *simple = mField.getInterface<Simple>();
190  auto *pipeline_mng = mField.getInterface<PipelineManager>();
191 
192  auto add_domain_ops_lhs = [&](auto &pip) {
195  CHKERR HenckyOps::opFactoryDomainLhs<SPACE_DIM, PETSC, GAUSS, DomainEleOp>(
196  mField, pip, "U", "MAT_ELASTIC", Sev::inform);
198  };
199 
200  auto add_domain_ops_rhs = [&](auto &pip) {
203  CHKERR HenckyOps::opFactoryDomainRhs<SPACE_DIM, PETSC, GAUSS, DomainEleOp>(
204  mField, pip, "U", "MAT_ELASTIC", Sev::inform);
206  };
207 
208  CHKERR add_domain_ops_lhs(pipeline_mng->getOpDomainLhsPipeline());
209  CHKERR add_domain_ops_rhs(pipeline_mng->getOpDomainRhsPipeline());
210 
212 }
213 //! [Push operators to pipeline]
214 
215 /**
216  * @brief Monitor solution
217  *
218  * This functions is called by TS solver at the end of each step. It is used
219  * to output results to the hard drive.
220  */
221 struct Monitor : public FEMethod {
222  Monitor(SmartPetscObj<DM> dm, boost::shared_ptr<PostProcEleBdy> post_proc)
223  : dM(dm), postProc(post_proc){};
226  constexpr int save_every_nth_step = 1;
227  if (ts_step % save_every_nth_step == 0) {
228  CHKERR DMoFEMLoopFiniteElements(dM, "bFE", postProc);
230  postProc->mField, postProc->getPostProcMesh(), {"U"});
231 
232  CHKERR postProc->writeFile(
233  "out_step_" + boost::lexical_cast<std::string>(ts_step) + ".h5m");
234  }
236  }
237 
238 private:
240  boost::shared_ptr<PostProcEleBdy> postProc;
241 };
242 
243 //! [Solve]
246  auto *simple = mField.getInterface<Simple>();
247  auto *pipeline_mng = mField.getInterface<PipelineManager>();
248 
249  auto dm = simple->getDM();
250  auto ts = pipeline_mng->createTSIM();
251 
252  // Setup postprocessing
253  auto create_post_proc_fe = [dm, this, simple]() {
254  auto post_proc_ele_domain = [dm, this](auto &pip_domain) {
256  auto common_ptr = commonDataFactory<SPACE_DIM, GAUSS, DomainEleOp>(
257  mField, pip_domain, "U", "MAT_ELASTIC", Sev::inform);
258  return common_ptr;
259  };
260 
261  auto post_proc_fe_bdy = boost::make_shared<PostProcEleBdy>(mField);
262  auto u_ptr = boost::make_shared<MatrixDouble>();
263  post_proc_fe_bdy->getOpPtrVector().push_back(
265  auto op_loop_side =
266  new OpLoopSide<SideEle>(mField, simple->getDomainFEName(), SPACE_DIM);
267  auto common_ptr = post_proc_ele_domain(op_loop_side->getOpPtrVector());
268  post_proc_fe_bdy->getOpPtrVector().push_back(op_loop_side);
269 
271 
272  post_proc_fe_bdy->getOpPtrVector().push_back(
273 
274  new OpPPMap(
275 
276  post_proc_fe_bdy->getPostProcMesh(),
277  post_proc_fe_bdy->getMapGaussPts(),
278 
279  {},
280 
281  {{"U", u_ptr}},
282 
283  {{"GRAD", common_ptr->matGradPtr},
284  {"FIRST_PIOLA", common_ptr->getMatFirstPiolaStress()}},
285 
286  {}
287 
288  )
289 
290  );
291  return post_proc_fe_bdy;
292  };
293 
294  auto add_extra_finite_elements_to_ksp_solver_pipelines = [&]() {
296 
297  auto pre_proc_ptr = boost::make_shared<FEMethod>();
298  auto post_proc_rhs_ptr = boost::make_shared<FEMethod>();
299  auto post_proc_lhs_ptr = boost::make_shared<FEMethod>();
300 
301  auto time_scale = boost::make_shared<TimeScale>();
302 
303  auto get_bc_hook_rhs = [this, pre_proc_ptr, time_scale]() {
306  {time_scale}, false)();
308  };
309 
310  pre_proc_ptr->preProcessHook = get_bc_hook_rhs;
311 
312  auto get_post_proc_hook_rhs = [this, post_proc_rhs_ptr]() {
315  mField, post_proc_rhs_ptr, nullptr, Sev::verbose)();
317  mField, post_proc_rhs_ptr, 1.)();
318  CHKERR EssentialPostProcRhs<MPCsType>(mField, post_proc_rhs_ptr)();
320  };
321  auto get_post_proc_hook_lhs = [this, post_proc_lhs_ptr]() {
324  mField, post_proc_lhs_ptr, 1.)();
325  CHKERR EssentialPostProcLhs<MPCsType>(mField, post_proc_lhs_ptr)();
327  };
328  post_proc_rhs_ptr->postProcessHook = get_post_proc_hook_rhs;
329  post_proc_lhs_ptr->postProcessHook = get_post_proc_hook_lhs;
330 
331  // This is low level pushing finite elements (pipelines) to solver
332  auto ts_ctx_ptr = getDMTsCtx(simple->getDM());
333  ts_ctx_ptr->getPreProcessIFunction().push_front(pre_proc_ptr);
334  ts_ctx_ptr->getPreProcessIJacobian().push_front(pre_proc_ptr);
335  ts_ctx_ptr->getPostProcessIFunction().push_back(post_proc_rhs_ptr);
336  ts_ctx_ptr->getPostProcessIJacobian().push_back(post_proc_lhs_ptr);
338  };
339 
340  // Add extra finite elements to SNES solver pipelines to resolve essential
341  // boundary conditions
342  CHKERR add_extra_finite_elements_to_ksp_solver_pipelines();
343 
344  auto create_monitor_fe = [dm](auto &&post_proc_fe) {
345  return boost::make_shared<Monitor>(dm, post_proc_fe);
346  };
347 
348  // Set monitor which postprocessing results and saves them to the hard drive
349  boost::shared_ptr<FEMethod> null_fe;
350  auto monitor_ptr = create_monitor_fe(create_post_proc_fe());
351  CHKERR DMMoFEMTSSetMonitor(dm, ts, simple->getDomainFEName(), null_fe,
352  null_fe, monitor_ptr);
353 
354  // Set time solver
355  double ftime = 1;
356  CHKERR TSSetDuration(ts, PETSC_DEFAULT, ftime);
357  CHKERR TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);
358 
359  auto D = createDMVector(simple->getDM());
360 
361  CHKERR TSSetSolution(ts, D);
362  CHKERR TSSetFromOptions(ts);
363 
364  CHKERR TSSolve(ts, NULL);
365  CHKERR TSGetTime(ts, &ftime);
366 
367  PetscInt steps, snesfails, rejects, nonlinits, linits;
368  CHKERR TSGetStepNumber(ts, &steps);
369  CHKERR TSGetSNESFailures(ts, &snesfails);
370  CHKERR TSGetStepRejections(ts, &rejects);
371  CHKERR TSGetSNESIterations(ts, &nonlinits);
372  CHKERR TSGetKSPIterations(ts, &linits);
373  MOFEM_LOG_C("EXAMPLE", Sev::inform,
374  "steps %d (%d rejected, %d SNES fails), ftime %g, nonlinits "
375  "%d, linits %d",
376  steps, rejects, snesfails, ftime, nonlinits, linits);
377 
379 }
380 //! [Solve]
381 
382 //! [Getting norms]
385 
386  auto simple = mField.getInterface<Simple>();
387  auto dm = simple->getDM();
388 
389  auto T = createDMVector(simple->getDM());
390  CHKERR DMoFEMMeshToLocalVector(simple->getDM(), T, INSERT_VALUES,
391  SCATTER_FORWARD);
392  double nrm2;
393  CHKERR VecNorm(T, NORM_2, &nrm2);
394  MOFEM_LOG("EXAMPLE", Sev::inform) << "Solution norm " << nrm2;
395 
396  auto post_proc_norm_fe = boost::make_shared<DomainEle>(mField);
397 
398  auto post_proc_norm_rule_hook = [](int, int, int p) -> int { return 2 * p; };
399  post_proc_norm_fe->getRuleHook = post_proc_norm_rule_hook;
400 
402  post_proc_norm_fe->getOpPtrVector(), {H1});
403 
404  enum NORMS { U_NORM_L2 = 0, PIOLA_NORM, LAST_NORM };
405  auto norms_vec =
406  createVectorMPI(mField.get_comm(),
407  (mField.get_comm_rank() == 0) ? LAST_NORM : 0, LAST_NORM);
408  CHKERR VecZeroEntries(norms_vec);
409 
410  auto u_ptr = boost::make_shared<MatrixDouble>();
411  post_proc_norm_fe->getOpPtrVector().push_back(
413 
414  post_proc_norm_fe->getOpPtrVector().push_back(
415  new OpCalcNormL2Tensor1<SPACE_DIM>(u_ptr, norms_vec, U_NORM_L2));
416 
417  auto common_ptr = commonDataFactory<SPACE_DIM, GAUSS, DomainEleOp>(
418  mField, post_proc_norm_fe->getOpPtrVector(), "U", "MAT_ELASTIC",
419  Sev::inform);
420 
421  post_proc_norm_fe->getOpPtrVector().push_back(
423  common_ptr->getMatFirstPiolaStress(), norms_vec, PIOLA_NORM));
424 
425  CHKERR DMoFEMLoopFiniteElements(dm, simple->getDomainFEName(),
426  post_proc_norm_fe);
427 
428  CHKERR VecAssemblyBegin(norms_vec);
429  CHKERR VecAssemblyEnd(norms_vec);
430 
431  MOFEM_LOG_CHANNEL("SELF"); // Clear channel from old tags
432  if (mField.get_comm_rank() == 0) {
433  const double *norms;
434  CHKERR VecGetArrayRead(norms_vec, &norms);
435  MOFEM_TAG_AND_LOG("SELF", Sev::inform, "example")
436  << "norm_u: " << std::scientific << std::sqrt(norms[U_NORM_L2]);
437  MOFEM_TAG_AND_LOG("SELF", Sev::inform, "example")
438  << "norm_piola: " << std::scientific << std::sqrt(norms[PIOLA_NORM]);
439  CHKERR VecRestoreArrayRead(norms_vec, &norms);
440  }
441 
443 }
444 //! [Getting norms]
445 
446 //! [Postprocessing results]
449  PetscInt test_nb = 0;
450  PetscBool test_flg = PETSC_FALSE;
451  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-test", &test_nb, &test_flg);
452 
453  if (test_flg) {
454  auto simple = mField.getInterface<Simple>();
455  auto T = createDMVector(simple->getDM());
456  CHKERR DMoFEMMeshToLocalVector(simple->getDM(), T, INSERT_VALUES,
457  SCATTER_FORWARD);
458  double nrm2;
459  CHKERR VecNorm(T, NORM_2, &nrm2);
460  MOFEM_LOG("EXAMPLE", Sev::verbose) << "Regression norm " << nrm2;
461  double regression_value = 0;
462  switch (test_nb) {
463  case 1:
464  regression_value = 1.02789;
465  break;
466  case 2:
467  regression_value = 1.8841e+00;
468  break;
469  case 3:
470  regression_value = 1.8841e+00;
471  break;
472  case 4: // just links
473  regression_value = 4.9625e+00;
474  break;
475  case 5: // link master
476  regression_value = 6.6394e+00;
477  break;
478  case 6: // link master swap
479  regression_value = 4.98764e+00;
480  break;
481  case 7: // link Y
482  regression_value = 4.9473e+00;
483  break;
484  case 8: // link_3D_repr
485  regression_value = 2.5749e-01;
486  break;
487 
488  default:
489  SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Wrong test number.");
490  break;
491  }
492  if (fabs(nrm2 - regression_value) > 1e-2)
493  SETERRQ2(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
494  "Regression test field; wrong norm value. %6.4e != %6.4e", nrm2,
495  regression_value);
496  }
498 }
499 //! [Postprocessing results]
500 
501 //! [Check]
505 }
506 //! [Check]
507 
508 static char help[] = "...\n\n";
509 
510 int main(int argc, char *argv[]) {
511 
512  // Initialisation of MoFEM/PETSc and MOAB data structures
513  const char param_file[] = "param_file.petsc";
514  MoFEM::Core::Initialize(&argc, &argv, param_file, help);
515 
516  // Add logging channel for example
517  auto core_log = logging::core::get();
518  core_log->add_sink(
520  LogManager::setLog("EXAMPLE");
521  MOFEM_LOG_TAG("EXAMPLE", "example");
522 
523  try {
524 
525  //! [Register MoFEM discrete manager in PETSc]
526  DMType dm_name = "DMMOFEM";
527  CHKERR DMRegister_MoFEM(dm_name);
528  //! [Register MoFEM discrete manager in PETSc
529 
530  //! [Create MoAB]
531  moab::Core mb_instance; ///< mesh database
532  moab::Interface &moab = mb_instance; ///< mesh database interface
533  //! [Create MoAB]
534 
535  //! [Create MoFEM]
536  MoFEM::Core core(moab); ///< finite element database
537  MoFEM::Interface &m_field = core; ///< finite element database interface
538  //! [Create MoFEM]
539 
540  //! [Example]
541  Example ex(m_field);
542  CHKERR ex.runProblem();
543  //! [Example]
544  }
545  CATCH_ERRORS;
546 
548 }
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:447
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce 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:127
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:151
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:59
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:130
Monitor::Monitor
Monitor(SmartPetscObj< DM > dm, boost::shared_ptr< PostProcEleBdy > post_proc)
Definition: nonlinear_elastic.cpp:222
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: adolc_plasticity.cpp:88
main
int main(int argc, char *argv[])
Definition: nonlinear_elastic.cpp:22
MoFEM::EssentialPostProcRhs< DisplacementCubitBcData >
Specialization for DisplacementCubitBcData.
Definition: EssentialDisplacementCubitBcData.hpp:111
MoFEM::PETSC
@ PETSC
Definition: FormsIntegrators.hpp:104
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:527
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
SideEle
PostProcEleByDim< SPACE_DIM >::SideEle SideEle
Definition: adolc_plasticity.cpp:98
MoFEM::EssentialPostProcRhs< MPCsType >
Specialization for DisplacementCubitBcData.
Definition: EssentialMPCsData.hpp:80
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:224
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:1975
EntData
EntitiesFieldData::EntData EntData
Definition: nonlinear_elastic.cpp:17
approx_order
static constexpr int approx_order
Definition: prism_polynomial_approximation.cpp:14
PostProcEleBdy
PostProcEleByDim< SPACE_DIM >::PostProcEleBdy PostProcEleBdy
Definition: adolc_plasticity.cpp:99
Example::readMesh
MoFEMErrorCode readMesh()
[Run problem]
Definition: dynamic_first_order_con_law.cpp:463
PostProcEleByDim
Definition: adolc_plasticity.cpp:82
MoFEM::OpCalcNormL2Tensor2
Get norm of input MatrixDouble for Tensor2.
Definition: NormsOperators.hpp:69
save_every_nth_step
int save_every_nth_step
Definition: photon_diffusion.cpp:67
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::createDMVector
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1018
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
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator
default operator for TRI element
Definition: FaceElementForcesAndSourcesCore.hpp:94
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:228
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
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:47
MoFEM::EssentialPostProcLhs< MPCsType >
Specialization for MPCsType.
Definition: EssentialMPCsData.hpp:99
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:42
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
MoFEM::AddFluxToRhsPipelineImpl
Definition: Natural.hpp:46
MoFEM::PostProcBrokenMeshInMoabBase
Definition: PostProcBrokenMeshInMoabBase.hpp:94
MoFEM::AddHOOps
Add operators pushing bases from local to physical configuration.
Definition: HODataOperators.hpp:503
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:372
DEMKOWICZ_JACOBI_BASE
@ DEMKOWICZ_JACOBI_BASE
Definition: definitions.h:66
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
Example::setupProblem
MoFEMErrorCode setupProblem()
[Run problem]
Definition: plastic.cpp:276
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: adolc_plasticity.cpp:333
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:1060
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
young_modulus
constexpr double young_modulus
Definition: nonlinear_elastic.cpp:51
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:508
MoFEM::createVectorMPI
auto createVectorMPI(MPI_Comm comm, PetscInt n, PetscInt N)
Create MPI Vector.
Definition: PetscSmartObj.hpp:198
Example::gettingNorms
MoFEMErrorCode gettingNorms()
[Solve]
Definition: nonlinear_elastic.cpp:383
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
Example::runProblem
MoFEMErrorCode runProblem()
[Run problem]
Definition: plastic.cpp:264
Monitor::postProc
boost::shared_ptr< PostProcEleBdy > postProc
Definition: nonlinear_elastic.cpp:240
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:440
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:590
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:416
MoFEM::getDMTsCtx
auto getDMTsCtx(DM dm)
Get TS context data structure used by DM.
Definition: DMMoFEM.hpp:1060
poisson_ratio
constexpr double poisson_ratio
Definition: nonlinear_elastic.cpp:52
PostProcEleByDim< 3 >::SideEle
PipelineManager::ElementsAndOpsByDim< 3 >::FaceSideEle SideEle
Definition: adolc_plasticity.cpp:94
MoFEM::OpLoopSide
Element used to execute operators on side of the element.
Definition: ForcesAndSourcesCore.hpp:1289
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
MoFEM::OpPostProcMapInMoab
Post post-proc data at points from hash maps.
Definition: PostProcBrokenMeshInMoabBase.hpp:698