Reference documentation for deal.II version 9.4.1
\(\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// Copyright (C) 2009 - 2020 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_distributed_grid_refinement_h
17#define dealii_distributed_grid_refinement_h
18
19
20#include <deal.II/base/config.h>
21
23
25
27
28#include <limits>
29#include <vector>
30
32
33namespace internal
34{
35 namespace parallel
36 {
37 namespace distributed
38 {
40 {
46 template <typename number>
47 std::pair<number, number>
49 const ::Vector<number> &criteria,
50 const MPI_Comm & mpi_communicator);
51
52 namespace RefineAndCoarsenFixedNumber
53 {
58 template <typename number>
59 number
60 compute_threshold(const ::Vector<number> & criteria,
61 const std::pair<double, double> &global_min_and_max,
62 const types::global_cell_index n_target_cells,
63 const MPI_Comm & mpi_communicator);
64 } // namespace RefineAndCoarsenFixedNumber
65
66 namespace RefineAndCoarsenFixedFraction
67 {
74 template <typename number>
75 number
76 compute_threshold(const ::Vector<number> & criteria,
77 const std::pair<double, double> &global_min_and_max,
78 const double target_error,
79 const MPI_Comm & mpi_communicator);
80 } // namespace RefineAndCoarsenFixedFraction
81 } // namespace GridRefinement
82 } // namespace distributed
83 } // namespace parallel
84} // namespace internal
85
86
87
88namespace parallel
89{
90 namespace distributed
91 {
107 {
152 template <int dim, typename Number, int spacedim>
153 void
156 const ::Vector<Number> & criteria,
157 const double top_fraction_of_cells,
158 const double bottom_fraction_of_cells,
159 const types::global_cell_index max_n_cells =
160 std::numeric_limits<types::global_cell_index>::max());
161
205 template <int dim, typename Number, int spacedim>
206 void
209 const ::Vector<Number> & criteria,
210 const double top_fraction_of_error,
211 const double bottom_fraction_of_error,
213 } // namespace GridRefinement
214 } // namespace distributed
215} // namespace parallel
216
217
219
220#endif // dealii_distributed_grid_refinement_h
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
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_number(parallel::distributed::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())
void refine_and_coarsen_fixed_fraction(parallel::distributed::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::NormType::L1_norm)
unsigned int global_cell_index
Definition: types.h:105
const ::Triangulation< dim, spacedim > & tria