Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
fe_data.cc
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2001 - 2021 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
17
18#include <deal.II/fe/fe.h>
19
21
22namespace internal
23{
25 expand(const unsigned int dim,
26 const std::vector<unsigned int> &dofs_per_object,
27 const ::ReferenceCell cell_type)
28 {
30
31 const unsigned int face_no = 0;
32
33 result.dofs_per_object_exclusive.resize(4, std::vector<unsigned int>(1));
34 result.dofs_per_object_inclusive.resize(4, std::vector<unsigned int>(1));
35 result.object_index.resize(4, std::vector<unsigned int>(1));
36 result.first_object_index_on_face.resize(3, std::vector<unsigned int>(1));
37
38 // dofs_per_vertex
39 const unsigned int dofs_per_vertex = dofs_per_object[0];
40 result.dofs_per_object_exclusive[0][0] = dofs_per_vertex;
41
42 // dofs_per_line
43 const unsigned int dofs_per_line = dofs_per_object[1];
44 result.dofs_per_object_exclusive[1][0] = dofs_per_line;
45
46 // dofs_per_quad
47 const unsigned int dofs_per_quad = dim > 1 ? dofs_per_object[2] : 0;
48 result.dofs_per_object_exclusive[2][0] = dofs_per_quad;
49
50 // dofs_per_hex
51 const unsigned int dofs_per_hex = dim > 2 ? dofs_per_object[3] : 0;
52 result.dofs_per_object_exclusive[3][0] = dofs_per_hex;
53
54
55 // first_line_index
56 const unsigned int first_line_index =
57 (cell_type.n_vertices() * dofs_per_vertex);
58 result.object_index[1][0] = first_line_index;
59
60 // first_quad_index
61 const unsigned int first_quad_index =
62 (first_line_index + cell_type.n_lines() * dofs_per_line);
63 result.object_index[2][0] = first_quad_index;
64
65 // first_hex_index
66 result.object_index[3][0] =
67 (first_quad_index +
68 (dim == 2 ? 1 : (dim == 3 ? cell_type.n_faces() : 0)) * dofs_per_quad);
69
70 // first_face_line_index
71 result.first_object_index_on_face[1][0] =
72 (cell_type.face_reference_cell(face_no).n_vertices() * dofs_per_vertex);
73
74 // first_face_quad_index
75 result.first_object_index_on_face[2][0] =
76 ((dim == 3 ? cell_type.face_reference_cell(face_no).n_vertices() *
77 dofs_per_vertex :
78 cell_type.n_vertices() * dofs_per_vertex) +
79 cell_type.face_reference_cell(face_no).n_lines() * dofs_per_line);
80
81 // dofs_per_face
82 result.dofs_per_object_inclusive[dim - 1][0] =
83 (cell_type.face_reference_cell(face_no).n_vertices() * dofs_per_vertex +
84 cell_type.face_reference_cell(face_no).n_lines() * dofs_per_line +
85 (dim == 3 ? 1 : 0) * dofs_per_quad);
86
87
88 // dofs_per_cell
89 result.dofs_per_object_inclusive[dim][0] =
90 (cell_type.n_vertices() * dofs_per_vertex +
91 cell_type.n_lines() * dofs_per_line +
92 (dim == 2 ? 1 : (dim == 3 ? cell_type.n_faces() : 0)) * dofs_per_quad +
93 (dim == 3 ? 1 : 0) * dofs_per_hex);
94
95 return result;
96 }
97} // namespace internal
98
99
100
101template <int dim>
103 const std::vector<unsigned int> &dofs_per_object,
104 const unsigned int n_components,
105 const unsigned int degree,
106 const Conformity conformity,
107 const BlockIndices & block_indices)
108 : FiniteElementData(dofs_per_object,
109 dim == 0 ?
111 (dim == 1 ? ReferenceCells::Line :
112 (dim == 2 ? ReferenceCells::Quadrilateral :
114 n_components,
115 degree,
116 conformity,
117 block_indices)
118{}
119
120template <int dim>
122 const std::vector<unsigned int> &dofs_per_object,
123 const ReferenceCell cell_type,
124 const unsigned int n_components,
125 const unsigned int degree,
126 const Conformity conformity,
127 const BlockIndices & block_indices)
128 : FiniteElementData(internal::expand(dim, dofs_per_object, cell_type),
129 cell_type,
130 n_components,
131 degree,
132 conformity,
133 block_indices)
134{}
135
136template <int dim>
140 const unsigned int n_components,
141 const unsigned int degree,
142 const Conformity conformity,
143 const BlockIndices & block_indices)
144 : reference_cell_kind(reference_cell)
145 , number_unique_quads(data.dofs_per_object_inclusive[2].size())
146 , number_unique_faces(data.dofs_per_object_inclusive[dim - 1].size())
147 , dofs_per_vertex(data.dofs_per_object_exclusive[0][0])
148 , dofs_per_line(data.dofs_per_object_exclusive[1][0])
149 , n_dofs_on_quad(dim > 1 ? data.dofs_per_object_exclusive[2] :
150 std::vector<unsigned int>{0})
151 , dofs_per_quad(n_dofs_on_quad[0])
152 , dofs_per_quad_max(
153 *max_element(n_dofs_on_quad.begin(), n_dofs_on_quad.end()))
154 , dofs_per_hex(dim > 2 ? data.dofs_per_object_exclusive[3][0] : 0)
155 , first_line_index(data.object_index[1][0])
156 , first_index_of_quads(data.object_index[2])
157 , first_quad_index(first_index_of_quads[0])
158 , first_hex_index(data.object_index[3][0])
159 , first_line_index_of_faces(data.first_object_index_on_face[1])
160 , first_face_line_index(first_line_index_of_faces[0])
161 , first_quad_index_of_faces(data.first_object_index_on_face[2])
162 , first_face_quad_index(first_quad_index_of_faces[0])
163 , n_dofs_on_face(data.dofs_per_object_inclusive[dim - 1])
164 , dofs_per_face(n_dofs_on_face[0])
165 , dofs_per_face_max(
166 *max_element(n_dofs_on_face.begin(), n_dofs_on_face.end()))
167 , dofs_per_cell(data.dofs_per_object_inclusive[dim][0])
168 , components(n_components)
169 , degree(degree)
170 , conforming_space(conformity)
171 , block_indices_data(block_indices.size() == 0 ?
172 BlockIndices(1, dofs_per_cell) :
173 block_indices)
174{}
175
176
177
178template <int dim>
179bool
181{
182 return ((dofs_per_vertex == f.dofs_per_vertex) &&
183 (dofs_per_line == f.dofs_per_line) &&
184 (dofs_per_quad == f.dofs_per_quad) &&
185 (dofs_per_hex == f.dofs_per_hex) && (components == f.components) &&
186 (degree == f.degree) && (conforming_space == f.conforming_space));
187}
188
189
190template class FiniteElementData<1>;
191template class FiniteElementData<2>;
192template class FiniteElementData<3>;
193
const unsigned int dofs_per_quad
Definition: fe_base.h:329
const unsigned int components
Definition: fe_base.h:429
const unsigned int degree
Definition: fe_base.h:435
bool operator==(const FiniteElementData &) const
Definition: fe_data.cc:180
const unsigned int dofs_per_line
Definition: fe_base.h:315
const unsigned int dofs_per_hex
Definition: fe_base.h:342
const Conformity conforming_space
Definition: fe_base.h:440
const unsigned int dofs_per_vertex
Definition: fe_base.h:309
FiniteElementData(const std::vector< unsigned int > &dofs_per_object, const unsigned int n_components, const unsigned int degree, const Conformity conformity=unknown, const BlockIndices &block_indices=BlockIndices())
Definition: fe_data.cc:102
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
graph_traits< Graph >::vertex_descriptor Vertex
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
constexpr const ReferenceCell Quadrilateral
constexpr const ReferenceCell Hexahedron
constexpr const ReferenceCell Line
internal::GenericDoFsPerObject expand(const unsigned int dim, const std::vector< unsigned int > &dofs_per_object, const ::ReferenceCell reference_cell)
Definition: fe_data.cc:25
STL namespace.
std::vector< std::vector< unsigned int > > object_index
Definition: fe_base.h:188
std::vector< std::vector< unsigned int > > first_object_index_on_face
Definition: fe_base.h:193
std::vector< std::vector< unsigned int > > dofs_per_object_inclusive
Definition: fe_base.h:183
std::vector< std::vector< unsigned int > > dofs_per_object_exclusive
Definition: fe_base.h:178