![]() |
Reference documentation for deal.II version 8.1.0
|
Functions | |
template<int dim, class Vector , int spacedim> | |
void | refine_and_coarsen_fixed_number (parallel::distributed::Triangulation< dim, spacedim > &tria, const Vector &criteria, const double top_fraction_of_cells, const double bottom_fraction_of_cells) |
template<int dim, class Vector , int spacedim> | |
void | refine_and_coarsen_fixed_fraction (parallel::distributed::Triangulation< dim, spacedim > &tria, const Vector &criteria, const double top_fraction_of_error, const double bottom_fraction_of_error) |
Collection of functions controlling refinement and coarsening of parallel::distributed::Triangulation objects. This namespace provides similar functionality to the GridRefinement namespace, except that it works for meshes that are parallel and distributed.
void parallel::distributed::GridRefinement::refine_and_coarsen_fixed_number | ( | parallel::distributed::Triangulation< dim, spacedim > & | tria, |
const Vector & | criteria, | ||
const double | top_fraction_of_cells, | ||
const double | bottom_fraction_of_cells | ||
) |
Like GridRefinement::refine_and_coarsen_fixed_number, but for parallel distributed triangulation.
The vector of criteria needs to be a vector of refinement criteria for all cells active on the current triangulation, i.e. tria.n_active_cells()
(and not tria.n_locally_owned_active_cells()
). However, the function will only look at the indicators that correspond to those cells that are actually locally owned, and ignore the indicators for all other cells. The function will then coordinate among all processors that store part of the triangulation so that at the end top_fraction_of_cells
are refined, where the fraction is enforced as a fraction of Triangulation::n_global_active_cells, not Triangulation::n_locally_active_cells on each processor individually. In other words, it may be that on some processors, no cells are refined at all.
The same is true for the fraction of cells that is coarsened.
void parallel::distributed::GridRefinement::refine_and_coarsen_fixed_fraction | ( | parallel::distributed::Triangulation< dim, spacedim > & | tria, |
const Vector & | criteria, | ||
const double | top_fraction_of_error, | ||
const double | bottom_fraction_of_error | ||
) |
Like GridRefinement::refine_and_coarsen_fixed_fraction, but for parallel distributed triangulation.
The vector of criteria needs to be a vector of refinement criteria for all cells active on the current triangulation, tria.n_active_cells()
(and not tria.n_locally_owned_active_cells()
). However, the function will only look at the indicators that correspond to those cells that are actually locally owned, and ignore the indicators for all other cells. The function will then coordinate among all processors that store part of the triangulation so that at the end the smallest fraction of Triangulation::n_global_active_cells (not Triangulation::n_locally_active_cells on each processor individually) is refined that together make up a total of top_fraction_of_error
of the total error. In other words, it may be that on some processors, no cells are refined at all.
The same is true for the fraction of cells that is coarsened.