Reference documentation for deal.II version 9.2.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\}}\)
grid_refinement.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 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_grid_refinement_h
17 #define dealii_grid_refinement_h
18 
19 
20 #include <deal.II/base/config.h>
21 
23 
24 #include <limits>
25 
27 
28 // forward declarations
29 #ifndef DOXYGEN
30 template <int dim, int spacedim>
31 class Triangulation;
32 template <typename Number>
33 class Vector;
34 #endif
35 
54 namespace GridRefinement
55 {
87  template <int dim>
88  std::pair<double, double>
90  const types::global_cell_index current_n_cells,
91  const types::global_cell_index max_n_cells,
92  const double top_fraction_of_cells,
93  const double bottom_fraction_of_cells);
94 
159  template <int dim, typename Number, int spacedim>
160  void
163  const Vector<Number> & criteria,
164  const double top_fraction_of_cells,
165  const double bottom_fraction_of_cells,
166  const unsigned int max_n_cells = std::numeric_limits<unsigned int>::max());
167 
224  template <int dim, typename Number, int spacedim>
225  void
228  const Vector<Number> & criteria,
229  const double top_fraction,
230  const double bottom_fraction,
231  const unsigned int max_n_cells = std::numeric_limits<unsigned int>::max());
232 
233 
234 
308  template <int dim, typename Number, int spacedim>
309  void
311  const Vector<Number> & criteria,
312  const unsigned int order = 2);
313 
328  template <int dim, typename Number, int spacedim>
329  void
331  const Vector<Number> & criteria,
332  const double threshold,
333  const unsigned int max_to_mark = numbers::invalid_unsigned_int);
334 
349  template <int dim, typename Number, int spacedim>
350  void
352  const Vector<Number> & criteria,
353  const double threshold);
354 
360 
366 } // namespace GridRefinement
367 
368 
369 
371 
372 #endif // dealii_grid_refinement_h
LinearAlgebraDealII::Vector
Vector< double > Vector
Definition: generic_linear_algebra.h:43
GridRefinement
Definition: grid_refinement.h:54
Triangulation
Definition: tria.h:1109
GridRefinement::refine_and_coarsen_fixed_number
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())
Definition: grid_refinement.cc:189
GridRefinement::adjust_refine_and_coarsen_number_fraction
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)
Definition: grid_refinement.cc:106
GridRefinement::coarsen
void coarsen(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold)
Definition: grid_refinement.cc:88
GridRefinement::refine_and_coarsen_fixed_fraction
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())
Definition: grid_refinement.cc:257
GridRefinement::ExcNegativeCriteria
static ::ExceptionBase & ExcNegativeCriteria()
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
GridRefinement::refine_and_coarsen_optimize
void refine_and_coarsen_optimize(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const unsigned int order=2)
Definition: grid_refinement.cc:383
exceptions.h
unsigned int
DeclException0
#define DeclException0(Exception0)
Definition: exceptions.h:473
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
config.h
triangulation
const typename ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
Definition: p4est_wrappers.cc:69
GridRefinement::ExcInvalidParameterValue
static ::ExceptionBase & ExcInvalidParameterValue()
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
GridRefinement::refine
void refine(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold, const unsigned int max_to_mark=numbers::invalid_unsigned_int)
Definition: grid_refinement.cc:41
Utilities::MPI::max
T max(const T &t, const MPI_Comm &mpi_communicator)