Reference documentation for deal.II version 9.1.0-pre
fe_q.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2017 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 #include <deal.II/base/quadrature_lib.h>
18 #include <deal.II/base/std_cxx14/memory.h>
19 
20 #include <deal.II/fe/fe_q.h>
21 
22 #include <deal.II/lac/vector.h>
23 
24 #include <sstream>
25 #include <vector>
26 
27 DEAL_II_NAMESPACE_OPEN
28 
29 
30 namespace internal
31 {
32  namespace FE_Q
33  {
34  namespace
35  {
36  std::vector<Point<1>>
37  get_QGaussLobatto_points(const unsigned int degree)
38  {
39  if (degree > 0)
40  return QGaussLobatto<1>(degree + 1).get_points();
41  else
42  {
43  using FEQ = ::FE_Q_Base<TensorProductPolynomials<1>, 1, 1>;
44  AssertThrow(false, FEQ::ExcFEQCannotHaveDegree0());
45  }
46  return std::vector<Point<1>>();
47  }
48  } // namespace
49  } // namespace FE_Q
50 } // namespace internal
51 
52 
53 
54 template <int dim, int spacedim>
55 FE_Q<dim, spacedim>::FE_Q(const unsigned int degree)
56  : FE_Q_Base<TensorProductPolynomials<dim>, dim, spacedim>(
58  Polynomials::generate_complete_Lagrange_basis(
59  internal::FE_Q::get_QGaussLobatto_points(degree))),
60  FiniteElementData<dim>(this->get_dpo_vector(degree),
61  1,
62  degree,
63  FiniteElementData<dim>::H1),
64  std::vector<bool>(1, false))
65 {
66  this->initialize(internal::FE_Q::get_QGaussLobatto_points(degree));
67 }
68 
69 
70 
71 template <int dim, int spacedim>
73  : FE_Q_Base<TensorProductPolynomials<dim>, dim, spacedim>(
75  Polynomials::generate_complete_Lagrange_basis(points.get_points())),
76  FiniteElementData<dim>(this->get_dpo_vector(points.size() - 1),
77  1,
78  points.size() - 1,
79  FiniteElementData<dim>::H1),
80  std::vector<bool>(1, false))
81 {
82  this->initialize(points.get_points());
83 }
84 
85 
86 
87 template <int dim, int spacedim>
88 std::string
90 {
91  // note that the FETools::get_fe_by_name function depends on the
92  // particular format of the string this function returns, so they have to be
93  // kept in synch
94 
95  std::ostringstream namebuf;
96  bool equidistant = true;
97  std::vector<double> points(this->degree + 1);
98 
99  // Decode the support points in one coordinate direction.
100  std::vector<unsigned int> lexicographic =
102  for (unsigned int j = 0; j <= this->degree; j++)
103  points[j] = this->unit_support_points[lexicographic[j]][0];
104 
105  // Check whether the support points are equidistant.
106  for (unsigned int j = 0; j <= this->degree; j++)
107  if (std::fabs(points[j] - (double)j / this->degree) > 1e-15)
108  {
109  equidistant = false;
110  break;
111  }
112 
113  if (equidistant == true)
114  {
115  if (this->degree > 2)
116  namebuf << "FE_Q<" << Utilities::dim_string(dim, spacedim)
117  << ">(QIterated(QTrapez()," << this->degree << "))";
118  else
119  namebuf << "FE_Q<" << Utilities::dim_string(dim, spacedim) << ">("
120  << this->degree << ")";
121  }
122  else
123  {
124  // Check whether the support points come from QGaussLobatto.
125  const QGaussLobatto<1> points_gl(this->degree + 1);
126  bool gauss_lobatto = true;
127  for (unsigned int j = 0; j <= this->degree; j++)
128  if (points[j] != points_gl.point(j)(0))
129  {
130  gauss_lobatto = false;
131  break;
132  }
133  if (gauss_lobatto == true)
134  namebuf << "FE_Q<" << Utilities::dim_string(dim, spacedim) << ">("
135  << this->degree << ")";
136  else
137  namebuf << "FE_Q<" << Utilities::dim_string(dim, spacedim)
138  << ">(QUnknownNodes(" << this->degree << "))";
139  }
140  return namebuf.str();
141 }
142 
143 
144 
145 template <int dim, int spacedim>
146 void
148  const std::vector<Vector<double>> &support_point_values,
149  std::vector<double> & nodal_values) const
150 {
151  AssertDimension(support_point_values.size(),
152  this->get_unit_support_points().size());
153  AssertDimension(support_point_values.size(), nodal_values.size());
154  AssertDimension(this->dofs_per_cell, nodal_values.size());
155 
156  for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
157  {
158  AssertDimension(support_point_values[i].size(), 1);
159 
160  nodal_values[i] = support_point_values[i](0);
161  }
162 }
163 
164 
165 
166 template <int dim, int spacedim>
167 std::unique_ptr<FiniteElement<dim, spacedim>>
169 {
170  return std_cxx14::make_unique<FE_Q<dim, spacedim>>(*this);
171 }
172 
173 
174 // explicit instantiations
175 #include "fe_q.inst"
176 
177 DEAL_II_NAMESPACE_CLOSE
const std::vector< Point< dim > > & get_unit_support_points() const
Definition: fe.cc:1000
virtual std::string get_name() const override
Definition: fe_q.cc:89
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1366
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
Definition: fe_q.cc:168
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
FE_Q(const unsigned int p)
Definition: fe_q.cc:55
const unsigned int degree
Definition: fe_base.h:313
STL namespace.
#define AssertThrow(cond, exc)
Definition: exceptions.h:1329
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_q.cc:147
const std::vector< unsigned int > & get_numbering_inverse() const
const Point< dim > & point(const unsigned int i) const
std::vector< Point< dim > > unit_support_points
Definition: fe.h:2457
const std::vector< Point< dim > > & get_points() const
Definition: fe_q.h:554
std::string dim_string(const int dim, const int spacedim)
Definition: utilities.cc:170
const unsigned int dofs_per_cell
Definition: fe_base.h:297
void initialize(const std::vector< Point< 1 >> &support_points_1d)