v0.14.0
Classes | Macros | Typedefs | Functions | Variables
HookeElement.hpp File Reference
#include <BasicFiniteElements.hpp>

Go to the source code of this file.

Classes

struct  NonlinearElasticElement
 structure grouping operators and data used for calculation of nonlinear elastic element More...
 
struct  NonlinearElasticElement::BlockData
 data for calculation heat conductivity and heat capacity elements More...
 
struct  ConvectiveMassElement
 structure grouping operators and data used for calculation of mass (convective) element More...
 
struct  ConvectiveMassElement::BlockData
 data for calculation inertia forces More...
 
struct  DataAtIntegrationPts
 
struct  OpCalculateStrain< D >
 
struct  OpCalculateStrainAle
 
struct  OpCalculateStress< S >
 
struct  OpCalculateEnergy
 
struct  OpCalculateEshelbyStress
 
struct  OpCalculateHomogeneousStiffness< S >
 
struct  OpCalculateMassMatrix
 Assemble mass matrix for elastic element TODO: CHANGE FORMULA *. More...
 
struct  OpCalculateStiffnessScaledByDensityField
 
struct  OpAssemble
 
struct  OpRhs_dx
 
struct  OpLhs_dx_dx< S >
 
struct  OpAleRhs_dx
 
struct  OpAleLhs_dx_dx< S >
 
struct  OpAleLhs_dx_dX< S >
 
struct  OpAleLhsWithDensity_dx_dX
 
struct  OpAleLhsWithDensity_dX_dX
 
struct  OpAleRhs_dX
 
struct  OpAleLhs_dX_dX< S >
 
struct  OpAleLhsPre_dX_dx< S >
 
struct  OpAleLhs_dX_dx
 
struct  OpAnalyticalInternalStrain_dx< S >
 
struct  OpAnalyticalInternalAleStrain_dX< S >
 
struct  OpAnalyticalInternalAleStrain_dx< S >
 
struct  OpPostProcHookeElement< ELEMENT >
 

Macros

#define MAT_TO_DDG(SM)
 

Typedefs

using MassBlockData = ConvectiveMassElement::BlockData
 
using EntData = EntitiesFieldData::EntData
 
using UserDataOperator = ForcesAndSourcesCore::UserDataOperator
 
using VolUserDataOperator = VolumeElementForcesAndSourcesCore::UserDataOperator
 

Functions

static MoFEMErrorCode setBlocks (MoFEM::Interface &m_field, boost::shared_ptr< map< int, BlockData >> &block_sets_ptr)
 
static MoFEMErrorCode addElasticElement (MoFEM::Interface &m_field, boost::shared_ptr< map< int, BlockData >> &block_sets_ptr, const std::string element_name, const std::string x_field, const std::string X_field, const bool ale)
 
static MoFEMErrorCode setOperators (boost::shared_ptr< ForcesAndSourcesCore > fe_lhs_ptr, boost::shared_ptr< ForcesAndSourcesCore > fe_rhs_ptr, boost::shared_ptr< map< int, BlockData >> block_sets_ptr, const std::string x_field, const std::string X_field, const bool ale, const bool field_disp, const EntityType type=MBTET, boost::shared_ptr< DataAtIntegrationPts > data_at_pts=nullptr)
 
static MoFEMErrorCode calculateEnergy (DM dm, boost::shared_ptr< map< int, BlockData >> block_sets_ptr, const std::string x_field, const std::string X_field, const bool ale, const bool field_disp, SmartPetscObj< Vec > &v_energy_ptr)
 

Variables

struct ConvectiveMassElement BlockData = NonlinearElasticElement::BlockData
 
MatrixDouble invJac
 

Macro Definition Documentation

◆ MAT_TO_DDG

#define MAT_TO_DDG (   SM)
Value:
&(*SM)(0, 0), &(*SM)(1, 0), &(*SM)(2, 0), &(*SM)(3, 0), &(*SM)(4, 0), \
&(*SM)(5, 0), &(*SM)(6, 0), &(*SM)(7, 0), &(*SM)(8, 0), &(*SM)(9, 0), \
&(*SM)(10, 0), &(*SM)(11, 0), &(*SM)(12, 0), &(*SM)(13, 0), \
&(*SM)(14, 0), &(*SM)(15, 0), &(*SM)(16, 0), &(*SM)(17, 0), \
&(*SM)(18, 0), &(*SM)(19, 0), &(*SM)(20, 0), &(*SM)(21, 0), \
&(*SM)(22, 0), &(*SM)(23, 0), &(*SM)(24, 0), &(*SM)(25, 0), \
&(*SM)(26, 0), &(*SM)(27, 0), &(*SM)(28, 0), &(*SM)(29, 0), \
&(*SM)(30, 0), &(*SM)(31, 0), &(*SM)(32, 0), &(*SM)(33, 0), \
&(*SM)(34, 0), &(*SM)(35, 0)
Examples
HookeElement.cpp, HookeElement.hpp, and HookeInternalStressElement.hpp.

Definition at line 142 of file HookeElement.hpp.

Typedef Documentation

◆ EntData

Examples
HookeElement.hpp.

Definition at line 74 of file HookeElement.hpp.

◆ MassBlockData

Definition at line 72 of file HookeElement.hpp.

◆ UserDataOperator

◆ VolUserDataOperator

Definition at line 77 of file HookeElement.hpp.

Function Documentation

◆ addElasticElement()

MoFEMErrorCode HookeElement::addElasticElement ( MoFEM::Interface m_field,
boost::shared_ptr< map< int, BlockData >> &  block_sets_ptr,
const std::string  element_name,
const std::string  x_field,
const std::string  X_field,
const bool  ale 
)
static
Examples
elasticity.cpp, HookeElement.cpp, HookeElement.hpp, mortar_contact_thermal.cpp, and simple_contact_thermal.cpp.

Definition at line 533 of file HookeElement.cpp.

537  {
539 
540  if (!block_sets_ptr)
541  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
542  "Pointer to block of sets is null");
543 
544  CHKERR m_field.add_finite_element(element_name, MF_ZERO);
545  CHKERR m_field.modify_finite_element_add_field_row(element_name, x_field);
546  CHKERR m_field.modify_finite_element_add_field_col(element_name, x_field);
547  CHKERR m_field.modify_finite_element_add_field_data(element_name, x_field);
548  if (m_field.check_field(X_field)) {
549  if (ale) {
550  CHKERR m_field.modify_finite_element_add_field_row(element_name, X_field);
551  CHKERR m_field.modify_finite_element_add_field_col(element_name, X_field);
552  }
553  CHKERR m_field.modify_finite_element_add_field_data(element_name, X_field);
554  }
555 
556  for (auto &m : (*block_sets_ptr)) {
557  CHKERR m_field.add_ents_to_finite_element_by_dim(m.second.tEts, 3,
558  element_name);
559  }
560 
562 }

◆ calculateEnergy()

static MoFEMErrorCode calculateEnergy ( DM  dm,
boost::shared_ptr< map< int, BlockData >>  block_sets_ptr,
const std::string  x_field,
const std::string  X_field,
const bool  ale,
const bool  field_disp,
SmartPetscObj< Vec > &  v_energy_ptr 
)
static

◆ setBlocks()

static MoFEMErrorCode setBlocks ( MoFEM::Interface m_field,
boost::shared_ptr< map< int, BlockData >> &  block_sets_ptr 
)
static

◆ setOperators()

static MoFEMErrorCode setOperators ( boost::shared_ptr< ForcesAndSourcesCore fe_lhs_ptr,
boost::shared_ptr< ForcesAndSourcesCore fe_rhs_ptr,
boost::shared_ptr< map< int, BlockData >>  block_sets_ptr,
const std::string  x_field,
const std::string  X_field,
const bool  ale,
const bool  field_disp,
const EntityType  type = MBTET,
boost::shared_ptr< DataAtIntegrationPts data_at_pts = nullptr 
)
static

Variable Documentation

◆ invJac

MatrixDouble invJac
private
Examples
HookeElement.hpp.

Definition at line 683 of file HookeElement.hpp.

MoFEM::CoreInterface::modify_finite_element_add_field_row
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::CoreInterface::add_finite_element
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
MoFEM::CoreInterface::modify_finite_element_add_field_col
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
MoFEM::CoreInterface::check_field
virtual bool check_field(const std::string &name) const =0
check if field is in database
MoFEM::CoreInterface::add_ents_to_finite_element_by_dim
virtual MoFEMErrorCode add_ents_to_finite_element_by_dim(const EntityHandle entities, const int dim, const std::string &name, const bool recursive=true)=0
add entities to finite element
MF_ZERO
@ MF_ZERO
Definition: definitions.h:98
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::CoreInterface::modify_finite_element_add_field_data
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346