fe_clough.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 
27 namespace libMesh
28 {
29 
30 // ------------------------------------------------------------
31 // Clough-specific implementations
32 
33 // Anonymous namespace for local helper functions
34 namespace {
35 
36 void clough_nodal_soln(const Elem * elem,
37  const Order order,
38  const std::vector<Number> & elem_soln,
39  std::vector<Number> & nodal_soln,
40  unsigned Dim)
41 {
42  const unsigned int n_nodes = elem->n_nodes();
43 
44  const ElemType elem_type = elem->type();
45 
46  nodal_soln.resize(n_nodes);
47 
48  const Order totalorder = static_cast<Order>(order + elem->p_level());
49 
50  // FEType object to be passed to various FEInterface functions below.
51  FEType fe_type(totalorder, CLOUGH);
52 
53  switch (totalorder)
54  {
55  // Piecewise cubic shape functions with linear flux edges
56  case SECOND:
57  // Piecewise cubic shape functions
58  case THIRD:
59  {
60 
61  const unsigned int n_sf =
62  // FE<Dim,T>::n_shape_functions(elem_type, totalorder);
63  FEInterface::n_shape_functions(Dim, fe_type, elem_type);
64 
65  std::vector<Point> refspace_nodes;
66  FEBase::get_refspace_nodes(elem_type,refspace_nodes);
67  libmesh_assert_equal_to (refspace_nodes.size(), n_nodes);
68 
69  for (unsigned int n=0; n<n_nodes; n++)
70  {
71  libmesh_assert_equal_to (elem_soln.size(), n_sf);
72 
73  // Zero before summation
74  nodal_soln[n] = 0;
75 
76  // u_i = Sum (alpha_i phi_i)
77  for (unsigned int i=0; i<n_sf; i++)
78  nodal_soln[n] += elem_soln[i] *
79  // FE<Dim,T>::shape(elem, order, i, mapped_point);
80  FEInterface::shape(Dim, fe_type, elem, i, refspace_nodes[n]);
81  }
82 
83  return;
84  }
85 
86  default:
87  libmesh_error_msg("ERROR: Invalid total order " << totalorder);
88  }
89 } // clough_nodal_soln()
90 
91 
92 
93 
94 unsigned int clough_n_dofs(const ElemType t, const Order o)
95 {
96  if (t == INVALID_ELEM)
97  return 0;
98 
99  switch (o)
100  {
101  // Piecewise cubic shape functions with linear flux edges
102  case SECOND:
103  {
104  switch (t)
105  {
106  case TRI6:
107  return 9;
108 
109  default:
110  libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
111  }
112  }
113  // Piecewise cubic Clough-Tocher element
114  case THIRD:
115  {
116  switch (t)
117  {
118  case TRI6:
119  return 12;
120 
121  default:
122  libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
123  }
124  }
125 
126  default:
127  libmesh_error_msg("ERROR: Invalid Order " << Utility::enum_to_string(o) << " selected for CLOUGH FE family!");
128  }
129 } // clough_n_dofs()
130 
131 
132 
133 
134 
135 unsigned int clough_n_dofs_at_node(const ElemType t,
136  const Order o,
137  const unsigned int n)
138 {
139  switch (o)
140  {
141  // Piecewise cubic shape functions with linear flux edges
142  case SECOND:
143  {
144  switch (t)
145  {
146  // The 2D Clough-Tocher defined on a 6-noded triangle
147  case TRI6:
148  {
149  switch (n)
150  {
151  case 0:
152  case 1:
153  case 2:
154  return 3;
155 
156  case 3:
157  case 4:
158  case 5:
159  return 0;
160 
161  default:
162  libmesh_error_msg("ERROR: Invalid node ID " << n << " selected for TRI6!");
163  }
164  }
165 
166  default:
167  libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
168  }
169  }
170  // The third-order clough shape functions
171  case THIRD:
172  {
173  switch (t)
174  {
175  // The 2D Clough-Tocher defined on a 6-noded triangle
176  case TRI6:
177  {
178  switch (n)
179  {
180  case 0:
181  case 1:
182  case 2:
183  return 3;
184 
185  case 3:
186  case 4:
187  case 5:
188  return 1;
189 
190  default:
191  libmesh_error_msg("ERROR: Invalid node ID " << n << " selected for TRI6!");
192  }
193  }
194 
195  default:
196  libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
197  }
198  }
199  default:
200  libmesh_error_msg("ERROR: Invalid Order " << Utility::enum_to_string(o) << " selected for CLOUGH FE family!");
201  }
202 } // clough_n_dofs_at_node()
203 
204 
205 
206 
207 unsigned int clough_n_dofs_per_elem(const ElemType t, const Order o)
208 {
209  switch (o)
210  {
211  // Piecewise cubic shape functions with linear flux edges
212  case SECOND:
213  // The third-order Clough-Tocher shape functions
214  case THIRD:
215  {
216  switch (t)
217  {
218  // The 2D clough defined on a 6-noded triangle
219  case TRI6:
220  return 0;
221 
222  default:
223  libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
224  }
225  }
226 
227  default:
228  libmesh_error_msg("ERROR: Invalid Order " << Utility::enum_to_string(o) << " selected for CLOUGH FE family!");
229  }
230 } // clough_n_dofs_per_elem()
231 
232 } // anonymous
233 
234 
235 
236 
237  // Do full-specialization of nodal_soln() function for every
238  // dimension, instead of explicit instantiation at the end of this
239  // file.
240  // This could be macro-ified so that it fits on one line...
241 template <>
242 void FE<0,CLOUGH>::nodal_soln(const Elem * elem,
243  const Order order,
244  const std::vector<Number> & elem_soln,
245  std::vector<Number> & nodal_soln)
246 { clough_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/0); }
247 
248 template <>
249 void FE<1,CLOUGH>::nodal_soln(const Elem * elem,
250  const Order order,
251  const std::vector<Number> & elem_soln,
252  std::vector<Number> & nodal_soln)
253 { clough_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/1); }
254 
255 template <>
256 void FE<2,CLOUGH>::nodal_soln(const Elem * elem,
257  const Order order,
258  const std::vector<Number> & elem_soln,
259  std::vector<Number> & nodal_soln)
260 { clough_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/2); }
261 
262 template <>
263 void FE<3,CLOUGH>::nodal_soln(const Elem * elem,
264  const Order order,
265  const std::vector<Number> & elem_soln,
266  std::vector<Number> & nodal_soln)
267 { clough_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/3); }
268 
269 
270 // Full specialization of n_dofs() function for every dimension
271 template <> unsigned int FE<0,CLOUGH>::n_dofs(const ElemType t, const Order o) { return clough_n_dofs(t, o); }
272 template <> unsigned int FE<1,CLOUGH>::n_dofs(const ElemType t, const Order o) { return clough_n_dofs(t, o); }
273 template <> unsigned int FE<2,CLOUGH>::n_dofs(const ElemType t, const Order o) { return clough_n_dofs(t, o); }
274 template <> unsigned int FE<3,CLOUGH>::n_dofs(const ElemType t, const Order o) { return clough_n_dofs(t, o); }
275 
276 
277 // Full specialization of n_dofs_at_node() function for every dimension.
278 template <> unsigned int FE<0,CLOUGH>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return clough_n_dofs_at_node(t, o, n); }
279 template <> unsigned int FE<1,CLOUGH>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return clough_n_dofs_at_node(t, o, n); }
280 template <> unsigned int FE<2,CLOUGH>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return clough_n_dofs_at_node(t, o, n); }
281 template <> unsigned int FE<3,CLOUGH>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return clough_n_dofs_at_node(t, o, n); }
282 
283 // Full specialization of n_dofs_per_elem() function for every dimension.
284 template <> unsigned int FE<0,CLOUGH>::n_dofs_per_elem(const ElemType t, const Order o) { return clough_n_dofs_per_elem(t, o); }
285 template <> unsigned int FE<1,CLOUGH>::n_dofs_per_elem(const ElemType t, const Order o) { return clough_n_dofs_per_elem(t, o); }
286 template <> unsigned int FE<2,CLOUGH>::n_dofs_per_elem(const ElemType t, const Order o) { return clough_n_dofs_per_elem(t, o); }
287 template <> unsigned int FE<3,CLOUGH>::n_dofs_per_elem(const ElemType t, const Order o) { return clough_n_dofs_per_elem(t, o); }
288 
289 // Clough FEMs are C^1 continuous
290 template <> FEContinuity FE<0,CLOUGH>::get_continuity() const { return C_ONE; }
291 template <> FEContinuity FE<1,CLOUGH>::get_continuity() const { return C_ONE; }
292 template <> FEContinuity FE<2,CLOUGH>::get_continuity() const { return C_ONE; }
293 template <> FEContinuity FE<3,CLOUGH>::get_continuity() const { return C_ONE; }
294 
295 // Clough FEMs are not (currently) hierarchic
296 template <> bool FE<0,CLOUGH>::is_hierarchic() const { return false; } // FIXME - this will be changed
297 template <> bool FE<1,CLOUGH>::is_hierarchic() const { return false; } // FIXME - this will be changed
298 template <> bool FE<2,CLOUGH>::is_hierarchic() const { return false; } // FIXME - this will be changed
299 template <> bool FE<3,CLOUGH>::is_hierarchic() const { return false; } // FIXME - this will be changed
300 
301 #ifdef LIBMESH_ENABLE_AMR
302 // compute_constraints() specializations are only needed for 2 and 3D
303 template <>
305  DofMap & dof_map,
306  const unsigned int variable_number,
307  const Elem * elem)
308 { compute_proj_constraints(constraints, dof_map, variable_number, elem); }
309 
310 template <>
312  DofMap & dof_map,
313  const unsigned int variable_number,
314  const Elem * elem)
315 { compute_proj_constraints(constraints, dof_map, variable_number, elem); }
316 #endif // #ifdef LIBMESH_ENABLE_AMR
317 
318 // Clough FEM shapes need reinit
319 template <> bool FE<0,CLOUGH>::shapes_need_reinit() const { return true; }
320 template <> bool FE<1,CLOUGH>::shapes_need_reinit() const { return true; }
321 template <> bool FE<2,CLOUGH>::shapes_need_reinit() const { return true; }
322 template <> bool FE<3,CLOUGH>::shapes_need_reinit() const { return true; }
323 
324 } // 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)