16 #ifndef _d2k_parsed_grid_refinement_h 17 #define _d2k_parsed_grid_refinement_h 22 #ifdef DEAL_II_WITH_MPI 23 #ifdef DEAL_II_WITH_P4EST 30 #include <deal2lkit/config.h> 47 const std::string &strategy=
"fraction",
48 const double &top_parameter=.3,
49 const double &bottom_parameter=.1,
50 const unsigned int &max_cells=0,
51 const unsigned int &order=2);
68 template<
int dim,
class Vector ,
int spacedim>
69 void mark_cells(
const Vector &criteria,
72 #ifdef DEAL_II_WITH_MPI 73 #ifdef DEAL_II_WITH_P4EST 85 template<
int dim,
class Vector ,
int spacedim>
86 void mark_cells(
const Vector &criteria,
104 #ifdef DEAL_II_WITH_MPI 105 #ifdef DEAL_II_WITH_P4EST 106 template<
int dim,
class Vector ,
int spacedim>
110 if (strategy ==
"number")
113 top_parameter, bottom_parameter,
114 max_cells ? max_cells :
115 std::numeric_limits< unsigned int >::max());
116 else if (strategy ==
"fraction")
119 top_parameter, bottom_parameter);
127 template<
int dim,
class Vector ,
int spacedim>
131 if (strategy ==
"number")
134 top_parameter, bottom_parameter,
135 max_cells ? max_cells :
136 std::numeric_limits< unsigned int >::max());
137 else if (strategy ==
"fraction")
140 top_parameter, bottom_parameter,
141 max_cells ? max_cells :
142 std::numeric_limits< unsigned int >::max());
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())
A parameter acceptor base class.
A wrapper for refinement strategies.
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())
#define Assert(cond, exc)
std::string strategy
Default expression of this function.
void refine_and_coarsen_fixed_fraction(parallel::distributed::Triangulation< dim, spacedim > &tria, const VectorType &criteria, const double top_fraction_of_error, const double bottom_fraction_of_error)
void refine_and_coarsen_fixed_number(parallel::distributed::Triangulation< dim, spacedim > &tria, 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::ExceptionBase & ExcInternalError()
void mark_cells(const Vector &criteria, Triangulation< dim, spacedim > &tria) const
Mark cells a the triangulation for refinement or coarsening, according to the given strategy applied ...