libMesh
fe_scalar.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/dof_map.h"
22 #include "libmesh/elem.h"
23 #include "libmesh/fe.h"
24 #include "libmesh/fe_macro.h"
25 
26 namespace libMesh
27 {
28 
29 // ------------------------------------------------------------
30 // SCALAR-specific implementations
31 
32 // Anonymous namespace for local helper functions
33 namespace {
34 
35 void scalar_nodal_soln(const Elem * elem,
36  const Order order,
37  const std::vector<Number> & elem_soln,
38  std::vector<Number> & nodal_soln,
39  const bool /*add_p_level*/)
40 {
41  const unsigned int n_nodes = elem->n_nodes();
42  nodal_soln.resize(n_nodes);
43 
44  // If the SCALAR order is CONSTANT, just set the nodal values
45  // to zero, otherwise, set to the value of the first SCALAR dof
46  for (unsigned int i=0; i<n_nodes; i++)
47  nodal_soln[i] = (order == CONSTANT) ? 0. : elem_soln[0];
48 } // scalar_nodal_soln()
49 
50 } // anonymous namespace
51 
52 
53 // Instantiate (side_) nodal_soln() function for every dimension
54 LIBMESH_FE_NODAL_SOLN(SCALAR, scalar_nodal_soln)
56 
57 
58 // Full specialization of n_dofs() function for every dimension
59 // The Order indicates the number of SCALAR dofs
60 template <> unsigned int FE<0,SCALAR>::n_dofs(const ElemType, const Order o) { return o; }
61 template <> unsigned int FE<1,SCALAR>::n_dofs(const ElemType, const Order o) { return o; }
62 template <> unsigned int FE<2,SCALAR>::n_dofs(const ElemType, const Order o) { return o; }
63 template <> unsigned int FE<3,SCALAR>::n_dofs(const ElemType, const Order o) { return o; }
64 
65 template <> unsigned int FE<0,SCALAR>::n_dofs(const Elem *, const Order o) { return o; }
66 template <> unsigned int FE<1,SCALAR>::n_dofs(const Elem *, const Order o) { return o; }
67 template <> unsigned int FE<2,SCALAR>::n_dofs(const Elem *, const Order o) { return o; }
68 template <> unsigned int FE<3,SCALAR>::n_dofs(const Elem *, const Order o) { return o; }
69 
70 // Full specialization of n_dofs_at_node() function for every dimension.
71 // SCALARs have no dofs at nodes
72 template <> unsigned int FE<0,SCALAR>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
73 template <> unsigned int FE<1,SCALAR>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
74 template <> unsigned int FE<2,SCALAR>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
75 template <> unsigned int FE<3,SCALAR>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
76 
77 template <> unsigned int FE<0,SCALAR>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
78 template <> unsigned int FE<1,SCALAR>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
79 template <> unsigned int FE<2,SCALAR>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
80 template <> unsigned int FE<3,SCALAR>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
81 
82 // Full specialization of n_dofs_per_elem() function for every dimension.
83 // SCALARs have no dofs per element
84 template <> unsigned int FE<0,SCALAR>::n_dofs_per_elem(const ElemType, const Order) { return 0; }
85 template <> unsigned int FE<1,SCALAR>::n_dofs_per_elem(const ElemType, const Order) { return 0; }
86 template <> unsigned int FE<2,SCALAR>::n_dofs_per_elem(const ElemType, const Order) { return 0; }
87 template <> unsigned int FE<3,SCALAR>::n_dofs_per_elem(const ElemType, const Order) { return 0; }
88 
89 template <> unsigned int FE<0,SCALAR>::n_dofs_per_elem(const Elem &, const Order) { return 0; }
90 template <> unsigned int FE<1,SCALAR>::n_dofs_per_elem(const Elem &, const Order) { return 0; }
91 template <> unsigned int FE<2,SCALAR>::n_dofs_per_elem(const Elem &, const Order) { return 0; }
92 template <> unsigned int FE<3,SCALAR>::n_dofs_per_elem(const Elem &, const Order) { return 0; }
93 
94 // Scalar FEMs are discontinuous
99 
100 // Scalar FEMs are not hierarchic
101 template <> bool FE<0,SCALAR>::is_hierarchic() const { return false; }
102 template <> bool FE<1,SCALAR>::is_hierarchic() const { return false; }
103 template <> bool FE<2,SCALAR>::is_hierarchic() const { return false; }
104 template <> bool FE<3,SCALAR>::is_hierarchic() const { return false; }
105 
106 
107 #ifdef LIBMESH_ENABLE_AMR
108 // compute_constraints() just returns for SCALAR FEMs
109 template <>
111  DofMap &,
112  const unsigned int,
113  const Elem *)
114 { }
115 
116 template <>
118  DofMap &,
119  const unsigned int,
120  const Elem *)
121 { }
122 #endif // #ifdef LIBMESH_ENABLE_AMR
123 
124 // Scalar FEM shapes do not need reinit
125 template <> bool FE<0,SCALAR>::shapes_need_reinit() const { return false; }
126 template <> bool FE<1,SCALAR>::shapes_need_reinit() const { return false; }
127 template <> bool FE<2,SCALAR>::shapes_need_reinit() const { return false; }
128 template <> bool FE<3,SCALAR>::shapes_need_reinit() const { return false; }
129 
130 } // namespace libMesh
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)
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_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...
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