Reference documentation for deal.II version 9.4.1
\(\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\}}\)
Loading...
Searching...
No Matches
shape_info.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2011 - 2022 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
16
17#ifndef dealii_matrix_free_shape_info_h
18#define dealii_matrix_free_shape_info_h
19
20
21#include <deal.II/base/config.h>
22
26#include <deal.II/base/table.h>
28
29
31
32
33// forward declaration
34template <int dim, int spacedim>
35class FiniteElement;
36
37
38namespace internal
39{
40 namespace MatrixFreeFunctions
41 {
59 {
69
76
83
88
94
101
108
112 tensor_none = 7
113
114
115 };
116
117
118
127 template <typename Number>
129 {
134
138 std::size_t
140
147
155
163
171
177
183
190
197
204
212
220
236
241
247 std::array<AlignedVector<Number>, 2> shape_data_on_face;
248
259 std::array<AlignedVector<Number>, 2> quadrature_data_on_face;
260
265 std::array<AlignedVector<Number>, 2> values_within_subface;
266
271 std::array<AlignedVector<Number>, 2> gradients_within_subface;
272
277 std::array<AlignedVector<Number>, 2> hessians_within_subface;
278
283 std::array<AlignedVector<Number>, 2> subface_interpolation_matrices;
284
289 std::array<AlignedVector<typename ::internal::VectorizedArrayTrait<
290 Number>::value_type>,
291 2>
293
299
303 unsigned int fe_degree;
304
308 unsigned int n_q_points_1d;
309
315
322
329 };
330
331
332
343 template <typename Number>
345 {
352
357
361 template <int dim, int spacedim, int dim_q>
364 const unsigned int base_element = 0);
365
374 template <int dim, int spacedim, int dim_q>
375 void
377 const FiniteElement<dim, spacedim> &fe_dim,
378 const unsigned int base_element = 0);
379
383 template <int dim, int spacedim>
384 static bool
386
392 compute_orientation_table(const unsigned int n_points_per_dim);
393
401 get_shape_data(const unsigned int dimension = 0,
402 const unsigned int component = 0) const;
403
407 std::size_t
409
417 std::vector<unsigned int> lexicographic_numbering;
418
423 std::vector<UnivariateShapeData<Number>> data;
424
431
435 unsigned int n_dimensions;
436
441 unsigned int n_components;
442
447 unsigned int n_q_points;
448
454
458 unsigned int n_q_points_face;
459
464 std::vector<unsigned int> n_q_points_faces;
465
470
507
547
555
563
564 private:
569 bool
571 UnivariateShapeData<Number> &univariate_shape_data);
572
579 bool
581 const UnivariateShapeData<Number> &univariate_shape_data) const;
582 };
583
584
585
586 // ------------------------------------------ inline functions
587
588 template <typename Number>
589 inline const UnivariateShapeData<Number> &
590 ShapeInfo<Number>::get_shape_data(const unsigned int dimension,
591 const unsigned int component) const
592 {
593 AssertDimension(n_dimensions, data_access.size(0));
594 AssertDimension(n_components, data_access.size(1));
595 AssertIndexRange(dimension, n_dimensions);
596 AssertIndexRange(component, n_components);
597 return *(data_access(dimension, component));
598 }
599
600 } // end of namespace MatrixFreeFunctions
601
602} // end of namespace internal
603
605
606#endif
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1667
#define AssertIndexRange(index, range)
Definition: exceptions.h:1732
static Table< 2, unsigned int > compute_orientation_table(const unsigned int n_points_per_dim)
ShapeInfo(const Quadrature< dim_q > &quad, const FiniteElement< dim, spacedim > &fe, const unsigned int base_element=0)
::Table< 2, unsigned int > face_to_cell_index_nodal
Definition: shape_info.h:506
void reinit(const Quadrature< dim_q > &quad, const FiniteElement< dim, spacedim > &fe_dim, const unsigned int base_element=0)
static bool is_supported(const FiniteElement< dim, spacedim > &fe)
::Table< 2, unsigned int > face_orientations_dofs
Definition: shape_info.h:554
bool check_1d_shapes_collocation(const UnivariateShapeData< Number > &univariate_shape_data) const
std::vector< UnivariateShapeData< Number > > data
Definition: shape_info.h:423
::Table< 2, unsigned int > face_orientations_quad
Definition: shape_info.h:562
const UnivariateShapeData< Number > & get_shape_data(const unsigned int dimension=0, const unsigned int component=0) const
Definition: shape_info.h:590
std::vector< unsigned int > lexicographic_numbering
Definition: shape_info.h:417
::Table< 2, unsigned int > face_to_cell_index_hermite
Definition: shape_info.h:546
std::vector< unsigned int > n_q_points_faces
Definition: shape_info.h:464
::Table< 2, UnivariateShapeData< Number > * > data_access
Definition: shape_info.h:430
bool check_1d_shapes_symmetric(UnivariateShapeData< Number > &univariate_shape_data)
std::array< AlignedVector< Number >, 2 > shape_data_on_face
Definition: shape_info.h:247
std::array< AlignedVector< typename ::internal::VectorizedArrayTrait< Number >::value_type >, 2 > subface_interpolation_matrices_scalar
Definition: shape_info.h:292
AlignedVector< Number > shape_gradients_collocation_eo
Definition: shape_info.h:211
std::array< AlignedVector< Number >, 2 > subface_interpolation_matrices
Definition: shape_info.h:283
std::array< AlignedVector< Number >, 2 > hessians_within_subface
Definition: shape_info.h:277
std::array< AlignedVector< Number >, 2 > values_within_subface
Definition: shape_info.h:265
std::array< AlignedVector< Number >, 2 > quadrature_data_on_face
Definition: shape_info.h:259
AlignedVector< Number > shape_hessians_collocation_eo
Definition: shape_info.h:219
std::array< AlignedVector< Number >, 2 > gradients_within_subface
Definition: shape_info.h:271