libMesh
inf_fe_base_radial.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 // Local includes
21 #include "libmesh/libmesh_config.h"
22 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
23 #include "libmesh/inf_fe.h"
24 #include "libmesh/inf_fe_macro.h"
25 #include "libmesh/fe.h"
26 #include "libmesh/elem.h"
27 
28 namespace libMesh
29 {
30 
31 
32 
33 // ------------------------------------------------------------
34 // InfFEBase class members
35 Elem * InfFEBase::build_elem (const Elem * inf_elem)
36 {
37  std::unique_ptr<const Elem> ape(inf_elem->build_side_ptr(0));
38 
39  // The incoming inf_elem is const, but this function is required to
40  // return a non-const Elem * so that it can be used by
41  // update_base_elem(). Therefore a const_cast seems to be
42  // unavoidable here.
43  return const_cast<Elem *>(ape.release());
44 }
45 
46 
47 
48 
50 {
51  switch (type)
52  {
53  // 3D infinite elements:
54  // with Dim=3 -> infinite elements on their own
55  case INFHEX8:
56  return QUAD4;
57 
58  case INFHEX16:
59  return QUAD8;
60 
61  case INFHEX18:
62  return QUAD9;
63 
64  case INFPRISM6:
65  return TRI3;
66 
67  case INFPRISM12:
68  return TRI6;
69 
70  // 2D infinite elements:
71  // with Dim=3 -> used as boundary condition,
72  // with Dim=2 -> infinite elements on their own
73  case INFQUAD4:
74  return EDGE2;
75 
76  case INFQUAD6:
77  return EDGE3;
78 
79  // 1D infinite elements:
80  // with Dim=2 -> used as boundary condition,
81  // with Dim=1 -> infinite elements on their own,
82  // but no base element!
83  case INFEDGE2:
84  return INVALID_ELEM;
85 
86  default:
87  libmesh_error_msg("ERROR: Unsupported element type!: " << type);
88  }
89 }
90 
91 
92 
93 
94 
95 unsigned int InfFEBase::n_base_mapping_sf (const Elem & base_elem,
96  const Order base_mapping_order)
97 {
98  switch (base_elem.dim())
99  {
100  case 0:
101  return 1;
102  case 1:
103  return FE<1,LAGRANGE>::n_shape_functions (base_elem.type(),
104  base_mapping_order);
105  case 2:
106  return FE<2,LAGRANGE>::n_shape_functions (base_elem.type(),
107  base_mapping_order);
108  default:
109  libmesh_error_msg("Unsupported base_elem dim = " << base_elem.dim());
110  }
111 }
112 
113 
114 
115 
116 
117 // ------------------------------------------------------------
118 // InfFERadial class members
119 unsigned int InfFERadial::n_dofs_at_node (const Order o_radial,
120  const unsigned int n_onion)
121 {
122  libmesh_assert_less (n_onion, 2);
123 
124  if (n_onion == 0)
125  /*
126  * in the base, no matter what, we have 1 node associated
127  * with radial direction
128  */
129  return 1;
130  else
131  /*
132  * this works, since for Order o_radial=CONST=0, we still
133  * have the (1-v)/2 mode, associated to the base
134  */
135  return static_cast<unsigned int>(o_radial);
136 }
137 
138 
139 } // namespace libMesh
140 
141 #endif //ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
libMesh::Elem::build_side_ptr
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i, bool proxy=true)=0
libMesh::FE::n_shape_functions
virtual unsigned int n_shape_functions() const override
Definition: fe.C:50
libMesh::INFHEX8
Definition: enum_elem_type.h:60
libMesh::INFQUAD4
Definition: enum_elem_type.h:58
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::Elem::dim
virtual unsigned short dim() const =0
libMesh::Order
Order
Definition: enum_order.h:40
libMesh::INFEDGE2
Definition: enum_elem_type.h:57
libMesh::INFPRISM6
Definition: enum_elem_type.h:63
libMesh::INFHEX18
Definition: enum_elem_type.h:62
libMesh::QUAD4
Definition: enum_elem_type.h:41
libMesh::TRI3
Definition: enum_elem_type.h:39
libMesh::INFHEX16
Definition: enum_elem_type.h:61
libMesh::INFPRISM12
Definition: enum_elem_type.h:64
libMesh::TRI6
Definition: enum_elem_type.h:40
libMesh::INVALID_ELEM
Definition: enum_elem_type.h:75
libMesh::InfFEBase::build_elem
static Elem * build_elem(const Elem *inf_elem)
Build the base element of an infinite element.
Definition: inf_fe_base_radial.C:35
libMesh::INFQUAD6
Definition: enum_elem_type.h:59
libMesh::InfFERadial::n_dofs_at_node
static unsigned int n_dofs_at_node(const Order o_radial, const unsigned int n_onion)
Definition: inf_fe_base_radial.C:119
libMesh::EDGE3
Definition: enum_elem_type.h:36
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::QUAD9
Definition: enum_elem_type.h:43
libMesh::Elem::type
virtual ElemType type() const =0
libMesh::InfFEBase::get_elem_type
static ElemType get_elem_type(const ElemType type)
Definition: inf_fe_base_radial.C:49
libMesh::EDGE2
Definition: enum_elem_type.h:35
libMesh::QUAD8
Definition: enum_elem_type.h:42
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33
libMesh::InfFEBase::n_base_mapping_sf
static unsigned int n_base_mapping_sf(const Elem &base_elem, const Order base_mapping_order)
Definition: inf_fe_base_radial.C:95