1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2020 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/fe.h"
23 #include "libmesh/elem.h"
24 
25 namespace libMesh
26 {
27 
28 // ------------------------------------------------------------
29 // SCALAR-specific implementations
30 
31 // Anonymous namespace for local helper functions
32 namespace {
33 
scalar_nodal_soln(const Elem * elem,const Order order,const std::vector<Number> & elem_soln,std::vector<Number> & nodal_soln)34 void scalar_nodal_soln(const Elem * elem,
35                        const Order order,
36                        const std::vector<Number> & elem_soln,
37                        std::vector<Number> &       nodal_soln)
38 {
39   const unsigned int n_nodes = elem->n_nodes();
40   nodal_soln.resize(n_nodes);
41 
42   // If the SCALAR order is CONSTANT, just set the nodal values
43   // to zero, otherwise, set to the value of the first SCALAR dof
44   for (unsigned int i=0; i<n_nodes; i++)
45     nodal_soln[i] = (order == CONSTANT) ? 0. : elem_soln[0];
46 } // scalar_nodal_soln()
47 
48 } // anonymous namespace
49 
50 
51 
52 
53   // Do full-specialization of nodal_soln() function for every
54   // dimension, instead of explicit instantiation at the end of this
55   // file.
56   // This could be macro-ified so that it fits on one line...
57 template <>
nodal_soln(const Elem * elem,const Order order,const std::vector<Number> & elem_soln,std::vector<Number> & nodal_soln)58 void FE<0,SCALAR>::nodal_soln(const Elem * elem,
59                               const Order order,
60                               const std::vector<Number> & elem_soln,
61                               std::vector<Number> & nodal_soln)
62 { scalar_nodal_soln(elem, order, elem_soln, nodal_soln); }
63 
64 template <>
nodal_soln(const Elem * elem,const Order order,const std::vector<Number> & elem_soln,std::vector<Number> & nodal_soln)65 void FE<1,SCALAR>::nodal_soln(const Elem * elem,
66                               const Order order,
67                               const std::vector<Number> & elem_soln,
68                               std::vector<Number> & nodal_soln)
69 { scalar_nodal_soln(elem, order, elem_soln, nodal_soln); }
70 
71 template <>
nodal_soln(const Elem * elem,const Order order,const std::vector<Number> & elem_soln,std::vector<Number> & nodal_soln)72 void FE<2,SCALAR>::nodal_soln(const Elem * elem,
73                               const Order order,
74                               const std::vector<Number> & elem_soln,
75                               std::vector<Number> & nodal_soln)
76 { scalar_nodal_soln(elem, order, elem_soln, nodal_soln); }
77 
78 template <>
nodal_soln(const Elem * elem,const Order order,const std::vector<Number> & elem_soln,std::vector<Number> & nodal_soln)79 void FE<3,SCALAR>::nodal_soln(const Elem * elem,
80                               const Order order,
81                               const std::vector<Number> & elem_soln,
82                               std::vector<Number> & nodal_soln)
83 { scalar_nodal_soln(elem, order, elem_soln, nodal_soln); }
84 
85 // Full specialization of n_dofs() function for every dimension
86 // The Order indicates the number of SCALAR dofs
n_dofs(const ElemType,const Order o)87 template <> unsigned int FE<0,SCALAR>::n_dofs(const ElemType, const Order o) { return o; }
n_dofs(const ElemType,const Order o)88 template <> unsigned int FE<1,SCALAR>::n_dofs(const ElemType, const Order o) { return o; }
n_dofs(const ElemType,const Order o)89 template <> unsigned int FE<2,SCALAR>::n_dofs(const ElemType, const Order o) { return o; }
n_dofs(const ElemType,const Order o)90 template <> unsigned int FE<3,SCALAR>::n_dofs(const ElemType, const Order o) { return o; }
91 
92 // Full specialization of n_dofs_at_node() function for every dimension.
93 // SCALARs have no dofs at nodes
n_dofs_at_node(const ElemType,const Order,const unsigned int)94 template <> unsigned int FE<0,SCALAR>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
n_dofs_at_node(const ElemType,const Order,const unsigned int)95 template <> unsigned int FE<1,SCALAR>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
n_dofs_at_node(const ElemType,const Order,const unsigned int)96 template <> unsigned int FE<2,SCALAR>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
n_dofs_at_node(const ElemType,const Order,const unsigned int)97 template <> unsigned int FE<3,SCALAR>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
98 
99 // Full specialization of n_dofs_per_elem() function for every dimension.
100 // SCALARs have no dofs per element
n_dofs_per_elem(const ElemType,const Order)101 template <> unsigned int FE<0,SCALAR>::n_dofs_per_elem(const ElemType, const Order) { return 0; }
n_dofs_per_elem(const ElemType,const Order)102 template <> unsigned int FE<1,SCALAR>::n_dofs_per_elem(const ElemType, const Order) { return 0; }
n_dofs_per_elem(const ElemType,const Order)103 template <> unsigned int FE<2,SCALAR>::n_dofs_per_elem(const ElemType, const Order) { return 0; }
n_dofs_per_elem(const ElemType,const Order)104 template <> unsigned int FE<3,SCALAR>::n_dofs_per_elem(const ElemType, const Order) { return 0; }
105 
106 // Scalar FEMs are discontinuous
get_continuity()107 template <> FEContinuity FE<0,SCALAR>::get_continuity() const { return DISCONTINUOUS; }
get_continuity()108 template <> FEContinuity FE<1,SCALAR>::get_continuity() const { return DISCONTINUOUS; }
get_continuity()109 template <> FEContinuity FE<2,SCALAR>::get_continuity() const { return DISCONTINUOUS; }
get_continuity()110 template <> FEContinuity FE<3,SCALAR>::get_continuity() const { return DISCONTINUOUS; }
111 
112 // Scalar FEMs are not hierarchic
is_hierarchic()113 template <> bool FE<0,SCALAR>::is_hierarchic() const { return false; }
is_hierarchic()114 template <> bool FE<1,SCALAR>::is_hierarchic() const { return false; }
is_hierarchic()115 template <> bool FE<2,SCALAR>::is_hierarchic() const { return false; }
is_hierarchic()116 template <> bool FE<3,SCALAR>::is_hierarchic() const { return false; }
117 
118 
119 #ifdef LIBMESH_ENABLE_AMR
120 // compute_constraints() just returns for SCALAR FEMs
121 template <>
compute_constraints(DofConstraints &,DofMap &,const unsigned int,const Elem *)122 void FE<2,SCALAR>::compute_constraints (DofConstraints &,
123                                         DofMap &,
124                                         const unsigned int,
125                                         const Elem *)
126 { }
127 
128 template <>
compute_constraints(DofConstraints &,DofMap &,const unsigned int,const Elem *)129 void FE<3,SCALAR>::compute_constraints (DofConstraints &,
130                                         DofMap &,
131                                         const unsigned int,
132                                         const Elem *)
133 { }
134 #endif // #ifdef LIBMESH_ENABLE_AMR
135 
136 // Scalar FEM shapes do not need reinit
shapes_need_reinit()137 template <> bool FE<0,SCALAR>::shapes_need_reinit() const { return false; }
shapes_need_reinit()138 template <> bool FE<1,SCALAR>::shapes_need_reinit() const { return false; }
shapes_need_reinit()139 template <> bool FE<2,SCALAR>::shapes_need_reinit() const { return false; }
shapes_need_reinit()140 template <> bool FE<3,SCALAR>::shapes_need_reinit() const { return false; }
141 
142 } // namespace libMesh
143