Reference documentation for deal.II version 8.4.2
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 #include <deal.II/grid/tria.h>
23 
24 #include <vector>
25 #include <limits>
26 
27 DEAL_II_NAMESPACE_OPEN
28 
29 // forward declarations
30 template <int dim, int spacedim> class Triangulation;
31 template <class T> class Vector;
32 
33 
48 namespace GridRefinement
49 {
81  template <int dim>
82  std::pair<double, double>
83  adjust_refine_and_coarsen_number_fraction (const unsigned int current_n_cells,
84  const unsigned int max_n_cells,
85  const double top_fraction_of_cells,
86  const double bottom_fraction_of_cells);
87 
153  template <int dim, class VectorType, int spacedim>
154  void
156  (Triangulation<dim,spacedim> &triangulation,
157  const VectorType &criteria,
158  const double top_fraction_of_cells,
159  const double bottom_fraction_of_cells,
160  const unsigned int max_n_cells = std::numeric_limits<unsigned int>::max());
161 
217  template <int dim, class VectorType, int spacedim>
218  void
221  const VectorType &criteria,
222  const double top_fraction,
223  const double bottom_fraction,
224  const unsigned int max_n_cells = std::numeric_limits<unsigned int>::max());
225 
226 
227 
300  template <int dim, class VectorType, int spacedim>
301  void
303  const VectorType &criteria,
304  const unsigned int order=2);
305 
320  template <int dim, class VectorType, int spacedim>
322  const VectorType &criteria,
323  const double threshold,
324  const unsigned int max_to_mark = numbers::invalid_unsigned_int);
325 
340  template <int dim, class VectorType, int spacedim>
342  const VectorType &criteria,
343  const double threshold);
344 
349  DeclException0(ExcNegativeCriteria);
350 
355  DeclException0 (ExcInvalidParameterValue);
356 }
357 
358 
359 
360 DEAL_II_NAMESPACE_CLOSE
361 
362 #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:164
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())
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)
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)
DeclException0(ExcNegativeCriteria)