deal.II version GIT relicensing-2169-gec1b43f35b 2024-11-22 07:10:00+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\}}\)
Loading...
Searching...
No Matches
fe_abf.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2003 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_fe_abf_h
16#define dealii_fe_abf_h
17
18#include <deal.II/base/config.h>
19
23#include <deal.II/base/table.h>
25
26#include <deal.II/fe/fe.h>
28
29#include <vector>
30
32
33
99template <int dim>
100class FE_ABF : public FE_PolyTensor<dim>
101{
102public:
106 FE_ABF(const unsigned int p);
107
113 virtual std::string
114 get_name() const override;
115
123 virtual bool
124 has_support_on_face(const unsigned int shape_index,
125 const unsigned int face_index) const override;
126
127 // documentation inherited from the base class
128 virtual void
130 const std::vector<Vector<double>> &support_point_values,
131 std::vector<double> &nodal_values) const override;
132
133 virtual std::size_t
134 memory_consumption() const override;
135
136 virtual std::unique_ptr<FiniteElement<dim, dim>>
137 clone() const override;
138
139private:
145 const unsigned int rt_order;
146
153 static std::vector<unsigned int>
154 get_dpo_vector(const unsigned int degree);
155
165 void
166 initialize_support_points(const unsigned int rt_degree);
167
174 void
176
183 class InternalData : public FiniteElement<dim>::InternalDataBase
184 {
185 public:
197 std::vector<std::vector<Tensor<1, dim>>> shape_values;
198
208 std::vector<std::vector<Tensor<2, dim>>> shape_gradients;
209 };
210
225
226
227
242
250 void
252
253 // Allow access from other dimensions.
254 template <int dim1>
255 friend class FE_ABF;
256};
257
258
259
264
265#endif
std::vector< std::vector< Tensor< 2, dim > > > shape_gradients
Definition fe_abf.h:208
std::vector< std::vector< Tensor< 1, dim > > > shape_values
Definition fe_abf.h:197
void initialize_quad_dof_index_permutation_and_sign_change()
Definition fe_abf.cc:122
Table< 2, double > boundary_weights
Definition fe_abf.h:217
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const override
Definition fe_abf.cc:155
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:563
virtual std::size_t memory_consumption() const override
Definition fe_abf.cc:648
void initialize_restriction()
Definition fe_abf.cc:341
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
Definition fe_abf.cc:521
virtual std::string get_name() const override
Definition fe_abf.cc:135
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition fe_abf.cc:488
Table< 2, double > boundary_weights_abf
Definition fe_abf.h:234
const unsigned int rt_order
Definition fe_abf.h:145
void initialize_support_points(const unsigned int rt_degree)
Definition fe_abf.cc:170
Table< 3, double > interior_weights
Definition fe_abf.h:224
friend class FE_ABF
Definition fe_abf.h:255
Table< 3, double > interior_weights_abf
Definition fe_abf.h:241
const unsigned int degree
Definition fe_data.h:452
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499