v0.14.0
Functions | Variables
nonlinear_elastic.cpp File Reference

Atom test for linear elastic dynamics. More...

#include <BasicFiniteElements.hpp>

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Variables

static char help [] = "...\n\n"
 

Detailed Description

Atom test for linear elastic dynamics.

This is not exactly procedure for linear elastic dynamics, since jacobian is evaluated at every time step and snes procedure is involved. However it is implemented like that, to test methodology for general nonlinear problem.

Definition in file nonlinear_elastic.cpp.

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 22 of file nonlinear_elastic.cpp.

22  {
23 
24  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
25 
26  try {
27 
28  moab::Core mb_instance;
29  moab::Interface &moab = mb_instance;
30  int rank;
31  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
32 
33  PetscBool flg = PETSC_TRUE;
34  char mesh_file_name[255];
35  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
36  mesh_file_name, 255, &flg);
37  if (flg != PETSC_TRUE) {
38  SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
39  }
40 
41  const char *option;
42  option = ""; //"PARALLEL=BCAST;";//;DEBUG_IO";
43  CHKERR moab.load_file(mesh_file_name, 0, option);
44 
45  MoFEM::Core core(moab);
46  MoFEM::Interface &m_field = core;
47 
48  // ref meshset ref level 0
49  BitRefLevel bit_level0;
50  bit_level0.set(0);
51  EntityHandle meshset_level0;
52  CHKERR moab.create_meshset(MESHSET_SET, meshset_level0);
53  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
54  0, 3, bit_level0);
55  CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByRefLevel(
56  bit_level0, BitRefLevel().set(), meshset_level0);
57 
58  // Fields
59  CHKERR m_field.add_field("SPATIAL_POSITION", H1, AINSWORTH_LEGENDRE_BASE,
60  3);
61  // add entitities (by tets) to the field
62  CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "SPATIAL_POSITION");
63  // set app. order
64  PetscInt order;
65  CHKERR PetscOptionsGetInt(PETSC_NULL, PETSC_NULL, "-my_order", &order,
66  &flg);
67  if (flg != PETSC_TRUE) {
68  order = 1;
69  }
70  CHKERR m_field.set_field_order(0, MBTET, "SPATIAL_POSITION", order);
71  CHKERR m_field.set_field_order(0, MBTRI, "SPATIAL_POSITION", order);
72  CHKERR m_field.set_field_order(0, MBEDGE, "SPATIAL_POSITION", order);
73  CHKERR m_field.set_field_order(0, MBVERTEX, "SPATIAL_POSITION", 1);
74 
75  NonlinearElasticElement elastic(m_field, 1);
76  boost::shared_ptr<
78  double_kirchhoff_material_ptr(
80  double>());
81  boost::shared_ptr<
83  adouble_kirchhoff_material_ptr(
85  adouble>());
86  CHKERR elastic.setBlocks(double_kirchhoff_material_ptr,
87  adouble_kirchhoff_material_ptr);
88  CHKERR elastic.addElement("ELASTIC", "SPATIAL_POSITION");
89  CHKERR elastic.setOperators("SPATIAL_POSITION");
90 
91  // define problems
92  CHKERR m_field.add_problem("ELASTIC_MECHANICS");
93  // set refinement level for problem
94  CHKERR m_field.modify_problem_ref_level_add_bit("ELASTIC_MECHANICS",
95  bit_level0);
96  // set finite elements for problems
97  CHKERR m_field.modify_problem_add_finite_element("ELASTIC_MECHANICS",
98  "ELASTIC");
99 
100  // build field
101  CHKERR m_field.build_fields();
102 
103  // use this to apply some strain field to the body (testing only)
104  double scale_positions = 2;
105  {
106  EntityHandle node = 0;
107  double coords[3];
108  for (_IT_GET_DOFS_FIELD_BY_NAME_FOR_LOOP_(m_field, "SPATIAL_POSITION",
109  dof_ptr)) {
110  if (dof_ptr->get()->getEntType() != MBVERTEX)
111  continue;
112  EntityHandle ent = dof_ptr->get()->getEnt();
113  int dof_rank = dof_ptr->get()->getDofCoeffIdx();
114  double &fval = dof_ptr->get()->getFieldData();
115  if (node != ent) {
116  CHKERR moab.get_coords(&ent, 1, coords);
117  node = ent;
118  }
119  fval = scale_positions * coords[dof_rank];
120  }
121  }
122 
123  // build finite elemnts
124  CHKERR m_field.build_finite_elements();
125  // build adjacencies
126  CHKERR m_field.build_adjacencies(bit_level0);
127 
128  ProblemsManager *prb_mng_ptr;
129  CHKERR m_field.getInterface(prb_mng_ptr);
130  // build problem
131  CHKERR prb_mng_ptr->buildProblem("ELASTIC_MECHANICS", true);
132  // partition
133  CHKERR prb_mng_ptr->partitionProblem("ELASTIC_MECHANICS");
134  CHKERR prb_mng_ptr->partitionFiniteElements("ELASTIC_MECHANICS");
135  CHKERR prb_mng_ptr->partitionGhostDofs("ELASTIC_MECHANICS");
136 
137  // create matrices
138  Vec F;
139  CHKERR m_field.getInterface<VecManager>()->vecCreateGhost(
140  "ELASTIC_MECHANICS", COL, &F);
141  Mat Aij;
143  ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>("ELASTIC_MECHANICS",
144  &Aij);
145 
146  elastic.getLoopFeRhs().snes_f = F;
147  elastic.getLoopFeLhs().snes_B = Aij;
148 
149  CHKERR m_field.loop_finite_elements("ELASTIC_MECHANICS", "ELASTIC",
150  elastic.getLoopFeRhs());
151  CHKERR VecGhostUpdateBegin(F, ADD_VALUES, SCATTER_REVERSE);
152  CHKERR VecGhostUpdateEnd(F, ADD_VALUES, SCATTER_REVERSE);
153  CHKERR VecAssemblyBegin(F);
154  CHKERR VecAssemblyEnd(F);
155 
156  CHKERR m_field.loop_finite_elements("ELASTIC_MECHANICS", "ELASTIC",
157  elastic.getLoopFeLhs());
158  CHKERR MatAssemblyBegin(Aij, MAT_FINAL_ASSEMBLY);
159  CHKERR MatAssemblyEnd(Aij, MAT_FINAL_ASSEMBLY);
160 
161  double sum = 0;
162  CHKERR VecSum(F, &sum);
163  CHKERR PetscPrintf(PETSC_COMM_WORLD, "sum = %4.3e\n", sum);
164  double fnorm;
165  CHKERR VecNorm(F, NORM_2, &fnorm);
166  CHKERR PetscPrintf(PETSC_COMM_WORLD, "fnorm = %9.8e\n", fnorm);
167 
168  double mnorm;
169  CHKERR MatNorm(Aij, NORM_1, &mnorm);
170  CHKERR PetscPrintf(PETSC_COMM_WORLD, "mnorm = %9.8e\n", mnorm);
171 
172  if (fabs(sum) > 1e-8) {
173  SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Failed to pass test");
174  }
175  if (fabs(fnorm - 5.12196914e+00) > 1e-6) {
176  SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Failed to pass test");
177  }
178  if (fabs(mnorm - 5.48280139e+01) > 1e-6) {
179  SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Failed to pass test");
180  }
181 
182  CHKERR VecDestroy(&F);
183  CHKERR MatDestroy(&Aij);
184  }
185  CATCH_ERRORS;
186 
188 
189  return 0;
190 }

Variable Documentation

◆ help

char help[] = "...\n\n"
static

Definition at line 20 of file nonlinear_elastic.cpp.

MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::CoreInterface::loop_finite_elements
virtual MoFEMErrorCode loop_finite_elements(const std::string problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
MoFEM::ProblemsManager::buildProblem
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
Definition: ProblemsManager.cpp:87
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
H1
@ H1
continuous field
Definition: definitions.h:85
EntityHandle
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
NonlinearElasticElement::FunctionsToCalculatePiolaKirchhoffI< double >
NonlinearElasticElement
structure grouping operators and data used for calculation of nonlinear elastic element
Definition: HookeElement.hpp:27
help
static char help[]
Definition: nonlinear_elastic.cpp:20
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::CoreInterface::add_ents_to_field_by_type
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
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
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
COL
@ COL
Definition: definitions.h:136
MoFEM::ProblemsManager::partitionFiniteElements
MoFEMErrorCode partitionFiniteElements(const std::string name, bool part_from_moab=false, int low_proc=-1, int hi_proc=-1, int verb=VERBOSE)
partition finite elements
Definition: ProblemsManager.cpp:2167
MoFEM::CoreInterface::modify_problem_ref_level_add_bit
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
MoFEM::VecManager
Vector manager is used to create vectors \mofem_vectors.
Definition: VecManager.hpp:23
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:22
MoFEM::MatrixManager
Matrix manager is used to build and partition problems.
Definition: MatrixManager.hpp:21
adouble
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
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
_IT_GET_DOFS_FIELD_BY_NAME_FOR_LOOP_
#define _IT_GET_DOFS_FIELD_BY_NAME_FOR_LOOP_(MFIELD, NAME, IT)
Definition: Interface.hpp:1878
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::CoreInterface::build_fields
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::CoreInterface::modify_problem_add_finite_element
virtual MoFEMErrorCode modify_problem_add_finite_element(const std::string name_problem, const std::string &fe_name)=0
add finite element to problem, this add entities assigned to finite element to a particular problem
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
MoFEM::CoreInterface::build_adjacencies
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
MoFEM::CoreInterface::set_field_order
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
MoFEM::ProblemsManager::partitionGhostDofs
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
Definition: ProblemsManager.cpp:2339
MoFEM::CoreInterface::add_problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
MoFEM::CoreInterface::add_field
virtual MoFEMErrorCode add_field(const std::string name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
MoFEM::ProblemsManager::partitionProblem
MoFEMErrorCode partitionProblem(const std::string name, int verb=VERBOSE)
partition problem dofs (collective)
Definition: ProblemsManager.cpp:1683
F
@ F
Definition: free_surface.cpp:394