Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Functions
parallel::distributed::GridRefinement Namespace Reference

Functions

template<int dim, typename Number , int spacedim>
void refine_and_coarsen_fixed_number (parallel::distributed::Triangulation< dim, spacedim > &tria, 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())
 
template<int dim, typename Number , int spacedim>
void refine_and_coarsen_fixed_fraction (parallel::distributed::Triangulation< dim, spacedim > &tria, const ::Vector< Number > &criteria, const double top_fraction_of_error, const double bottom_fraction_of_error)
 

Detailed Description

This namespace provides a collection of functions that aid in refinement and coarsening of triangulations. Despite the name of the namespace, the functions do not actually refine the triangulation, but only mark cells for refinement or coarsening. In other words, they perform the "mark" part of the typical "solve-estimate-mark-refine" cycle of the adaptive finite element loop.

In contrast to the functions in namespace GridRefinement, the functions in the current namespace are intended for distributed meshes, i.e., objects of type parallel::distributed::Triangulation.

Author
Wolfgang Bangerth, 2009

Function Documentation

◆ refine_and_coarsen_fixed_number()

template<int dim, typename Number , int spacedim>
void parallel::distributed::GridRefinement::refine_and_coarsen_fixed_number ( parallel::distributed::Triangulation< dim, spacedim > &  tria,
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() 
)

Like GridRefinement::refine_and_coarsen_fixed_number, but for parallel distributed triangulations.

The vector of criteria needs to be a vector of refinement criteria for all cells active on the current triangulation, i.e., it needs to be of length tria.n_active_cells() (and not tria.n_locally_owned_active_cells()). In other words, the vector needs to include entries for ghost and artificial cells. However, the current 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 a fraction top_fraction_of_cells of all Triangulation::n_global_active_cells() active cells are refined, rather than a fraction of the 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.

Definition at line 429 of file grid_refinement.cc.

◆ refine_and_coarsen_fixed_fraction()

template<int dim, typename Number , int spacedim>
void parallel::distributed::GridRefinement::refine_and_coarsen_fixed_fraction ( parallel::distributed::Triangulation< dim, spacedim > &  tria,
const ::Vector< Number > &  criteria,
const double  top_fraction_of_error,
const double  bottom_fraction_of_error 
)

Like GridRefinement::refine_and_coarsen_fixed_fraction, but for parallel distributed triangulations.

The vector of criteria needs to be a vector of refinement criteria for all cells active on the current triangulation, i.e., it needs to be of length tria.n_active_cells() (and not tria.n_locally_owned_active_cells()). In other words, the vector needs to include entries for ghost and artificial cells. However, the current 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_owned_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.

Definition at line 501 of file grid_refinement.cc.