v0.14.0
nonlinear_elastic.cpp

Plane stress elastic dynamic problem

/** \file nonlinear_elastic.cpp
\brief 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.
*/
using namespace MoFEM;
namespace bio = boost::iostreams;
using bio::stream;
using bio::tee_device;
static char help[] = "...\n\n";
int main(int argc, char *argv[]) {
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
try {
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
int rank;
MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
PetscBool flg = PETSC_TRUE;
char mesh_file_name[255];
CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
mesh_file_name, 255, &flg);
if (flg != PETSC_TRUE) {
SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
}
const char *option;
option = ""; //"PARALLEL=BCAST;";//;DEBUG_IO";
CHKERR moab.load_file(mesh_file_name, 0, option);
MoFEM::Core core(moab);
MoFEM::Interface &m_field = core;
// ref meshset ref level 0
BitRefLevel bit_level0;
bit_level0.set(0);
EntityHandle meshset_level0;
CHKERR moab.create_meshset(MESHSET_SET, meshset_level0);
CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
0, 3, bit_level0);
CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByRefLevel(
bit_level0, BitRefLevel().set(), meshset_level0);
// Fields
CHKERR m_field.add_field("SPATIAL_POSITION", H1, AINSWORTH_LEGENDRE_BASE,
3);
// add entitities (by tets) to the field
CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "SPATIAL_POSITION");
// set app. order
PetscInt order;
CHKERR PetscOptionsGetInt(PETSC_NULL, PETSC_NULL, "-my_order", &order,
&flg);
if (flg != PETSC_TRUE) {
order = 1;
}
CHKERR m_field.set_field_order(0, MBTET, "SPATIAL_POSITION", order);
CHKERR m_field.set_field_order(0, MBTRI, "SPATIAL_POSITION", order);
CHKERR m_field.set_field_order(0, MBEDGE, "SPATIAL_POSITION", order);
CHKERR m_field.set_field_order(0, MBVERTEX, "SPATIAL_POSITION", 1);
NonlinearElasticElement elastic(m_field, 1);
boost::shared_ptr<
double_kirchhoff_material_ptr(
double>());
boost::shared_ptr<
adouble_kirchhoff_material_ptr(
adouble>());
CHKERR elastic.setBlocks(double_kirchhoff_material_ptr,
adouble_kirchhoff_material_ptr);
CHKERR elastic.addElement("ELASTIC", "SPATIAL_POSITION");
CHKERR elastic.setOperators("SPATIAL_POSITION");
// define problems
CHKERR m_field.add_problem("ELASTIC_MECHANICS");
// set refinement level for problem
CHKERR m_field.modify_problem_ref_level_add_bit("ELASTIC_MECHANICS",
bit_level0);
// set finite elements for problems
CHKERR m_field.modify_problem_add_finite_element("ELASTIC_MECHANICS",
"ELASTIC");
// build field
CHKERR m_field.build_fields();
// use this to apply some strain field to the body (testing only)
double scale_positions = 2;
{
EntityHandle node = 0;
double coords[3];
for (_IT_GET_DOFS_FIELD_BY_NAME_FOR_LOOP_(m_field, "SPATIAL_POSITION",
dof_ptr)) {
if (dof_ptr->get()->getEntType() != MBVERTEX)
continue;
EntityHandle ent = dof_ptr->get()->getEnt();
int dof_rank = dof_ptr->get()->getDofCoeffIdx();
double &fval = dof_ptr->get()->getFieldData();
if (node != ent) {
CHKERR moab.get_coords(&ent, 1, coords);
node = ent;
}
fval = scale_positions * coords[dof_rank];
}
}
// build finite elemnts
// build adjacencies
CHKERR m_field.build_adjacencies(bit_level0);
ProblemsManager *prb_mng_ptr;
CHKERR m_field.getInterface(prb_mng_ptr);
// build problem
CHKERR prb_mng_ptr->buildProblem("ELASTIC_MECHANICS", true);
// partition
CHKERR prb_mng_ptr->partitionProblem("ELASTIC_MECHANICS");
CHKERR prb_mng_ptr->partitionFiniteElements("ELASTIC_MECHANICS");
CHKERR prb_mng_ptr->partitionGhostDofs("ELASTIC_MECHANICS");
// create matrices
Vec F;
CHKERR m_field.getInterface<VecManager>()->vecCreateGhost(
"ELASTIC_MECHANICS", COL, &F);
Mat Aij;
->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>("ELASTIC_MECHANICS",
&Aij);
elastic.getLoopFeRhs().snes_f = F;
elastic.getLoopFeLhs().snes_B = Aij;
CHKERR m_field.loop_finite_elements("ELASTIC_MECHANICS", "ELASTIC",
elastic.getLoopFeRhs());
CHKERR VecGhostUpdateBegin(F, ADD_VALUES, SCATTER_REVERSE);
CHKERR VecGhostUpdateEnd(F, ADD_VALUES, SCATTER_REVERSE);
CHKERR VecAssemblyBegin(F);
CHKERR VecAssemblyEnd(F);
CHKERR m_field.loop_finite_elements("ELASTIC_MECHANICS", "ELASTIC",
elastic.getLoopFeLhs());
CHKERR MatAssemblyBegin(Aij, MAT_FINAL_ASSEMBLY);
CHKERR MatAssemblyEnd(Aij, MAT_FINAL_ASSEMBLY);
double sum = 0;
CHKERR VecSum(F, &sum);
CHKERR PetscPrintf(PETSC_COMM_WORLD, "sum = %4.3e\n", sum);
double fnorm;
CHKERR VecNorm(F, NORM_2, &fnorm);
CHKERR PetscPrintf(PETSC_COMM_WORLD, "fnorm = %9.8e\n", fnorm);
double mnorm;
CHKERR MatNorm(Aij, NORM_1, &mnorm);
CHKERR PetscPrintf(PETSC_COMM_WORLD, "mnorm = %9.8e\n", mnorm);
if (fabs(sum) > 1e-8) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Failed to pass test");
}
if (fabs(fnorm - 5.12196914e+00) > 1e-6) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Failed to pass test");
}
if (fabs(mnorm - 5.48280139e+01) > 1e-6) {
SETERRQ(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID, "Failed to pass test");
}
CHKERR VecDestroy(&F);
CHKERR MatDestroy(&Aij);
}
return 0;
}
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
BasicFiniteElements.hpp
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
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
main
int main(int argc, char *argv[])
Definition: nonlinear_elastic.cpp:508
COL
@ COL
Definition: definitions.h:136
NonlinearElasticElement::addElement
MoFEMErrorCode addElement(const std::string element_name, const std::string spatial_position_field_name, const std::string material_position_field_name="MESH_NODE_POSITIONS", const bool ale=false)
Definition: NonLinearElasticElement.cpp:1120
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
NonlinearElasticElement::setOperators
MoFEMErrorCode setOperators(const std::string spatial_position_field_name, const std::string material_position_field_name="MESH_NODE_POSITIONS", const bool ale=false, const bool field_disp=false)
Set operators to calculate left hand tangent matrix and right hand residual.
Definition: NonLinearElasticElement.cpp:1153
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
help
static char help[]
[Check]
Definition: nonlinear_elastic.cpp:506
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.
NonlinearElasticElement::setBlocks
MoFEMErrorCode setBlocks(boost::shared_ptr< FunctionsToCalculatePiolaKirchhoffI< double >> materialDoublePtr, boost::shared_ptr< FunctionsToCalculatePiolaKirchhoffI< adouble >> materialAdoublePtr)
Definition: NonLinearElasticElement.cpp:1086
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