Reference documentation for deal.II version 8.5.1
mapping_q_generic.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 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 #ifndef dealii__mapping_q_generic_h
17 #define dealii__mapping_q_generic_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/derivative_form.h>
22 #include <deal.II/base/table.h>
23 #include <deal.II/base/quadrature_lib.h>
24 #include <deal.II/grid/tria_iterator.h>
25 #include <deal.II/fe/mapping.h>
26 #include <deal.II/fe/fe_q.h>
27 
28 #include <deal.II/base/std_cxx11/array.h>
29 
30 #include <cmath>
31 
32 DEAL_II_NAMESPACE_OPEN
33 
34 template <int,int> class MappingQ;
35 
36 
39 
40 
136 template <int dim, int spacedim=dim>
137 class MappingQGeneric : public Mapping<dim,spacedim>
138 {
139 public:
145  MappingQGeneric (const unsigned int polynomial_degree);
146 
151 
152  // for documentation, see the Mapping base class
153  virtual
154  Mapping<dim,spacedim> *clone () const;
155 
160  unsigned int get_degree () const;
161 
166  virtual
167  bool preserves_vertex_locations () const;
168 
174  // for documentation, see the Mapping base class
175  virtual
178  const Point<dim> &p) const;
179 
180  // for documentation, see the Mapping base class
181  virtual
182  Point<dim>
184  const Point<spacedim> &p) const;
185 
195  // for documentation, see the Mapping base class
196  virtual
197  void
198  transform (const ArrayView<const Tensor<1,dim> > &input,
199  const MappingType type,
200  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
201  const ArrayView<Tensor<1,spacedim> > &output) const;
202 
203  // for documentation, see the Mapping base class
204  virtual
205  void
207  const MappingType type,
208  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
209  const ArrayView<Tensor<2,spacedim> > &output) const;
210 
211  // for documentation, see the Mapping base class
212  virtual
213  void
214  transform (const ArrayView<const Tensor<2, dim> > &input,
215  const MappingType type,
216  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
217  const ArrayView<Tensor<2,spacedim> > &output) const;
218 
219  // for documentation, see the Mapping base class
220  virtual
221  void
223  const MappingType type,
224  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
225  const ArrayView<Tensor<3,spacedim> > &output) const;
226 
227  // for documentation, see the Mapping base class
228  virtual
229  void
230  transform (const ArrayView<const Tensor<3, dim> > &input,
231  const MappingType type,
232  const typename Mapping<dim,spacedim>::InternalDataBase &internal,
233  const ArrayView<Tensor<3,spacedim> > &output) const;
234 
244 public:
256  class InternalData : public Mapping<dim,spacedim>::InternalDataBase
257  {
258  public:
263  InternalData(const unsigned int polynomial_degree);
264 
273  void
274  initialize (const UpdateFlags update_flags,
275  const Quadrature<dim> &quadrature,
276  const unsigned int n_original_q_points);
277 
283  void
284  initialize_face (const UpdateFlags update_flags,
285  const Quadrature<dim> &quadrature,
286  const unsigned int n_original_q_points);
287 
303  void compute_shape_function_values (const std::vector<Point<dim> > &unit_points);
304 
305 
310  const double &shape (const unsigned int qpoint,
311  const unsigned int shape_nr) const;
312 
316  double &shape (const unsigned int qpoint,
317  const unsigned int shape_nr);
318 
322  const Tensor<1,dim> &derivative (const unsigned int qpoint,
323  const unsigned int shape_nr) const;
324 
328  Tensor<1,dim> &derivative (const unsigned int qpoint,
329  const unsigned int shape_nr);
330 
334  const Tensor<2,dim> &second_derivative (const unsigned int qpoint,
335  const unsigned int shape_nr) const;
336 
340  Tensor<2,dim> &second_derivative (const unsigned int qpoint,
341  const unsigned int shape_nr);
342 
346  const Tensor<3,dim> &third_derivative (const unsigned int qpoint,
347  const unsigned int shape_nr) const;
348 
352  Tensor<3,dim> &third_derivative (const unsigned int qpoint,
353  const unsigned int shape_nr);
354 
358  const Tensor<4,dim> &fourth_derivative (const unsigned int qpoint,
359  const unsigned int shape_nr) const;
360 
364  Tensor<4,dim> &fourth_derivative (const unsigned int qpoint,
365  const unsigned int shape_nr);
366 
370  virtual std::size_t memory_consumption () const;
371 
377  std::vector<double> shape_values;
378 
384  std::vector<Tensor<1,dim> > shape_derivatives;
385 
392  std::vector<Tensor<2,dim> > shape_second_derivatives;
393 
400  std::vector<Tensor<3,dim> > shape_third_derivatives;
401 
408  std::vector<Tensor<4,dim> > shape_fourth_derivatives;
409 
423  std_cxx11::array<std::vector<Tensor<1,dim> >, GeometryInfo<dim>::faces_per_cell *(dim-1)> unit_tangentials;
424 
429  const unsigned int polynomial_degree;
430 
440  const unsigned int n_shape_functions;
441 
442  /*
443  * The default line support points. Is used in when the shape function
444  * values are computed.
445  *
446  * The number of quadrature points depends on the degree of this
447  * class, and it matches the number of degrees of freedom of an
448  * FE_Q<1>(this->degree).
449  */
450  QGaussLobatto<1> line_support_points;
451 
461  mutable std::vector<DerivativeForm<1,dim, spacedim > > covariant;
462 
470  mutable std::vector< DerivativeForm<1,dim,spacedim> > contravariant;
471 
475  mutable std::vector<std::vector<Tensor<1,spacedim> > > aux;
476 
481  mutable std::vector<Point<spacedim> > mapping_support_points;
482 
487 
492  mutable std::vector<double> volume_elements;
493  };
494 
495 
496  // documentation can be found in Mapping::requires_update_flags()
497  virtual
499  requires_update_flags (const UpdateFlags update_flags) const;
500 
501  // documentation can be found in Mapping::get_data()
502  virtual
503  InternalData *
504  get_data (const UpdateFlags,
505  const Quadrature<dim> &quadrature) const;
506 
507  // documentation can be found in Mapping::get_face_data()
508  virtual
509  InternalData *
510  get_face_data (const UpdateFlags flags,
511  const Quadrature<dim-1>& quadrature) const;
512 
513  // documentation can be found in Mapping::get_subface_data()
514  virtual
515  InternalData *
516  get_subface_data (const UpdateFlags flags,
517  const Quadrature<dim-1>& quadrature) const;
518 
519  // documentation can be found in Mapping::fill_fe_values()
520  virtual
523  const CellSimilarity::Similarity cell_similarity,
524  const Quadrature<dim> &quadrature,
525  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
527 
528  // documentation can be found in Mapping::fill_fe_face_values()
529  virtual void
531  const unsigned int face_no,
532  const Quadrature<dim-1> &quadrature,
533  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
535 
536  // documentation can be found in Mapping::fill_fe_subface_values()
537  virtual void
539  const unsigned int face_no,
540  const unsigned int subface_no,
541  const Quadrature<dim-1> &quadrature,
542  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
544 
549 protected:
550 
555  const unsigned int polynomial_degree;
556 
557  /*
558  * The default line support points. These are used when computing
559  * the location in real space of the support points on lines and
560  * quads, which are asked to the Manifold<dim,spacedim> class.
561  *
562  * The number of quadrature points depends on the degree of this
563  * class, and it matches the number of degrees of freedom of an
564  * FE_Q<1>(this->degree).
565  */
566  QGaussLobatto<1> line_support_points;
567 
573  const std_cxx11::unique_ptr<FE_Q<dim> > fe_q;
574 
594  std::vector<Table<2,double> > support_point_weights_perimeter_to_interior;
595 
609 
638  virtual
639  std::vector<Point<spacedim> >
641 
646  Point<dim>
648  const Point<spacedim> &p,
649  const Point<dim> &initial_p_unit) const;
650 
665  virtual
666  void
668  std::vector<Point<spacedim> > &a) const;
669 
684  virtual
685  void
687  std::vector<Point<spacedim> > &a) const;
688 
693  template <int, int> friend class MappingQ;
694 };
695 
696 
697 
700 /*----------------------------------------------------------------------*/
701 
702 #ifndef DOXYGEN
703 
704 template<int dim, int spacedim>
705 inline
706 const double &
707 MappingQGeneric<dim,spacedim>::InternalData::shape (const unsigned int qpoint,
708  const unsigned int shape_nr) const
709 {
710  Assert(qpoint*n_shape_functions + shape_nr < shape_values.size(),
711  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
712  shape_values.size()));
713  return shape_values [qpoint*n_shape_functions + shape_nr];
714 }
715 
716 
717 
718 template<int dim, int spacedim>
719 inline
720 double &
721 MappingQGeneric<dim,spacedim>::InternalData::shape (const unsigned int qpoint,
722  const unsigned int shape_nr)
723 {
724  Assert(qpoint*n_shape_functions + shape_nr < shape_values.size(),
725  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
726  shape_values.size()));
727  return shape_values [qpoint*n_shape_functions + shape_nr];
728 }
729 
730 
731 template<int dim, int spacedim>
732 inline
733 const Tensor<1,dim> &
735  const unsigned int shape_nr) const
736 {
737  Assert(qpoint*n_shape_functions + shape_nr < shape_derivatives.size(),
738  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
739  shape_derivatives.size()));
740  return shape_derivatives [qpoint*n_shape_functions + shape_nr];
741 }
742 
743 
744 
745 template<int dim, int spacedim>
746 inline
749  const unsigned int shape_nr)
750 {
751  Assert(qpoint*n_shape_functions + shape_nr < shape_derivatives.size(),
752  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
753  shape_derivatives.size()));
754  return shape_derivatives [qpoint*n_shape_functions + shape_nr];
755 }
756 
757 
758 template <int dim, int spacedim>
759 inline
760 const Tensor<2,dim> &
762  const unsigned int shape_nr) const
763 {
764  Assert(qpoint*n_shape_functions + shape_nr < shape_second_derivatives.size(),
765  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
766  shape_second_derivatives.size()));
767  return shape_second_derivatives [qpoint*n_shape_functions + shape_nr];
768 }
769 
770 
771 template <int dim, int spacedim>
772 inline
775  const unsigned int shape_nr)
776 {
777  Assert(qpoint*n_shape_functions + shape_nr < shape_second_derivatives.size(),
778  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
779  shape_second_derivatives.size()));
780  return shape_second_derivatives [qpoint*n_shape_functions + shape_nr];
781 }
782 
783 template <int dim, int spacedim>
784 inline
785 const Tensor<3,dim> &
787  const unsigned int shape_nr) const
788 {
789  Assert(qpoint*n_shape_functions + shape_nr < shape_third_derivatives.size(),
790  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
791  shape_third_derivatives.size()));
792  return shape_third_derivatives [qpoint*n_shape_functions + shape_nr];
793 }
794 
795 
796 template <int dim, int spacedim>
797 inline
800  const unsigned int shape_nr)
801 {
802  Assert(qpoint*n_shape_functions + shape_nr < shape_third_derivatives.size(),
803  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
804  shape_third_derivatives.size()));
805  return shape_third_derivatives [qpoint*n_shape_functions + shape_nr];
806 }
807 
808 
809 template <int dim, int spacedim>
810 inline
811 const Tensor<4,dim> &
813  const unsigned int shape_nr) const
814 {
815  Assert(qpoint*n_shape_functions + shape_nr < shape_fourth_derivatives.size(),
816  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
817  shape_fourth_derivatives.size()));
818  return shape_fourth_derivatives [qpoint*n_shape_functions + shape_nr];
819 }
820 
821 
822 template <int dim, int spacedim>
823 inline
826  const unsigned int shape_nr)
827 {
828  Assert(qpoint*n_shape_functions + shape_nr < shape_fourth_derivatives.size(),
829  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
830  shape_fourth_derivatives.size()));
831  return shape_fourth_derivatives [qpoint*n_shape_functions + shape_nr];
832 }
833 
834 
835 
836 template <int dim, int spacedim>
837 inline
838 bool
840 {
841  return true;
842 }
843 
844 #endif // DOXYGEN
845 
846 /* -------------- declaration of explicit specializations ------------- */
847 
848 
849 DEAL_II_NAMESPACE_CLOSE
850 
851 #endif
std_cxx11::array< std::vector< Tensor< 1, dim > >, GeometryInfo< dim >::faces_per_cell *(dim-1)> unit_tangentials
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
const std_cxx11::unique_ptr< FE_Q< dim > > fe_q
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< Table< 2, double > > support_point_weights_perimeter_to_interior
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
Table< 2, double > support_point_weights_cell
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)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
const Tensor< 4, dim > & fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
InternalData(const unsigned int polynomial_degree)
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:313
UpdateFlags
Abstract base class for mapping classes.
Definition: dof_tools.h:46
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
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< double > shape_values