Reference documentation for deal.II version Git f15f581df6 2020-07-10 15:30:09 -0400
\(\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 
75 template <int dim>
77 {
78 public:
83  using SubQuadrature = Quadrature<dim - 1>;
84 
90  static void
91  project_to_face(const SubQuadrature & quadrature,
92  const unsigned int face_no,
93  std::vector<Point<dim>> &q_points);
94 
100  static Quadrature<dim>
101  project_to_face(const SubQuadrature &quadrature, const unsigned int face_no);
102 
112  static void
113  project_to_subface(const SubQuadrature & quadrature,
114  const unsigned int face_no,
115  const unsigned int subface_no,
116  std::vector<Point<dim>> & q_points,
117  const RefinementCase<dim - 1> &ref_case =
119 
129  static Quadrature<dim>
130  project_to_subface(const SubQuadrature & quadrature,
131  const unsigned int face_no,
132  const unsigned int subface_no,
133  const RefinementCase<dim - 1> &ref_case =
135 
152  static Quadrature<dim>
153  project_to_all_faces(const SubQuadrature &quadrature);
154 
168  static Quadrature<dim>
169  project_to_all_subfaces(const SubQuadrature &quadrature);
170 
181  static Quadrature<dim>
182  project_to_child(const Quadrature<dim> &quadrature,
183  const unsigned int child_no);
184 
194  static Quadrature<dim>
195  project_to_all_children(const Quadrature<dim> &quadrature);
196 
201  static Quadrature<dim>
202  project_to_line(const Quadrature<1> &quadrature,
203  const Point<dim> & p1,
204  const Point<dim> & p2);
205 
217  {
218  public:
225 
232  static DataSetDescriptor
233  cell();
234 
245  static DataSetDescriptor
246  face(const unsigned int face_no,
247  const bool face_orientation,
248  const bool face_flip,
249  const bool face_rotation,
250  const unsigned int n_quadrature_points);
251 
264  static DataSetDescriptor
265  subface(const unsigned int face_no,
266  const unsigned int subface_no,
267  const bool face_orientation,
268  const bool face_flip,
269  const bool face_rotation,
270  const unsigned int n_quadrature_points,
271  const internal::SubfaceCase<dim> ref_case =
273 
280  operator unsigned int() const;
281 
282  private:
286  const unsigned int dataset_offset;
287 
292  DataSetDescriptor(const unsigned int dataset_offset);
293  };
294 
295 private:
303  static Quadrature<2>
304  reflect(const Quadrature<2> &q);
305 
315  static Quadrature<2>
316  rotate(const Quadrature<2> &q, const unsigned int n_times);
317 };
318 
322 // ------------------- inline and template functions ----------------
323 
324 
325 
326 template <int dim>
328  const unsigned int dataset_offset)
329  : dataset_offset(dataset_offset)
330 {}
331 
332 
333 template <int dim>
336 {}
337 
338 
339 
340 template <int dim>
343 {
344  return 0;
345 }
346 
347 
348 
349 template <int dim>
351 {
352  return dataset_offset;
353 }
354 
355 
356 /* -------------- declaration of explicit specializations ------------- */
357 
358 #ifndef DOXYGEN
359 
360 
361 template <>
362 void
364  const unsigned int,
365  std::vector<Point<1>> &);
366 template <>
367 void
369  const unsigned int face_no,
370  std::vector<Point<2>> &q_points);
371 template <>
372 void
374  const unsigned int face_no,
375  std::vector<Point<3>> &q_points);
376 
377 template <>
380 
381 
382 template <>
383 void
385  const unsigned int,
386  const unsigned int,
387  std::vector<Point<1>> &,
388  const RefinementCase<0> &);
389 template <>
390 void
392  const unsigned int face_no,
393  const unsigned int subface_no,
394  std::vector<Point<2>> &q_points,
395  const RefinementCase<1> &);
396 template <>
397 void
399  const unsigned int face_no,
400  const unsigned int subface_no,
401  std::vector<Point<3>> & q_points,
402  const RefinementCase<2> &face_ref_case);
403 
404 template <>
407 
408 template <>
409 QIterated<1>::QIterated(const Quadrature<1> &base_quadrature,
410  const unsigned int n_copies);
411 
412 
413 #endif // DOXYGEN
415 
416 #endif
static const unsigned int invalid_unsigned_int
Definition: types.h:191
static Quadrature< dim > project_to_child(const Quadrature< dim > &quadrature, const unsigned int child_no)
Definition: quadrature.cc:1064
const unsigned int dataset_offset
Definition: qprojector.h:286
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)
static Quadrature< 2 > rotate(const Quadrature< 2 > &q, const unsigned int n_times)
Definition: quadrature.cc:451
static DataSetDescriptor cell()
Definition: qprojector.h:342
static Quadrature< dim > project_to_all_subfaces(const SubQuadrature &quadrature)
static Quadrature< dim > project_to_line(const Quadrature< 1 > &quadrature, const Point< dim > &p1, const Point< dim > &p2)
Definition: quadrature.cc:1116
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:363
static Quadrature< 2 > reflect(const Quadrature< 2 > &q)
Definition: quadrature.cc:433
static Quadrature< dim > project_to_all_faces(const SubQuadrature &quadrature)
static void project_to_face(const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim >> &q_points)
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:362
QIterated(const Quadrature< 1 > &base_quadrature, const unsigned int n_copies)
Definition: quadrature.cc:1776
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)
static Quadrature< dim > project_to_all_children(const Quadrature< dim > &quadrature)
Definition: quadrature.cc:1090
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:1139