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\}}\)
immersed_surface_quadrature.cc
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 1998 - 2019 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
17
19namespace NonMatching
20{
21 template <int dim>
23 const std::vector<Point<dim>> & points,
24 const std::vector<double> & weights,
25 const std::vector<Tensor<1, dim>> &normals)
26 : Quadrature<dim>(points, weights)
27 , normals(normals)
28 {
29 AssertDimension(weights.size(), points.size());
30 AssertDimension(normals.size(), points.size());
31 for (const auto &normal : normals)
32 {
33 (void)normal;
34 Assert(std::abs(normal.norm() - 1.0) < 1e-9,
35 ExcMessage("Normal is not normalized."));
36 }
37 }
38
39
40
41 template <int dim>
42 void
44 const double weight,
45 const Tensor<1, dim> &normal)
46 {
47 this->quadrature_points.push_back(point);
48 this->weights.push_back(weight);
49 this->normals.push_back(normal);
50 Assert(std::abs(normal.norm() - 1.0) < 1e-9,
51 ExcMessage("Normal is not normalized."));
52 }
53
54
55
56 template <int dim>
57 const Tensor<1, dim> &
59 {
60 AssertIndexRange(i, this->size());
61 return normals[i];
62 }
63
64
65
66 template <int dim>
67 const std::vector<Tensor<1, dim>> &
69 {
70 return normals;
71 }
72
73
74
75 template class ImmersedSurfaceQuadrature<1>;
76 template class ImmersedSurfaceQuadrature<2>;
77 template class ImmersedSurfaceQuadrature<3>;
78
79} // namespace NonMatching
const Tensor< 1, dim > & normal_vector(const unsigned int i) const
void push_back(const Point< dim > &point, const double weight, const Tensor< 1, dim > &normal)
const std::vector< Tensor< 1, dim > > & get_normal_vectors() const
Definition: point.h:111
std::vector< double > weights
Definition: quadrature.h:289
numbers::NumberTraits< Number >::real_type norm() const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define Assert(cond, exc)
Definition: exceptions.h:1465
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1622
#define AssertIndexRange(index, range)
Definition: exceptions.h:1690
static ::ExceptionBase & ExcMessage(std::string arg1)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition: utilities.cc:188
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
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)