16 #ifndef dealii_quadrature_h
17 #define dealii_quadrature_h
140 explicit Quadrature(
const unsigned int n_quadrature_points = 0);
190 const std::vector<
double> &
weights);
245 const std::vector<
double> &
weights);
268 const std::vector<
Point<dim>> &
280 const std::vector<
double> &
295 template <class Archive>
326 typename std::conditional<dim == 1,
327 std::array<Quadrature<1>, dim>,
328 const std::array<Quadrature<1>, dim> &>::type
330 const std::array<Quadrature<1>, dim> &
443 const std::vector<
Point<1>> &intervals);
449 "The quadrature formula you provided cannot be used "
450 "as the basis for iteration.");
466 return weights.size();
475 return weights.empty();
501 inline const std::vector<Point<dim>> &
510 inline const std::vector<double> &
522 return is_tensor_product_flag;
528 template <
class Archive>
561 const unsigned int n_copies);
QAnisotropic(const Quadrature< 1 > &qx)
QIterated(const Quadrature< 1 > &base_quadrature, const unsigned int n_copies)
std::vector< Point< dim > > quadrature_points
std::unique_ptr< std::array< Quadrature< 1 >, dim > > tensor_basis
const std::vector< Point< dim > > & get_points() const
bool is_tensor_product() const
std::size_t memory_consumption() const
void initialize(const std::vector< Point< dim >> &points, const std::vector< double > &weights)
bool is_tensor_product_flag
Quadrature(const unsigned int n_quadrature_points=0)
const std::vector< double > & get_weights() const
double weight(const unsigned int i) const
const std::array< Quadrature< 1 >, dim > & get_tensor_basis() const
const Point< dim > & point(const unsigned int i) const
Quadrature(Quadrature< dim > &&) noexcept=default
std::vector< double > weights
unsigned int size() const
void serialize(Archive &ar, const unsigned int version)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInvalidQuadratureFormula()
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={})