Reference documentation for deal.II version 8.4.2
fe_q_dg0.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2012 - 2015 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 
29 DEAL_II_NAMESPACE_OPEN
30 
31 
32 template <int dim, int spacedim>
33 FE_Q_DG0<dim,spacedim>::FE_Q_DG0 (const unsigned int degree)
34  :
35  FE_Q_Base<TensorProductPolynomialsConst<dim>, dim, spacedim> (
36  TensorProductPolynomialsConst<dim>(Polynomials::LagrangeEquidistant::generate_complete_basis(degree)),
37  FiniteElementData<dim>(get_dpo_vector(degree),
38  1, degree,
39  FiniteElementData<dim>::L2),
40  get_riaf_vector(degree))
41 {
42  Assert (degree > 0,
43  ExcMessage ("This element can only be used for polynomial degrees "
44  "greater than zero"));
45 
46  std::vector<Point<1> > support_points_1d(degree+1);
47  for (unsigned int i=0; i<=degree; ++i)
48  support_points_1d[i][0] = static_cast<double>(i)/degree;
49 
50  this->initialize(support_points_1d);
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  :
65  FE_Q_Base<TensorProductPolynomialsConst<dim>, dim, spacedim> (
66  TensorProductPolynomialsConst<dim>(Polynomials::generate_complete_Lagrange_basis(points.get_points())),
67  FiniteElementData<dim>(get_dpo_vector(points.size()-1),
68  1, points.size()-1,
69  FiniteElementData<dim>::L2),
70  get_riaf_vector(points.size()-1))
71 {
72  const int degree = points.size()-1;
73  (void)degree;
74 
75  Assert (degree > 0,
76  ExcMessage ("This element can only be used for polynomial degrees "
77  "at least zero"));
78 
79  this->initialize(points.get_points());
80 
81  // adjust unit support point for discontinuous node
82  Point<dim> point;
83  for (unsigned int d=0; d<dim; ++d)
84  point[d] = 0.5;
85  this->unit_support_points.push_back(point);
87 }
88 
89 
90 
91 template <int dim, int spacedim>
92 std::string
94 {
95  // note that the FETools::get_fe_from_name function depends on the
96  // particular format of the string this function returns, so they have to be
97  // kept in synch
98 
99  std::ostringstream namebuf;
100  bool type = true;
101  const unsigned int n_points = this->degree +1;
102  std::vector<double> points(n_points);
103  const unsigned int dofs_per_cell = this->dofs_per_cell;
104  const std::vector<Point<dim> > &unit_support_points = this->unit_support_points;
105  unsigned int index = 0;
106 
107  // Decode the support points in one coordinate direction.
108  for (unsigned int j=0; j<dofs_per_cell; j++)
109  {
110  if ((dim>1) ? (unit_support_points[j](1)==0 &&
111  ((dim>2) ? unit_support_points[j](2)==0: true)) : true)
112  {
113  if (index == 0)
114  points[index] = unit_support_points[j](0);
115  else if (index == 1)
116  points[n_points-1] = unit_support_points[j](0);
117  else
118  points[index-1] = unit_support_points[j](0);
119 
120  index++;
121  }
122  }
123  // Do not consider the discontinuous node for dimension 1
124  Assert (index == n_points || (dim==1 && index == n_points+1),
125  ExcMessage ("Could not decode support points in one coordinate direction."));
126 
127  // Check whether the support points are equidistant.
128  for (unsigned int j=0; j<n_points; j++)
129  if (std::fabs(points[j] - (double)j/this->degree) > 1e-15)
130  {
131  type = false;
132  break;
133  }
134 
135  if (type == true)
136  namebuf << "FE_Q_DG0<"
137  << Utilities::dim_string(dim,spacedim)
138  << ">(" << this->degree << ")";
139  else
140  {
141 
142  // Check whether the support points come from QGaussLobatto.
143  const QGaussLobatto<1> points_gl(n_points);
144  type = true;
145  for (unsigned int j=0; j<n_points; j++)
146  if (points[j] != points_gl.point(j)(0))
147  {
148  type = false;
149  break;
150  }
151  if (type == true)
152  namebuf << "FE_Q_DG0<"
153  << Utilities::dim_string(dim,spacedim)
154  << ">(QGaussLobatto(" << this->degree+1 << "))";
155  else
156  namebuf << "FE_Q_DG0<"
157  << Utilities::dim_string(dim,spacedim)
158  << ">(QUnknownNodes(" << this->degree << "))";
159  }
160  return namebuf.str();
161 }
162 
163 
164 
165 template <int dim, int spacedim>
168 {
169  return new FE_Q_DG0<dim,spacedim>(*this);
170 }
171 
172 
173 
174 template <int dim, int spacedim>
175 void
176 FE_Q_DG0<dim,spacedim>::interpolate(std::vector<double> &local_dofs,
177  const std::vector<double> &values) const
178 {
179  Assert (values.size() == this->unit_support_points.size(),
180  ExcDimensionMismatch(values.size(),
181  this->unit_support_points.size()));
182  Assert (local_dofs.size() == this->dofs_per_cell,
183  ExcDimensionMismatch(local_dofs.size(),this->dofs_per_cell));
184  Assert (this->n_components() == 1,
185  ExcDimensionMismatch(this->n_components(), 1));
186 
187  std::copy(values.begin(), values.end(), local_dofs.begin());
188 
189  // We don't need the discontinuous function for local interpolation
190  local_dofs[local_dofs.size()-1] = 0.;
191 }
192 
193 
194 
195 template <int dim, int spacedim>
196 void
197 FE_Q_DG0<dim,spacedim>::interpolate(std::vector<double> &local_dofs,
198  const std::vector<Vector<double> > &values,
199  unsigned int offset) const
200 {
201  Assert (values.size() == this->unit_support_points.size(),
202  ExcDimensionMismatch(values.size(),
203  this->unit_support_points.size()));
204  Assert (local_dofs.size() == this->dofs_per_cell,
205  ExcDimensionMismatch(local_dofs.size(),this->dofs_per_cell));
206  Assert (values[0].size() >= offset+this->n_components(),
207  ExcDimensionMismatch(values[0].size(),offset+this->n_components()));
208 
209  for (unsigned int i=0; i<this->dofs_per_cell-1; ++i)
210  {
211  const std::pair<unsigned int, unsigned int> index
212  = this->system_to_component_index(i);
213  local_dofs[i] = values[i](offset+index.first);
214  }
215 
216  // We don't need the discontinuous function for local interpolation
217  local_dofs[local_dofs.size()-1] = 0.;
218 }
219 
220 
221 
222 template <int dim, int spacedim>
223 void
225  std::vector<double> &local_dofs,
226  const VectorSlice<const std::vector<std::vector<double> > > &values) const
227 {
228  Assert (values[0].size() == this->unit_support_points.size(),
229  ExcDimensionMismatch(values.size(),
230  this->unit_support_points.size()));
231  Assert (local_dofs.size() == this->dofs_per_cell,
232  ExcDimensionMismatch(local_dofs.size(),this->dofs_per_cell));
233  Assert (values.size() == this->n_components(),
234  ExcDimensionMismatch(values.size(), this->n_components()));
235 
236  for (unsigned int i=0; i<this->dofs_per_cell-1; ++i)
237  {
238  const std::pair<unsigned int, unsigned int> index
239  = this->system_to_component_index(i);
240  local_dofs[i] = values[index.first][i];
241  }
242 
243  // We don't need the discontinuous function for local interpolation
244  local_dofs[local_dofs.size()-1] = 0.;
245 }
246 
247 
248 
249 template <int dim, int spacedim>
250 void
253  FullMatrix<double> &interpolation_matrix) const
254 {
255  // this is only implemented, if the source FE is also a Q_DG0 element
256  typedef FE_Q_DG0<dim,spacedim> FEQDG0;
257  typedef FiniteElement<dim,spacedim> FEL;
258 
259  AssertThrow ((x_source_fe.get_name().find ("FE_Q_DG0<") == 0)
260  ||
261  (dynamic_cast<const FEQDG0 *>(&x_source_fe) != 0),
262  typename FEL::
263  ExcInterpolationNotImplemented());
264 
265  Assert (interpolation_matrix.m() == this->dofs_per_cell,
266  ExcDimensionMismatch (interpolation_matrix.m(),
267  this->dofs_per_cell));
268  Assert (interpolation_matrix.n() == x_source_fe.dofs_per_cell,
269  ExcDimensionMismatch (interpolation_matrix.m(),
270  x_source_fe.dofs_per_cell));
271 
273  get_interpolation_matrix(x_source_fe, interpolation_matrix);
274 }
275 
276 
277 
278 template <int dim, int spacedim>
279 std::vector<bool>
281 {
282  std::vector<bool> riaf(Utilities::fixed_power<dim>(deg+1)+1,false);
283  riaf[riaf.size()-1]=true;
284  return riaf;
285 }
286 
287 
288 
289 template <int dim, int spacedim>
290 std::vector<unsigned int>
292 {
293  std::vector<unsigned int> dpo(dim+1, 1U);
294  for (unsigned int i=1; i<dpo.size(); ++i)
295  dpo[i]=dpo[i-1]*(deg-1);
296 
297  dpo[dim]++;//we need an additional DG0-node for a dim-dimensional object
298  return dpo;
299 }
300 
301 
302 
303 template <int dim, int spacedim>
304 bool
305 FE_Q_DG0<dim,spacedim>::has_support_on_face (const unsigned int shape_index,
306  const unsigned int face_index) const
307 {
308  // discontinuous function has support on all faces
309  if (shape_index == this->dofs_per_cell-1)
310  return true;
311  else
312  return FE_Q_Base<TensorProductPolynomialsConst<dim>,dim,spacedim>::has_support_on_face(shape_index, face_index);
313 }
314 
315 
316 
317 template <int dim, int spacedim>
318 std::pair<Table<2,bool>, std::vector<unsigned int> >
320 {
321  Table<2,bool> constant_modes(2, this->dofs_per_cell);
322 
323  // 1 represented by FE_Q part
324  for (unsigned int i=0; i<this->dofs_per_cell-1; ++i)
325  constant_modes(0, i) = true;
326 
327  // 1 represented by DG0 part
328  constant_modes(1, this->dofs_per_cell-1) = true;
329 
330  return std::pair<Table<2,bool>, std::vector<unsigned int> >
331  (constant_modes, std::vector<unsigned int> (2, 0));
332 }
333 
334 
335 
336 // explicit instantiations
337 #include "fe_q_dg0.inst"
338 
339 DEAL_II_NAMESPACE_CLOSE
size_type m() const
virtual void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) const
Definition: fe_q_dg0.cc:176
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1052
const std::vector< Point< dim > > & get_points() const
::ExceptionBase & ExcMessage(std::string arg1)
const unsigned int degree
Definition: fe_base.h:299
const Point< dim > & point(const unsigned int i) const
#define AssertThrow(cond, exc)
Definition: exceptions.h:358
std::vector< Point< dim > > unit_support_points
Definition: fe.h:2144
virtual std::string get_name() const
Definition: fe_q_dg0.cc:93
size_type n() const
#define Assert(cond, exc)
Definition: exceptions.h:294
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Definition: fe_q_dg0.cc:305
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
Definition: fe_q_dg0.cc:319
virtual std::string get_name() const =0
std::string dim_string(const int dim, const int spacedim)
Definition: utilities.cc:158
unsigned int size() const
const unsigned int dofs_per_cell
Definition: fe_base.h:283
static std::vector< bool > get_riaf_vector(const unsigned int degree)
Definition: fe_q_dg0.cc:280
void initialize(const std::vector< Point< 1 > > &support_points_1d)
Definition: fe_q_base.cc:428
FE_Q_DG0(const unsigned int p)
Definition: fe_q_dg0.cc:33
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
Definition: fe.h:2742
unsigned int n_components() const
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe_q_dg0.cc:252
virtual FiniteElement< dim, spacedim > * clone() const
Definition: fe_q_dg0.cc:167
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition: fe_q_dg0.cc:291