Reference documentation for deal.II version 9.1.0-pre
shape_info.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2011 - 2018 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/aligned_vector.h>
22 #include <deal.II/base/exceptions.h>
23 #include <deal.II/base/quadrature_lib.h>
24 
25 #include <deal.II/fe/fe.h>
26 
27 
28 DEAL_II_NAMESPACE_OPEN
29 
30 namespace internal
31 {
32  namespace MatrixFreeFunctions
33  {
39  enum ElementType
40  {
49  tensor_symmetric_collocation = 0,
55  tensor_symmetric_hermite = 1,
61  tensor_symmetric = 2,
65  tensor_general = 3,
70  truncated_tensor = 4,
76  tensor_symmetric_plus_dg0 = 5
77  };
78 
89  template <typename Number>
90  struct ShapeInfo
91  {
95  ShapeInfo();
96 
100  template <int dim>
101  ShapeInfo(const Quadrature<1> & quad,
102  const FiniteElement<dim> &fe,
103  const unsigned int base_element = 0);
104 
113  template <int dim>
114  void
115  reinit(const Quadrature<1> & quad,
116  const FiniteElement<dim> &fe_dim,
117  const unsigned int base_element = 0);
118 
122  std::size_t
123  memory_consumption() const;
124 
130  ElementType element_type;
131 
140 
149 
158 
165 
172 
179 
186 
193 
200 
206 
212 
218 
226  std::vector<unsigned int> lexicographic_numbering;
227 
231  unsigned int fe_degree;
232 
236  unsigned int n_q_points_1d;
237 
242  unsigned int n_q_points;
243 
249 
253  unsigned int n_q_points_face;
254 
259 
265 
302 
342 
347  bool
348  check_1d_shapes_symmetric(const unsigned int n_q_points_1d);
349 
356  bool
358  };
359 
360 
361 
362  // ------------------------------------------ inline functions
363 
364  template <typename Number>
365  template <int dim>
367  const FiniteElement<dim> &fe_in,
368  const unsigned int base_element_number)
369  : element_type(tensor_general)
370  , fe_degree(0)
371  , n_q_points_1d(0)
372  , n_q_points(0)
374  , n_q_points_face(0)
376  , nodal_at_cell_boundaries(false)
377  {
378  reinit(quad, fe_in, base_element_number);
379  }
380 
381  } // end of namespace MatrixFreeFunctions
382 
383 } // end of namespace internal
384 
385 DEAL_II_NAMESPACE_CLOSE
386 
387 #endif
void reinit(const Quadrature< 1 > &quad, const FiniteElement< dim > &fe_dim, const unsigned int base_element=0)
AlignedVector< Number > shape_values_eo
Definition: shape_info.h:164
AlignedVector< Number > shape_hessians
Definition: shape_info.h:157
AlignedVector< Number > shape_values
Definition: shape_info.h:139
AlignedVector< Number > hessians_within_subface[2]
Definition: shape_info.h:217
bool check_1d_shapes_symmetric(const unsigned int n_q_points_1d)
AlignedVector< Number > shape_gradients_collocation_eo
Definition: shape_info.h:185
AlignedVector< Number > shape_data_on_face[2]
Definition: shape_info.h:199
AlignedVector< Number > shape_hessians_eo
Definition: shape_info.h:178
::Table< 2, unsigned int > face_to_cell_index_hermite
Definition: shape_info.h:341
AlignedVector< Number > gradients_within_subface[2]
Definition: shape_info.h:211
::Table< 2, unsigned int > face_to_cell_index_nodal
Definition: shape_info.h:301
AlignedVector< Number > shape_gradients_eo
Definition: shape_info.h:171
AlignedVector< Number > shape_gradients
Definition: shape_info.h:148
AlignedVector< Number > values_within_subface[2]
Definition: shape_info.h:205
std::vector< unsigned int > lexicographic_numbering
Definition: shape_info.h:226
AlignedVector< Number > shape_hessians_collocation_eo
Definition: shape_info.h:192