Reference documentation for deal.II version GIT bef9f326eb 2022-07-03 20:25:02+00:00
\(\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 - 2021 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 
28 
30 
31 
34 
35 
79 template <int dim>
81 {
82 public:
87  using SubQuadrature = Quadrature<dim - 1>;
88 
94  static void
96  const SubQuadrature & quadrature,
97  const unsigned int face_no,
98  std::vector<Point<dim>> &q_points);
99 
105  static Quadrature<dim>
107  const SubQuadrature &quadrature,
108  const unsigned int face_no);
109 
119  static void
121  const SubQuadrature & quadrature,
122  const unsigned int face_no,
123  const unsigned int subface_no,
124  std::vector<Point<dim>> & q_points,
125  const RefinementCase<dim - 1> &ref_case =
127 
141  static Quadrature<dim>
143  const SubQuadrature & quadrature,
144  const unsigned int face_no,
145  const unsigned int subface_no,
146  const RefinementCase<dim - 1> &ref_case =
148 
165  static Quadrature<dim>
167  const hp::QCollection<dim - 1> &quadrature);
168 
173  static Quadrature<dim>
175  const Quadrature<dim - 1> &quadrature);
176 
190  static Quadrature<dim>
192  const SubQuadrature &quadrature);
193 
204  static Quadrature<dim>
206  const Quadrature<dim> &quadrature,
207  const unsigned int child_no);
208 
218  static Quadrature<dim>
220  const Quadrature<dim> &quadrature);
221 
226  static Quadrature<dim>
228  const Quadrature<1> &quadrature,
229  const Point<dim> & p1,
230  const Point<dim> & p2);
231 
243  {
244  public:
251 
258  static DataSetDescriptor
259  cell();
260 
271  static DataSetDescriptor
273  const unsigned int face_no,
274  const bool face_orientation,
275  const bool face_flip,
276  const bool face_rotation,
277  const unsigned int n_quadrature_points);
278 
283  static DataSetDescriptor
285  const unsigned int face_no,
286  const bool face_orientation,
287  const bool face_flip,
288  const bool face_rotation,
289  const hp::QCollection<dim - 1> &quadrature);
290 
303  static DataSetDescriptor
305  const unsigned int face_no,
306  const unsigned int subface_no,
307  const bool face_orientation,
308  const bool face_flip,
309  const bool face_rotation,
310  const unsigned int n_quadrature_points,
311  const internal::SubfaceCase<dim> ref_case =
313 
320  operator unsigned int() const;
321 
322  private:
326  const unsigned int dataset_offset;
327 
332  DataSetDescriptor(const unsigned int dataset_offset);
333  };
334 };
335 
339 // ------------------- inline and template functions ----------------
340 
341 
342 
343 template <int dim>
345  const unsigned int dataset_offset)
346  : dataset_offset(dataset_offset)
347 {}
348 
349 
350 template <int dim>
352  : dataset_offset(numbers::invalid_unsigned_int)
353 {}
354 
355 
356 
357 template <int dim>
360 {
361  return 0;
362 }
363 
364 
365 
366 template <int dim>
368 {
369  return dataset_offset;
370 }
371 
372 
373 
374 template <int dim>
377  const Quadrature<dim - 1> &quadrature)
378 {
380  hp::QCollection<dim - 1>(quadrature));
381 }
382 
383 
384 /* -------------- declaration of explicit specializations ------------- */
385 
386 #ifndef DOXYGEN
387 
388 
389 template <>
390 void
392  const Quadrature<0> &,
393  const unsigned int,
394  std::vector<Point<1>> &);
395 template <>
396 void
398  const Quadrature<1> & quadrature,
399  const unsigned int face_no,
400  std::vector<Point<2>> &q_points);
401 template <>
402 void
404  const Quadrature<2> & quadrature,
405  const unsigned int face_no,
406  std::vector<Point<3>> &q_points);
407 
408 template <>
411  const hp::QCollection<0> &quadrature);
412 
413 
414 template <>
415 void
417  const Quadrature<0> &,
418  const unsigned int,
419  const unsigned int,
420  std::vector<Point<1>> &,
421  const RefinementCase<0> &);
422 template <>
423 void
425  const Quadrature<1> & quadrature,
426  const unsigned int face_no,
427  const unsigned int subface_no,
428  std::vector<Point<2>> &q_points,
429  const RefinementCase<1> &);
430 template <>
431 void
433  const Quadrature<2> & quadrature,
434  const unsigned int face_no,
435  const unsigned int subface_no,
436  std::vector<Point<3>> & q_points,
437  const RefinementCase<2> &face_ref_case);
438 
439 template <>
442  const Quadrature<0> &quadrature);
443 
444 
445 #endif // DOXYGEN
447 
448 #endif
Definition: point.h:111
const unsigned int dataset_offset
Definition: qprojector.h:326
static DataSetDescriptor subface(const ReferenceCell reference_cell, 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 DataSetDescriptor cell()
Definition: qprojector.h:359
static DataSetDescriptor face(const ReferenceCell reference_cell, const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points)
Definition: qprojector.cc:1222
static Quadrature< dim > project_to_line(const ReferenceCell reference_cell, const Quadrature< 1 > &quadrature, const Point< dim > &p1, const Point< dim > &p2)
Definition: qprojector.cc:1194
static Quadrature< dim > project_to_child(const ReferenceCell reference_cell, const Quadrature< dim > &quadrature, const unsigned int child_no)
Definition: qprojector.cc:1130
static void project_to_subface(const ReferenceCell reference_cell, 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< dim > project_to_all_subfaces(const ReferenceCell reference_cell, const SubQuadrature &quadrature)
static void project_to_face(const ReferenceCell reference_cell, const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim >> &q_points)
static Quadrature< dim > project_to_all_faces(const ReferenceCell reference_cell, const hp::QCollection< dim - 1 > &quadrature)
static Quadrature< dim > project_to_all_children(const ReferenceCell reference_cell, const Quadrature< dim > &quadrature)
Definition: qprojector.cc:1162
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
static const unsigned int invalid_unsigned_int
Definition: types.h:201