Reference documentation for deal.II version 9.1.0-pre
fe_dgp_nonparametric.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2002 - 2018 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.h>
18 #include <deal.II/base/std_cxx14/memory.h>
19 
20 #include <deal.II/dofs/dof_accessor.h>
21 
22 #include <deal.II/fe/fe.h>
23 #include <deal.II/fe/fe_dgp_nonparametric.h>
24 #include <deal.II/fe/fe_values.h>
25 #include <deal.II/fe/mapping.h>
26 
27 #include <deal.II/grid/tria.h>
28 #include <deal.II/grid/tria_iterator.h>
29 
30 #include <sstream>
31 
32 
33 DEAL_II_NAMESPACE_OPEN
34 
35 template <int dim, int spacedim>
37  const unsigned int degree)
38  : FiniteElement<dim, spacedim>(
39  FiniteElementData<dim>(get_dpo_vector(degree),
40  1,
41  degree,
42  FiniteElementData<dim>::L2),
43  std::vector<bool>(
44  FiniteElementData<dim>(get_dpo_vector(degree), 1, degree).dofs_per_cell,
45  true),
46  std::vector<ComponentMask>(
47  FiniteElementData<dim>(get_dpo_vector(degree), 1, degree).dofs_per_cell,
48  std::vector<bool>(1, true)))
49  , polynomial_space(Polynomials::Legendre::generate_complete_basis(degree))
50 {
51  const unsigned int n_dofs = this->dofs_per_cell;
52  for (unsigned int ref_case = RefinementCase<dim>::cut_x;
53  ref_case < RefinementCase<dim>::isotropic_refinement + 1;
54  ++ref_case)
55  {
56  if (dim != 2 && ref_case != RefinementCase<dim>::isotropic_refinement)
57  // do nothing, as anisotropic
58  // refinement is not
59  // implemented so far
60  continue;
61 
62  const unsigned int nc =
64  for (unsigned int i = 0; i < nc; ++i)
65  {
66  this->prolongation[ref_case - 1][i].reinit(n_dofs, n_dofs);
67  // Fill prolongation matrices with
68  // embedding operators
69  for (unsigned int j = 0; j < n_dofs; ++j)
70  this->prolongation[ref_case - 1][i](j, j) = 1.;
71  }
72  }
73 
74  // restriction can be defined
75  // through projection for
76  // discontinuous elements, but is
77  // presently not implemented for DGPNonparametric
78  // elements.
79  //
80  // if it were, then the following
81  // snippet would be the right code
82  // if ((degree < Matrices::n_projection_matrices) &&
83  // (Matrices::projection_matrices[degree] != 0))
84  // {
85  // restriction[0].fill (Matrices::projection_matrices[degree]);
86  // }
87  // else
88  // // matrix undefined, set size to zero
89  // for (unsigned int i=0;i<GeometryInfo<dim>::max_children_per_cell;++i)
90  // restriction[i].reinit(0, 0);
91  // since not implemented, set to
92  // "empty". however, that is done in the
93  // default constructor already, so do nothing
94  // for (unsigned int i=0;i<GeometryInfo<dim>::max_children_per_cell;++i)
95  // this->restriction[i].reinit(0, 0);
96 
97  // note further, that these
98  // elements have neither support
99  // nor face-support points, so
100  // leave these fields empty
101 }
102 
103 
104 
105 template <int dim, int spacedim>
106 std::string
108 {
109  // note that the
110  // FETools::get_fe_by_name
111  // function depends on the
112  // particular format of the string
113  // this function returns, so they
114  // have to be kept in synch
115 
116  std::ostringstream namebuf;
117  namebuf << "FE_DGPNonparametric<" << Utilities::dim_string(dim, spacedim)
118  << ">(" << this->degree << ")";
119 
120  return namebuf.str();
121 }
122 
123 
124 
125 template <int dim, int spacedim>
126 std::unique_ptr<FiniteElement<dim, spacedim>>
128 {
129  return std_cxx14::make_unique<FE_DGPNonparametric<dim, spacedim>>(*this);
130 }
131 
132 
133 
134 template <int dim, int spacedim>
135 double
137  const Point<dim> & p) const
138 {
139  (void)i;
140  (void)p;
141  Assert(i < this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
142  AssertThrow(false,
144  return 0;
145 }
146 
147 
148 
149 template <int dim, int spacedim>
150 double
152  const unsigned int i,
153  const Point<dim> & p,
154  const unsigned int component) const
155 {
156  (void)i;
157  (void)p;
158  (void)component;
159  Assert(i < this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
160  Assert(component == 0, ExcIndexRange(component, 0, 1));
161  AssertThrow(false,
163  return 0;
164 }
165 
166 
167 
168 template <int dim, int spacedim>
171  const Point<dim> & p) const
172 {
173  (void)i;
174  (void)p;
175  Assert(i < this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
176  AssertThrow(false,
178  return Tensor<1, dim>();
179 }
180 
181 
182 template <int dim, int spacedim>
185  const unsigned int i,
186  const Point<dim> & p,
187  const unsigned int component) const
188 {
189  (void)i;
190  (void)p;
191  (void)component;
192  Assert(i < this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
193  Assert(component == 0, ExcIndexRange(component, 0, 1));
194  AssertThrow(false,
196  return Tensor<1, dim>();
197 }
198 
199 
200 
201 template <int dim, int spacedim>
204  const Point<dim> & p) const
205 {
206  (void)i;
207  (void)p;
208  Assert(i < this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
209  AssertThrow(false,
211  return Tensor<2, dim>();
212 }
213 
214 
215 
216 template <int dim, int spacedim>
219  const unsigned int i,
220  const Point<dim> & p,
221  const unsigned int component) const
222 {
223  (void)i;
224  (void)p;
225  (void)component;
226  Assert(i < this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
227  Assert(component == 0, ExcIndexRange(component, 0, 1));
228  AssertThrow(false,
230  return Tensor<2, dim>();
231 }
232 
233 
234 //---------------------------------------------------------------------------
235 // Auxiliary functions
236 //---------------------------------------------------------------------------
237 
238 
239 template <int dim, int spacedim>
240 std::vector<unsigned int>
242 {
243  std::vector<unsigned int> dpo(dim + 1, static_cast<unsigned int>(0));
244  dpo[dim] = deg + 1;
245  for (unsigned int i = 1; i < dim; ++i)
246  {
247  dpo[dim] *= deg + 1 + i;
248  dpo[dim] /= i + 1;
249  }
250  return dpo;
251 }
252 
253 
254 
255 template <int dim, int spacedim>
258  const UpdateFlags flags) const
259 {
260  UpdateFlags out = flags;
261 
264 
265  return out;
266 }
267 
268 
269 
270 //---------------------------------------------------------------------------
271 // Data field initialization
272 //---------------------------------------------------------------------------
273 
274 template <int dim, int spacedim>
275 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
277  const UpdateFlags update_flags,
278  const Mapping<dim, spacedim> &,
279  const Quadrature<dim> &,
281  spacedim>
282  & /*output_data*/) const
283 {
284  // generate a new data object
285  auto data = std_cxx14::make_unique<
287  data->update_each = requires_update_flags(update_flags);
288 
289  // other than that, there is nothing we can add here as discussed
290  // in the general documentation of this class
291 
292  return std::move(data);
293 }
294 
295 
296 
297 //---------------------------------------------------------------------------
298 // Fill data of FEValues
299 //---------------------------------------------------------------------------
300 
301 template <int dim, int spacedim>
302 void
306  const Quadrature<dim> &,
307  const Mapping<dim, spacedim> &,
309  const ::internal::FEValuesImplementation::MappingRelatedData<dim,
310  spacedim>
311  & mapping_data,
312  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
314  spacedim>
315  &output_data) const
316 {
318  ExcInternalError());
319 
320  const unsigned int n_q_points = mapping_data.quadrature_points.size();
321 
322  std::vector<double> values(
323  (fe_internal.update_each & update_values) ? this->dofs_per_cell : 0);
324  std::vector<Tensor<1, dim>> grads(
325  (fe_internal.update_each & update_gradients) ? this->dofs_per_cell : 0);
326  std::vector<Tensor<2, dim>> grad_grads(
327  (fe_internal.update_each & update_hessians) ? this->dofs_per_cell : 0);
328  std::vector<Tensor<3, dim>> empty_vector_of_3rd_order_tensors;
329  std::vector<Tensor<4, dim>> empty_vector_of_4th_order_tensors;
330 
331  if (fe_internal.update_each & (update_values | update_gradients))
332  for (unsigned int i = 0; i < n_q_points; ++i)
333  {
334  polynomial_space.compute(mapping_data.quadrature_points[i],
335  values,
336  grads,
337  grad_grads,
338  empty_vector_of_3rd_order_tensors,
339  empty_vector_of_4th_order_tensors);
340 
341  if (fe_internal.update_each & update_values)
342  for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
343  output_data.shape_values[k][i] = values[k];
344 
345  if (fe_internal.update_each & update_gradients)
346  for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
347  output_data.shape_gradients[k][i] = grads[k];
348 
349  if (fe_internal.update_each & update_hessians)
350  for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
351  output_data.shape_hessians[k][i] = grad_grads[k];
352  }
353 }
354 
355 
356 
357 template <int dim, int spacedim>
358 void
361  const unsigned int,
362  const Quadrature<dim - 1> &,
363  const Mapping<dim, spacedim> &,
365  const ::internal::FEValuesImplementation::MappingRelatedData<dim,
366  spacedim>
367  & mapping_data,
368  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
370  spacedim>
371  &output_data) const
372 {
374  ExcInternalError());
375 
376  const unsigned int n_q_points = mapping_data.quadrature_points.size();
377 
378  std::vector<double> values(
379  (fe_internal.update_each & update_values) ? this->dofs_per_cell : 0);
380  std::vector<Tensor<1, dim>> grads(
381  (fe_internal.update_each & update_gradients) ? this->dofs_per_cell : 0);
382  std::vector<Tensor<2, dim>> grad_grads(
383  (fe_internal.update_each & update_hessians) ? this->dofs_per_cell : 0);
384  std::vector<Tensor<3, dim>> empty_vector_of_3rd_order_tensors;
385  std::vector<Tensor<4, dim>> empty_vector_of_4th_order_tensors;
386 
387  if (fe_internal.update_each & (update_values | update_gradients))
388  for (unsigned int i = 0; i < n_q_points; ++i)
389  {
390  polynomial_space.compute(mapping_data.quadrature_points[i],
391  values,
392  grads,
393  grad_grads,
394  empty_vector_of_3rd_order_tensors,
395  empty_vector_of_4th_order_tensors);
396 
397  if (fe_internal.update_each & update_values)
398  for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
399  output_data.shape_values[k][i] = values[k];
400 
401  if (fe_internal.update_each & update_gradients)
402  for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
403  output_data.shape_gradients[k][i] = grads[k];
404 
405  if (fe_internal.update_each & update_hessians)
406  for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
407  output_data.shape_hessians[k][i] = grad_grads[k];
408  }
409 }
410 
411 
412 
413 template <int dim, int spacedim>
414 void
417  const unsigned int,
418  const unsigned int,
419  const Quadrature<dim - 1> &,
420  const Mapping<dim, spacedim> &,
422  const ::internal::FEValuesImplementation::MappingRelatedData<dim,
423  spacedim>
424  & mapping_data,
425  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
427  spacedim>
428  &output_data) const
429 {
431  ExcInternalError());
432 
433  const unsigned int n_q_points = mapping_data.quadrature_points.size();
434 
435  std::vector<double> values(
436  (fe_internal.update_each & update_values) ? this->dofs_per_cell : 0);
437  std::vector<Tensor<1, dim>> grads(
438  (fe_internal.update_each & update_gradients) ? this->dofs_per_cell : 0);
439  std::vector<Tensor<2, dim>> grad_grads(
440  (fe_internal.update_each & update_hessians) ? this->dofs_per_cell : 0);
441  std::vector<Tensor<3, dim>> empty_vector_of_3rd_order_tensors;
442  std::vector<Tensor<4, dim>> empty_vector_of_4th_order_tensors;
443 
444  if (fe_internal.update_each & (update_values | update_gradients))
445  for (unsigned int i = 0; i < n_q_points; ++i)
446  {
447  polynomial_space.compute(mapping_data.quadrature_points[i],
448  values,
449  grads,
450  grad_grads,
451  empty_vector_of_3rd_order_tensors,
452  empty_vector_of_4th_order_tensors);
453 
454  if (fe_internal.update_each & update_values)
455  for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
456  output_data.shape_values[k][i] = values[k];
457 
458  if (fe_internal.update_each & update_gradients)
459  for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
460  output_data.shape_gradients[k][i] = grads[k];
461 
462  if (fe_internal.update_each & update_hessians)
463  for (unsigned int k = 0; k < this->dofs_per_cell; ++k)
464  output_data.shape_hessians[k][i] = grad_grads[k];
465  }
466 }
467 
468 
469 
470 template <int dim, int spacedim>
471 void
473  const FiniteElement<dim, spacedim> &x_source_fe,
474  FullMatrix<double> & interpolation_matrix) const
475 {
476  // this is only implemented, if the source
477  // FE is also a DGPNonparametric element. in that case,
478  // both elements have no dofs on their
479  // faces and the face interpolation matrix
480  // is necessarily empty -- i.e. there isn't
481  // much we need to do here.
482  (void)interpolation_matrix;
483  AssertThrow(
484  (x_source_fe.get_name().find("FE_DGPNonparametric<") == 0) ||
485  (dynamic_cast<const FE_DGPNonparametric<dim, spacedim> *>(&x_source_fe) !=
486  nullptr),
488 
489  Assert(interpolation_matrix.m() == 0,
490  ExcDimensionMismatch(interpolation_matrix.m(), 0));
491  Assert(interpolation_matrix.n() == 0,
492  ExcDimensionMismatch(interpolation_matrix.n(), 0));
493 }
494 
495 
496 
497 template <int dim, int spacedim>
498 void
500  const FiniteElement<dim, spacedim> &x_source_fe,
501  const unsigned int,
502  FullMatrix<double> &interpolation_matrix) const
503 {
504  // this is only implemented, if the source
505  // FE is also a DGPNonparametric element. in that case,
506  // both elements have no dofs on their
507  // faces and the face interpolation matrix
508  // is necessarily empty -- i.e. there isn't
509  // much we need to do here.
510  (void)interpolation_matrix;
511  AssertThrow(
512  (x_source_fe.get_name().find("FE_DGPNonparametric<") == 0) ||
513  (dynamic_cast<const FE_DGPNonparametric<dim, spacedim> *>(&x_source_fe) !=
514  nullptr),
516 
517  Assert(interpolation_matrix.m() == 0,
518  ExcDimensionMismatch(interpolation_matrix.m(), 0));
519  Assert(interpolation_matrix.n() == 0,
520  ExcDimensionMismatch(interpolation_matrix.n(), 0));
521 }
522 
523 
524 
525 template <int dim, int spacedim>
526 bool
528 {
529  return true;
530 }
531 
532 
533 
534 template <int dim, int spacedim>
535 std::vector<std::pair<unsigned int, unsigned int>>
537  const FiniteElement<dim, spacedim> &fe_other) const
538 {
539  // there are no such constraints for DGPNonparametric
540  // elements at all
541  if (dynamic_cast<const FE_DGPNonparametric<dim, spacedim> *>(&fe_other) !=
542  nullptr)
543  return std::vector<std::pair<unsigned int, unsigned int>>();
544  else
545  {
546  Assert(false, ExcNotImplemented());
547  return std::vector<std::pair<unsigned int, unsigned int>>();
548  }
549 }
550 
551 
552 
553 template <int dim, int spacedim>
554 std::vector<std::pair<unsigned int, unsigned int>>
556  const FiniteElement<dim, spacedim> &fe_other) const
557 {
558  // there are no such constraints for DGPNonparametric
559  // elements at all
560  if (dynamic_cast<const FE_DGPNonparametric<dim, spacedim> *>(&fe_other) !=
561  nullptr)
562  return std::vector<std::pair<unsigned int, unsigned int>>();
563  else
564  {
565  Assert(false, ExcNotImplemented());
566  return std::vector<std::pair<unsigned int, unsigned int>>();
567  }
568 }
569 
570 
571 
572 template <int dim, int spacedim>
573 std::vector<std::pair<unsigned int, unsigned int>>
575  const FiniteElement<dim, spacedim> &fe_other) const
576 {
577  // there are no such constraints for DGPNonparametric
578  // elements at all
579  if (dynamic_cast<const FE_DGPNonparametric<dim, spacedim> *>(&fe_other) !=
580  nullptr)
581  return std::vector<std::pair<unsigned int, unsigned int>>();
582  else
583  {
584  Assert(false, ExcNotImplemented());
585  return std::vector<std::pair<unsigned int, unsigned int>>();
586  }
587 }
588 
589 
590 
591 template <int dim, int spacedim>
594  const FiniteElement<dim, spacedim> &fe_other) const
595 {
596  // check whether both are discontinuous
597  // elements, see the description of
598  // FiniteElementDomination::Domination
599  if (dynamic_cast<const FE_DGPNonparametric<dim, spacedim> *>(&fe_other) !=
600  nullptr)
602 
603  Assert(false, ExcNotImplemented());
605 }
606 
607 
608 
609 template <int dim, int spacedim>
610 bool
612  const unsigned int,
613  const unsigned int) const
614 {
615  return true;
616 }
617 
618 
619 
620 template <int dim, int spacedim>
621 std::size_t
623 {
624  Assert(false, ExcNotImplemented());
625  return 0;
626 }
627 
628 
629 
630 template <int dim, int spacedim>
631 unsigned int
633 {
634  return this->degree;
635 }
636 
637 
638 
639 // explicit instantiations
640 #include "fe_dgp_nonparametric.inst"
641 
642 
643 DEAL_II_NAMESPACE_CLOSE
Shape function values.
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const override
friend class FE_DGPNonparametric
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
virtual std::size_t memory_consumption() const override
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const override
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
const unsigned int degree
Definition: fe_base.h:313
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
const PolynomialSpace< dim > polynomial_space
STL namespace.
Transformed quadrature points.
#define AssertThrow(cond, exc)
Definition: exceptions.h:1329
static::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
virtual bool hp_constraints_are_implemented() const override
unsigned int get_degree() const
size_type n() const
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition: fe.h:2433
#define Assert(cond, exc)
Definition: exceptions.h:1227
UpdateFlags
static::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
Definition: dof_tools.h:57
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature,::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual std::string get_name() const =0
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
size_type m() const
Second derivatives of shape functions.
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 unsigned int n_children(const RefinementCase< dim > &refinement_case)
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const override
Shape function gradients.
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const override
virtual std::string get_name() const override
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const override
static::ExceptionBase & ExcNotImplemented()
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
UpdateFlags update_each
Definition: fe.h:713
static::ExceptionBase & ExcInternalError()