libMesh
fe_rational.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/libmesh_config.h"
22 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
23 
24 #include "libmesh/fe.h"
25 #include "libmesh/elem.h"
26 #include "libmesh/fe_interface.h"
27 #include "libmesh/enum_to_string.h"
28 
29 // Anonymous namespace for local helper functions
30 namespace {
31 
32 using namespace libMesh;
33 
34 static const FEFamily _underlying_fe_family = BERNSTEIN;
35 
36 void rational_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, _underlying_fe_family);
52 
53  switch (totalorder)
54  {
55  // Constant shape functions
56  case CONSTANT:
57  {
58  libmesh_assert_equal_to (elem_soln.size(), 1);
59 
60  const Number val = elem_soln[0];
61 
62  for (unsigned int n=0; n<n_nodes; n++)
63  nodal_soln[n] = val;
64 
65  return;
66  }
67 
68 
69  // For other bases do interpolation at the nodes
70  // explicitly.
71  default:
72  {
73  const unsigned int n_sf =
74  FEInterface::n_shape_functions(Dim, fe_type, elem_type);
75 
76  std::vector<Point> refspace_nodes;
77  FEBase::get_refspace_nodes(elem_type,refspace_nodes);
78  libmesh_assert_equal_to (refspace_nodes.size(), n_nodes);
79  libmesh_assert_equal_to (n_sf, n_nodes);
80  libmesh_assert_equal_to (elem_soln.size(), n_sf);
81 
82  std::vector<Real> node_weights(n_nodes);
83 
84  for (unsigned int n=0; n<n_nodes; n++)
85  node_weights[n] = elem->node_ref(n).get_extra_integer(0);
86 
87  for (unsigned int n=0; n<n_nodes; n++)
88  {
89  std::vector<Real> weighted_shape(n_sf);
90  Real weighted_sum = 0;
91 
92  for (unsigned int i=0; i<n_sf; i++)
93  {
94  weighted_shape[i] = node_weights[i] *
95  FEInterface::shape(Dim, fe_type, elem, i, refspace_nodes[n]);
96  weighted_sum += weighted_shape[i];
97  }
98 
99  // Zero before summation
100  nodal_soln[n] = 0;
101 
102  // u_i = Sum (alpha_i w_i phi_i) / Sum (w_j phi_j)
103  for (unsigned int i=0; i<n_sf; i++)
104  nodal_soln[n] += elem_soln[i] * weighted_shape[i];
105  nodal_soln[n] /= weighted_sum;
106  }
107 
108  return;
109  }
110  }
111 } // rational_nodal_soln()
112 
113 } // anon namespace
114 
115 
116 namespace libMesh {
117 
118 template <>
120  const Order order,
121  const std::vector<Number> & elem_soln,
122  std::vector<Number> & nodal_soln)
123 { rational_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/0); }
124 
125 template <>
127  const Order order,
128  const std::vector<Number> & elem_soln,
129  std::vector<Number> & nodal_soln)
130 { rational_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/1); }
131 
132 template <>
134  const Order order,
135  const std::vector<Number> & elem_soln,
136  std::vector<Number> & nodal_soln)
137 { rational_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/2); }
138 
139 template <>
141  const Order order,
142  const std::vector<Number> & elem_soln,
143  std::vector<Number> & nodal_soln)
144 { rational_nodal_soln(elem, order, elem_soln, nodal_soln, /*Dim=*/3); }
145 
146 
147 // Full specialization of n_dofs() function for every dimension
148 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); }
149 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); }
150 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); }
151 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); }
152 
153 // Full specialization of n_dofs_at_node() function for every dimension.
154 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); }
155 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); }
156 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); }
157 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); }
158 
159 // Full specialization of n_dofs_per_elem() function for every dimension.
164 
165 // Our current rational FEMs are C^0 continuous
170 
171 // Rational FEMs are in general not hierarchic
172 template <> bool FE<0,RATIONAL_BERNSTEIN>::is_hierarchic() const { return false; }
173 template <> bool FE<1,RATIONAL_BERNSTEIN>::is_hierarchic() const { return false; }
174 template <> bool FE<2,RATIONAL_BERNSTEIN>::is_hierarchic() const { return false; }
175 template <> bool FE<3,RATIONAL_BERNSTEIN>::is_hierarchic() const { return false; }
176 
177 #ifdef LIBMESH_ENABLE_AMR
178 // compute_constraints() specializations are only needed for 2 and 3D
179 template <>
181  DofMap & dof_map,
182  const unsigned int variable_number,
183  const Elem * elem)
184 { compute_proj_constraints(constraints, dof_map, variable_number, elem); }
185 
186 template <>
188  DofMap & dof_map,
189  const unsigned int variable_number,
190  const Elem * elem)
191 { compute_proj_constraints(constraints, dof_map, variable_number, elem); }
192 #endif // #ifdef LIBMESH_ENABLE_AMR
193 
194 // Bernstein shapes need reinit only for approximation orders >= 3,
195 // and at the moment we're only supporting equal-order shapes due to
196 // the need for nodal weights.
197 template <> bool FE<0,RATIONAL_BERNSTEIN>::shapes_need_reinit() const { return false; }
198 template <> bool FE<1,RATIONAL_BERNSTEIN>::shapes_need_reinit() const { return false; }
199 template <> bool FE<2,RATIONAL_BERNSTEIN>::shapes_need_reinit() const { return false; }
200 template <> bool FE<3,RATIONAL_BERNSTEIN>::shapes_need_reinit() const { return false; }
201 
202 } // namespace libMesh
203 
204 #endif //LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
libMesh::DofConstraints
The constraint matrix storage format.
Definition: dof_map.h:105
libMesh::Number
Real Number
Definition: libmesh_common.h:195
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::Elem::n_nodes
virtual unsigned int n_nodes() const =0
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::Elem::p_level
unsigned int p_level() const
Definition: elem.h:2512
libMesh::C_ZERO
Definition: enum_fe_family.h:79
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...
libMesh::CONSTANT
Definition: enum_order.h:41
libMesh::BERNSTEIN
Definition: enum_fe_family.h:43
n_nodes
const dof_id_type n_nodes
Definition: tecplot_io.C:68
libMesh::FE::shapes_need_reinit
virtual bool shapes_need_reinit() const override
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
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::FEFamily
FEFamily
Definition: enum_fe_family.h:34
libMesh::FE::get_continuity
virtual FEContinuity get_continuity() const override
libMesh::DofObject::get_extra_integer
dof_id_type get_extra_integer(const unsigned int index) const
Gets the value on this object of the extra integer associated with index, which should have been obta...
Definition: dof_object.h:1026
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::FE::n_dofs
static unsigned int n_dofs(const ElemType t, const Order o)
libMesh::Elem::node_ref
const Node & node_ref(const unsigned int i) const
Definition: elem.h:2031
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::FEContinuity
FEContinuity
Definition: enum_fe_family.h:77
libMesh::Elem::type
virtual ElemType type() const =0
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