libMesh
fe_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 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, 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 } // hierarchic_nodal_soln()
99 
100 
101 
102 
103 unsigned int 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  libmesh_assert_less (o, 2);
116  libmesh_fallthrough();
117  case QUAD8:
118  case QUADSHELL8:
119  case QUAD9:
120  return ((o+1)*(o+1));
121  case HEX8:
122  libmesh_assert_less (o, 2);
123  libmesh_fallthrough();
124  case HEX20:
125  libmesh_assert_less (o, 2);
126  libmesh_fallthrough();
127  case HEX27:
128  return ((o+1)*(o+1)*(o+1));
129  case TRI3:
130  libmesh_assert_less (o, 2);
131  libmesh_fallthrough();
132  case TRI6:
133  return ((o+1)*(o+2)/2);
134  case INVALID_ELEM:
135  return 0;
136  default:
137  libmesh_error_msg("ERROR: Invalid ElemType " << Utility::enum_to_string(t) << " selected for HIERARCHIC FE family!");
138  }
139 } // hierarchic_n_dofs()
140 
141 
142 
143 
144 unsigned int hierarchic_n_dofs_at_node(const ElemType t,
145  const Order o,
146  const unsigned int n)
147 {
148  libmesh_assert_greater (o, 0);
149  switch (t)
150  {
151  case NODEELEM:
152  return 1;
153  case EDGE2:
154  case EDGE3:
155  switch (n)
156  {
157  case 0:
158  case 1:
159  return 1;
160  // Internal DoFs are associated with the elem, not its nodes
161  case 2:
162  libmesh_assert_equal_to(t, EDGE3);
163  return 0;
164  default:
165  libmesh_error_msg("ERROR: Invalid node ID " << n << " selected for EDGE2/3!");
166  }
167  case TRI3:
168  libmesh_assert_less (n, 3);
169  libmesh_assert_less (o, 2);
170  libmesh_fallthrough();
171  case TRI6:
172  switch (n)
173  {
174  case 0:
175  case 1:
176  case 2:
177  return 1;
178 
179  case 3:
180  case 4:
181  case 5:
182  return (o-1);
183 
184  // Internal DoFs are associated with the elem, not its nodes
185  default:
186  libmesh_error_msg("ERROR: Invalid node ID " << n << " selected for TRI6!");
187  }
188  case QUAD4:
189  case QUADSHELL4:
190  libmesh_assert_less (n, 4);
191  libmesh_assert_less (o, 2);
192  libmesh_fallthrough();
193  case QUAD8:
194  case QUADSHELL8:
195  case QUAD9:
196  switch (n)
197  {
198  case 0:
199  case 1:
200  case 2:
201  case 3:
202  return 1;
203 
204  case 4:
205  case 5:
206  case 6:
207  case 7:
208  return (o-1);
209 
210  // Internal DoFs are associated with the elem, not its nodes
211  case 8:
212  return 0;
213 
214  default:
215  libmesh_error_msg("ERROR: Invalid node ID " << n << " selected for QUAD4/8/9!");
216  }
217  case HEX8:
218  libmesh_assert_less (n, 8);
219  libmesh_assert_less (o, 2);
220  libmesh_fallthrough();
221  case HEX20:
222  libmesh_assert_less (n, 20);
223  libmesh_assert_less (o, 2);
224  libmesh_fallthrough();
225  case HEX27:
226  switch (n)
227  {
228  case 0:
229  case 1:
230  case 2:
231  case 3:
232  case 4:
233  case 5:
234  case 6:
235  case 7:
236  return 1;
237 
238  case 8:
239  case 9:
240  case 10:
241  case 11:
242  case 12:
243  case 13:
244  case 14:
245  case 15:
246  case 16:
247  case 17:
248  case 18:
249  case 19:
250  return (o-1);
251 
252  case 20:
253  case 21:
254  case 22:
255  case 23:
256  case 24:
257  case 25:
258  return ((o-1)*(o-1));
259 
260  // Internal DoFs are associated with the elem, not its nodes
261  case 26:
262  return 0;
263  default:
264  libmesh_error_msg("ERROR: Invalid node ID " << n << " selected for HEX8/20/27!");
265  }
266 
267  case INVALID_ELEM:
268  return 0;
269 
270  default:
271  libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
272  }
273 } // hierarchic_n_dofs_at_node()
274 
275 
276 
277 
278 unsigned int hierarchic_n_dofs_per_elem(const ElemType t,
279  const Order o)
280 {
281  libmesh_assert_greater (o, 0);
282  switch (t)
283  {
284  case NODEELEM:
285  return 0;
286  case EDGE2:
287  case EDGE3:
288  return (o-1);
289  case TRI3:
290  case TRISHELL3:
291  case QUAD4:
292  case QUADSHELL4:
293  return 0;
294  case TRI6:
295  return ((o-1)*(o-2)/2);
296  case QUAD8:
297  case QUADSHELL8:
298  case QUAD9:
299  return ((o-1)*(o-1));
300  case HEX8:
301  case HEX20:
302  libmesh_assert_less (o, 2);
303  return 0;
304  case HEX27:
305  return ((o-1)*(o-1)*(o-1));
306  case INVALID_ELEM:
307  return 0;
308  default:
309  libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
310  }
311 } // hierarchic_n_dofs_per_elem()
312 
313 } // anonymous namespace
314 
315 
316 
317 
318  // Do full-specialization of nodal_soln() function for every
319  // dimension, instead of explicit instantiation at the end of this
320  // file.
321  // This could be macro-ified so that it fits on one line...
322 template <>
324  const Order order,
325  const std::vector<Number> & elem_soln,
326  std::vector<Number> & nodal_soln)
327 { hierarchic_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/0); }
328 
329 template <>
331  const Order order,
332  const std::vector<Number> & elem_soln,
333  std::vector<Number> & nodal_soln)
334 { hierarchic_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/1); }
335 
336 template <>
338  const Order order,
339  const std::vector<Number> & elem_soln,
340  std::vector<Number> & nodal_soln)
341 { hierarchic_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/2); }
342 
343 template <>
345  const Order order,
346  const std::vector<Number> & elem_soln,
347  std::vector<Number> & nodal_soln)
348 { hierarchic_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/3); }
349 
350 
351 // Full specialization of n_dofs() function for every dimension
352 template <> unsigned int FE<0,HIERARCHIC>::n_dofs(const ElemType t, const Order o) { return hierarchic_n_dofs(t, o); }
353 template <> unsigned int FE<1,HIERARCHIC>::n_dofs(const ElemType t, const Order o) { return hierarchic_n_dofs(t, o); }
354 template <> unsigned int FE<2,HIERARCHIC>::n_dofs(const ElemType t, const Order o) { return hierarchic_n_dofs(t, o); }
355 template <> unsigned int FE<3,HIERARCHIC>::n_dofs(const ElemType t, const Order o) { return hierarchic_n_dofs(t, o); }
356 
357 // Full specialization of n_dofs_at_node() function for every dimension.
358 template <> unsigned int FE<0,HIERARCHIC>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return hierarchic_n_dofs_at_node(t, o, n); }
359 template <> unsigned int FE<1,HIERARCHIC>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return hierarchic_n_dofs_at_node(t, o, n); }
360 template <> unsigned int FE<2,HIERARCHIC>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return hierarchic_n_dofs_at_node(t, o, n); }
361 template <> unsigned int FE<3,HIERARCHIC>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return hierarchic_n_dofs_at_node(t, o, n); }
362 
363 // Full specialization of n_dofs_per_elem() function for every dimension.
364 template <> unsigned int FE<0,HIERARCHIC>::n_dofs_per_elem(const ElemType t, const Order o) { return hierarchic_n_dofs_per_elem(t, o); }
365 template <> unsigned int FE<1,HIERARCHIC>::n_dofs_per_elem(const ElemType t, const Order o) { return hierarchic_n_dofs_per_elem(t, o); }
366 template <> unsigned int FE<2,HIERARCHIC>::n_dofs_per_elem(const ElemType t, const Order o) { return hierarchic_n_dofs_per_elem(t, o); }
367 template <> unsigned int FE<3,HIERARCHIC>::n_dofs_per_elem(const ElemType t, const Order o) { return hierarchic_n_dofs_per_elem(t, o); }
368 
369 // Hierarchic FEMs are C^0 continuous
370 template <> FEContinuity FE<0,HIERARCHIC>::get_continuity() const { return C_ZERO; }
371 template <> FEContinuity FE<1,HIERARCHIC>::get_continuity() const { return C_ZERO; }
372 template <> FEContinuity FE<2,HIERARCHIC>::get_continuity() const { return C_ZERO; }
373 template <> FEContinuity FE<3,HIERARCHIC>::get_continuity() const { return C_ZERO; }
374 
375 // Hierarchic FEMs are hierarchic (duh!)
376 template <> bool FE<0,HIERARCHIC>::is_hierarchic() const { return true; }
377 template <> bool FE<1,HIERARCHIC>::is_hierarchic() const { return true; }
378 template <> bool FE<2,HIERARCHIC>::is_hierarchic() const { return true; }
379 template <> bool FE<3,HIERARCHIC>::is_hierarchic() const { return true; }
380 
381 #ifdef LIBMESH_ENABLE_AMR
382 // compute_constraints() specializations are only needed for 2 and 3D
383 template <>
385  DofMap & dof_map,
386  const unsigned int variable_number,
387  const Elem * elem)
388 { compute_proj_constraints(constraints, dof_map, variable_number, elem); }
389 
390 template <>
392  DofMap & dof_map,
393  const unsigned int variable_number,
394  const Elem * elem)
395 { compute_proj_constraints(constraints, dof_map, variable_number, elem); }
396 #endif // #ifdef LIBMESH_ENABLE_AMR
397 
398 // Hierarchic FEM shapes need reinit
399 template <> bool FE<0,HIERARCHIC>::shapes_need_reinit() const { return true; }
400 template <> bool FE<1,HIERARCHIC>::shapes_need_reinit() const { return true; }
401 template <> bool FE<2,HIERARCHIC>::shapes_need_reinit() const { return true; }
402 template <> bool FE<3,HIERARCHIC>::shapes_need_reinit() const { return true; }
403 
404 } // 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::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::TRI3
Definition: enum_elem_type.h:39
libMesh::C_ZERO
Definition: enum_fe_family.h:79
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::HIERARCHIC
Definition: enum_fe_family.h:37
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::TRISHELL3
Definition: enum_elem_type.h:71
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