Reference documentation for deal.II version 9.3.3
\(\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\}}\)
quadrature.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 1998 - 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_quadrature_h
17#define dealii_quadrature_h
18
19
20#include <deal.II/base/config.h>
21
22#include <deal.II/base/point.h>
24
25#include <array>
26#include <memory>
27#include <vector>
28
30
33
82template <int dim>
83class Quadrature : public Subscriptor
84{
85public:
90 using SubQuadrature = Quadrature<dim - 1>;
91
100 explicit Quadrature(const unsigned int n_quadrature_points = 0);
101
110 Quadrature(const SubQuadrature &, const Quadrature<1> &);
111
128 explicit Quadrature(const Quadrature<dim != 1 ? 1 : 0> &quadrature_1d);
129
133 Quadrature(const Quadrature<dim> &q);
134
139 Quadrature(Quadrature<dim> &&) noexcept = default;
140
146 Quadrature(const std::vector<Point<dim>> &points,
147 const std::vector<double> & weights);
148
156 Quadrature(const std::vector<Point<dim>> &points);
157
162 Quadrature(const Point<dim> &point);
163
167 virtual ~Quadrature() override = default;
168
173 Quadrature &
174 operator=(const Quadrature<dim> &);
175
180 Quadrature &
181 operator=(Quadrature<dim> &&) = default; // NOLINT
182
186 bool
187 operator==(const Quadrature<dim> &p) const;
188
193 void
194 initialize(const std::vector<Point<dim>> &points,
195 const std::vector<double> & weights);
196
200 unsigned int
201 size() const;
202
206 const Point<dim> &
207 point(const unsigned int i) const;
208
212 const std::vector<Point<dim>> &
213 get_points() const;
214
218 double
219 weight(const unsigned int i) const;
220
224 const std::vector<double> &
225 get_weights() const;
226
231 std::size_t
232 memory_consumption() const;
233
239 template <class Archive>
240 void
241 serialize(Archive &ar, const unsigned int version);
242
248 bool
250
269#ifndef DOXYGEN
270 typename std::conditional<dim == 1,
271 std::array<Quadrature<1>, dim>,
272 const std::array<Quadrature<1>, dim> &>::type
273#else
274 const std::array<Quadrature<1>, dim> &
275#endif
276 get_tensor_basis() const;
277
278protected:
283 std::vector<Point<dim>> quadrature_points;
284
289 std::vector<double> weights;
290
299
304 std::unique_ptr<std::array<Quadrature<1>, dim>> tensor_basis;
305};
306
307
316template <int dim>
317class QAnisotropic : public Quadrature<dim>
318{
319public:
324 QAnisotropic(const Quadrature<1> &qx);
325
329 QAnisotropic(const Quadrature<1> &qx, const Quadrature<1> &qy);
330
334 QAnisotropic(const Quadrature<1> &qx,
335 const Quadrature<1> &qy,
336 const Quadrature<1> &qz);
337};
338
339
363template <int dim>
364class QIterated : public Quadrature<dim>
365{
366public:
371 QIterated(const Quadrature<1> &base_quadrature, const unsigned int n_copies);
372
377 "The quadrature formula you provided cannot be used "
378 "as the basis for iteration.");
379};
380
381
382
385#ifndef DOXYGEN
386
387// ------------------- inline and template functions ----------------
388
389
390template <int dim>
391inline unsigned int
393{
394 return weights.size();
395}
396
397
398template <int dim>
399inline const Point<dim> &
400Quadrature<dim>::point(const unsigned int i) const
401{
402 AssertIndexRange(i, size());
403 return quadrature_points[i];
404}
405
406
407
408template <int dim>
409double
410Quadrature<dim>::weight(const unsigned int i) const
411{
412 AssertIndexRange(i, size());
413 return weights[i];
414}
415
416
417
418template <int dim>
419inline const std::vector<Point<dim>> &
421{
422 return quadrature_points;
423}
424
425
426
427template <int dim>
428inline const std::vector<double> &
430{
431 return weights;
432}
433
434
435
436template <int dim>
437inline bool
439{
440 return is_tensor_product_flag;
441}
442
443
444
445template <int dim>
446template <class Archive>
447inline void
448Quadrature<dim>::serialize(Archive &ar, const unsigned int)
449{
450 // forward to serialization
451 // function in the base class.
452 ar &static_cast<Subscriptor &>(*this);
453
454 ar &quadrature_points &weights;
455}
456
457
458
459/* -------------- declaration of explicit specializations ------------- */
460
461template <>
462Quadrature<0>::Quadrature(const unsigned int);
463template <>
465template <>
467template <>
469
470template <>
472
473template <>
475
476template <>
477QIterated<1>::QIterated(const Quadrature<1> &base_quadrature,
478 const unsigned int n_copies);
479
480#endif // DOXYGEN
482
483#endif
Definition: point.h:111
QAnisotropic(const Quadrature< 1 > &qx)
Definition: quadrature.cc:351
QIterated(const Quadrature< 1 > &base_quadrature, const unsigned int n_copies)
Definition: quadrature.cc:561
std::vector< Point< dim > > quadrature_points
Definition: quadrature.h:283
std::unique_ptr< std::array< Quadrature< 1 >, dim > > tensor_basis
Definition: quadrature.h:304
bool is_tensor_product() const
std::size_t memory_consumption() const
Definition: quadrature.cc:311
const Point< dim > & point(const unsigned int i) const
bool is_tensor_product_flag
Definition: quadrature.h:298
Quadrature(const unsigned int n_quadrature_points=0)
Definition: quadrature.cc:40
const std::array< Quadrature< 1 >, dim > & get_tensor_basis() const
Definition: quadrature.cc:323
double weight(const unsigned int i) const
const std::vector< double > & get_weights() const
void initialize(const std::vector< Point< dim > > &points, const std::vector< double > &weights)
Definition: quadrature.cc:50
Quadrature(Quadrature< dim > &&) noexcept=default
std::vector< double > weights
Definition: quadrature.h:289
const std::vector< Point< dim > > & get_points() const
unsigned int size() const
void serialize(Archive &ar, const unsigned int version)
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
static ::ExceptionBase & ExcInvalidQuadratureFormula()
#define AssertIndexRange(index, range)
Definition: exceptions.h:1690
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:493
void quadrature_points(const Triangulation< dim, spacedim > &triangulation, const Quadrature< dim > &quadrature, const std::vector< std::vector< BoundingBox< spacedim > > > &global_bounding_boxes, ParticleHandler< dim, spacedim > &particle_handler, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()), const std::vector< std::vector< double > > &properties={})
Definition: generators.cc:451
STL namespace.