Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
qprojector.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2005 - 2019 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 #ifndef dealii_qprojector_h
17 #define dealii_qprojector_h
18 
19 
20 #include <deal.II/base/config.h>
21 
24 
26 
27 
30 
31 
77 template <int dim>
79 {
80 public:
85  using SubQuadrature = Quadrature<dim - 1>;
86 
92  static void
93  project_to_face(const SubQuadrature & quadrature,
94  const unsigned int face_no,
95  std::vector<Point<dim>> &q_points);
96 
102  static Quadrature<dim>
103  project_to_face(const SubQuadrature &quadrature, const unsigned int face_no);
104 
114  static void
115  project_to_subface(const SubQuadrature & quadrature,
116  const unsigned int face_no,
117  const unsigned int subface_no,
118  std::vector<Point<dim>> & q_points,
119  const RefinementCase<dim - 1> &ref_case =
121 
131  static Quadrature<dim>
132  project_to_subface(const SubQuadrature & quadrature,
133  const unsigned int face_no,
134  const unsigned int subface_no,
135  const RefinementCase<dim - 1> &ref_case =
137 
154  static Quadrature<dim>
155  project_to_all_faces(const SubQuadrature &quadrature);
156 
170  static Quadrature<dim>
171  project_to_all_subfaces(const SubQuadrature &quadrature);
172 
183  static Quadrature<dim>
184  project_to_child(const Quadrature<dim> &quadrature,
185  const unsigned int child_no);
186 
196  static Quadrature<dim>
197  project_to_all_children(const Quadrature<dim> &quadrature);
198 
203  static Quadrature<dim>
204  project_to_line(const Quadrature<1> &quadrature,
205  const Point<dim> & p1,
206  const Point<dim> & p2);
207 
221  {
222  public:
229 
236  static DataSetDescriptor
237  cell();
238 
249  static DataSetDescriptor
250  face(const unsigned int face_no,
251  const bool face_orientation,
252  const bool face_flip,
253  const bool face_rotation,
254  const unsigned int n_quadrature_points);
255 
268  static DataSetDescriptor
269  subface(const unsigned int face_no,
270  const unsigned int subface_no,
271  const bool face_orientation,
272  const bool face_flip,
273  const bool face_rotation,
274  const unsigned int n_quadrature_points,
275  const internal::SubfaceCase<dim> ref_case =
277 
284  operator unsigned int() const;
285 
286  private:
290  const unsigned int dataset_offset;
291 
296  DataSetDescriptor(const unsigned int dataset_offset);
297  };
298 
299 private:
307  static Quadrature<2>
308  reflect(const Quadrature<2> &q);
309 
319  static Quadrature<2>
320  rotate(const Quadrature<2> &q, const unsigned int n_times);
321 };
322 
326 // ------------------- inline and template functions ----------------
327 
328 
329 
330 template <int dim>
332  const unsigned int dataset_offset)
333  : dataset_offset(dataset_offset)
334 {}
335 
336 
337 template <int dim>
339  : dataset_offset(numbers::invalid_unsigned_int)
340 {}
341 
342 
343 
344 template <int dim>
347 {
348  return 0;
349 }
350 
351 
352 
353 template <int dim>
355 {
356  return dataset_offset;
357 }
358 
359 
360 /* -------------- declaration of explicit specializations ------------- */
361 
362 #ifndef DOXYGEN
363 
364 
365 template <>
366 void
368  const unsigned int,
369  std::vector<Point<1>> &);
370 template <>
371 void
373  const unsigned int face_no,
374  std::vector<Point<2>> &q_points);
375 template <>
376 void
378  const unsigned int face_no,
379  std::vector<Point<3>> &q_points);
380 
381 template <>
384 
385 
386 template <>
387 void
389  const unsigned int,
390  const unsigned int,
391  std::vector<Point<1>> &,
392  const RefinementCase<0> &);
393 template <>
394 void
396  const unsigned int face_no,
397  const unsigned int subface_no,
398  std::vector<Point<2>> &q_points,
399  const RefinementCase<1> &);
400 template <>
401 void
403  const unsigned int face_no,
404  const unsigned int subface_no,
405  std::vector<Point<3>> & q_points,
406  const RefinementCase<2> &face_ref_case);
407 
408 template <>
411 
412 template <>
413 QIterated<1>::QIterated(const Quadrature<1> &base_quadrature,
414  const unsigned int n_copies);
415 
416 
417 #endif // DOXYGEN
419 
420 #endif
QProjector::DataSetDescriptor::DataSetDescriptor
DataSetDescriptor()
Definition: qprojector.h:338
QProjector::reflect
static Quadrature< 2 > reflect(const Quadrature< 2 > &q)
Definition: quadrature.cc:435
QProjector::project_to_all_children
static Quadrature< dim > project_to_all_children(const Quadrature< dim > &quadrature)
Definition: quadrature.cc:1092
QProjector::project_to_all_faces
static Quadrature< dim > project_to_all_faces(const SubQuadrature &quadrature)
QProjector::project_to_face
static void project_to_face(const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim >> &q_points)
QProjector
Definition: qprojector.h:78
geometry_info.h
QProjector::project_to_subface
static void project_to_subface(const SubQuadrature &quadrature, const unsigned int face_no, const unsigned int subface_no, std::vector< Point< dim >> &q_points, const RefinementCase< dim - 1 > &ref_case=RefinementCase< dim - 1 >::isotropic_refinement)
internal::SubfaceCase
Definition: geometry_info.h:1172
QProjector::DataSetDescriptor
Definition: qprojector.h:220
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
QProjector::project_to_child
static Quadrature< dim > project_to_child(const Quadrature< dim > &quadrature, const unsigned int child_no)
Definition: quadrature.cc:1066
QProjector::DataSetDescriptor::subface
static DataSetDescriptor subface(const unsigned int face_no, const unsigned int subface_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points, const internal::SubfaceCase< dim > ref_case=internal::SubfaceCase< dim >::case_isotropic)
RefinementCase
Definition: geometry_info.h:795
QIterated::QIterated
QIterated(const Quadrature< 1 > &base_quadrature, const unsigned int n_copies)
Definition: quadrature.cc:1778
QProjector::rotate
static Quadrature< 2 > rotate(const Quadrature< 2 > &q, const unsigned int n_times)
Definition: quadrature.cc:453
numbers
Definition: numbers.h:207
QProjector::DataSetDescriptor::face
static DataSetDescriptor face(const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points)
Definition: quadrature.cc:1141
QProjector::DataSetDescriptor::dataset_offset
const unsigned int dataset_offset
Definition: qprojector.h:290
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
QProjector::DataSetDescriptor::cell
static DataSetDescriptor cell()
Definition: qprojector.h:346
Point< dim >
QProjector::project_to_all_subfaces
static Quadrature< dim > project_to_all_subfaces(const SubQuadrature &quadrature)
quadrature.h
config.h
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
Quadrature< dim - 1 >
QProjector::project_to_line
static Quadrature< dim > project_to_line(const Quadrature< 1 > &quadrature, const Point< dim > &p1, const Point< dim > &p2)
Definition: quadrature.cc:1118
int