Line data Source code
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 13458 : 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 13458 : const unsigned int n_nodes = elem->n_nodes();
43 :
44 13458 : const ElemType elem_type = elem->type();
45 :
46 13458 : nodal_soln.resize(n_nodes);
47 :
48 14587 : const Order totalorder = order + add_p_level*elem->p_level();
49 :
50 : // FEType object to be passed to various FEInterface functions below.
51 1129 : FEType fe_type(order, CLOUGH);
52 :
53 13458 : switch (totalorder)
54 : {
55 : // Piecewise cubic shape functions with linear flux edges
56 13458 : case SECOND:
57 : // Piecewise cubic shape functions
58 : case THIRD:
59 : {
60 : const unsigned int n_sf =
61 13458 : FEInterface::n_shape_functions(fe_type, elem);
62 :
63 2258 : std::vector<Point> refspace_nodes;
64 13458 : FEBase::get_refspace_nodes(elem_type,refspace_nodes);
65 1129 : libmesh_assert_equal_to (refspace_nodes.size(), n_nodes);
66 1129 : libmesh_assert_equal_to (elem_soln.size(), n_sf);
67 :
68 : // Zero before summation
69 1129 : std::fill(nodal_soln.begin(), nodal_soln.end(), 0);
70 :
71 94206 : for (unsigned int n=0; n<n_nodes; n++)
72 : // u_i = Sum (alpha_i phi_i)
73 1049724 : for (unsigned int i=0; i<n_sf; i++)
74 1050264 : nodal_soln[n] += elem_soln[i] *
75 1050264 : FEInterface::shape(fe_type, elem, i, refspace_nodes[n]);
76 :
77 14587 : return;
78 : }
79 :
80 0 : default:
81 0 : libmesh_error_msg("ERROR: Invalid total order " << totalorder);
82 : }
83 : } // clough_nodal_soln()
84 :
85 :
86 :
87 574420 : unsigned int CLOUGH_n_dofs(const ElemType t, const Order o)
88 : {
89 574420 : if (t == INVALID_ELEM)
90 0 : return 0;
91 :
92 574420 : switch (o)
93 : {
94 : // Piecewise cubic shape functions with linear flux edges
95 0 : case SECOND:
96 : {
97 0 : switch (t)
98 : {
99 0 : case TRI6:
100 : case TRI7:
101 0 : return 9;
102 :
103 0 : default:
104 0 : 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 574420 : case THIRD:
109 : {
110 574420 : switch (t)
111 : {
112 83488 : case TRI6:
113 : case TRI7:
114 83488 : return 12;
115 :
116 0 : default:
117 0 : libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
118 : }
119 : }
120 :
121 0 : default:
122 0 : 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 532608 : unsigned int CLOUGH_n_dofs(const Elem * e, const Order o)
129 : {
130 83488 : libmesh_assert(e);
131 574420 : return CLOUGH_n_dofs(e->type(), o);
132 : }
133 :
134 :
135 :
136 5165194 : unsigned int CLOUGH_n_dofs_at_node(const ElemType t,
137 : const Order o,
138 : const unsigned int n)
139 : {
140 5165194 : switch (o)
141 : {
142 : // Piecewise cubic shape functions with linear flux edges
143 0 : case SECOND:
144 : {
145 0 : switch (t)
146 : {
147 : // The 2D Clough-Tocher defined on a 6-noded triangle
148 0 : case TRI6:
149 : case TRI7:
150 : {
151 0 : switch (n)
152 : {
153 0 : case 0:
154 : case 1:
155 : case 2:
156 0 : return 3;
157 :
158 0 : case 3:
159 : case 4:
160 : case 5:
161 : case 6:
162 0 : return 0;
163 :
164 0 : default:
165 0 : libmesh_error_msg("ERROR: Invalid node ID " << n << " selected for TRI!");
166 : }
167 : }
168 :
169 0 : default:
170 0 : 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 5165194 : case THIRD:
175 : {
176 5165194 : switch (t)
177 : {
178 : // The 2D Clough-Tocher defined on a 6-noded triangle
179 5165194 : case TRI6:
180 : case TRI7:
181 : {
182 4744829 : switch (n)
183 : {
184 217299 : case 0:
185 : case 1:
186 : case 2:
187 217299 : return 3;
188 :
189 203459 : case 3:
190 : case 4:
191 : case 5:
192 203459 : return 1;
193 :
194 405 : case 6:
195 405 : return 0;
196 :
197 0 : default:
198 0 : libmesh_error_msg("ERROR: Invalid node ID " << n << " selected for TRI6!");
199 : }
200 : }
201 :
202 0 : default:
203 0 : libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
204 : }
205 : }
206 0 : default:
207 0 : 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 901652 : unsigned int CLOUGH_n_dofs_at_node(const Elem & e,
214 : const Order o,
215 : const unsigned int n)
216 : {
217 984082 : return CLOUGH_n_dofs_at_node(e.type(), o, n);
218 : }
219 :
220 :
221 :
222 792496 : unsigned int CLOUGH_n_dofs_per_elem(const ElemType t, const Order o)
223 : {
224 792496 : switch (o)
225 : {
226 : // Piecewise cubic shape functions with linear flux edges
227 792496 : case SECOND:
228 : // The third-order Clough-Tocher shape functions
229 : case THIRD:
230 : {
231 792496 : switch (t)
232 : {
233 : // The 2D clough defined on a 6-noded triangle
234 792496 : case TRI6:
235 : case TRI7:
236 792496 : return 0;
237 :
238 0 : default:
239 0 : libmesh_error_msg("ERROR: Bad ElemType = " << Utility::enum_to_string(t) << " for " << Utility::enum_to_string(o) << " order approximation!");
240 : }
241 : }
242 :
243 0 : default:
244 0 : 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 727943 : unsigned int CLOUGH_n_dofs_per_elem(const Elem & e, const Order o)
251 : {
252 792496 : 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 13458 : LIBMESH_FE_NODAL_SOLN(CLOUGH, clough_nodal_soln)
262 0 : LIBMESH_FE_SIDE_NODAL_SOLN(CLOUGH)
263 :
264 :
265 : // Instantiate n_dofs*() functions for every dimension
266 6532110 : LIBMESH_DEFAULT_NDOFS(CLOUGH)
267 :
268 :
269 : // Clough FEMs are C^1 continuous
270 0 : template <> FEContinuity FE<0,CLOUGH>::get_continuity() const { return C_ONE; }
271 0 : template <> FEContinuity FE<1,CLOUGH>::get_continuity() const { return C_ONE; }
272 75634 : template <> FEContinuity FE<2,CLOUGH>::get_continuity() const { return C_ONE; }
273 0 : template <> FEContinuity FE<3,CLOUGH>::get_continuity() const { return C_ONE; }
274 :
275 : // Clough FEMs are not (currently) hierarchic
276 0 : template <> bool FE<0,CLOUGH>::is_hierarchic() const { return false; } // FIXME - this will be changed
277 0 : template <> bool FE<1,CLOUGH>::is_hierarchic() const { return false; } // FIXME - this will be changed
278 44 : template <> bool FE<2,CLOUGH>::is_hierarchic() const { return false; } // FIXME - this will be changed
279 0 : 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 <>
284 38496 : void FE<2,CLOUGH>::compute_constraints (DofConstraints & constraints,
285 : DofMap & dof_map,
286 : const unsigned int variable_number,
287 : const Elem * elem)
288 38496 : { compute_proj_constraints(constraints, dof_map, variable_number, elem); }
289 :
290 : template <>
291 0 : void FE<3,CLOUGH>::compute_constraints (DofConstraints & constraints,
292 : DofMap & dof_map,
293 : const unsigned int variable_number,
294 : const Elem * elem)
295 0 : { compute_proj_constraints(constraints, dof_map, variable_number, elem); }
296 : #endif // #ifdef LIBMESH_ENABLE_AMR
297 :
298 : // Clough FEM shapes need reinit
299 0 : template <> bool FE<0,CLOUGH>::shapes_need_reinit() const { return true; }
300 0 : template <> bool FE<1,CLOUGH>::shapes_need_reinit() const { return true; }
301 386408 : template <> bool FE<2,CLOUGH>::shapes_need_reinit() const { return true; }
302 0 : template <> bool FE<3,CLOUGH>::shapes_need_reinit() const { return true; }
303 :
304 : } // namespace libMesh
|