inf_fe_base_radial.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 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 // InfFE::Base class members
35 template <unsigned int Dim, FEFamily T_radial, InfMapType T_base>
37 {
38  std::unique_ptr<const Elem> ape(inf_elem->build_side_ptr(0));
39 
40  // The incoming inf_elem is const, but this function is required to
41  // return a non-const Elem * so that it can be used by
42  // update_base_elem(). Therefore a const_cast seems to be
43  // unavoidable here.
44  return const_cast<Elem *>(ape.release());
45 }
46 
47 
48 
49 
50 template <unsigned int Dim, FEFamily T_radial, InfMapType T_base>
52 {
53  switch (type)
54  {
55  // 3D infinite elements:
56  // with Dim=3 -> infinite elements on their own
57  case INFHEX8:
58  return QUAD4;
59 
60  case INFHEX16:
61  return QUAD8;
62 
63  case INFHEX18:
64  return QUAD9;
65 
66  case INFPRISM6:
67  return TRI3;
68 
69  case INFPRISM12:
70  return TRI6;
71 
72  // 2D infinite elements:
73  // with Dim=3 -> used as boundary condition,
74  // with Dim=2 -> infinite elements on their own
75  case INFQUAD4:
76  return EDGE2;
77 
78  case INFQUAD6:
79  return EDGE3;
80 
81  // 1D infinite elements:
82  // with Dim=2 -> used as boundary condition,
83  // with Dim=1 -> infinite elements on their own,
84  // but no base element!
85  case INFEDGE2:
86  return INVALID_ELEM;
87 
88  default:
89  libmesh_error_msg("ERROR: Unsupported element type!: " << type);
90  }
91 }
92 
93 
94 
95 
96 
97 template <unsigned int Dim, FEFamily T_radial, InfMapType T_base>
99  const Order base_mapping_order)
100 {
101  if (Dim == 1)
102  return 1;
103 
104  else if (Dim == 2)
105  return FE<1,LAGRANGE>::n_shape_functions (base_elem_type,
106  base_mapping_order);
107  else if (Dim == 3)
108  return FE<2,LAGRANGE>::n_shape_functions (base_elem_type,
109  base_mapping_order);
110  else
111  {
112  libmesh_error_msg("Unsupported Dim = " << Dim);
113  return 0;
114  }
115 }
116 
117 
118 
119 
120 
121 // ------------------------------------------------------------
122 // InfFE::Radial class members
123 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
125  const unsigned int n_onion)
126 {
127  libmesh_assert_less (n_onion, 2);
128 
129  if (n_onion == 0)
130  /*
131  * in the base, no matter what, we have 1 node associated
132  * with radial direction
133  */
134  return 1;
135  else
136  /*
137  * this works, since for Order o_radial=CONST=0, we still
138  * have the (1-v)/2 mode, associated to the base
139  */
140  return static_cast<unsigned int>(o_radial);
141 }
142 
143 
144 
145 
146 
147 
148 //--------------------------------------------------------------
149 // Explicit instantiations
159 INSTANTIATE_INF_FE_MBRF(1,CARTESIAN,unsigned int,Radial::n_dofs_at_node (const Order,const unsigned int));
160 INSTANTIATE_INF_FE_MBRF(2,CARTESIAN,unsigned int,Radial::n_dofs_at_node (const Order,const unsigned int));
161 INSTANTIATE_INF_FE_MBRF(3,CARTESIAN,unsigned int,Radial::n_dofs_at_node (const Order,const unsigned int));
162 
163 } // namespace libMesh
164 
165 #endif //ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
static Elem * build_elem(const Elem *inf_elem)
The base class for all geometric element types.
Definition: elem.h:100
INSTANTIATE_INF_FE_MBRF(1, CARTESIAN, Elem *, Base::build_elem(const Elem *))
static unsigned int n_dofs_at_node(const Order o_radial, const unsigned int n_onion)
static ElemType get_elem_type(const ElemType type)
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i, bool proxy=true)=0
virtual unsigned int n_shape_functions() const override
Definition: fe.C:50
static unsigned int n_base_mapping_sf(const ElemType base_elem_type, const Order base_mapping_order)