Reference documentation for deal.II version 9.6.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// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2022 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_cgal_intersections_h
16#define dealii_cgal_intersections_h
17
18#include <deal.II/base/config.h>
19
21#include <deal.II/base/point.h>
22
23#include <deal.II/fe/mapping.h>
24
25#include <deal.II/grid/tria.h>
26
27
29namespace CGALWrappers
30{
52 template <int structdim0, int structdim1, int spacedim>
53 std::vector<std::array<Point<spacedim>, structdim1 + 1>>
57 const Mapping<structdim0, spacedim> &mapping0,
58 const Mapping<structdim1, spacedim> &mapping1,
59 const double tol = 1e-9);
60
61
70 template <int structdim0, int structdim1, int spacedim>
71 std::vector<std::array<Point<spacedim>, structdim1 + 1>>
73 const ArrayView<const Point<spacedim>> &vertices0,
74 const ArrayView<const Point<spacedim>> &vertices1,
75 const double tol = 1e-9);
76
77} // namespace CGALWrappers
78
80
81#endif
Abstract base class for mapping classes.
Definition mapping.h:318
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
std::vector< std::array< Point< spacedim >, structdim1+1 > > compute_intersection_of_cells(const typename Triangulation< structdim0, spacedim >::cell_iterator &cell0, const typename Triangulation< structdim1, spacedim >::cell_iterator &cell1, const Mapping< structdim0, spacedim > &mapping0, const Mapping< structdim1, spacedim > &mapping1, const double tol=1e-9)