libMesh
fe_l2_hierarchic.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/elem.h"
22 #include "libmesh/fe.h"
23 #include "libmesh/fe_interface.h"
24 #include "libmesh/enum_to_string.h"
25 
26 namespace libMesh
27 {
28 
29 // ------------------------------------------------------------
30 // Hierarchic-specific implementations
31 
32 // Anonymous namespace for local helper functions
33 namespace {
34 
35 void l2_hierarchic_nodal_soln(const Elem * elem,
36  const Order order,
37  const std::vector<Number> & elem_soln,
38  std::vector<Number> & nodal_soln,
39  unsigned Dim)
40 {
41  const unsigned int n_nodes = elem->n_nodes();
42 
43  const ElemType elem_type = elem->type();
44 
45  nodal_soln.resize(n_nodes);
46 
47  const Order totalorder = static_cast<Order>(order + elem->p_level());
48 
49  // FEType object to be passed to various FEInterface functions below.
50  FEType fe_type(totalorder, L2_HIERARCHIC);
51 
52  switch (totalorder)
53  {
54  // Constant shape functions
55  case CONSTANT:
56  {
57  libmesh_assert_equal_to (elem_soln.size(), 1);
58 
59  const Number val = elem_soln[0];
60 
61  for (unsigned int n=0; n<n_nodes; n++)
62  nodal_soln[n] = val;
63 
64  return;
65  }
66 
67 
68  // For other orders do interpolation at the nodes
69  // explicitly.
70  default:
71  {
72 
73  const unsigned int n_sf =
74  // FE<Dim,T>::n_shape_functions(elem_type, totalorder);
75  FEInterface::n_shape_functions(Dim, fe_type, elem_type);
76 
77  std::vector<Point> refspace_nodes;
78  FEBase::get_refspace_nodes(elem_type,refspace_nodes);
79  libmesh_assert_equal_to (refspace_nodes.size(), n_nodes);
80 
81  for (unsigned int n=0; n<n_nodes; n++)
82  {
83  libmesh_assert_equal_to (elem_soln.size(), n_sf);
84 
85  // Zero before summation
86  nodal_soln[n] = 0;
87 
88  // u_i = Sum (alpha_i phi_i)
89  for (unsigned int i=0; i<n_sf; i++)
90  nodal_soln[n] += elem_soln[i] *
91  // FE<Dim,T>::shape(elem, order, i, mapped_point);
92  FEInterface::shape(Dim, fe_type, elem, i, refspace_nodes[n]);
93  }
94 
95  return;
96  }
97  }
98 } // l2_hierarchic_nodal_soln()
99 
100 
101 
102 
103 unsigned int l2_hierarchic_n_dofs(const ElemType t, const Order o)
104 {
105  libmesh_assert_greater (o, 0);
106  switch (t)
107  {
108  case NODEELEM:
109  return 1;
110  case EDGE2:
111  case EDGE3:
112  return (o+1);
113  case QUAD4:
114  case QUADSHELL4:
115  case QUAD8:
116  case QUADSHELL8:
117  case QUAD9:
118  return ((o+1)*(o+1));
119  case HEX8:
120  case HEX20:
121  case HEX27:
122  return ((o+1)*(o+1)*(o+1));
123  case TRI3:
124  libmesh_assert_less (o, 2);
125  libmesh_fallthrough();
126  case TRI6:
127  return ((o+1)*(o+2)/2);
128  case INVALID_ELEM:
129  return 0;
130  default:
131  libmesh_error_msg("ERROR: Invalid ElemType " << Utility::enum_to_string(t) << " selected for L2_HIERARCHIC FE family!");
132  }
133 } // l2_hierarchic_n_dofs()
134 
135 
136 } // anonymous namespace
137 
138 
139 
140 
141  // Do full-specialization of nodal_soln() function for every
142  // dimension, instead of explicit instantiation at the end of this
143  // file.
144  // This could be macro-ified so that it fits on one line...
145 template <>
147  const Order order,
148  const std::vector<Number> & elem_soln,
149  std::vector<Number> & nodal_soln)
150 { l2_hierarchic_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/0); }
151 
152 template <>
154  const Order order,
155  const std::vector<Number> & elem_soln,
156  std::vector<Number> & nodal_soln)
157 { l2_hierarchic_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/1); }
158 
159 template <>
161  const Order order,
162  const std::vector<Number> & elem_soln,
163  std::vector<Number> & nodal_soln)
164 { l2_hierarchic_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/2); }
165 
166 template <>
168  const Order order,
169  const std::vector<Number> & elem_soln,
170  std::vector<Number> & nodal_soln)
171 { l2_hierarchic_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/3); }
172 
173 // Full specialization of n_dofs() function for every dimension
174 template <> unsigned int FE<0,L2_HIERARCHIC>::n_dofs(const ElemType t, const Order o) { return l2_hierarchic_n_dofs(t, o); }
175 template <> unsigned int FE<1,L2_HIERARCHIC>::n_dofs(const ElemType t, const Order o) { return l2_hierarchic_n_dofs(t, o); }
176 template <> unsigned int FE<2,L2_HIERARCHIC>::n_dofs(const ElemType t, const Order o) { return l2_hierarchic_n_dofs(t, o); }
177 template <> unsigned int FE<3,L2_HIERARCHIC>::n_dofs(const ElemType t, const Order o) { return l2_hierarchic_n_dofs(t, o); }
178 
179 // Full specialization of n_dofs_at_node() function for every dimension.
180 // Discontinuous L2 elements only have interior nodes
181 template <> unsigned int FE<0,L2_HIERARCHIC>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
182 template <> unsigned int FE<1,L2_HIERARCHIC>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
183 template <> unsigned int FE<2,L2_HIERARCHIC>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
184 template <> unsigned int FE<3,L2_HIERARCHIC>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
185 
186 // Full specialization of n_dofs_per_elem() function for every dimension.
187 template <> unsigned int FE<0,L2_HIERARCHIC>::n_dofs_per_elem(const ElemType t, const Order o) { return l2_hierarchic_n_dofs(t, o); }
188 template <> unsigned int FE<1,L2_HIERARCHIC>::n_dofs_per_elem(const ElemType t, const Order o) { return l2_hierarchic_n_dofs(t, o); }
189 template <> unsigned int FE<2,L2_HIERARCHIC>::n_dofs_per_elem(const ElemType t, const Order o) { return l2_hierarchic_n_dofs(t, o); }
190 template <> unsigned int FE<3,L2_HIERARCHIC>::n_dofs_per_elem(const ElemType t, const Order o) { return l2_hierarchic_n_dofs(t, o); }
191 
192 // L2 Hierarchic FEMs are C^0 continuous
197 
198 // L2 Hierarchic FEMs are hierarchic (duh!)
199 template <> bool FE<0,L2_HIERARCHIC>::is_hierarchic() const { return true; }
200 template <> bool FE<1,L2_HIERARCHIC>::is_hierarchic() const { return true; }
201 template <> bool FE<2,L2_HIERARCHIC>::is_hierarchic() const { return true; }
202 template <> bool FE<3,L2_HIERARCHIC>::is_hierarchic() const { return true; }
203 
204 #ifdef LIBMESH_ENABLE_AMR
205 // compute_constraints() is a NOOP for DISCONTINUOUS FE's
206 template <>
208  DofMap &,
209  const unsigned int,
210  const Elem *)
211 { }
212 
213 template <>
215  DofMap &,
216  const unsigned int,
217  const Elem *)
218 { }
219 #endif // #ifdef LIBMESH_ENABLE_AMR
220 
221 // L2-Hierarchic FEM shapes need reinit
222 template <> bool FE<0,L2_HIERARCHIC>::shapes_need_reinit() const { return true; }
223 template <> bool FE<1,L2_HIERARCHIC>::shapes_need_reinit() const { return true; }
224 template <> bool FE<2,L2_HIERARCHIC>::shapes_need_reinit() const { return true; }
225 template <> bool FE<3,L2_HIERARCHIC>::shapes_need_reinit() const { return true; }
226 
227 } // namespace libMesh
libMesh::DofConstraints
The constraint matrix storage format.
Definition: dof_map.h:105
libMesh::HEX20
Definition: enum_elem_type.h:48
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::FEInterface::n_shape_functions
static unsigned int n_shape_functions(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:446
libMesh::L2_HIERARCHIC
Definition: enum_fe_family.h:40
libMesh::HEX8
Definition: enum_elem_type.h:47
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::Order
Order
Definition: enum_order.h:40
libMesh::HEX27
Definition: enum_elem_type.h:49
libMesh::QUAD4
Definition: enum_elem_type.h:41
libMesh::DISCONTINUOUS
Definition: enum_fe_family.h:78
libMesh::TRI3
Definition: enum_elem_type.h:39
libMesh::FE::compute_constraints
static void compute_constraints(DofConstraints &constraints, DofMap &dof_map, const unsigned int variable_number, const Elem *elem)
Computes the constraint matrix contributions (for non-conforming adapted meshes) corresponding to var...
libMesh::CONSTANT
Definition: enum_order.h:41
n_nodes
const dof_id_type n_nodes
Definition: tecplot_io.C:68
libMesh::Utility::enum_to_string
std::string enum_to_string(const T e)
libMesh::TRI6
Definition: enum_elem_type.h:40
libMesh::FE::shapes_need_reinit
virtual bool shapes_need_reinit() const override
libMesh::INVALID_ELEM
Definition: enum_elem_type.h:75
libMesh::FE::n_dofs_at_node
static unsigned int n_dofs_at_node(const ElemType t, const Order o, const unsigned int n)
libMesh::FE::nodal_soln
static void nodal_soln(const Elem *elem, const Order o, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln)
Build the nodal soln from the element soln.
libMesh::FE::get_continuity
virtual FEContinuity get_continuity() const override
libMesh::QUADSHELL8
Definition: enum_elem_type.h:73
libMesh::DofMap
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:176
libMesh::FE::n_dofs_per_elem
static unsigned int n_dofs_per_elem(const ElemType t, const Order o)
libMesh::EDGE3
Definition: enum_elem_type.h:36
libMesh::QUADSHELL4
Definition: enum_elem_type.h:72
libMesh::FE::is_hierarchic
virtual bool is_hierarchic() const override
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::NODEELEM
Definition: enum_elem_type.h:66
libMesh::FEAbstract::get_refspace_nodes
static void get_refspace_nodes(const ElemType t, std::vector< Point > &nodes)
Definition: fe_abstract.C:308
libMesh::FE::n_dofs
static unsigned int n_dofs(const ElemType t, const Order o)
libMesh::QUAD9
Definition: enum_elem_type.h:43
libMesh::FEContinuity
FEContinuity
Definition: enum_fe_family.h:77
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::FEInterface::shape
static Real shape(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p)
Definition: fe_interface.C:687