28 #include "libmesh/mesh_base.h" 29 #include "libmesh/boundary_info.h" 34 _use_local_elem (false),
36 _local_elem (nullptr) {
121 libmesh_assert_equal_to(y_vec.size(), 3);
141 std::unique_ptr<MAST::FEBase>
143 bool init_second_order_derivative,
144 int extra_quadrature_order)
const {
149 fe->set_extra_quadrature_order(extra_quadrature_order);
150 fe->set_evaluate_second_order_derivatives(init_second_order_derivative);
152 fe->init(*
this, init_grads);
158 std::unique_ptr<MAST::FEBase>
161 int extra_quadrature_order)
const {
164 fe->set_extra_quadrature_order(extra_quadrature_order);
166 fe->init_for_side(*
this, s, init_grads);
175 (std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc,
176 std::map<
unsigned int, std::vector<MAST::BoundaryConditionBase*>>& loads)
const {
180 typedef std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*> maptype;
181 std::pair<maptype::const_iterator, maptype::const_iterator> range;
185 const libMesh::BoundaryInfo&
188 std::vector<libMesh::boundary_id_type> bids;
190 for (
unsigned int i=0; i<
_ref_elem->n_sides(); i++) {
195 for (
unsigned int j=0; j<bids.size(); j++) {
197 range = bc.equal_range(bids[j]);
199 maptype::const_iterator it = range.first;
200 for ( ; it != range.second; it++) loads[i].push_back(it->second);
209 (
unsigned int s, std::vector<libMesh::boundary_id_type>& bc_ids)
const {
239 libMesh::Point& global_pt)
const {
246 for (
unsigned int j=0; j<3; j++)
247 for (
unsigned int k=0; k<3; k++)
248 global_pt(j) +=
_T_mat(j,k)*local_pt(k);
259 libMesh::Point& global_vec)
const {
266 for (
unsigned int j=0; j<3; j++)
267 for (
unsigned int k=0; k<3; k++)
268 global_vec(j) +=
_T_mat(j,k)*local_vec(k);
276 libMesh::Point& local_vec)
const {
283 for (
unsigned int j=0; j<3; j++)
284 for (
unsigned int k=0; k<3; k++)
285 local_vec(j) +=
_T_mat(k,j) * global_vec(k);
297 global_vec =
_T_mat * local_vec;
308 local_vec =
_T_mat.transpose() * global_vec;
339 for (
unsigned int i=1; i<
_ref_elem->n_nodes(); i++) {
340 if (std::fabs(y-
_ref_elem->point(i)(1)) != 0. ||
341 std::fabs(z-
_ref_elem->point(i)(2)) != 0.)
358 for (
unsigned int i=1; i<
_ref_elem->n_nodes(); i++) {
359 if (std::fabs(z-
_ref_elem->point(i)(2)) != 0.)
386 libmesh_assert_equal_to(
_ref_elem->dim(), 1);
391 libMesh::Point v1, v2, v3, p;
393 for (
unsigned int i=0; i<3; i++) v2(i) =
_local_y(i);
395 libmesh_assert_greater(v3.norm(), 0.);
397 v2 = v3.cross(v1); v2 /= v2.norm();
399 _T_mat = RealMatrixX::Zero(3,3);
403 for (
unsigned int i=0; i<
_ref_elem->n_nodes(); i++) {
414 for (
unsigned int i=0; i<3; i++) {
424 for (
unsigned int j=0; j<3; j++)
425 for (
unsigned int k=0; k<3; k++)
437 libmesh_assert_equal_to(
_ref_elem->dim(), 2);
441 libMesh::Point v1, v2, v3, p;
444 v3 = v1.cross(v2); v3 /= v3.norm();
445 v2 = v3.cross(v1); v2 /= v2.norm();
449 for (
unsigned int i=0; i<3; i++)
452 _T_mat = RealMatrixX::Zero(3,3);
456 for (
unsigned int i=0; i<
_ref_elem->n_nodes(); i++) {
467 for (
unsigned int i=0; i<3; i++) {
477 for (
unsigned int j=0; j<3; j++)
478 for (
unsigned int k=0; k<3; k++)
bool use_local_elem() const
Vector and matrix quantities defined on one- and two-dimensional elements that are oriented in two or...
MAST::NonlinearSystem & system()
const MAST::SystemInitialization * _sys_init
system initialization object for this element
void _init_local_elem_1d()
initializes the local element
virtual void clear() libmesh_override
Clear all the data structures associated with the system.
const libMesh::Elem * _ref_elem
reference element in the mesh for which the data structure is initialized
RealVectorX _local_y
the y-axis that is used to define the coordinate system for a 1-D element.
void _init_local_elem()
initializes the local element
RealMatrixX _T_mat
Transformation matrix defines T_ij = V_i^t .
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 al...
virtual const libMesh::Elem & get_reference_elem() const
libMesh::Elem * _local_elem
a local element is created if
const libMesh::FEType & fetype(unsigned int i) const
unsigned int n_sides_quadrature_elem() const
number of sides on quadrature element.
virtual const libMesh::Elem & get_reference_local_elem() const
libMesh::FEType get_fe_type(unsigned int i) const
Matrix< Real, Dynamic, Dynamic > RealMatrixX
const RealMatrixX & T_matrix() const
O.
RealVectorX _domain_surface_normal
surface normal of the 1D/2D element.
Matrix< Real, Dynamic, 1 > RealVectorX
std::vector< libMesh::Node * > _local_nodes
nodes for local element
virtual void get_boundary_ids_on_quadrature_elem_side(unsigned int s, std::vector< libMesh::boundary_id_type > &bc_ids) const
void transform_point_to_global_coordinate(const libMesh::Point &local_pt, libMesh::Point &global_pt) const
virtual std::unique_ptr< MAST::FEBase > init_side_fe(unsigned int s, bool init_grads, int extra_quadrature_order=0) const
initializes the finite element shape function and quadrature object for the side with the order of qu...
const RealVectorX & domain_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 lo...
virtual const libMesh::Elem & get_quadrature_elem() const
void _init_local_elem_2d()
initializes the local element
void transform_vector_to_local_coordinate(const libMesh::Point &global_vec, libMesh::Point &local_vec) const
virtual void init(const libMesh::Elem &elem, const MAST::SystemInitialization &sys_init)
initialize the object for the specified reference elem.
virtual const libMesh::Elem & get_quadrature_local_elem() const
virtual std::unique_ptr< MAST::FEBase > init_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...
void transform_vector_to_global_coordinate(const libMesh::Point &local_vec, libMesh::Point &global_vec) const