Reference documentation for deal.II version GIT 6a72d26406 2023-06-07 13:05:01+00:00
\(\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 
22 namespace 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 
98  unsigned int
99  number_unique_entries(const std::vector<unsigned int> &vector)
100  {
101  if (std::all_of(vector.begin(), vector.end(), [&](const auto &e) {
102  return e == vector.front();
103  }))
104  {
105  return 1;
106  }
107  else
108  return vector.size();
109  }
110 } // namespace internal
111 
112 
113 
114 template <int dim>
116  const std::vector<unsigned int> &dofs_per_object,
117  const unsigned int n_components,
118  const unsigned int degree,
119  const Conformity conformity,
120  const BlockIndices & block_indices)
121  : FiniteElementData(dofs_per_object,
122  dim == 0 ?
124  (dim == 1 ? ReferenceCells::Line :
125  (dim == 2 ? ReferenceCells::Quadrilateral :
127  n_components,
128  degree,
129  conformity,
130  block_indices)
131 {}
132 
133 template <int dim>
135  const std::vector<unsigned int> &dofs_per_object,
136  const ReferenceCell cell_type,
137  const unsigned int n_components,
138  const unsigned int degree,
139  const Conformity conformity,
140  const BlockIndices & block_indices)
141  : FiniteElementData(internal::expand(dim, dofs_per_object, cell_type),
142  cell_type,
143  n_components,
144  degree,
145  conformity,
146  block_indices)
147 {}
148 
149 
150 
151 template <int dim>
153  const internal::GenericDoFsPerObject &data,
155  const unsigned int n_components,
156  const unsigned int degree,
157  const Conformity conformity,
158  const BlockIndices & block_indices)
159  : reference_cell_kind(reference_cell)
160  , number_unique_quads(
161  internal::number_unique_entries(data.dofs_per_object_inclusive[2]))
162  , number_unique_faces(
163  internal::number_unique_entries(data.dofs_per_object_inclusive[dim - 1]))
164  , dofs_per_vertex(data.dofs_per_object_exclusive[0][0])
165  , dofs_per_line(data.dofs_per_object_exclusive[1][0])
166  , n_dofs_on_quad(data.dofs_per_object_exclusive[2])
167  , dofs_per_quad(n_dofs_on_quad[0])
168  , dofs_per_quad_max(
169  *max_element(n_dofs_on_quad.begin(), n_dofs_on_quad.end()))
170  , dofs_per_hex(data.dofs_per_object_exclusive[3][0])
171  , first_line_index(data.object_index[1][0])
172  , first_index_of_quads(data.object_index[2])
173  , first_quad_index(first_index_of_quads[0])
174  , first_hex_index(data.object_index[3][0])
175  , first_line_index_of_faces(data.first_object_index_on_face[1])
176  , first_face_line_index(first_line_index_of_faces[0])
177  , first_quad_index_of_faces(data.first_object_index_on_face[2])
178  , first_face_quad_index(first_quad_index_of_faces[0])
179  , n_dofs_on_face(data.dofs_per_object_inclusive[dim - 1])
180  , dofs_per_face(n_dofs_on_face[0])
181  , dofs_per_face_max(
182  *max_element(n_dofs_on_face.begin(), n_dofs_on_face.end()))
183  , dofs_per_cell(data.dofs_per_object_inclusive[dim][0])
184  , components(n_components)
185  , degree(degree)
186  , conforming_space(conformity)
187  , block_indices_data(block_indices.size() == 0 ?
188  BlockIndices(1, dofs_per_cell) :
189  block_indices)
190 {}
191 
192 
193 
194 template <int dim>
195 bool
197 {
198  return ((dofs_per_vertex == f.dofs_per_vertex) &&
199  (dofs_per_line == f.dofs_per_line) &&
200  (dofs_per_quad == f.dofs_per_quad) &&
201  (dofs_per_hex == f.dofs_per_hex) && (components == f.components) &&
202  (degree == f.degree) && (conforming_space == f.conforming_space));
203 }
204 
205 
206 template class FiniteElementData<1>;
207 template class FiniteElementData<2>;
208 template class FiniteElementData<3>;
209 
const unsigned int dofs_per_quad
Definition: fe_data.h:343
const unsigned int components
Definition: fe_data.h:443
const unsigned int degree
Definition: fe_data.h:449
bool operator==(const FiniteElementData &) const
Definition: fe_data.cc:196
const unsigned int dofs_per_line
Definition: fe_data.h:329
const unsigned int dofs_per_hex
Definition: fe_data.h:356
const Conformity conforming_space
Definition: fe_data.h:454
const unsigned int dofs_per_vertex
Definition: fe_data.h:323
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:115
unsigned int n_vertices() const
unsigned int n_faces() const
unsigned int n_lines() const
ReferenceCell face_reference_cell(const unsigned int face_no) const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:474
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:475
graph_traits< Graph >::vertex_descriptor Vertex
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
constexpr const ReferenceCell Quadrilateral
constexpr const ReferenceCell Hexahedron
constexpr const ReferenceCell Line
VectorType::value_type * begin(VectorType &V)
VectorType::value_type * end(VectorType &V)
unsigned int number_unique_entries(const std::vector< unsigned int > &vector)
Definition: fe_data.cc:99
internal::GenericDoFsPerObject expand(const unsigned int dim, const std::vector< unsigned int > &dofs_per_object, const ReferenceCell reference_cell)
Definition: fe_data.cc:25
std::vector< std::vector< unsigned int > > object_index
Definition: fe_data.h:195
std::vector< std::vector< unsigned int > > first_object_index_on_face
Definition: fe_data.h:200
std::vector< std::vector< unsigned int > > dofs_per_object_inclusive
Definition: fe_data.h:190
std::vector< std::vector< unsigned int > > dofs_per_object_exclusive
Definition: fe_data.h:185