Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
immersed_surface_quadrature.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2018 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_non_matching_immersed_surface_quadrature
17 #define dealii_non_matching_immersed_surface_quadrature
18 
19 #include <deal.II/base/config.h>
20 
21 #include <deal.II/base/point.h>
22 #include <deal.II/base/quadrature.h>
23 #include <deal.II/base/subscriptor.h>
24 #include <deal.II/base/tensor.h>
25 
26 #include <vector>
27 
28 DEAL_II_NAMESPACE_OPEN
29 namespace NonMatching
30 {
75  template <int dim>
77  {
78  public:
83  ImmersedSurfaceQuadrature() = default;
84 
90  ImmersedSurfaceQuadrature(const std::vector<Point<dim>> & points,
91  const std::vector<double> & weights,
92  const std::vector<Tensor<1, dim>> &normals);
93 
109  void
110  push_back(const Point<dim> & point,
111  const double weight,
112  const Tensor<1, dim> &normal);
113 
117  const Tensor<1, dim> &
118  normal_vector(const unsigned int i) const;
119 
123  const std::vector<Tensor<1, dim>> &
124  get_normal_vectors() const;
125 
126  protected:
130  std::vector<Tensor<1, dim>> normals;
131  };
132 
133 } // namespace NonMatching
134 DEAL_II_NAMESPACE_CLOSE
135 
136 #endif
std::vector< double > weights
Definition: quadrature.h:290
const Point< dim > & point(const unsigned int i) const
const Tensor< 1, dim > & normal_vector(const unsigned int i) const
const std::vector< Tensor< 1, dim > > & get_normal_vectors() const
void push_back(const Point< dim > &point, const double weight, const Tensor< 1, dim > &normal)
double weight(const unsigned int i) const