libMesh
fe_clough.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 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/enum_to_string.h"
23 #include "libmesh/fe.h"
24 #include "libmesh/fe_interface.h"
25 #include "libmesh/fe_macro.h"
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  const bool add_p_level)
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 = order + add_p_level*elem->p_level();
49 
50  // FEType object to be passed to various FEInterface functions below.
51  FEType fe_type(order, 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  const unsigned int n_sf =
61  FEInterface::n_shape_functions(fe_type, elem);
62 
63  std::vector<Point> refspace_nodes;
64  FEBase::get_refspace_nodes(elem_type,refspace_nodes);
65  libmesh_assert_equal_to (refspace_nodes.size(), n_nodes);
66  libmesh_assert_equal_to (elem_soln.size(), n_sf);
67 
68  // Zero before summation
69  std::fill(nodal_soln.begin(), nodal_soln.end(), 0);
70 
71  for (unsigned int n=0; n<n_nodes; n++)
72  // u_i = Sum (alpha_i phi_i)
73  for (unsigned int i=0; i<n_sf; i++)
74  nodal_soln[n] += elem_soln[i] *
75  FEInterface::shape(fe_type, elem, i, refspace_nodes[n]);
76 
77  return;
78  }
79 
80  default:
81  libmesh_error_msg("ERROR: Invalid total order " << totalorder);
82  }
83 } // clough_nodal_soln()
84 
85 
86 
87 unsigned int CLOUGH_n_dofs(const ElemType t, const Order o)
88 {
89  if (t == INVALID_ELEM)
90  return 0;
91 
92  switch (o)
93  {
94  // Piecewise cubic shape functions with linear flux edges
95  case SECOND:
96  {
97  switch (t)
98  {
99  case TRI6:
100  case TRI7:
101  return 9;
102 
103  default:
104  libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
105  }
106  }
107  // Piecewise cubic Clough-Tocher element
108  case THIRD:
109  {
110  switch (t)
111  {
112  case TRI6:
113  case TRI7:
114  return 12;
115 
116  default:
117  libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
118  }
119  }
120 
121  default:
122  libmesh_error_msg("ERROR: Invalid Order " << Utility::enum_to_string(o) << " selected for CLOUGH FE family!");
123  }
124 } // CLOUGH_n_dofs()
125 
126 
127 
128 unsigned int CLOUGH_n_dofs(const Elem * e, const Order o)
129 {
130  libmesh_assert(e);
131  return CLOUGH_n_dofs(e->type(), o);
132 }
133 
134 
135 
136 unsigned int CLOUGH_n_dofs_at_node(const ElemType t,
137  const Order o,
138  const unsigned int n)
139 {
140  switch (o)
141  {
142  // Piecewise cubic shape functions with linear flux edges
143  case SECOND:
144  {
145  switch (t)
146  {
147  // The 2D Clough-Tocher defined on a 6-noded triangle
148  case TRI6:
149  case TRI7:
150  {
151  switch (n)
152  {
153  case 0:
154  case 1:
155  case 2:
156  return 3;
157 
158  case 3:
159  case 4:
160  case 5:
161  case 6:
162  return 0;
163 
164  default:
165  libmesh_error_msg("ERROR: Invalid node ID " << n << " selected for TRI!");
166  }
167  }
168 
169  default:
170  libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
171  }
172  }
173  // The third-order clough shape functions
174  case THIRD:
175  {
176  switch (t)
177  {
178  // The 2D Clough-Tocher defined on a 6-noded triangle
179  case TRI6:
180  case TRI7:
181  {
182  switch (n)
183  {
184  case 0:
185  case 1:
186  case 2:
187  return 3;
188 
189  case 3:
190  case 4:
191  case 5:
192  return 1;
193 
194  case 6:
195  return 0;
196 
197  default:
198  libmesh_error_msg("ERROR: Invalid node ID " << n << " selected for TRI6!");
199  }
200  }
201 
202  default:
203  libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
204  }
205  }
206  default:
207  libmesh_error_msg("ERROR: Invalid Order " << Utility::enum_to_string(o) << " selected for CLOUGH FE family!");
208  }
209 } // CLOUGH_n_dofs_at_node()
210 
211 
212 
213 unsigned int CLOUGH_n_dofs_at_node(const Elem & e,
214  const Order o,
215  const unsigned int n)
216 {
217  return CLOUGH_n_dofs_at_node(e.type(), o, n);
218 }
219 
220 
221 
222 unsigned int CLOUGH_n_dofs_per_elem(const ElemType t, const Order o)
223 {
224  switch (o)
225  {
226  // Piecewise cubic shape functions with linear flux edges
227  case SECOND:
228  // The third-order Clough-Tocher shape functions
229  case THIRD:
230  {
231  switch (t)
232  {
233  // The 2D clough defined on a 6-noded triangle
234  case TRI6:
235  case TRI7:
236  return 0;
237 
238  default:
239  libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
240  }
241  }
242 
243  default:
244  libmesh_error_msg("ERROR: Invalid Order " << Utility::enum_to_string(o) << " selected for CLOUGH FE family!");
245  }
246 } // CLOUGH_n_dofs_per_elem()
247 
248 
249 
250 unsigned int CLOUGH_n_dofs_per_elem(const Elem & e, const Order o)
251 {
252  return CLOUGH_n_dofs_per_elem(e.type(), o);
253 }
254 
255 } // anonymous
256 
257 
258 
259 
260 // Instantiate (side_) nodal_soln() function for every dimension
261 LIBMESH_FE_NODAL_SOLN(CLOUGH, clough_nodal_soln)
263 
264 
265 // Instantiate n_dofs*() functions for every dimension
266 LIBMESH_DEFAULT_NDOFS(CLOUGH)
267 
268 
269 // Clough FEMs are C^1 continuous
270 template <> FEContinuity FE<0,CLOUGH>::get_continuity() const { return C_ONE; }
271 template <> FEContinuity FE<1,CLOUGH>::get_continuity() const { return C_ONE; }
272 template <> FEContinuity FE<2,CLOUGH>::get_continuity() const { return C_ONE; }
273 template <> FEContinuity FE<3,CLOUGH>::get_continuity() const { return C_ONE; }
274 
275 // Clough FEMs are not (currently) hierarchic
276 template <> bool FE<0,CLOUGH>::is_hierarchic() const { return false; } // FIXME - this will be changed
277 template <> bool FE<1,CLOUGH>::is_hierarchic() const { return false; } // FIXME - this will be changed
278 template <> bool FE<2,CLOUGH>::is_hierarchic() const { return false; } // FIXME - this will be changed
279 template <> bool FE<3,CLOUGH>::is_hierarchic() const { return false; } // FIXME - this will be changed
280 
281 #ifdef LIBMESH_ENABLE_AMR
282 // compute_constraints() specializations are only needed for 2 and 3D
283 template <>
285  DofMap & dof_map,
286  const unsigned int variable_number,
287  const Elem * elem)
288 { compute_proj_constraints(constraints, dof_map, variable_number, elem); }
289 
290 template <>
292  DofMap & dof_map,
293  const unsigned int variable_number,
294  const Elem * elem)
295 { compute_proj_constraints(constraints, dof_map, variable_number, elem); }
296 #endif // #ifdef LIBMESH_ENABLE_AMR
297 
298 // Clough FEM shapes need reinit
299 template <> bool FE<0,CLOUGH>::shapes_need_reinit() const { return true; }
300 template <> bool FE<1,CLOUGH>::shapes_need_reinit() const { return true; }
301 template <> bool FE<2,CLOUGH>::shapes_need_reinit() const { return true; }
302 template <> bool FE<3,CLOUGH>::shapes_need_reinit() const { return true; }
303 
304 } // namespace libMesh
ElemType
Defines an enum for geometric element types.
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
virtual bool shapes_need_reinit() const override
The libMesh namespace provides an interface to certain functionality in the library.
LIBMESH_FE_NODAL_SOLN(LIBMESH_FE_SIDE_NODAL_SOLN() LIBMESH_DEFAULT_NDOFS(BERNSTEIN) template<> FEContinuity FE< 0 BERNSTEIN, bernstein_nodal_soln)
Definition: fe_bernstein.C:415
virtual bool is_hierarchic() const override
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:179
LIBMESH_FE_SIDE_NODAL_SOLN(HIERARCHIC_VEC)
const dof_id_type n_nodes
Definition: tecplot_io.C:67
static unsigned int n_shape_functions(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:272
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:760
static void get_refspace_nodes(const ElemType t, std::vector< Point > &nodes)
Definition: fe_abstract.C:400
libmesh_assert(ctx)
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)
Computes the constraint matrix contributions (for non-conforming adapted meshes) corresponding to var...
FEContinuity
defines an enum for finite element types to libmesh_assert a certain level (or type? Hcurl?) of continuity.
The constraint matrix storage format.
Definition: dof_map.h:108