Reference documentation for deal.II version GIT relicensing-422-gb369f187d8 2024-04-17 18:10:02+00:00
\(\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
grid_refinement.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2010 - 2023 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_distributed_grid_refinement_h
16#define dealii_distributed_grid_refinement_h
17
18
19#include <deal.II/base/config.h>
20
22
24
26
27#include <limits>
28#include <vector>
29
31
32namespace internal
33{
34 namespace parallel
35 {
36 namespace distributed
37 {
39 {
45 template <typename number>
46 std::pair<number, number>
48 const ::Vector<number> &criteria,
49 const MPI_Comm mpi_communicator);
50
51 namespace RefineAndCoarsenFixedNumber
52 {
57 template <typename number>
58 number
59 compute_threshold(const ::Vector<number> &criteria,
60 const std::pair<double, double> &global_min_and_max,
61 const types::global_cell_index n_target_cells,
62 const MPI_Comm mpi_communicator);
63 } // namespace RefineAndCoarsenFixedNumber
64
65 namespace RefineAndCoarsenFixedFraction
66 {
75 template <typename number>
76 number
77 compute_threshold(const ::Vector<number> &criteria,
78 const std::pair<double, double> &global_min_and_max,
79 const double target_error,
80 const MPI_Comm mpi_communicator);
81 } // namespace RefineAndCoarsenFixedFraction
82 } // namespace GridRefinement
83 } // namespace distributed
84 } // namespace parallel
85} // namespace internal
86
87
88
89namespace parallel
90{
91 namespace distributed
92 {
109 {
155 template <int dim, typename Number, int spacedim>
156 void
159 const ::Vector<Number> &criteria,
160 const double top_fraction_of_cells,
161 const double bottom_fraction_of_cells,
162 const types::global_cell_index max_n_cells =
163 std::numeric_limits<types::global_cell_index>::max());
164
209 template <int dim, typename Number, int spacedim>
210 void
213 const ::Vector<Number> &criteria,
214 const double top_fraction_of_error,
215 const double bottom_fraction_of_error,
217 } // namespace GridRefinement
218 } // namespace distributed
219} // namespace parallel
220
221
223
224#endif // dealii_distributed_grid_refinement_h
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
number compute_threshold(const ::Vector< number > &criteria, const std::pair< double, double > &global_min_and_max, const double target_error, const MPI_Comm mpi_communicator)
number compute_threshold(const ::Vector< number > &criteria, const std::pair< double, double > &global_min_and_max, const types::global_cell_index n_target_cells, const MPI_Comm mpi_communicator)
std::pair< number, number > compute_global_min_and_max_at_root(const ::Vector< number > &criteria, const MPI_Comm mpi_communicator)
void refine_and_coarsen_fixed_fraction(::Triangulation< dim, spacedim > &tria, const ::Vector< Number > &criteria, const double top_fraction_of_error, const double bottom_fraction_of_error, const VectorTools::NormType norm_type=VectorTools::L1_norm)
void refine_and_coarsen_fixed_number(::Triangulation< dim, spacedim > &tria, const ::Vector< Number > &criteria, const double top_fraction_of_cells, const double bottom_fraction_of_cells, const types::global_cell_index max_n_cells=std::numeric_limits< types::global_cell_index >::max())