Reference documentation for deal.II version 9.1.0-pre
fe_q_dg0.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2012 - 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/qprojector.h>
18 #include <deal.II/base/quadrature.h>
19 #include <deal.II/base/quadrature_lib.h>
20 #include <deal.II/base/std_cxx14/memory.h>
21 #include <deal.II/base/template_constraints.h>
22 
23 #include <deal.II/dofs/dof_accessor.h>
24 
25 #include <deal.II/fe/fe_nothing.h>
26 #include <deal.II/fe/fe_q_dg0.h>
27 
28 #include <sstream>
29 #include <vector>
30 
31 DEAL_II_NAMESPACE_OPEN
32 
33 
34 template <int dim, int spacedim>
35 FE_Q_DG0<dim, spacedim>::FE_Q_DG0(const unsigned int degree)
36  : FE_Q_Base<TensorProductPolynomialsConst<dim>, dim, spacedim>(
38  Polynomials::generate_complete_Lagrange_basis(
39  QGaussLobatto<1>(degree + 1).get_points())),
40  FiniteElementData<dim>(get_dpo_vector(degree),
41  1,
42  degree,
43  FiniteElementData<dim>::L2),
44  get_riaf_vector(degree))
45 {
46  Assert(degree > 0,
47  ExcMessage("This element can only be used for polynomial degrees "
48  "greater than zero"));
49 
50  this->initialize(QGaussLobatto<1>(degree + 1).get_points());
51 
52  // adjust unit support point for discontinuous node
53  Point<dim> point;
54  for (unsigned int d = 0; d < dim; ++d)
55  point[d] = 0.5;
56  this->unit_support_points.push_back(point);
58 }
59 
60 
61 
62 template <int dim, int spacedim>
64  : FE_Q_Base<TensorProductPolynomialsConst<dim>, dim, spacedim>(
66  Polynomials::generate_complete_Lagrange_basis(points.get_points())),
67  FiniteElementData<dim>(get_dpo_vector(points.size() - 1),
68  1,
69  points.size() - 1,
70  FiniteElementData<dim>::L2),
71  get_riaf_vector(points.size() - 1))
72 {
73  const int degree = points.size() - 1;
74  (void)degree;
75 
76  Assert(degree > 0,
77  ExcMessage("This element can only be used for polynomial degrees "
78  "at least zero"));
79 
80  this->initialize(points.get_points());
81 
82  // adjust unit support point for discontinuous node
83  Point<dim> point;
84  for (unsigned int d = 0; d < dim; ++d)
85  point[d] = 0.5;
86  this->unit_support_points.push_back(point);
88 }
89 
90 
91 
92 template <int dim, int spacedim>
93 std::string
95 {
96  // note that the FETools::get_fe_by_name function depends on the
97  // particular format of the string this function returns, so they have to be
98  // kept in synch
99 
100  std::ostringstream namebuf;
101  bool type = true;
102  const unsigned int n_points = this->degree + 1;
103  std::vector<double> points(n_points);
104  const unsigned int dofs_per_cell = this->dofs_per_cell;
105  const std::vector<Point<dim>> &unit_support_points =
106  this->unit_support_points;
107  unsigned int index = 0;
108 
109  // Decode the support points in one coordinate direction.
110  for (unsigned int j = 0; j < dofs_per_cell; j++)
111  {
112  if ((dim > 1) ? (unit_support_points[j](1) == 0 &&
113  ((dim > 2) ? unit_support_points[j](2) == 0 : true)) :
114  true)
115  {
116  if (index == 0)
117  points[index] = unit_support_points[j](0);
118  else if (index == 1)
119  points[n_points - 1] = unit_support_points[j](0);
120  else
121  points[index - 1] = unit_support_points[j](0);
122 
123  index++;
124  }
125  }
126  // Do not consider the discontinuous node for dimension 1
127  Assert(index == n_points || (dim == 1 && index == n_points + 1),
128  ExcMessage(
129  "Could not decode support points in one coordinate direction."));
130 
131  // Check whether the support points are equidistant.
132  for (unsigned int j = 0; j < n_points; j++)
133  if (std::fabs(points[j] - (double)j / this->degree) > 1e-15)
134  {
135  type = false;
136  break;
137  }
138 
139  if (type == true)
140  {
141  if (this->degree > 2)
142  namebuf << "FE_Q_DG0<" << Utilities::dim_string(dim, spacedim)
143  << ">(QIterated(QTrapez()," << this->degree << "))";
144  else
145  namebuf << "FE_Q_DG0<" << Utilities::dim_string(dim, spacedim) << ">("
146  << this->degree << ")";
147  }
148  else
149  {
150  // Check whether the support points come from QGaussLobatto.
151  const QGaussLobatto<1> points_gl(n_points);
152  type = true;
153  for (unsigned int j = 0; j < n_points; j++)
154  if (points[j] != points_gl.point(j)(0))
155  {
156  type = false;
157  break;
158  }
159  if (type == true)
160  namebuf << "FE_Q_DG0<" << Utilities::dim_string(dim, spacedim) << ">("
161  << this->degree << ")";
162  else
163  namebuf << "FE_Q_DG0<" << Utilities::dim_string(dim, spacedim)
164  << ">(QUnknownNodes(" << this->degree << "))";
165  }
166  return namebuf.str();
167 }
168 
169 
170 
171 template <int dim, int spacedim>
172 std::unique_ptr<FiniteElement<dim, spacedim>>
174 {
175  return std_cxx14::make_unique<FE_Q_DG0<dim, spacedim>>(*this);
176 }
177 
178 
179 
180 template <int dim, int spacedim>
181 void
183  const std::vector<Vector<double>> &support_point_values,
184  std::vector<double> & nodal_dofs) const
185 {
186  Assert(support_point_values.size() == this->unit_support_points.size(),
187  ExcDimensionMismatch(support_point_values.size(),
188  this->unit_support_points.size()));
189  Assert(nodal_dofs.size() == this->dofs_per_cell,
190  ExcDimensionMismatch(nodal_dofs.size(), this->dofs_per_cell));
191  Assert(support_point_values[0].size() == this->n_components(),
192  ExcDimensionMismatch(support_point_values[0].size(),
193  this->n_components()));
194 
195  for (unsigned int i = 0; i < this->dofs_per_cell - 1; ++i)
196  {
197  const std::pair<unsigned int, unsigned int> index =
198  this->system_to_component_index(i);
199  nodal_dofs[i] = support_point_values[i](index.first);
200  }
201 
202  // We don't need the discontinuous function for local interpolation
203  nodal_dofs[nodal_dofs.size() - 1] = 0.;
204 }
205 
206 
207 
208 template <int dim, int spacedim>
209 void
211  const FiniteElement<dim, spacedim> &x_source_fe,
212  FullMatrix<double> & interpolation_matrix) const
213 {
214  // this is only implemented, if the source FE is also a Q_DG0 element
215  using FEQDG0 = FE_Q_DG0<dim, spacedim>;
216 
217  AssertThrow(
218  (x_source_fe.get_name().find("FE_Q_DG0<") == 0) ||
219  (dynamic_cast<const FEQDG0 *>(&x_source_fe) != nullptr),
221 
222  Assert(interpolation_matrix.m() == this->dofs_per_cell,
223  ExcDimensionMismatch(interpolation_matrix.m(), this->dofs_per_cell));
224  Assert(interpolation_matrix.n() == x_source_fe.dofs_per_cell,
225  ExcDimensionMismatch(interpolation_matrix.m(),
226  x_source_fe.dofs_per_cell));
227 
229  get_interpolation_matrix(x_source_fe, interpolation_matrix);
230 }
231 
232 
233 
234 template <int dim, int spacedim>
235 std::vector<bool>
237 {
238  std::vector<bool> riaf(Utilities::fixed_power<dim>(deg + 1) + 1, false);
239  riaf[riaf.size() - 1] = true;
240  return riaf;
241 }
242 
243 
244 
245 template <int dim, int spacedim>
246 std::vector<unsigned int>
248 {
249  std::vector<unsigned int> dpo(dim + 1, 1U);
250  for (unsigned int i = 1; i < dpo.size(); ++i)
251  dpo[i] = dpo[i - 1] * (deg - 1);
252 
253  dpo[dim]++; // we need an additional DG0-node for a dim-dimensional object
254  return dpo;
255 }
256 
257 
258 
259 template <int dim, int spacedim>
260 bool
262  const unsigned int shape_index,
263  const unsigned int face_index) const
264 {
265  // discontinuous function has support on all faces
266  if (shape_index == this->dofs_per_cell - 1)
267  return true;
268  else
270  has_support_on_face(shape_index, face_index);
271 }
272 
273 
274 
275 template <int dim, int spacedim>
276 std::pair<Table<2, bool>, std::vector<unsigned int>>
278 {
279  Table<2, bool> constant_modes(2, this->dofs_per_cell);
280 
281  // 1 represented by FE_Q part
282  for (unsigned int i = 0; i < this->dofs_per_cell - 1; ++i)
283  constant_modes(0, i) = true;
284 
285  // 1 represented by DG0 part
286  constant_modes(1, this->dofs_per_cell - 1) = true;
287 
288  return std::pair<Table<2, bool>, std::vector<unsigned int>>(
289  constant_modes, std::vector<unsigned int>(2, 0));
290 }
291 
292 
293 
294 // explicit instantiations
295 #include "fe_q_dg0.inst"
296 
297 DEAL_II_NAMESPACE_CLOSE
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1366
const unsigned int degree
Definition: fe_base.h:313
#define AssertThrow(cond, exc)
Definition: exceptions.h:1329
const Point< dim > & point(const unsigned int i) const
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
Definition: fe_q_dg0.cc:277
std::vector< Point< dim > > unit_support_points
Definition: fe.h:2457
const std::vector< Point< dim > > & get_points() const
unsigned int size() const
static::ExceptionBase & ExcMessage(std::string arg1)
size_type n() const
#define Assert(cond, exc)
Definition: exceptions.h:1227
static::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
virtual std::string get_name() const =0
unsigned int n_components() const
size_type m() const
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
static std::vector< bool > get_riaf_vector(const unsigned int degree)
Definition: fe_q_dg0.cc:236
void initialize(const std::vector< Point< 1 >> &support_points_1d)
FE_Q_DG0(const unsigned int p)
Definition: fe_q_dg0.cc:35
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
Definition: fe.h:3087
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
Definition: fe_q_dg0.cc:173
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
Definition: fe_q_dg0.cc:261
Definition: table.h:37
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const override
Definition: fe_q_dg0.cc:210
virtual std::string get_name() const override
Definition: fe_q_dg0.cc:94
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_dg0.cc:182
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition: fe_q_dg0.cc:247