Reference documentation for deal.II version 8.4.2
mapping_q_generic.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2016 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 #ifndef dealii__mapping_q_generic_h
17 #define dealii__mapping_q_generic_h
18 
19 
20 #include <deal.II/base/derivative_form.h>
21 #include <deal.II/base/config.h>
22 #include <deal.II/base/table.h>
23 #include <deal.II/base/quadrature_lib.h>
24 #include <deal.II/base/qprojector.h>
25 #include <deal.II/grid/tria_iterator.h>
26 #include <deal.II/dofs/dof_accessor.h>
27 #include <deal.II/fe/mapping.h>
28 #include <deal.II/fe/fe_q.h>
29 
30 #include <cmath>
31 
32 DEAL_II_NAMESPACE_OPEN
33 
34 template <int,int> class MappingQ;
35 
36 
39 
40 
87 template <int dim, int spacedim=dim>
88 class MappingQGeneric : public Mapping<dim,spacedim>
89 {
90 public:
96  MappingQGeneric (const unsigned int polynomial_degree);
97 
102 
103  // for documentation, see the Mapping base class
104  virtual
105  Mapping<dim,spacedim> *clone () const;
106 
111  unsigned int get_degree () const;
112 
117  virtual
118  bool preserves_vertex_locations () const;
119 
125  // for documentation, see the Mapping base class
126  virtual
129  const Point<dim> &p) const;
130 
131  // for documentation, see the Mapping base class
132  virtual
133  Point<dim>
135  const Point<spacedim> &p) const;
136 
146  // for documentation, see the Mapping base class
147  virtual
148  void
149  transform (const ArrayView<const Tensor<1,dim> > &input,
150  const MappingType type,
151  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
152  const ArrayView<Tensor<1,spacedim> > &output) const;
153 
154  // for documentation, see the Mapping base class
155  virtual
156  void
158  const MappingType type,
159  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
160  const ArrayView<Tensor<2,spacedim> > &output) const;
161 
162  // for documentation, see the Mapping base class
163  virtual
164  void
165  transform (const ArrayView<const Tensor<2, dim> > &input,
166  const MappingType type,
167  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
168  const ArrayView<Tensor<2,spacedim> > &output) const;
169 
170  // for documentation, see the Mapping base class
171  virtual
172  void
174  const MappingType type,
175  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
176  const ArrayView<Tensor<3,spacedim> > &output) const;
177 
178  // for documentation, see the Mapping base class
179  virtual
180  void
181  transform (const ArrayView<const Tensor<3, dim> > &input,
182  const MappingType type,
183  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
184  const ArrayView<Tensor<3,spacedim> > &output) const;
185 
195 public:
207  class InternalData : public Mapping<dim,spacedim>::InternalDataBase
208  {
209  public:
214  InternalData(const unsigned int polynomial_degree);
215 
224  void
225  initialize (const UpdateFlags update_flags,
226  const Quadrature<dim> &quadrature,
227  const unsigned int n_original_q_points);
228 
234  void
235  initialize_face (const UpdateFlags update_flags,
236  const Quadrature<dim> &quadrature,
237  const unsigned int n_original_q_points);
238 
254  void compute_shape_function_values (const std::vector<Point<dim> > &unit_points);
255 
256 
261  const double &shape (const unsigned int qpoint,
262  const unsigned int shape_nr) const;
263 
267  double &shape (const unsigned int qpoint,
268  const unsigned int shape_nr);
269 
273  const Tensor<1,dim> &derivative (const unsigned int qpoint,
274  const unsigned int shape_nr) const;
275 
279  Tensor<1,dim> &derivative (const unsigned int qpoint,
280  const unsigned int shape_nr);
281 
285  const Tensor<2,dim> &second_derivative (const unsigned int qpoint,
286  const unsigned int shape_nr) const;
287 
291  Tensor<2,dim> &second_derivative (const unsigned int qpoint,
292  const unsigned int shape_nr);
293 
297  const Tensor<3,dim> &third_derivative (const unsigned int qpoint,
298  const unsigned int shape_nr) const;
299 
303  Tensor<3,dim> &third_derivative (const unsigned int qpoint,
304  const unsigned int shape_nr);
305 
309  const Tensor<4,dim> &fourth_derivative (const unsigned int qpoint,
310  const unsigned int shape_nr) const;
311 
315  Tensor<4,dim> &fourth_derivative (const unsigned int qpoint,
316  const unsigned int shape_nr);
317 
321  virtual std::size_t memory_consumption () const;
322 
328  std::vector<double> shape_values;
329 
335  std::vector<Tensor<1,dim> > shape_derivatives;
336 
343  std::vector<Tensor<2,dim> > shape_second_derivatives;
344 
351  std::vector<Tensor<3,dim> > shape_third_derivatives;
352 
359  std::vector<Tensor<4,dim> > shape_fourth_derivatives;
360 
374  std::vector<std::vector<Tensor<1,dim> > > unit_tangentials;
375 
380  unsigned int polynomial_degree;
381 
391  const unsigned int n_shape_functions;
392 
402  mutable std::vector<DerivativeForm<1,dim, spacedim > > covariant;
403 
411  mutable std::vector< DerivativeForm<1,dim,spacedim> > contravariant;
412 
416  mutable std::vector<std::vector<Tensor<1,spacedim> > > aux;
417 
422  mutable std::vector<Point<spacedim> > mapping_support_points;
423 
428 
433  mutable std::vector<double> volume_elements;
434  };
435 
436 
437  // documentation can be found in Mapping::requires_update_flags()
438  virtual
440  requires_update_flags (const UpdateFlags update_flags) const;
441 
442  // documentation can be found in Mapping::get_data()
443  virtual
444  InternalData *
445  get_data (const UpdateFlags,
446  const Quadrature<dim> &quadrature) const;
447 
448  // documentation can be found in Mapping::get_face_data()
449  virtual
450  InternalData *
451  get_face_data (const UpdateFlags flags,
452  const Quadrature<dim-1>& quadrature) const;
453 
454  // documentation can be found in Mapping::get_subface_data()
455  virtual
456  InternalData *
457  get_subface_data (const UpdateFlags flags,
458  const Quadrature<dim-1>& quadrature) const;
459 
460  // documentation can be found in Mapping::fill_fe_values()
461  virtual
462  CellSimilarity::Similarity
464  const CellSimilarity::Similarity cell_similarity,
465  const Quadrature<dim> &quadrature,
466  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
468 
469  // documentation can be found in Mapping::fill_fe_face_values()
470  virtual void
472  const unsigned int face_no,
473  const Quadrature<dim-1> &quadrature,
474  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
476 
477  // documentation can be found in Mapping::fill_fe_subface_values()
478  virtual void
480  const unsigned int face_no,
481  const unsigned int subface_no,
482  const Quadrature<dim-1> &quadrature,
483  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
485 
490 protected:
491 
496  const unsigned int polynomial_degree;
497 
498  /*
499  * The default line support points. These are used when computing
500  * the location in real space of the support points on lines and
501  * quads, which are asked to the Manifold<dim,spacedim> class.
502  *
503  * The number of quadrature points depends on the degree of this
504  * class, and it matches the number of degrees of freedom of an
505  * FE_Q<1>(this->degree).
506  */
507  QGaussLobatto<1> line_support_points;
508 
514  const std_cxx11::unique_ptr<FE_Q<dim> > fe_q;
515 
530 
540 
569  virtual
570  std::vector<Point<spacedim> >
572 
577  Point<dim>
579  const Point<spacedim> &p,
580  const Point<dim> &initial_p_unit) const;
581 
596  virtual
597  void
599  std::vector<Point<spacedim> > &a) const;
600 
615  virtual
616  void
618  std::vector<Point<spacedim> > &a) const;
619 
624  template <int, int> friend class MappingQ;
625 };
626 
627 
628 
631 /*----------------------------------------------------------------------*/
632 
633 #ifndef DOXYGEN
634 
635 template<int dim, int spacedim>
636 inline
637 const double &
638 MappingQGeneric<dim,spacedim>::InternalData::shape (const unsigned int qpoint,
639  const unsigned int shape_nr) const
640 {
641  Assert(qpoint*n_shape_functions + shape_nr < shape_values.size(),
642  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
643  shape_values.size()));
644  return shape_values [qpoint*n_shape_functions + shape_nr];
645 }
646 
647 
648 
649 template<int dim, int spacedim>
650 inline
651 double &
652 MappingQGeneric<dim,spacedim>::InternalData::shape (const unsigned int qpoint,
653  const unsigned int shape_nr)
654 {
655  Assert(qpoint*n_shape_functions + shape_nr < shape_values.size(),
656  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
657  shape_values.size()));
658  return shape_values [qpoint*n_shape_functions + shape_nr];
659 }
660 
661 
662 template<int dim, int spacedim>
663 inline
664 const Tensor<1,dim> &
666  const unsigned int shape_nr) const
667 {
668  Assert(qpoint*n_shape_functions + shape_nr < shape_derivatives.size(),
669  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
670  shape_derivatives.size()));
671  return shape_derivatives [qpoint*n_shape_functions + shape_nr];
672 }
673 
674 
675 
676 template<int dim, int spacedim>
677 inline
680  const unsigned int shape_nr)
681 {
682  Assert(qpoint*n_shape_functions + shape_nr < shape_derivatives.size(),
683  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
684  shape_derivatives.size()));
685  return shape_derivatives [qpoint*n_shape_functions + shape_nr];
686 }
687 
688 
689 template <int dim, int spacedim>
690 inline
691 const Tensor<2,dim> &
693  const unsigned int shape_nr) const
694 {
695  Assert(qpoint*n_shape_functions + shape_nr < shape_second_derivatives.size(),
696  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
697  shape_second_derivatives.size()));
698  return shape_second_derivatives [qpoint*n_shape_functions + shape_nr];
699 }
700 
701 
702 template <int dim, int spacedim>
703 inline
706  const unsigned int shape_nr)
707 {
708  Assert(qpoint*n_shape_functions + shape_nr < shape_second_derivatives.size(),
709  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
710  shape_second_derivatives.size()));
711  return shape_second_derivatives [qpoint*n_shape_functions + shape_nr];
712 }
713 
714 template <int dim, int spacedim>
715 inline
716 const Tensor<3,dim> &
718  const unsigned int shape_nr) const
719 {
720  Assert(qpoint*n_shape_functions + shape_nr < shape_third_derivatives.size(),
721  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
722  shape_third_derivatives.size()));
723  return shape_third_derivatives [qpoint*n_shape_functions + shape_nr];
724 }
725 
726 
727 template <int dim, int spacedim>
728 inline
731  const unsigned int shape_nr)
732 {
733  Assert(qpoint*n_shape_functions + shape_nr < shape_third_derivatives.size(),
734  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
735  shape_third_derivatives.size()));
736  return shape_third_derivatives [qpoint*n_shape_functions + shape_nr];
737 }
738 
739 
740 template <int dim, int spacedim>
741 inline
742 const Tensor<4,dim> &
744  const unsigned int shape_nr) const
745 {
746  Assert(qpoint*n_shape_functions + shape_nr < shape_fourth_derivatives.size(),
747  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
748  shape_fourth_derivatives.size()));
749  return shape_fourth_derivatives [qpoint*n_shape_functions + shape_nr];
750 }
751 
752 
753 template <int dim, int spacedim>
754 inline
757  const unsigned int shape_nr)
758 {
759  Assert(qpoint*n_shape_functions + shape_nr < shape_fourth_derivatives.size(),
760  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
761  shape_fourth_derivatives.size()));
762  return shape_fourth_derivatives [qpoint*n_shape_functions + shape_nr];
763 }
764 
765 
766 
767 template <int dim, int spacedim>
768 inline
769 bool
771 {
772  return true;
773 }
774 
775 #endif // DOXYGEN
776 
777 /* -------------- declaration of explicit specializations ------------- */
778 
779 
780 DEAL_II_NAMESPACE_CLOSE
781 
782 #endif
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
const std_cxx11::unique_ptr< FE_Q< dim > > fe_q
Table< 2, double > support_point_weights_on_hex
std::vector< Tensor< 2, dim > > shape_second_derivatives
const unsigned int polynomial_degree
virtual CellSimilarity::Similarity fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValues::MappingRelatedData< dim, spacedim > &output_data) const
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim-1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValues::MappingRelatedData< dim, spacedim > &output_data) const
virtual InternalData * get_subface_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const
MappingType
Definition: mapping.h:50
virtual InternalData * get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const
std::vector< Tensor< 4, dim > > shape_fourth_derivatives
virtual Mapping< dim, spacedim > * clone() const
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const
Triangulation< dim, spacedim >::cell_iterator cell_of_current_support_points
const Tensor< 3, dim > & third_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
MappingQGeneric(const unsigned int polynomial_degree)
const Tensor< 4, dim > & fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
InternalData(const unsigned int polynomial_degree)
std::vector< Tensor< 1, dim > > shape_derivatives
virtual void add_quad_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, std::vector< Point< spacedim > > &a) const
void compute_shape_function_values(const std::vector< Point< dim > > &unit_points)
#define Assert(cond, exc)
Definition: exceptions.h:294
UpdateFlags
Abstract base class for mapping classes.
Definition: dof_tools.h:52
virtual void transform(const ArrayView< const Tensor< 1, dim > > &input, const MappingType type, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 1, spacedim > > &output) const
std::vector< Point< spacedim > > mapping_support_points
std::vector< DerivativeForm< 1, dim, spacedim > > contravariant
std::vector< std::vector< Tensor< 1, spacedim > > > aux
void initialize(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
std::vector< double > volume_elements
const unsigned int n_shape_functions
Point< dim > transform_real_to_unit_cell_internal(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p, const Point< dim > &initial_p_unit) const
virtual std::size_t memory_consumption() const
unsigned int get_degree() const
const double & shape(const unsigned int qpoint, const unsigned int shape_nr) const
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const
Table< 2, double > support_point_weights_on_quad
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const Quadrature< dim-1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValues::MappingRelatedData< dim, spacedim > &output_data) const
virtual std::vector< Point< spacedim > > compute_mapping_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
virtual void add_line_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, std::vector< Point< spacedim > > &a) const
virtual InternalData * get_face_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const
virtual bool preserves_vertex_locations() const
std::vector< Tensor< 3, dim > > shape_third_derivatives
const Tensor< 2, dim > & second_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
const Tensor< 1, dim > & derivative(const unsigned int qpoint, const unsigned int shape_nr) const
void initialize_face(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
std::vector< std::vector< Tensor< 1, dim > > > unit_tangentials
std::vector< double > shape_values