Reference documentation for deal.II version 8.5.1
qprojector.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2005 - 2015 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__qprojector_h
17 #define dealii__qprojector_h
18 
19 
20 #include <deal.II/base/quadrature.h>
21 #include <deal.II/base/geometry_info.h>
22 
23 DEAL_II_NAMESPACE_OPEN
24 
25 
28 
29 
75 template <int dim>
77 {
78 public:
83  typedef Quadrature<dim-1> SubQuadrature;
84 
90  static void project_to_face (const SubQuadrature &quadrature,
91  const unsigned int face_no,
92  std::vector<Point<dim> > &q_points);
93 
99  static Quadrature<dim>
100  project_to_face (const SubQuadrature &quadrature,
101  const unsigned int face_no);
102 
112  static void project_to_subface (const SubQuadrature &quadrature,
113  const unsigned int face_no,
114  const unsigned int subface_no,
115  std::vector<Point<dim> > &q_points,
117 
127  static Quadrature<dim>
128  project_to_subface (const SubQuadrature &quadrature,
129  const unsigned int face_no,
130  const unsigned int subface_no,
132 
149  static Quadrature<dim>
150  project_to_all_faces (const SubQuadrature &quadrature);
151 
165  static Quadrature<dim>
166  project_to_all_subfaces (const SubQuadrature &quadrature);
167 
178  static
180  project_to_child (const Quadrature<dim> &quadrature,
181  const unsigned int child_no);
182 
192  static
194  project_to_all_children (const Quadrature<dim> &quadrature);
195 
200  static
202  project_to_line(const Quadrature<1> &quadrature,
203  const Point<dim> &p1,
204  const Point<dim> &p2);
205 
219  {
220  public:
227 
234  static DataSetDescriptor cell ();
235 
246  static
248  face (const unsigned int face_no,
249  const bool face_orientation,
250  const bool face_flip,
251  const bool face_rotation,
252  const unsigned int n_quadrature_points);
253 
266  static
268  subface (const unsigned int face_no,
269  const unsigned int subface_no,
270  const bool face_orientation,
271  const bool face_flip,
272  const bool face_rotation,
273  const unsigned int n_quadrature_points,
275 
282  operator unsigned int () const;
283 
284  private:
288  const unsigned int dataset_offset;
289 
294  DataSetDescriptor (const unsigned int dataset_offset);
295  };
296 
297 private:
305  static Quadrature<2> reflect (const Quadrature<2> &q);
306 
316  static Quadrature<2> rotate (const Quadrature<2> &q,
317  const unsigned int n_times);
318 };
319 
323 // ------------------- inline and template functions ----------------
324 
325 
326 
327 template <int dim>
328 inline
330 DataSetDescriptor (const unsigned int dataset_offset)
331  :
332  dataset_offset (dataset_offset)
333 {}
334 
335 
336 template <int dim>
337 inline
340  :
341  dataset_offset (numbers::invalid_unsigned_int)
342 {}
343 
344 
345 
346 template <int dim>
349 {
350  return 0;
351 }
352 
353 
354 
355 template <int dim>
356 inline
358 {
359  return dataset_offset;
360 }
361 
362 
363 /* -------------- declaration of explicit specializations ------------- */
364 
365 #ifndef DOXYGEN
366 
367 
368 template <>
369 void
371  const unsigned int,
372  std::vector<Point<1> > &);
373 template <>
374 void
376  const unsigned int face_no,
377  std::vector<Point<2> > &q_points);
378 template <>
379 void
381  const unsigned int face_no,
382  std::vector<Point<3> > &q_points);
383 
384 template <>
387 
388 
389 template <>
390 void
392  const unsigned int,
393  const unsigned int,
394  std::vector<Point<1> > &,
395  const RefinementCase<0> &);
396 template <>
397 void
399  const unsigned int face_no,
400  const unsigned int subface_no,
401  std::vector<Point<2> > &q_points,
402  const RefinementCase<1> &);
403 template <>
404 void
406  const unsigned int face_no,
407  const unsigned int subface_no,
408  std::vector<Point<3> > &q_points,
409  const RefinementCase<2> &face_ref_case);
410 
411 template <>
414 
415 
416 template <>
417 bool
418 QIterated<1>::uses_both_endpoints (const Quadrature<1> &base_quadrature);
419 
420 template <>
421 QIterated<1>::QIterated (const Quadrature<1> &base_quadrature,
422  const unsigned int n_copies);
423 
424 
425 #endif // DOXYGEN
426 DEAL_II_NAMESPACE_CLOSE
427 
428 #endif
static Quadrature< dim > project_to_child(const Quadrature< dim > &quadrature, const unsigned int child_no)
Definition: quadrature.cc:1003
const unsigned int dataset_offset
Definition: qprojector.h:288
static Quadrature< 2 > rotate(const Quadrature< 2 > &q, const unsigned int n_times)
Definition: quadrature.cc:372
static bool uses_both_endpoints(const Quadrature< 1 > &base_quadrature)
Quadrature< dim-1 > SubQuadrature
Definition: qprojector.h:83
static DataSetDescriptor cell()
Definition: qprojector.h:348
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:1056
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 > reflect(const Quadrature< 2 > &q)
Definition: quadrature.cc:354
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)
QIterated(const Quadrature< 1 > &base_quadrature, const unsigned int n_copies)
Definition: quadrature.cc:1758
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:1029
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:1081