fe_hierarchic.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/elem.h"
22 #include "libmesh/fe.h"
23 #include "libmesh/fe_interface.h"
24 #include "libmesh/string_to_enum.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  return 0;
163  default:
164  libmesh_error_msg("ERROR: Invalid node ID " << n << " selected for EDGE2/3!");
165  }
166  case TRI3:
167  libmesh_assert_less (n, 3);
168  libmesh_assert_less (o, 2);
169  libmesh_fallthrough();
170  case TRI6:
171  switch (n)
172  {
173  case 0:
174  case 1:
175  case 2:
176  return 1;
177 
178  case 3:
179  case 4:
180  case 5:
181  return (o-1);
182 
183  // Internal DoFs are associated with the elem, not its nodes
184  default:
185  libmesh_error_msg("ERROR: Invalid node ID " << n << " selected for TRI6!");
186  }
187  case QUAD4:
188  case QUADSHELL4:
189  libmesh_assert_less (n, 4);
190  libmesh_assert_less (o, 2);
191  libmesh_fallthrough();
192  case QUAD8:
193  case QUADSHELL8:
194  case QUAD9:
195  switch (n)
196  {
197  case 0:
198  case 1:
199  case 2:
200  case 3:
201  return 1;
202 
203  case 4:
204  case 5:
205  case 6:
206  case 7:
207  return (o-1);
208 
209  // Internal DoFs are associated with the elem, not its nodes
210  case 8:
211  return 0;
212 
213  default:
214  libmesh_error_msg("ERROR: Invalid node ID " << n << " selected for QUAD4/8/9!");
215  }
216  case HEX8:
217  libmesh_assert_less (n, 8);
218  libmesh_assert_less (o, 2);
219  libmesh_fallthrough();
220  case HEX20:
221  libmesh_assert_less (n, 20);
222  libmesh_assert_less (o, 2);
223  libmesh_fallthrough();
224  case HEX27:
225  switch (n)
226  {
227  case 0:
228  case 1:
229  case 2:
230  case 3:
231  case 4:
232  case 5:
233  case 6:
234  case 7:
235  return 1;
236 
237  case 8:
238  case 9:
239  case 10:
240  case 11:
241  case 12:
242  case 13:
243  case 14:
244  case 15:
245  case 16:
246  case 17:
247  case 18:
248  case 19:
249  return (o-1);
250 
251  case 20:
252  case 21:
253  case 22:
254  case 23:
255  case 24:
256  case 25:
257  return ((o-1)*(o-1));
258 
259  // Internal DoFs are associated with the elem, not its nodes
260  case 26:
261  return 0;
262  default:
263  libmesh_error_msg("ERROR: Invalid node ID " << n << " selected for HEX8/20/27!");
264  }
265 
266  case INVALID_ELEM:
267  return 0;
268 
269  default:
270  libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
271  }
272 } // hierarchic_n_dofs_at_node()
273 
274 
275 
276 
277 unsigned int hierarchic_n_dofs_per_elem(const ElemType t,
278  const Order o)
279 {
280  libmesh_assert_greater (o, 0);
281  switch (t)
282  {
283  case NODEELEM:
284  return 0;
285  case EDGE2:
286  case EDGE3:
287  return (o-1);
288  case TRI3:
289  case TRISHELL3:
290  case QUAD4:
291  case QUADSHELL4:
292  return 0;
293  case TRI6:
294  return ((o-1)*(o-2)/2);
295  case QUAD8:
296  case QUADSHELL8:
297  case QUAD9:
298  return ((o-1)*(o-1));
299  case HEX8:
300  case HEX20:
301  libmesh_assert_less (o, 2);
302  return 0;
303  case HEX27:
304  return ((o-1)*(o-1)*(o-1));
305  case INVALID_ELEM:
306  return 0;
307  default:
308  libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
309  }
310 } // hierarchic_n_dofs_per_elem()
311 
312 } // anonymous namespace
313 
314 
315 
316 
317  // Do full-specialization of nodal_soln() function for every
318  // dimension, instead of explicit instantiation at the end of this
319  // file.
320  // This could be macro-ified so that it fits on one line...
321 template <>
323  const Order order,
324  const std::vector<Number> & elem_soln,
325  std::vector<Number> & nodal_soln)
326 { hierarchic_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/0); }
327 
328 template <>
330  const Order order,
331  const std::vector<Number> & elem_soln,
332  std::vector<Number> & nodal_soln)
333 { hierarchic_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/1); }
334 
335 template <>
337  const Order order,
338  const std::vector<Number> & elem_soln,
339  std::vector<Number> & nodal_soln)
340 { hierarchic_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/2); }
341 
342 template <>
344  const Order order,
345  const std::vector<Number> & elem_soln,
346  std::vector<Number> & nodal_soln)
347 { hierarchic_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/3); }
348 
349 
350 // Full specialization of n_dofs() function for every dimension
351 template <> unsigned int FE<0,HIERARCHIC>::n_dofs(const ElemType t, const Order o) { return hierarchic_n_dofs(t, o); }
352 template <> unsigned int FE<1,HIERARCHIC>::n_dofs(const ElemType t, const Order o) { return hierarchic_n_dofs(t, o); }
353 template <> unsigned int FE<2,HIERARCHIC>::n_dofs(const ElemType t, const Order o) { return hierarchic_n_dofs(t, o); }
354 template <> unsigned int FE<3,HIERARCHIC>::n_dofs(const ElemType t, const Order o) { return hierarchic_n_dofs(t, o); }
355 
356 // Full specialization of n_dofs_at_node() function for every dimension.
357 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); }
358 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); }
359 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); }
360 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); }
361 
362 // Full specialization of n_dofs_per_elem() function for every dimension.
363 template <> unsigned int FE<0,HIERARCHIC>::n_dofs_per_elem(const ElemType t, const Order o) { return hierarchic_n_dofs_per_elem(t, o); }
364 template <> unsigned int FE<1,HIERARCHIC>::n_dofs_per_elem(const ElemType t, const Order o) { return hierarchic_n_dofs_per_elem(t, o); }
365 template <> unsigned int FE<2,HIERARCHIC>::n_dofs_per_elem(const ElemType t, const Order o) { return hierarchic_n_dofs_per_elem(t, o); }
366 template <> unsigned int FE<3,HIERARCHIC>::n_dofs_per_elem(const ElemType t, const Order o) { return hierarchic_n_dofs_per_elem(t, o); }
367 
368 // Hierarchic FEMs are C^0 continuous
369 template <> FEContinuity FE<0,HIERARCHIC>::get_continuity() const { return C_ZERO; }
370 template <> FEContinuity FE<1,HIERARCHIC>::get_continuity() const { return C_ZERO; }
371 template <> FEContinuity FE<2,HIERARCHIC>::get_continuity() const { return C_ZERO; }
372 template <> FEContinuity FE<3,HIERARCHIC>::get_continuity() const { return C_ZERO; }
373 
374 // Hierarchic FEMs are hierarchic (duh!)
375 template <> bool FE<0,HIERARCHIC>::is_hierarchic() const { return true; }
376 template <> bool FE<1,HIERARCHIC>::is_hierarchic() const { return true; }
377 template <> bool FE<2,HIERARCHIC>::is_hierarchic() const { return true; }
378 template <> bool FE<3,HIERARCHIC>::is_hierarchic() const { return true; }
379 
380 #ifdef LIBMESH_ENABLE_AMR
381 // compute_constraints() specializations are only needed for 2 and 3D
382 template <>
384  DofMap & dof_map,
385  const unsigned int variable_number,
386  const Elem * elem)
387 { compute_proj_constraints(constraints, dof_map, variable_number, elem); }
388 
389 template <>
391  DofMap & dof_map,
392  const unsigned int variable_number,
393  const Elem * elem)
394 { compute_proj_constraints(constraints, dof_map, variable_number, elem); }
395 #endif // #ifdef LIBMESH_ENABLE_AMR
396 
397 // Hierarchic FEM shapes need reinit
398 template <> bool FE<0,HIERARCHIC>::shapes_need_reinit() const { return true; }
399 template <> bool FE<1,HIERARCHIC>::shapes_need_reinit() const { return true; }
400 template <> bool FE<2,HIERARCHIC>::shapes_need_reinit() const { return true; }
401 template <> bool FE<3,HIERARCHIC>::shapes_need_reinit() const { return true; }
402 
403 } // namespace libMesh
static unsigned int n_dofs(const ElemType t, const Order o)
The base class for all geometric element types.
Definition: elem.h:100
static unsigned int n_dofs_at_node(const ElemType t, const Order o, const unsigned int n)
virtual bool shapes_need_reinit() const override
virtual bool is_hierarchic() const override
Manages the degrees of freedom (DOFs) in a simulation.
Definition: dof_map.h:176
const dof_id_type n_nodes
Definition: tecplot_io.C:68
static unsigned int n_shape_functions(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:428
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:657
static void get_refspace_nodes(const ElemType t, std::vector< Point > &nodes)
Definition: fe_abstract.C:283
static unsigned int n_dofs_per_elem(const ElemType t, const Order o)
virtual FEContinuity get_continuity() const override
std::string enum_to_string(const T e)
static void compute_constraints(DofConstraints &constraints, DofMap &dof_map, const unsigned int variable_number, const Elem *elem)
static void nodal_soln(const Elem *elem, const Order o, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln)