deal.II version GIT relicensing-2165-gc91f007519 2024-11-20 01:40: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_q_iso_q1.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2013 - 2024 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
17
18#include <deal.II/fe/fe_dgq.h>
21
22#include <deal.II/lac/vector.h>
23
24#include <memory>
25#include <sstream>
26#include <vector>
27
29
30
31
32template <int dim, int spacedim>
33FE_Q_iso_Q1<dim, spacedim>::FE_Q_iso_Q1(const unsigned int subdivisions)
34 : FE_Q_Base<dim, spacedim>(
35 TensorProductPolynomials<dim, Polynomials::PiecewisePolynomial<double>>(
36 Polynomials::generate_complete_Lagrange_basis_on_subdivisions(
37 subdivisions,
38 1)),
39 FiniteElementData<dim>(this->get_dpo_vector(subdivisions),
40 1,
41 subdivisions,
42 FiniteElementData<dim>::H1),
43 std::vector<bool>(1, false))
44{
45 Assert(subdivisions > 0,
46 ExcMessage("This element can only be used with a positive number of "
47 "subelements"));
48
49 const QTrapezoid<1> trapez;
50 const QIterated<1> points(trapez, subdivisions);
51
52 this->initialize(points.get_points());
53}
54
55
56
57template <int dim, int spacedim>
59 const std::vector<Point<1>> &support_points)
60 : FE_Q_Base<dim, spacedim>(
61 TensorProductPolynomials<dim, Polynomials::PiecewisePolynomial<double>>(
62 Polynomials::generate_complete_linear_basis_on_subdivisions(
63 support_points)),
64 FiniteElementData<dim>(this->get_dpo_vector(support_points.size() - 1),
65 1,
66 support_points.size() - 1,
67 FiniteElementData<dim>::H1),
68 std::vector<bool>(1, false))
69{
70 Assert(support_points.size() > 1,
71 ExcMessage("This element can only be used with a positive number of "
72 "subelements"));
73
74 this->initialize(support_points);
75}
76
77
78
79template <int dim, int spacedim>
80std::string
82{
83 // note that the FETools::get_fe_by_name function depends on the
84 // particular format of the string this function returns, so they have to be
85 // kept in sync
86
87 std::ostringstream namebuf;
88 namebuf << "FE_Q_iso_Q1<" << Utilities::dim_string(dim, spacedim) << ">("
89 << this->degree << ")";
90 return namebuf.str();
91}
92
93
94
95template <int dim, int spacedim>
96void
99 const std::vector<Vector<double>> &support_point_values,
100 std::vector<double> &nodal_values) const
101{
102 AssertDimension(support_point_values.size(),
103 this->get_unit_support_points().size());
104 AssertDimension(support_point_values.size(), nodal_values.size());
105 AssertDimension(this->n_dofs_per_cell(), nodal_values.size());
106
107 for (unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
108 {
109 AssertDimension(support_point_values[i].size(), 1);
110
111 nodal_values[i] = support_point_values[i](0);
112 }
113}
114
115
116
117template <int dim, int spacedim>
118std::unique_ptr<FiniteElement<dim, spacedim>>
120{
121 return std::make_unique<FE_Q_iso_Q1<dim, spacedim>>(*this);
122}
123
124
125
126template <int dim, int spacedim>
129 const FiniteElement<dim, spacedim> &fe_other,
130 const unsigned int codim) const
131{
132 Assert(codim <= dim, ExcImpossibleInDim(dim));
133 (void)codim;
134
135 // vertex/line/face domination
136 // (if fe_other is derived from FE_DGQ)
137 // ------------------------------------
138 if (codim > 0)
139 if (dynamic_cast<const FE_DGQ<dim, spacedim> *>(&fe_other) != nullptr)
140 // there are no requirements between continuous and discontinuous elements
142
143 // vertex/line/face domination
144 // (if fe_other is not derived from FE_DGQ)
145 // & cell domination
146 // ----------------------------------------
147 if (const FE_Q_iso_Q1<dim, spacedim> *fe_q_iso_q1_other =
148 dynamic_cast<const FE_Q_iso_Q1<dim, spacedim> *>(&fe_other))
149 {
150 // different behavior as in FE_Q: as FE_Q_iso_Q1(2) is not a subspace of
151 // FE_Q_iso_Q1(3), need that the element degrees are multiples of each
152 // other
153 if (this->degree < fe_q_iso_q1_other->degree &&
154 fe_q_iso_q1_other->degree % this->degree == 0)
156 else if (this->degree == fe_q_iso_q1_other->degree)
158 else if (this->degree > fe_q_iso_q1_other->degree &&
159 this->degree % fe_q_iso_q1_other->degree == 0)
161 else
163 }
164 else if (const FE_Nothing<dim> *fe_nothing =
165 dynamic_cast<const FE_Nothing<dim> *>(&fe_other))
166 {
167 if (fe_nothing->is_dominating())
169 else
170 // the FE_Nothing has no degrees of freedom and it is typically used
171 // in a context where we don't require any continuity along the
172 // interface
174 }
175
178}
179
180
181// explicit instantiations
182#include "fe_q_iso_q1.inst"
183
void initialize(const std::vector< Point< 1 > > &support_points_1d)
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
virtual std::string get_name() const override
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
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
FE_Q_iso_Q1(const unsigned int n_subdivisions)
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define DEAL_II_NOT_IMPLEMENTED()
std::size_t size
Definition mpi.cc:734
std::string dim_string(const int dim, const int spacedim)
Definition utilities.cc:555
STL namespace.