Reference documentation for deal.II version 9.5.0
\(\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
fe_abf.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2003 - 2023 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#ifndef dealii_fe_abf_h
17#define dealii_fe_abf_h
18
19#include <deal.II/base/config.h>
20
24#include <deal.II/base/table.h>
26
27#include <deal.II/fe/fe.h>
29
30#include <vector>
31
33
34
100template <int dim>
101class FE_ABF : public FE_PolyTensor<dim>
102{
103public:
107 FE_ABF(const unsigned int p);
108
114 virtual std::string
115 get_name() const override;
116
124 virtual bool
125 has_support_on_face(const unsigned int shape_index,
126 const unsigned int face_index) const override;
127
128 // documentation inherited from the base class
129 virtual void
131 const std::vector<Vector<double>> &support_point_values,
132 std::vector<double> & nodal_values) const override;
133
134 virtual std::size_t
135 memory_consumption() const override;
136
137 virtual std::unique_ptr<FiniteElement<dim, dim>>
138 clone() const override;
139
140private:
146 const unsigned int rt_order;
147
154 static std::vector<unsigned int>
155 get_dpo_vector(const unsigned int degree);
156
166 void
167 initialize_support_points(const unsigned int rt_degree);
168
175 void
177
184 class InternalData : public FiniteElement<dim>::InternalDataBase
185 {
186 public:
198 std::vector<std::vector<Tensor<1, dim>>> shape_values;
199
209 std::vector<std::vector<Tensor<2, dim>>> shape_gradients;
210 };
211
226
227
228
243
251 void
253
254 // Allow access from other dimensions.
255 template <int dim1>
256 friend class FE_ABF;
257};
258
259
260
265
266#endif
std::vector< std::vector< Tensor< 2, dim > > > shape_gradients
Definition fe_abf.h:209
std::vector< std::vector< Tensor< 1, dim > > > shape_values
Definition fe_abf.h:198
void initialize_quad_dof_index_permutation_and_sign_change()
Definition fe_abf.cc:123
Table< 2, double > boundary_weights
Definition fe_abf.h:218
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const override
Definition fe_abf.cc:156
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const override
Definition fe_abf.cc:564
virtual std::size_t memory_consumption() const override
Definition fe_abf.cc:646
void initialize_restriction()
Definition fe_abf.cc:342
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
Definition fe_abf.cc:522
virtual std::string get_name() const override
Definition fe_abf.cc:136
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition fe_abf.cc:489
Table< 2, double > boundary_weights_abf
Definition fe_abf.h:235
const unsigned int rt_order
Definition fe_abf.h:146
void initialize_support_points(const unsigned int rt_degree)
Definition fe_abf.cc:171
Table< 3, double > interior_weights
Definition fe_abf.h:225
friend class FE_ABF
Definition fe_abf.h:256
Table< 3, double > interior_weights_abf
Definition fe_abf.h:242
const unsigned int degree
Definition fe_data.h:453
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:472
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473