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