MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
MAST::GeomElem Class Reference

Detailed Description

This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface for cases where:

  • two-dimensional elements may exist in this three-dimensional space, and
  • level-set intersection may create sub-elements inside the element on either side of the interface. A physics implementation should be averse to these variants. This class and it children shoudl honor the following definitions:
    • a reference element is the element in the finite element mesh. In a standard FE analysis a quadrature is performed on this element.
    • a quadrature element is the element that is a subset of the reference element on which the quadrature points are defined. In a standard FE analysis the reference and quadrature elements are the same.
    • a local reference element is the element transformed in a coordinate system that is defined in along the 1D element or in the plane of the 2D element.
    • a local quadrature element is similarly defined for the quadrature element.

Definition at line 59 of file geom_elem.h.

#include <geom_elem.h>

Inheritance diagram for MAST::GeomElem:
Collaboration diagram for MAST::GeomElem:

Public Member Functions

 GeomElem ()
 
virtual ~GeomElem ()
 
unsigned int dim () const
 
const RealVectorXdomain_surface_normal () const
 
virtual void external_side_loads_for_quadrature_elem (std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase *> &bc, std::map< unsigned int, std::vector< MAST::BoundaryConditionBase *>> &loads) const
 From the given list of boundary loads, this identifies the sides of the quadrature element and the loads in bc that are to be applied on it. More...
 
virtual void get_boundary_ids_on_quadrature_elem_side (unsigned int s, std::vector< libMesh::boundary_id_type > &bc_ids) const
 
libMesh::FEType get_fe_type (unsigned int i) const
 
virtual const libMesh::Elem & get_quadrature_elem () const
 
virtual const libMesh::Elem & get_quadrature_local_elem () const
 
virtual const libMesh::Elem & get_reference_elem () const
 
virtual const libMesh::Elem & get_reference_local_elem () const
 
virtual void init (const libMesh::Elem &elem, const MAST::SystemInitialization &sys_init)
 initialize the object for the specified reference elem. More...
 
virtual std::unique_ptr< MAST::FEBaseinit_fe (bool init_grads, bool init_second_order_derivative, int extra_quadrature_order=0) const
 initializes the finite element shape function and quadrature object with the order of quadrature rule changed based on the extra_quadrature_order. More...
 
virtual std::unique_ptr< MAST::FEBaseinit_side_fe (unsigned int s, bool init_grads, bool init_second_order_derivative, int extra_quadrature_order=0) const
 initializes the finite element shape function and quadrature object for the side with the order of quadrature rule changed based on the extra_quadrature_order More...
 
unsigned int n_sides_quadrature_elem () const
 number of sides on quadrature element. More...
 
void set_bending (bool onoff)
 This sets the 1D elements to extension/torsional stiffness only. More...
 
void set_local_y_vector (const RealVectorX &y_vec)
 for 1D elements the transformed coordinate system attached to the element defines the local x-axis along the length of the element. More...
 
const RealMatrixXT_matrix () const
 O. More...
 
void transform_point_to_global_coordinate (const libMesh::Point &local_pt, libMesh::Point &global_pt) const
 
void transform_vector_to_global_coordinate (const libMesh::Point &local_vec, libMesh::Point &global_vec) const
 
void transform_vector_to_global_coordinate (const RealVectorX &local_vec, RealVectorX &global_vec) const
 
void transform_vector_to_local_coordinate (const libMesh::Point &global_vec, libMesh::Point &local_vec) const
 
void transform_vector_to_local_coordinate (const RealVectorX &global_vec, RealVectorX &local_vec) const
 
bool use_local_elem () const
 Vector and matrix quantities defined on one- and two-dimensional elements that are oriented in two or three-dimensional spaces may need to be transformed to/from element coordinate system. More...
 

Protected Member Functions

void _init_local_elem ()
 initializes the local element More...
 
void _init_local_elem_1d ()
 initializes the local element More...
 
void _init_local_elem_2d ()
 initializes the local element More...
 

Protected Attributes

bool _bending = true
 Defines if bending is used in this element or not. More...
 
RealVectorX _domain_surface_normal
 surface normal of the 1D/2D element. More...
 
libMesh::Elem * _local_elem
 a local element is created if More...
 
std::vector< libMesh::Node * > _local_nodes
 nodes for local element More...
 
RealVectorX _local_y
 the y-axis that is used to define the coordinate system for a 1-D element. More...
 
const libMesh::Elem * _ref_elem
 reference element in the mesh for which the data structure is initialized More...
 
const MAST::SystemInitialization_sys_init
 system initialization object for this element More...
 
RealMatrixX _T_mat
 Transformation matrix defines T_ij = V_i^t . More...
 
bool _use_local_elem
 

Constructor & Destructor Documentation

◆ GeomElem()

MAST::GeomElem::GeomElem ( )

Definition at line 32 of file geom_elem.cpp.

◆ ~GeomElem()

MAST::GeomElem::~GeomElem ( )
virtual

Definition at line 41 of file geom_elem.cpp.

Member Function Documentation

◆ _init_local_elem()

void MAST::GeomElem::_init_local_elem ( )
protected

initializes the local element

Definition at line 331 of file geom_elem.cpp.

◆ _init_local_elem_1d()

void MAST::GeomElem::_init_local_elem_1d ( )
protected

initializes the local element

Definition at line 388 of file geom_elem.cpp.

◆ _init_local_elem_2d()

void MAST::GeomElem::_init_local_elem_2d ( )
protected

initializes the local element

Definition at line 461 of file geom_elem.cpp.

◆ dim()

unsigned int MAST::GeomElem::dim ( ) const
Returns
dimension of the element.

Definition at line 91 of file geom_elem.cpp.

◆ domain_surface_normal()

const RealVectorX & MAST::GeomElem::domain_surface_normal ( ) const

Definition at line 227 of file geom_elem.cpp.

◆ external_side_loads_for_quadrature_elem()

void MAST::GeomElem::external_side_loads_for_quadrature_elem ( std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase *> &  bc,
std::map< unsigned int, std::vector< MAST::BoundaryConditionBase *>> &  loads 
) const
virtual

From the given list of boundary loads, this identifies the sides of the quadrature element and the loads in bc that are to be applied on it.

The map of side and loads is returned in loads.

Definition at line 183 of file geom_elem.cpp.

◆ get_boundary_ids_on_quadrature_elem_side()

void MAST::GeomElem::get_boundary_ids_on_quadrature_elem_side ( unsigned int  s,
std::vector< libMesh::boundary_id_type > &  bc_ids 
) const
virtual

Definition at line 217 of file geom_elem.cpp.

◆ get_fe_type()

libMesh::FEType MAST::GeomElem::get_fe_type ( unsigned int  i) const
Returns
FEType for the underlying FE discretization for the i-th variable

Definition at line 110 of file geom_elem.cpp.

◆ get_quadrature_elem()

const libMesh::Elem & MAST::GeomElem::get_quadrature_elem ( ) const
virtual
Returns
a reference to quadrature element.

Reimplemented in MAST::LevelSetIntersectedElem.

Definition at line 72 of file geom_elem.cpp.

◆ get_quadrature_local_elem()

const libMesh::Elem & MAST::GeomElem::get_quadrature_local_elem ( ) const
virtual
Returns
a reference to quadrature element.

Reimplemented in MAST::LevelSetIntersectedElem.

Definition at line 81 of file geom_elem.cpp.

◆ get_reference_elem()

const libMesh::Elem & MAST::GeomElem::get_reference_elem ( ) const
virtual
Returns
a reference to element in the mesh.

Definition at line 54 of file geom_elem.cpp.

◆ get_reference_local_elem()

const libMesh::Elem & MAST::GeomElem::get_reference_local_elem ( ) const
virtual
Returns
a reference to element in the mesh.

Definition at line 63 of file geom_elem.cpp.

◆ init()

void MAST::GeomElem::init ( const libMesh::Elem &  elem,
const MAST::SystemInitialization sys_init 
)
virtual

initialize the object for the specified reference elem.

Reimplemented in MAST::LevelSetIntersectedElem.

Definition at line 134 of file geom_elem.cpp.

◆ init_fe()

std::unique_ptr< MAST::FEBase > MAST::GeomElem::init_fe ( bool  init_grads,
bool  init_second_order_derivative,
int  extra_quadrature_order = 0 
) const
virtual

initializes the finite element shape function and quadrature object with the order of quadrature rule changed based on the extra_quadrature_order.

Reimplemented in MAST::LevelSetIntersectedElem.

Definition at line 148 of file geom_elem.cpp.

◆ init_side_fe()

std::unique_ptr< MAST::FEBase > MAST::GeomElem::init_side_fe ( unsigned int  s,
bool  init_grads,
bool  init_second_order_derivative,
int  extra_quadrature_order = 0 
) const
virtual

initializes the finite element shape function and quadrature object for the side with the order of quadrature rule changed based on the extra_quadrature_order

Reimplemented in MAST::LevelSetIntersectedElem.

Definition at line 165 of file geom_elem.cpp.

◆ n_sides_quadrature_elem()

unsigned int MAST::GeomElem::n_sides_quadrature_elem ( ) const

number of sides on quadrature element.

Definition at line 101 of file geom_elem.cpp.

◆ set_bending()

void MAST::GeomElem::set_bending ( bool  onoff)

This sets the 1D elements to extension/torsional stiffness only.

This is useful when modeling truss structures (e.g. using CROD elements in Nastran) which do not require an orientation vector like beam elements do. By default, 1D elements include bending.

Definition at line 128 of file geom_elem.cpp.

◆ set_local_y_vector()

void MAST::GeomElem::set_local_y_vector ( const RealVectorX y_vec)

for 1D elements the transformed coordinate system attached to the element defines the local x-axis along the length of the element.

The local-y axis should be specified using this method. This then provides the complete information about setting up the local three-dimensional coordinate system.

Definition at line 119 of file geom_elem.cpp.

◆ T_matrix()

const RealMatrixX & MAST::GeomElem::T_matrix ( ) const

O.

Returns
the 3x3 transformation matrix to transform a vector from the element to global coordinate system: $ v_{global} = T v_{elem} $.

Definition at line 236 of file geom_elem.cpp.

◆ transform_point_to_global_coordinate()

void MAST::GeomElem::transform_point_to_global_coordinate ( const libMesh::Point &  local_pt,
libMesh::Point &  global_pt 
) const

Definition at line 246 of file geom_elem.cpp.

◆ transform_vector_to_global_coordinate() [1/2]

void MAST::GeomElem::transform_vector_to_global_coordinate ( const libMesh::Point &  local_vec,
libMesh::Point &  global_vec 
) const

Definition at line 266 of file geom_elem.cpp.

◆ transform_vector_to_global_coordinate() [2/2]

void MAST::GeomElem::transform_vector_to_global_coordinate ( const RealVectorX local_vec,
RealVectorX global_vec 
) const

Definition at line 300 of file geom_elem.cpp.

◆ transform_vector_to_local_coordinate() [1/2]

void MAST::GeomElem::transform_vector_to_local_coordinate ( const libMesh::Point &  global_vec,
libMesh::Point &  local_vec 
) const

Definition at line 283 of file geom_elem.cpp.

◆ transform_vector_to_local_coordinate() [2/2]

void MAST::GeomElem::transform_vector_to_local_coordinate ( const RealVectorX global_vec,
RealVectorX local_vec 
) const

Definition at line 311 of file geom_elem.cpp.

◆ use_local_elem()

bool MAST::GeomElem::use_local_elem ( ) const

Vector and matrix quantities defined on one- and two-dimensional elements that are oriented in two or three-dimensional spaces may need to be transformed to/from element coordinate system.

This is required for one-dimensional elements that are not along the x-axis or two-dimensional elements that are not in the x-y plane.

Returns
true if transformation is required, false otherwise.

Definition at line 322 of file geom_elem.cpp.

Member Data Documentation

◆ _bending

bool MAST::GeomElem::_bending = true
protected

Defines if bending is used in this element or not.

True by default Added for github issue #40

Definition at line 275 of file geom_elem.h.

◆ _domain_surface_normal

RealVectorX MAST::GeomElem::_domain_surface_normal
protected

surface normal of the 1D/2D element.

Definition at line 255 of file geom_elem.h.

◆ _local_elem

libMesh::Elem* MAST::GeomElem::_local_elem
protected

a local element is created if

Definition at line 244 of file geom_elem.h.

◆ _local_nodes

std::vector<libMesh::Node*> MAST::GeomElem::_local_nodes
protected

nodes for local element

Definition at line 260 of file geom_elem.h.

◆ _local_y

RealVectorX MAST::GeomElem::_local_y
protected

the y-axis that is used to define the coordinate system for a 1-D element.

This must be provided a local element is required.

Definition at line 250 of file geom_elem.h.

◆ _ref_elem

const libMesh::Elem* MAST::GeomElem::_ref_elem
protected

reference element in the mesh for which the data structure is initialized

Definition at line 239 of file geom_elem.h.

◆ _sys_init

const MAST::SystemInitialization* MAST::GeomElem::_sys_init
protected

system initialization object for this element

Definition at line 231 of file geom_elem.h.

◆ _T_mat

RealMatrixX MAST::GeomElem::_T_mat
protected

Transformation matrix defines T_ij = V_i^t .

Vn_j, where V_i are the unit vectors of the global cs, and Vn_j are the unit vectors of the local cs. To transform a vector from global to local cs, an_j = T^t a_i, and the reverse transformation is obtained as a_j = T an_i

Definition at line 269 of file geom_elem.h.

◆ _use_local_elem

bool MAST::GeomElem::_use_local_elem
protected

Definition at line 233 of file geom_elem.h.


The documentation for this class was generated from the following files: