Reference documentation for deal.II version 8.5.1
grid_refinement.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2016 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 at
12 // the top level of the deal.II distribution.
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 #include <deal.II/base/exceptions.h>
22 
23 #include <limits>
24 
25 DEAL_II_NAMESPACE_OPEN
26 
27 // forward declarations
28 template <int dim, int spacedim> class Triangulation;
29 
30 
45 namespace GridRefinement
46 {
78  template <int dim>
79  std::pair<double, double>
80  adjust_refine_and_coarsen_number_fraction (const unsigned int current_n_cells,
81  const unsigned int max_n_cells,
82  const double top_fraction_of_cells,
83  const double bottom_fraction_of_cells);
84 
148  template <int dim, class VectorType, int spacedim>
149  void
151  (Triangulation<dim,spacedim> &triangulation,
152  const VectorType &criteria,
153  const double top_fraction_of_cells,
154  const double bottom_fraction_of_cells,
155  const unsigned int max_n_cells = std::numeric_limits<unsigned int>::max());
156 
212  template <int dim, class VectorType, int spacedim>
213  void
216  const VectorType &criteria,
217  const double top_fraction,
218  const double bottom_fraction,
219  const unsigned int max_n_cells = std::numeric_limits<unsigned int>::max());
220 
221 
222 
295  template <int dim, class VectorType, int spacedim>
296  void
298  const VectorType &criteria,
299  const unsigned int order=2);
300 
315  template <int dim, class VectorType, int spacedim>
317  const VectorType &criteria,
318  const double threshold,
319  const unsigned int max_to_mark = numbers::invalid_unsigned_int);
320 
335  template <int dim, class VectorType, int spacedim>
337  const VectorType &criteria,
338  const double threshold);
339 
345 
351 }
352 
353 
354 
355 DEAL_II_NAMESPACE_CLOSE
356 
357 #endif //dealii__grid_refinement_h
void refine_and_coarsen_fixed_number(Triangulation< dim, spacedim > &triangulation, const VectorType &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())
static const unsigned int invalid_unsigned_int
Definition: types.h:170
void refine_and_coarsen_fixed_fraction(Triangulation< dim, spacedim > &tria, const VectorType &criteria, const double top_fraction, const double bottom_fraction, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max())
static ::ExceptionBase & ExcNegativeCriteria()
#define DeclException0(Exception0)
Definition: exceptions.h:541
void coarsen(Triangulation< dim, spacedim > &tria, const VectorType &criteria, const double threshold)
std::pair< double, double > adjust_refine_and_coarsen_number_fraction(const unsigned int current_n_cells, const unsigned int max_n_cells, const double top_fraction_of_cells, const double bottom_fraction_of_cells)
static ::ExceptionBase & ExcInvalidParameterValue()
void refine(Triangulation< dim, spacedim > &tria, const VectorType &criteria, const double threshold, const unsigned int max_to_mark=numbers::invalid_unsigned_int)
void refine_and_coarsen_optimize(Triangulation< dim, spacedim > &tria, const VectorType &criteria, const unsigned int order=2)