Reference documentation for deal.II version GIT c14369f203 2023-10-01 07:40: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\}}\)
grid_refinement.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 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_grid_refinement_h
17 #define dealii_grid_refinement_h
18 
19 
20 #include <deal.II/base/config.h>
21 
23 
25 
26 #include <limits>
27 
29 
30 // forward declarations
31 #ifndef DOXYGEN
32 template <int dim, int spacedim>
33 DEAL_II_CXX20_REQUIRES((concepts::is_valid_dim_spacedim<dim, spacedim>))
34 class Triangulation;
35 template <typename Number>
36 class Vector;
37 #endif
38 
56 namespace GridRefinement
57 {
89  template <int dim>
90  std::pair<double, double>
92  const types::global_cell_index current_n_cells,
93  const types::global_cell_index max_n_cells,
94  const double top_fraction_of_cells,
95  const double bottom_fraction_of_cells);
96 
161  template <int dim, typename Number, int spacedim>
162  void
165  const Vector<Number> &criteria,
166  const double top_fraction_of_cells,
167  const double bottom_fraction_of_cells,
168  const unsigned int max_n_cells = std::numeric_limits<unsigned int>::max());
169 
232  template <int dim, typename Number, int spacedim>
233  void
236  const Vector<Number> &criteria,
237  const double top_fraction,
238  const double bottom_fraction,
239  const unsigned int max_n_cells = std::numeric_limits<unsigned int>::max(),
241 
242 
243 
315  template <int dim, typename Number, int spacedim>
316  void
318  const Vector<Number> &criteria,
319  const unsigned int order = 2);
320 
335  template <int dim, typename Number, int spacedim>
336  void
338  const Vector<Number> &criteria,
339  const double threshold,
340  const unsigned int max_to_mark = numbers::invalid_unsigned_int);
341 
356  template <int dim, typename Number, int spacedim>
357  void
359  const Vector<Number> &criteria,
360  const double threshold);
361 
367 
373 } // namespace GridRefinement
374 
375 
376 
378 
379 #endif // dealii_grid_refinement_h
Definition: vector.h:110
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:477
#define DEAL_II_CXX20_REQUIRES(condition)
Definition: config.h:166
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:478
#define DeclException0(Exception0)
Definition: exceptions.h:467
static ::ExceptionBase & ExcNegativeCriteria()
static ::ExceptionBase & ExcInvalidParameterValue()
void refine(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold, const unsigned int max_to_mark=numbers::invalid_unsigned_int)
void refine_and_coarsen_fixed_number(Triangulation< dim, spacedim > &triangulation, const Vector< Number > &criteria, const double top_fraction_of_cells, const double bottom_fraction_of_cells, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max())
void refine_and_coarsen_optimize(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const unsigned int order=2)
void coarsen(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold)
std::pair< double, double > adjust_refine_and_coarsen_number_fraction(const types::global_cell_index current_n_cells, const types::global_cell_index max_n_cells, const double top_fraction_of_cells, const double bottom_fraction_of_cells)
void refine_and_coarsen_fixed_fraction(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double top_fraction, const double bottom_fraction, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max(), const VectorTools::NormType norm_type=VectorTools::NormType::L1_norm)
static const unsigned int invalid_unsigned_int
Definition: types.h:213
unsigned int global_cell_index
Definition: types.h:118
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
const ::Triangulation< dim, spacedim > & tria