DOLFINx
DOLFINx C++ interface
plaza.h
1// Copyright (C) 2014-2018 Chris Richardson
2//
3// This file is part of DOLFINx (https://www.fenicsproject.org)
4//
5// SPDX-License-Identifier: LGPL-3.0-or-later
6
7#include <cstdint>
8#include <dolfinx/graph/AdjacencyList.h>
9#include <utility>
10#include <vector>
11#include <xtl/xspan.hpp>
12
13#pragma once
14
15namespace dolfinx::mesh
16{
17class Mesh;
18template <typename T>
19class MeshTags;
20} // namespace dolfinx::mesh
21
22namespace dolfinx::refinement
23{
24
28namespace plaza
29{
30
36enum class RefinementOptions : int
37{
38 none = 0,
39 parent_cell = 1,
40 parent_facet = 2,
41 parent_cell_and_facet = 3
42};
43
54std::tuple<mesh::Mesh, std::vector<std::int32_t>, std::vector<std::int8_t>>
55refine(const mesh::Mesh& mesh, bool redistribute, RefinementOptions options);
56
69std::tuple<mesh::Mesh, std::vector<std::int32_t>, std::vector<std::int8_t>>
70refine(const mesh::Mesh& mesh, const xtl::span<const std::int32_t>& edges,
71 bool redistribute, RefinementOptions options);
72
81std::tuple<graph::AdjacencyList<std::int64_t>, xt::xtensor<double, 2>,
82 std::vector<std::int32_t>, std::vector<std::int8_t>>
84
95std::tuple<graph::AdjacencyList<std::int64_t>, xt::xtensor<double, 2>,
96 std::vector<std::int32_t>, std::vector<std::int8_t>>
98 const xtl::span<const std::int32_t>& edges,
99 RefinementOptions options);
100
101} // namespace plaza
102} // namespace dolfinx::refinement
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition: Mesh.h:33
Mesh data structures and algorithms on meshes.
Definition: DofMap.h:30
RefinementOptions
Selection of options when refining a Mesh. parent_cell will output a list containing the local parent...
Definition: plaza.h:37
std::tuple< graph::AdjacencyList< std::int64_t >, xt::xtensor< double, 2 >, std::vector< std::int32_t >, std::vector< std::int8_t > > compute_refinement_data(const mesh::Mesh &mesh, RefinementOptions options)
Refine mesh returning new mesh data.
Definition: plaza.cpp:689
std::tuple< mesh::Mesh, std::vector< std::int32_t >, std::vector< std::int8_t > > refine(const mesh::Mesh &mesh, bool redistribute, RefinementOptions options)
Uniform refine, optionally redistributing and optionally calculating the parent-child relationships,...
Definition: plaza.cpp:619
Mesh refinement algorithms.