libMesh
fe_rational.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/libmesh_config.h"
22 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
23 
24 #include "libmesh/elem.h"
25 #include "libmesh/enum_to_string.h"
26 #include "libmesh/fe.h"
27 #include "libmesh/fe_interface.h"
28 #include "libmesh/fe_macro.h"
29 
30 // Anonymous namespace for local helper functions
31 namespace {
32 
33 using namespace libMesh;
34 
35 static const FEFamily _underlying_fe_family = BERNSTEIN;
36 
37 void rational_nodal_soln(const Elem * elem,
38  const Order order,
39  const std::vector<Number> & elem_soln,
40  std::vector<Number> & nodal_soln,
41  const bool add_p_level)
42 {
43  const unsigned int n_nodes = elem->n_nodes();
44 
45  const ElemType elem_type = elem->type();
46 
47  nodal_soln.resize(n_nodes);
48 
49  const Order totalorder = order + add_p_level*elem->p_level();
50 
51  // FEType object to be passed to various FEInterface functions below.
52  FEType fe_type(order, _underlying_fe_family);
53  FEType p_refined_fe_type(totalorder, _underlying_fe_family);
54 
55  switch (totalorder)
56  {
57  // Constant shape functions
58  case CONSTANT:
59  {
60  libmesh_assert_equal_to (elem_soln.size(), 1);
61 
62  std::fill(nodal_soln.begin(), nodal_soln.end(), elem_soln[0]);
63 
64  return;
65  }
66 
67 
68  // For other bases do interpolation at the nodes
69  // explicitly.
70  default:
71  {
72  const unsigned int n_sf =
73  FEInterface::n_shape_functions(fe_type, elem);
74 
75  std::vector<Point> refspace_nodes;
76  FEBase::get_refspace_nodes(elem_type,refspace_nodes);
77  libmesh_assert_equal_to (refspace_nodes.size(), n_nodes);
78  libmesh_assert_equal_to (n_sf, n_nodes);
79  libmesh_assert_equal_to (elem_soln.size(), n_sf);
80 
81  std::vector<Real> node_weights(n_nodes);
82 
83  const unsigned char weight_index = elem->mapping_data();
84 
85  for (unsigned int n=0; n<n_nodes; n++)
86  node_weights[n] =
87  elem->node_ref(n).get_extra_datum<Real>(weight_index);
88 
89  // Zero before summation
90  std::fill(nodal_soln.begin(), nodal_soln.end(), 0);
91 
92  for (unsigned int n=0; n<n_nodes; n++)
93  {
94  std::vector<Real> weighted_shape(n_sf);
95  Real weighted_sum = 0;
96 
97  for (unsigned int i=0; i<n_sf; i++)
98  {
99  weighted_shape[i] = node_weights[i] *
100  FEInterface::shape(fe_type, elem, i, refspace_nodes[n]);
101  weighted_sum += weighted_shape[i];
102  }
103 
104  // u_i = Sum (alpha_i w_i phi_i) / Sum (w_j phi_j)
105  for (unsigned int i=0; i<n_sf; i++)
106  nodal_soln[n] += elem_soln[i] * weighted_shape[i];
107  nodal_soln[n] /= weighted_sum;
108  }
109 
110  return;
111  }
112  }
113 } // rational_nodal_soln()
114 
115 } // anon namespace
116 
117 
118 namespace libMesh {
119 
120 
121 // Instantiate (side_) nodal_soln() function for every dimension
124 
125 
126 // Full specialization of n_dofs() function for every dimension
127 template <> unsigned int FE<0,RATIONAL_BERNSTEIN>::n_dofs(const ElemType t, const Order o) { return FE<0,_underlying_fe_family>::n_dofs(t, o); }
128 template <> unsigned int FE<1,RATIONAL_BERNSTEIN>::n_dofs(const ElemType t, const Order o) { return FE<1,_underlying_fe_family>::n_dofs(t, o); }
129 template <> unsigned int FE<2,RATIONAL_BERNSTEIN>::n_dofs(const ElemType t, const Order o) { return FE<2,_underlying_fe_family>::n_dofs(t, o); }
130 template <> unsigned int FE<3,RATIONAL_BERNSTEIN>::n_dofs(const ElemType t, const Order o) { return FE<3,_underlying_fe_family>::n_dofs(t, o); }
131 
132 template <> unsigned int FE<0,RATIONAL_BERNSTEIN>::n_dofs(const Elem * e, const Order o) { return FE<0,_underlying_fe_family>::n_dofs(e, o); }
133 template <> unsigned int FE<1,RATIONAL_BERNSTEIN>::n_dofs(const Elem * e, const Order o) { return FE<1,_underlying_fe_family>::n_dofs(e, o); }
134 template <> unsigned int FE<2,RATIONAL_BERNSTEIN>::n_dofs(const Elem * e, const Order o) { return FE<2,_underlying_fe_family>::n_dofs(e, o); }
135 template <> unsigned int FE<3,RATIONAL_BERNSTEIN>::n_dofs(const Elem * e, const Order o) { return FE<3,_underlying_fe_family>::n_dofs(e, o); }
136 
137 // Full specialization of n_dofs_at_node() function for every dimension.
138 template <> unsigned int FE<0,RATIONAL_BERNSTEIN>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return FE<0,_underlying_fe_family>::n_dofs_at_node(t, o, n); }
139 template <> unsigned int FE<1,RATIONAL_BERNSTEIN>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return FE<1,_underlying_fe_family>::n_dofs_at_node(t, o, n); }
140 template <> unsigned int FE<2,RATIONAL_BERNSTEIN>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return FE<2,_underlying_fe_family>::n_dofs_at_node(t, o, n); }
141 template <> unsigned int FE<3,RATIONAL_BERNSTEIN>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return FE<3,_underlying_fe_family>::n_dofs_at_node(t, o, n); }
142 
143 template <> unsigned int FE<0,RATIONAL_BERNSTEIN>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return FE<0,_underlying_fe_family>::n_dofs_at_node(e, o, n); }
144 template <> unsigned int FE<1,RATIONAL_BERNSTEIN>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return FE<1,_underlying_fe_family>::n_dofs_at_node(e, o, n); }
145 template <> unsigned int FE<2,RATIONAL_BERNSTEIN>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return FE<2,_underlying_fe_family>::n_dofs_at_node(e, o, n); }
146 template <> unsigned int FE<3,RATIONAL_BERNSTEIN>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return FE<3,_underlying_fe_family>::n_dofs_at_node(e, o, n); }
147 
148 // Full specialization of n_dofs_per_elem() function for every dimension.
153 
154 template <> unsigned int FE<0,RATIONAL_BERNSTEIN>::n_dofs_per_elem(const Elem & e, const Order o) { return FE<0,_underlying_fe_family>::n_dofs_per_elem(e, o); }
155 template <> unsigned int FE<1,RATIONAL_BERNSTEIN>::n_dofs_per_elem(const Elem & e, const Order o) { return FE<1,_underlying_fe_family>::n_dofs_per_elem(e, o); }
156 template <> unsigned int FE<2,RATIONAL_BERNSTEIN>::n_dofs_per_elem(const Elem & e, const Order o) { return FE<2,_underlying_fe_family>::n_dofs_per_elem(e, o); }
157 template <> unsigned int FE<3,RATIONAL_BERNSTEIN>::n_dofs_per_elem(const Elem & e, const Order o) { return FE<3,_underlying_fe_family>::n_dofs_per_elem(e, o); }
158 
159 // Our current rational FEMs are C^0 continuous
164 
165 // Rational FEMs are in general not hierarchic
166 template <> bool FE<0,RATIONAL_BERNSTEIN>::is_hierarchic() const { return false; }
167 template <> bool FE<1,RATIONAL_BERNSTEIN>::is_hierarchic() const { return false; }
168 template <> bool FE<2,RATIONAL_BERNSTEIN>::is_hierarchic() const { return false; }
169 template <> bool FE<3,RATIONAL_BERNSTEIN>::is_hierarchic() const { return false; }
170 
171 #ifdef LIBMESH_ENABLE_AMR
172 // compute_constraints() specializations are only needed for 2 and 3D
173 template <>
175  DofMap & dof_map,
176  const unsigned int variable_number,
177  const Elem * elem)
178 { compute_proj_constraints(constraints, dof_map, variable_number, elem); }
179 
180 template <>
182  DofMap & dof_map,
183  const unsigned int variable_number,
184  const Elem * elem)
185 { compute_proj_constraints(constraints, dof_map, variable_number, elem); }
186 #endif // #ifdef LIBMESH_ENABLE_AMR
187 
188 // Bernstein shapes need reinit only for approximation orders >= 3,
189 // but for *RATIONAL* shapes we could need reinit anywhere because our
190 // nodal weights can change from element to element.
191 template <> bool FE<0,RATIONAL_BERNSTEIN>::shapes_need_reinit() const { return true; }
192 template <> bool FE<1,RATIONAL_BERNSTEIN>::shapes_need_reinit() const { return true; }
193 template <> bool FE<2,RATIONAL_BERNSTEIN>::shapes_need_reinit() const { return true; }
194 template <> bool FE<3,RATIONAL_BERNSTEIN>::shapes_need_reinit() const { return true; }
195 
196 } // namespace libMesh
197 
198 #endif //LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:196
unsigned char mapping_data() const
Definition: elem.h:3136
static unsigned int n_dofs(const ElemType t, const Order o)
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
static unsigned int n_dofs_at_node(const ElemType t, const Order o, const unsigned int n)
unsigned int p_level() const
Definition: elem.h:3108
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
const Node & node_ref(const unsigned int i) const
Definition: elem.h:2529
static unsigned int n_shape_functions(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:272
virtual unsigned int n_nodes() const =0
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
static unsigned int n_dofs_per_elem(const ElemType t, const Order o)
virtual FEContinuity get_continuity() const override
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...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
FEContinuity
defines an enum for finite element types to libmesh_assert a certain level (or type? Hcurl?) of continuity.
T get_extra_datum(const unsigned int index) const
Gets the value on this object of the extra datum associated with index, which should have been obtain...
Definition: dof_object.h:1146
FEFamily
defines an enum for finite element families.
virtual ElemType type() const =0
The constraint matrix storage format.
Definition: dof_map.h:108