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