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