Reference documentation for deal.II version 9.5.0
\(\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\}}\)
Loading...
Searching...
No Matches
intersections.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2022 - 2023 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_cgal_intersections_h
17#define dealii_cgal_intersections_h
18
19#include <deal.II/base/config.h>
20
22#include <deal.II/base/point.h>
23
24#include <deal.II/fe/mapping.h>
25
26#include <deal.II/grid/tria.h>
27
28
30namespace CGALWrappers
31{
53 template <int dim0, int dim1, int spacedim>
54 std::vector<std::array<Point<spacedim>, dim1 + 1>>
58 const Mapping<dim0, spacedim> & mapping0,
59 const Mapping<dim1, spacedim> & mapping1,
60 const double tol = 1e-9);
61
62
71 template <int dim0, int dim1, int spacedim>
72 std::vector<std::array<Point<spacedim>, dim1 + 1>>
74 const ArrayView<const Point<spacedim>> &vertices0,
75 const ArrayView<const Point<spacedim>> &vertices1,
76 const double tol = 1e-9);
77
78} // namespace CGALWrappers
79
81
82#endif
Abstract base class for mapping classes.
Definition mapping.h:317
Definition point.h:112
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:472
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473
std::vector< std::array< Point< spacedim >, dim1+1 > > compute_intersection_of_cells(const typename Triangulation< dim0, spacedim >::cell_iterator &cell0, const typename Triangulation< dim1, spacedim >::cell_iterator &cell1, const Mapping< dim0, spacedim > &mapping0, const Mapping< dim1, spacedim > &mapping1, const double tol=1e-9)