v0.14.0
eigen_elastic.cpp
Go to the documentation of this file.
1 /**
2  * \file eigen_elastic.cpp
3  * \example eigen_elastic.cpp
4  *
5  * Calculate natural frequencies in 2d and 3d problems.
6  *
7  */
8 
9 #include <MoFEM.hpp>
10 #undef EPS
11 #include <slepceps.h>
12 
13 using namespace MoFEM;
14 
15 template <int DIM> struct ElementsAndOps {};
16 
17 template <> struct ElementsAndOps<2> {
19 };
20 
21 template <> struct ElementsAndOps<3> {
23 };
24 
25 constexpr int SPACE_DIM =
26  EXECUTABLE_DIMENSION; //< Space dimension of problem, mesh
27 
32 
34  GAUSS>::OpGradSymTensorGrad<1, SPACE_DIM, SPACE_DIM, 0>;
37 
38 double rho = 7829e-12;
39 double young_modulus = 207e3;
40 double poisson_ratio = 0.33;
41 
42 double bulk_modulus_K = young_modulus / (3 * (1 - 2 * poisson_ratio));
44 
45 int order = 1;
46 
47 struct Example {
48 
49  Example(MoFEM::Interface &m_field) : mField(m_field) {}
50 
51  MoFEMErrorCode runProblem();
52 
53 private:
54  MoFEM::Interface &mField;
55 
56  MoFEMErrorCode readMesh();
57  MoFEMErrorCode setupProblem();
58  MoFEMErrorCode createCommonData();
59  MoFEMErrorCode boundaryCondition();
60  MoFEMErrorCode assembleSystem();
61  MoFEMErrorCode solveSystem();
62  MoFEMErrorCode outputResults();
63  MoFEMErrorCode checkResults();
64 
65  boost::shared_ptr<MatrixDouble> matDPtr;
66 
70 
71  std::array<SmartPetscObj<Vec>, 6> rigidBodyMotion;
72 };
73 
74 //! [Create common data]
77 
78  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-rho", &rho, PETSC_NULL);
79  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-young_modulus", &young_modulus,
80  PETSC_NULL);
81  CHKERR PetscOptionsGetScalar(PETSC_NULL, "", "-poisson_ratio", &poisson_ratio,
82  PETSC_NULL);
83 
84  auto set_matrial_stiffens = [&]() {
91  auto t_D = getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, 0>(*matDPtr);
92 
93  const double A = (SPACE_DIM == 2)
94  ? 2 * shear_modulus_G /
95  (bulk_modulus_K + (4. / 3.) * shear_modulus_G)
96  : 1;
97  t_D(i, j, k, l) = 2 * shear_modulus_G * ((t_kd(i, k) ^ t_kd(j, l)) / 4.) +
98  A * (bulk_modulus_K - (2. / 3.) * shear_modulus_G) *
99  t_kd(i, j) * t_kd(k, l);
100 
102  };
103 
104  matDPtr = boost::make_shared<MatrixDouble>();
105 
106  constexpr auto size_symm = (SPACE_DIM * (SPACE_DIM + 1)) / 2;
107  matDPtr->resize(size_symm * size_symm, 1);
108 
109  CHKERR set_matrial_stiffens();
110 
112 }
113 //! [Create common data]
114 
115 //! [Run problem]
118  CHKERR readMesh();
119  CHKERR setupProblem();
120  CHKERR createCommonData();
121  CHKERR boundaryCondition();
122  CHKERR assembleSystem();
123  CHKERR solveSystem();
124  CHKERR outputResults();
125  CHKERR checkResults();
127 }
128 //! [Run problem]
129 
130 //! [Read mesh]
133  auto simple = mField.getInterface<Simple>();
134 
135  MOFEM_LOG("EXAMPLE", Sev::inform)
136  << "Read mesh for problem in " << EXECUTABLE_DIMENSION;
137  CHKERR simple->getOptions();
138  CHKERR simple->loadFile();
140 }
141 //! [Read mesh]
142 
143 //! [Set up problem]
145  auto *simple = mField.getInterface<Simple>();
147  // Add field
148  CHKERR simple->addDomainField("U", H1, AINSWORTH_BERNSTEIN_BEZIER_BASE,
149  SPACE_DIM);
150  CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-order", &order, PETSC_NULL);
151  CHKERR simple->setFieldOrder("U", order);
152  CHKERR simple->setUp();
154 }
155 //! [Set up problem]
156 
157 //! [Boundary condition]
159  auto *simple = mField.getInterface<Simple>();
161 
162  rigidBodyMotion[0] = createDMVector(simple->getDM());
163  for (int n = 1; n != 6; ++n)
164  rigidBodyMotion[n] = vectorDuplicate(rigidBodyMotion[0]);
165 
166  // Create space of vectors or rigid motion
167  auto problem_ptr = mField.get_problem(simple->getProblemName());
168  auto dofs = problem_ptr->getNumeredRowDofsPtr();
169 
170  // Get all vertices
171  const auto bit_number = mField.get_field_bit_number("U");
172  auto lo_uid =
173  DofEntity::getUniqueIdCalculate(0, get_id_for_min_type<MBVERTEX>());
174  auto hi_uid =
175  DofEntity::getUniqueIdCalculate(2, get_id_for_max_type<MBVERTEX>());
176 
177  auto hi = dofs->upper_bound(lo_uid);
178  std::array<double, 3> coords;
179 
180  for (auto lo = dofs->lower_bound(lo_uid); lo != hi; ++lo) {
181 
182  if ((*lo)->getPart() == mField.get_comm_rank()) {
183 
184  auto ent = (*lo)->getEnt();
185  CHKERR mField.get_moab().get_coords(&ent, 1, coords.data());
186 
187  if ((*lo)->getDofCoeffIdx() == 0) {
188  CHKERR VecSetValue(rigidBodyMotion[0], (*lo)->getPetscGlobalDofIdx(), 1,
189  INSERT_VALUES);
190  CHKERR VecSetValue(rigidBodyMotion[3], (*lo)->getPetscGlobalDofIdx(),
191  -coords[1], INSERT_VALUES);
192  if (SPACE_DIM == 3)
193  CHKERR VecSetValue(rigidBodyMotion[4], (*lo)->getPetscGlobalDofIdx(),
194  -coords[2], INSERT_VALUES);
195 
196  } else if ((*lo)->getDofCoeffIdx() == 1) {
197  CHKERR VecSetValue(rigidBodyMotion[1], (*lo)->getPetscGlobalDofIdx(), 1,
198  INSERT_VALUES);
199  CHKERR VecSetValue(rigidBodyMotion[3], (*lo)->getPetscGlobalDofIdx(),
200  coords[0], INSERT_VALUES);
201  if (SPACE_DIM == 3)
202  CHKERR VecSetValue(rigidBodyMotion[5], (*lo)->getPetscGlobalDofIdx(),
203  -coords[2], INSERT_VALUES);
204 
205  } else if ((*lo)->getDofCoeffIdx() == 2) {
206  if (SPACE_DIM == 3) {
207  CHKERR VecSetValue(rigidBodyMotion[2], (*lo)->getPetscGlobalDofIdx(),
208  1, INSERT_VALUES);
209  CHKERR VecSetValue(rigidBodyMotion[4], (*lo)->getPetscGlobalDofIdx(),
210  coords[0], INSERT_VALUES);
211  CHKERR VecSetValue(rigidBodyMotion[5], (*lo)->getPetscGlobalDofIdx(),
212  coords[1], INSERT_VALUES);
213  }
214  }
215  }
216  }
217 
218  for (int n = 0; n != rigidBodyMotion.size(); ++n) {
219  CHKERR VecAssemblyBegin(rigidBodyMotion[n]);
220  CHKERR VecAssemblyEnd(rigidBodyMotion[n]);
221  CHKERR VecGhostUpdateBegin(rigidBodyMotion[n], INSERT_VALUES,
222  SCATTER_FORWARD);
223  CHKERR VecGhostUpdateEnd(rigidBodyMotion[n], INSERT_VALUES,
224  SCATTER_FORWARD);
225  }
226 
228 }
229 //! [Boundary condition]
230 
231 //! [Push operators to pipeline]
234  auto *simple = mField.getInterface<Simple>();
235  auto *pipeline_mng = mField.getInterface<PipelineManager>();
236 
237  auto dm = simple->getDM();
239  M = matDuplicate(K, MAT_SHARE_NONZERO_PATTERN);
240 
241  auto calculate_stiffness_matrix = [&]() {
243  pipeline_mng->getDomainLhsFE().reset();
244 
245  auto det_ptr = boost::make_shared<VectorDouble>();
246  auto jac_ptr = boost::make_shared<MatrixDouble>();
247  auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
248  pipeline_mng->getOpDomainLhsPipeline().push_back(
249  new OpCalculateHOJac<SPACE_DIM>(jac_ptr));
250  pipeline_mng->getOpDomainLhsPipeline().push_back(
251  new OpInvertMatrix<SPACE_DIM>(jac_ptr, det_ptr, inv_jac_ptr));
252  pipeline_mng->getOpDomainLhsPipeline().push_back(
253  new OpSetHOInvJacToScalarBases<SPACE_DIM>(H1, inv_jac_ptr));
254  pipeline_mng->getOpDomainLhsPipeline().push_back(
255  new OpSetHOWeights(det_ptr));
256 
257  pipeline_mng->getOpDomainLhsPipeline().push_back(
258  new OpK("U", "U", matDPtr));
259  auto integration_rule = [](int, int, int approx_order) {
260  return 2 * (approx_order - 1);
261  };
262 
263  CHKERR pipeline_mng->setDomainLhsIntegrationRule(integration_rule);
264  pipeline_mng->getDomainLhsFE()->B = K;
265  CHKERR MatZeroEntries(K);
266  CHKERR pipeline_mng->loopFiniteElements();
267  CHKERR MatAssemblyBegin(K, MAT_FINAL_ASSEMBLY);
268  CHKERR MatAssemblyEnd(K, MAT_FINAL_ASSEMBLY);
270  };
271 
272  auto calculate_mass_matrix = [&]() {
274  pipeline_mng->getDomainLhsFE().reset();
275 
276  auto det_ptr = boost::make_shared<VectorDouble>();
277  auto jac_ptr = boost::make_shared<MatrixDouble>();
278  auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
279  pipeline_mng->getOpDomainLhsPipeline().push_back(
280  new OpCalculateHOJac<SPACE_DIM>(jac_ptr));
281  pipeline_mng->getOpDomainLhsPipeline().push_back(
282  new OpInvertMatrix<SPACE_DIM>(jac_ptr, det_ptr, inv_jac_ptr));
283  pipeline_mng->getOpDomainLhsPipeline().push_back(
284  new OpSetHOWeights(det_ptr));
285 
286  auto get_rho = [](const double, const double, const double) { return rho; };
287  pipeline_mng->getOpDomainLhsPipeline().push_back(
288  new OpMass("U", "U", get_rho));
289  auto integration_rule = [](int, int, int approx_order) {
290  return 2 * approx_order;
291  };
292  CHKERR pipeline_mng->setDomainLhsIntegrationRule(integration_rule);
293  CHKERR MatZeroEntries(M);
294  pipeline_mng->getDomainLhsFE()->B = M;
295  CHKERR pipeline_mng->loopFiniteElements();
296  CHKERR MatAssemblyBegin(M, MAT_FINAL_ASSEMBLY);
297  CHKERR MatAssemblyEnd(M, MAT_FINAL_ASSEMBLY);
299  };
300 
301  CHKERR calculate_stiffness_matrix();
302  CHKERR calculate_mass_matrix();
303 
305 }
306 //! [Push operators to pipeline]
307 
308 //! [Solve]
310  auto simple = mField.getInterface<Simple>();
311  auto is_manager = mField.getInterface<ISManager>();
313 
314  auto create_eps = [](MPI_Comm comm) {
315  EPS eps;
316  CHKERR EPSCreate(comm, &eps);
317  return SmartPetscObj<EPS>(eps);
318  };
319 
320  auto deflate_vectors = [&]() {
322  // Deflate vectors
323  std::array<Vec, 6> deflate_vectors;
324  for (int n = 0; n != 6; ++n) {
325  deflate_vectors[n] = rigidBodyMotion[n];
326  }
327  CHKERR EPSSetDeflationSpace(ePS, 6, &deflate_vectors[0]);
329  };
330 
331  auto print_info = [&]() {
333  ST st;
334  EPSType type;
335  PetscReal tol;
336  PetscInt nev, maxit, its;
337  // Optional: Get some information from the solver and display it
338  CHKERR EPSGetIterationNumber(ePS, &its);
339  MOFEM_LOG_C("EXAMPLE", Sev::inform,
340  " Number of iterations of the method: %d", its);
341  CHKERR EPSGetST(ePS, &st);
342  CHKERR EPSGetType(ePS, &type);
343  MOFEM_LOG_C("EXAMPLE", Sev::inform, " Solution method: %s", type);
344  CHKERR EPSGetDimensions(ePS, &nev, NULL, NULL);
345  MOFEM_LOG_C("EXAMPLE", Sev::inform, " Number of requested eigenvalues: %d",
346  nev);
347  CHKERR EPSGetTolerances(ePS, &tol, &maxit);
348  MOFEM_LOG_C("EXAMPLE", Sev::inform,
349  " Stopping condition: tol=%.4g, maxit=%d", (double)tol, maxit);
350 
351  PetscScalar eigr, eigi;
352  for (int nn = 0; nn < nev; nn++) {
353  CHKERR EPSGetEigenpair(ePS, nn, &eigr, &eigi, PETSC_NULL, PETSC_NULL);
354  MOFEM_LOG_C("EXAMPLE", Sev::inform,
355  " ncov = %d eigr = %.4g eigi = %.4g (inv eigr = %.4g)", nn,
356  eigr, eigi, 1. / eigr);
357  }
358 
360  };
361 
362  auto setup_eps = [&]() {
364  CHKERR EPSSetProblemType(ePS, EPS_GHEP);
365  CHKERR EPSSetWhichEigenpairs(ePS, EPS_SMALLEST_MAGNITUDE);
366  CHKERR EPSSetFromOptions(ePS);
368  };
369 
370  // Create eigensolver context
371  ePS = create_eps(mField.get_comm());
372  CHKERR EPSSetOperators(ePS, K, M);
373 
374  // Setup eps
375  CHKERR setup_eps();
376 
377  // Deflate vectors
378  CHKERR deflate_vectors();
379 
380  // Solve problem
381  CHKERR EPSSolve(ePS);
382 
383  // Print info
384  CHKERR print_info();
385 
387 }
388 //! [Solve]
389 
390 //! [Postprocess results]
393  auto *pipeline_mng = mField.getInterface<PipelineManager>();
394  auto *simple = mField.getInterface<Simple>();
395 
396  pipeline_mng->getDomainLhsFE().reset();
397  auto post_proc_fe = boost::make_shared<PostProcEle>(mField);
398 
399  auto det_ptr = boost::make_shared<VectorDouble>();
400  auto jac_ptr = boost::make_shared<MatrixDouble>();
401  auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
402  post_proc_fe->getOpPtrVector().push_back(
403  new OpCalculateHOJac<SPACE_DIM>(jac_ptr));
404  post_proc_fe->getOpPtrVector().push_back(
405  new OpInvertMatrix<SPACE_DIM>(jac_ptr, det_ptr, inv_jac_ptr));
406  post_proc_fe->getOpPtrVector().push_back(
407  new OpSetHOInvJacToScalarBases<SPACE_DIM>(H1, inv_jac_ptr));
408 
409  auto u_ptr = boost::make_shared<MatrixDouble>();
410  auto grad_ptr = boost::make_shared<MatrixDouble>();
411  auto strain_ptr = boost::make_shared<MatrixDouble>();
412  auto stress_ptr = boost::make_shared<MatrixDouble>();
413 
414  post_proc_fe->getOpPtrVector().push_back(
416  post_proc_fe->getOpPtrVector().push_back(
418  post_proc_fe->getOpPtrVector().push_back(
419  new OpSymmetrizeTensor<SPACE_DIM>("U", grad_ptr, strain_ptr));
420  post_proc_fe->getOpPtrVector().push_back(
422  "U", strain_ptr, stress_ptr, matDPtr));
423 
425 
426  post_proc_fe->getOpPtrVector().push_back(
427 
428  new OpPPMap(
429  post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
430 
432 
433  OpPPMap::DataMapMat{{"U", u_ptr}},
434 
436 
437  OpPPMap::DataMapMat{{"STRAIN", strain_ptr}, {"STRESS", stress_ptr}}
438 
439  )
440 
441  );
442 
443  pipeline_mng->getDomainRhsFE() = post_proc_fe;
444 
445  auto dm = simple->getDM();
446  auto D = createDMVector(dm);
447 
448  PetscInt nev;
449  CHKERR EPSGetDimensions(ePS, &nev, NULL, NULL);
450  PetscScalar eigr, eigi, nrm2r;
451  for (int nn = 0; nn < nev; nn++) {
452  CHKERR EPSGetEigenpair(ePS, nn, &eigr, &eigi, D, PETSC_NULL);
453  CHKERR VecGhostUpdateBegin(D, INSERT_VALUES, SCATTER_FORWARD);
454  CHKERR VecGhostUpdateEnd(D, INSERT_VALUES, SCATTER_FORWARD);
455  CHKERR VecNorm(D, NORM_2, &nrm2r);
456  MOFEM_LOG_C("EXAMPLE", Sev::inform,
457  " ncov = %d omega2 = %.8g omega = %.8g frequency = %.8g", nn,
458  eigr, std::sqrt(std::abs(eigr)),
459  std::sqrt(std::abs(eigr)) / (2 * M_PI));
460  CHKERR DMoFEMMeshToLocalVector(dm, D, INSERT_VALUES, SCATTER_REVERSE);
461  CHKERR pipeline_mng->loopFiniteElements();
462  post_proc_fe->writeFile("out_eig_" + boost::lexical_cast<std::string>(nn) +
463  ".h5m");
464  }
465 
467 }
468 //! [Postprocess results]
469 
470 //! [Check]
473  PetscBool test_flg = PETSC_FALSE;
474  CHKERR PetscOptionsGetBool(PETSC_NULL, "", "-test", &test_flg, PETSC_NULL);
475  if (test_flg) {
476  PetscScalar eigr, eigi;
477  CHKERR EPSGetEigenpair(ePS, 0, &eigr, &eigi, PETSC_NULL, PETSC_NULL);
478  constexpr double regression_value = 12579658;
479  if (fabs(eigr - regression_value) > 1)
480  SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
481  "Regression test faileed; wrong eigen value.");
482  }
484 }
485 //! [Check]
486 
487 static char help[] = "...\n\n";
488 
489 int main(int argc, char *argv[]) {
490 
491  // Initialisation of MoFEM/PETSc and MOAB data structures
492  const char param_file[] = "param_file.petsc";
493  SlepcInitialize(&argc, &argv, param_file, help);
494  MoFEM::Core::Initialize(&argc, &argv, param_file, help);
495 
496  // Add logging channel for example
497  auto core_log = logging::core::get();
498  core_log->add_sink(
499  LogManager::createSink(LogManager::getStrmWorld(), "EXAMPLE"));
500  LogManager::setLog("EXAMPLE");
501  MOFEM_LOG_TAG("EXAMPLE", "example");
502 
503  try {
504 
505  //! [Register MoFEM discrete manager in PETSc]
506  DMType dm_name = "DMMOFEM";
507  CHKERR DMRegister_MoFEM(dm_name);
508  //! [Register MoFEM discrete manager in PETSc
509 
510  //! [Create MoAB]
511  moab::Core mb_instance; ///< mesh database
512  moab::Interface &moab = mb_instance; ///< mesh database interface
513  //! [Create MoAB]
514 
515  //! [Create MoFEM]
516  MoFEM::Core core(moab); ///< finite element database
517  MoFEM::Interface &m_field = core; ///< finite element database insterface
518  //! [Create MoFEM]
519 
520  //! [Example]
521  Example ex(m_field);
522  CHKERR ex.runProblem();
523  //! [Example]
524  }
525  CATCH_ERRORS;
526 
527  SlepcFinalize();
529 }
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::K
VectorDouble K
Definition: Projection10NodeCoordsOnField.cpp:125
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
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
OpK
Definition: simple_elasticity.cpp:16
Example::Example
Example(MoFEM::Interface &m_field)
Definition: eigen_elastic.cpp:49
MoFEM::OpPostProcMapInMoab::DataMapVec
std::map< std::string, boost::shared_ptr< VectorDouble > > DataMapVec
Definition: PostProcBrokenMeshInMoabBase.hpp:700
shear_modulus_G
double shear_modulus_G
Definition: eigen_elastic.cpp:43
MoFEM::DofEntity::getUniqueIdCalculate
static UId getUniqueIdCalculate(const DofIdx dof, UId ent_uid)
Definition: DofsMultiIndices.hpp:54
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
order
int order
Definition: eigen_elastic.cpp:45
MoFEM::PETSC
@ PETSC
Definition: FormsIntegrators.hpp:105
MoFEM::PipelineManager
PipelineManager interface.
Definition: PipelineManager.hpp:24
MoFEM.hpp
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
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::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
Example::M
SmartPetscObj< Mat > M
Definition: eigen_elastic.cpp:67
MoFEM::Simple
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::OpCalculateHOJac
Definition: HODataOperators.hpp:267
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
Example::readMesh
MoFEMErrorCode readMesh()
[Run problem]
Definition: dynamic_first_order_con_law.cpp:463
MoFEM::PostProcBrokenMeshInMoab
Definition: PostProcBrokenMeshInMoabBase.hpp:667
PlasticOps::M
FTensor::Index< 'M', 3 > M
Definition: PlasticOps.hpp:117
MoFEM::DMCreateMatrix_MoFEM
PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M)
Definition: DMMoFEM.cpp:1197
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
bulk_modulus_K
double bulk_modulus_K
Definition: eigen_elastic.cpp:42
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
OpMass
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpMass
[Only used with Hooke equation (linear material model)]
Definition: seepage.cpp:57
help
static char help[]
[Check]
Definition: eigen_elastic.cpp:487
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
MoFEM::OpSetHOInvJacToScalarBases
Set inverse jacobian to base functions.
Definition: HODataOperators.hpp:73
Example
[Example]
Definition: plastic.cpp:177
double
convert.type
type
Definition: convert.py:64
MoFEM::FormsIntegrators::Assembly
Assembly methods.
Definition: FormsIntegrators.hpp:317
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
Example::K
SmartPetscObj< Mat > K
Definition: eigen_elastic.cpp:68
Example::solveSystem
MoFEMErrorCode solveSystem()
[Solve]
Definition: dynamic_first_order_con_law.cpp:893
MoFEM::DMRegister_MoFEM
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition: DMMoFEM.cpp:43
MoFEM::GAUSS
@ GAUSS
Definition: FormsIntegrators.hpp:136
rho
double rho
Definition: eigen_elastic.cpp:38
MoFEM::FaceElementForcesAndSourcesCore
Face finite element.
Definition: FaceElementForcesAndSourcesCore.hpp:23
EntData
EntitiesFieldData::EntData EntData
Definition: eigen_elastic.cpp:28
main
int main(int argc, char *argv[])
Definition: eigen_elastic.cpp:489
OpMass
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpMass
Definition: eigen_elastic.cpp:36
Example::ePS
SmartPetscObj< EPS > ePS
Definition: eigen_elastic.cpp:69
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
size_symm
constexpr auto size_symm
Definition: plastic.cpp:42
MoFEM::matDuplicate
SmartPetscObj< Mat > matDuplicate(Mat mat, MatDuplicateOption op)
Definition: PetscSmartObj.hpp:234
MoFEM::VolumeElementForcesAndSourcesCore
Volume finite element base.
Definition: VolumeElementForcesAndSourcesCore.hpp:26
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
MoFEM::OpPostProcMapInMoab::DataMapMat
std::map< std::string, boost::shared_ptr< MatrixDouble > > DataMapMat
Definition: PostProcBrokenMeshInMoabBase.hpp:701
MoFEM::OpTensorTimesSymmetricTensor
Calculate .
Definition: UserDataOperators.hpp:1883
t_kd
constexpr auto t_kd
Definition: free_surface.cpp:137
BiLinearForm
FTensor::Index< 'i', SPACE_DIM >
AINSWORTH_BERNSTEIN_BEZIER_BASE
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
Definition: definitions.h:64
MoFEM::OpCalculateVectorFieldGradient
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Definition: UserDataOperators.hpp:1535
convert.n
n
Definition: convert.py:82
integration_rule
auto integration_rule
Definition: free_surface.cpp:185
ElementsAndOps
Definition: child_and_parent.cpp:18
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_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::vectorDuplicate
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
Definition: PetscSmartObj.hpp:221
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
Example::matDPtr
boost::shared_ptr< MatrixDouble > matDPtr
Definition: eigen_elastic.cpp:65
MoFEM::OpSymmetrizeTensor
Definition: UserDataOperators.hpp:1948
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
OpK
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradSymTensorGrad< 1, SPACE_DIM, SPACE_DIM, 0 > OpK
Definition: eigen_elastic.cpp:34
Example::setupProblem
MoFEMErrorCode setupProblem()
[Run problem]
Definition: plastic.cpp:225
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
Example::boundaryCondition
MoFEMErrorCode boundaryCondition()
[Set up problem]
Definition: dynamic_first_order_con_law.cpp:536
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
eps
static const double eps
Definition: check_base_functions_derivatives_on_tet.cpp:11
approx_order
int approx_order
Definition: test_broken_space.cpp:50
MoFEM::OpSetHOWeights
Set inverse jacobian to base functions.
Definition: HODataOperators.hpp:159
MoFEM::OpInvertMatrix
Definition: UserDataOperators.hpp:3249
Example::runProblem
MoFEMErrorCode runProblem()
[Run problem]
Definition: plastic.cpp:213
Example::rigidBodyMotion
std::array< SmartPetscObj< Vec >, 6 > rigidBodyMotion
Definition: eigen_elastic.cpp:71
ReactionDiffusionEquation::D
const double D
diffusivity
Definition: reaction_diffusion.cpp:20
young_modulus
double young_modulus
Definition: eigen_elastic.cpp:39
FTensor::Kronecker_Delta_symmetric
Kronecker Delta class symmetric.
Definition: Kronecker_Delta.hpp:49
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
MoFEM::PetscOptionsGetScalar
PetscErrorCode PetscOptionsGetScalar(PetscOptions *, const char pre[], const char name[], PetscScalar *dval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:162
Example::createCommonData
MoFEMErrorCode createCommonData()
[Set up problem]
Definition: plastic.cpp:396
DomainEle
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
Definition: child_and_parent.cpp:34
MoFEM::SmartPetscObj< Mat >
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
Example::outputResults
MoFEMErrorCode outputResults()
[Solve]
Definition: dynamic_first_order_con_law.cpp:1183
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
poisson_ratio
double poisson_ratio
Definition: eigen_elastic.cpp:40
SPACE_DIM
constexpr int SPACE_DIM
Definition: eigen_elastic.cpp:25
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
l
FTensor::Index< 'l', 3 > l
Definition: matrix_function.cpp:21
tol
double tol
Definition: mesh_smoothing.cpp:26
MoFEM::PetscOptionsGetBool
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:182
MoFEM::OpPostProcMapInMoab
Post post-proc data at points from hash maps.
Definition: PostProcBrokenMeshInMoabBase.hpp:698