Reference documentation for deal.II version GIT relicensing-136-gb80d0be4af 2024-03-18 08:20:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Functions
GridRefinement Namespace Reference

Functions

template<int dim>
std::pair< double, double > adjust_refine_and_coarsen_number_fraction (const types::global_cell_index current_n_cells, const types::global_cell_index max_n_cells, const double top_fraction_of_cells, const double bottom_fraction_of_cells)
 
template<int dim, typename Number , int spacedim>
void refine_and_coarsen_fixed_number (Triangulation< dim, spacedim > &triangulation, 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 (Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double top_fraction, const double bottom_fraction, const unsigned int max_n_cells=std::numeric_limits< unsigned int >::max(), const VectorTools::NormType norm_type=VectorTools::L1_norm)
 
template<int dim, typename Number , int spacedim>
void refine_and_coarsen_optimize (Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const unsigned int order=2)
 
template<int dim, typename Number , int spacedim>
void refine (Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold, const unsigned int max_to_mark=numbers::invalid_unsigned_int)
 
template<int dim, typename Number , int spacedim>
void coarsen (Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold)
 
static ::ExceptionBaseExcNegativeCriteria ()
 
static ::ExceptionBaseExcInvalidParameterValue ()
 

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.

The functions in this namespace form two categories. There are the auxiliary functions refine() and coarsen(). More important for users are the other functions, which implement refinement strategies, as being found in the literature on adaptive finite element methods. For mathematical discussion of these methods, consider works by Dörfler, Morin, Nochetto, Rannacher, Stevenson, and others.

Function Documentation

◆ adjust_refine_and_coarsen_number_fraction()

template<int dim>
std::pair< double, double > GridRefinement::adjust_refine_and_coarsen_number_fraction ( const types::global_cell_index  current_n_cells,
const types::global_cell_index  max_n_cells,
const double  top_fraction_of_cells,
const double  bottom_fraction_of_cells 
)

Return a pair of double values of which the first is adjusted refinement fraction of cells and the second is adjusted coarsening fraction of cells.

Parameters
[in]current_n_cellsThe current cell number.
[in]max_n_cellsThe maximal number of cells. If current cell number current_n_cells is already exceeded maximal cell number max_n_cells, refinement fraction of cells will be set to zero and coarsening fraction of cells will be adjusted to reduce cell number to max_n_cells. If cell number is going to be exceeded only upon refinement, then refinement and coarsening fractions are going to be adjusted with a same ratio in an attempt to reach the maximum number of cells. Be aware though that through proliferation of refinement due to Triangulation::MeshSmoothing, this number is only an indicator. The default value of this argument is to impose no limit on the number of cells.
[in]top_fraction_of_cellsThe requested fraction of active cells to be refined.
[in]bottom_fraction_of_cellsThe requested fraction of active cells to be coarsened.
Note
Usually you do not need to call this function explicitly. Pass max_n_cells to function refine_and_coarsen_fixed_number() or function refine_and_coarsen_fixed_fraction() and they will call this function if necessary.

Definition at line 235 of file grid_refinement.cc.

◆ refine_and_coarsen_fixed_number()

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

This function provides a strategy to mark cells for refinement and coarsening with the goal of providing predictable growth in the size of the mesh by refining a given fraction of all cells.

The function takes a vector of refinement criteria and two values between zero and one denoting the fractions of cells to be refined and coarsened. It flags cells for further processing by Triangulation::execute_coarsening_and_refinement() according to the following greedy algorithm:

  1. Sort the cells according to descending values of criteria.

  2. Mark the top_fraction_of_cells times Triangulation::n_active_cells() active cells with the largest refinement criteria for refinement.

  3. Mark the bottom_fraction_of_cells times Triangulation::n_active_cells() active cells with the smallest refinement criteria for coarsening.

As an example, with no coarsening, setting top_fraction_of_cells to 1/3 will result in approximately doubling the number of cells in two dimensions. That is because each of these 1/3 of cells will be replaced by its four children, resulting in \(4\times \frac 13 N\) cells, whereas the remaining 2/3 of cells remains untouched – thus yielding a total of \(4\times \frac 13 N + \frac 23 N = 2N\) cells. The same effect in three dimensions is achieved by refining 1/7th of the cells. These values are therefore frequently used because they ensure that the cost of computations on subsequent meshes become expensive sufficiently quickly that the fraction of time spent on the coarse meshes is not too large. On the other hand, the fractions are small enough that mesh adaptation does not refine too many cells in each step.

Note
This function only sets the coarsening and refinement flags. The mesh is not changed until you call Triangulation::execute_coarsening_and_refinement().
Parameters
[in,out]triangulationThe triangulation whose cells this function is supposed to mark for coarsening and refinement.
[in]criteriaThe refinement criterion for each mesh cell. Entries may not be negative.
[in]top_fraction_of_cellsThe fraction of cells to be refined. If this number is zero, no cells will be refined. If it equals one, the result will be flagging for global refinement.
[in]bottom_fraction_of_cellsThe fraction of cells to be coarsened. If this number is zero, no cells will be coarsened.
[in]max_n_cellsThis argument can be used to specify a maximal number of cells. If this number is going to be exceeded upon refinement, then refinement and coarsening fractions are going to be adjusted in an attempt to reach the maximum number of cells. Be aware though that through proliferation of refinement due to Triangulation::MeshSmoothing, this number is only an indicator. The default value of this argument is to impose no limit on the number of cells.

Definition at line 318 of file grid_refinement.cc.

◆ refine_and_coarsen_fixed_fraction()

template<int dim, typename Number , int spacedim>
void GridRefinement::refine_and_coarsen_fixed_fraction ( Triangulation< dim, spacedim > &  tria,
const Vector< Number > &  criteria,
const double  top_fraction,
const double  bottom_fraction,
const unsigned int  max_n_cells = std::numeric_limits<unsigned int>::max(),
const VectorTools::NormType  norm_type = VectorTools::L1_norm 
)

This function provides a strategy to mark cells for refinement and coarsening with the goal of controlling the reduction of the error estimate.

Also known as the bulk criterion or Dörfler marking, this function computes the thresholds for refinement and coarsening such that the criteria of cells getting flagged for refinement make up for a certain fraction of the total error. We explain its operation for refinement, coarsening works analogously.

Let cK be the criterion of cell K. Then the total error estimate is computed by the formula

\[ E = \sum_{K\in \cal T} c_K. \]

If 0 < a < 1 is top_fraction, then we refine the smallest subset \(\cal M\) of the Triangulation \(\cal T\) such that

\[ a E \le \sum_{K\in \cal M} c_K \]

The algorithm is performed by the greedy algorithm described in refine_and_coarsen_fixed_number().

Note
The often used formula with squares on the left and right is recovered by actually storing the square of cK in the vector criteria.

From the point of view of implementation, this time we really need to sort the array of criteria. Just like the other strategy described above, this function only computes the threshold values and then passes over to refine() and coarsen().

Parameters
[in,out]triaThe triangulation whose cells this function is supposed to mark for coarsening and refinement.
[in]criteriaThe refinement criterion computed on each mesh cell. Entries may not be negative.
[in]top_fractionThe fraction of the total estimate which should be refined. If this number is zero, no cells will be refined. If it equals one, the result will be flagging for global refinement.
[in]bottom_fractionThe fraction of the estimate coarsened. If this number is zero, no cells will be coarsened.
[in]max_n_cellsThis argument can be used to specify a maximal number of cells. If this number is going to be exceeded upon refinement, then refinement and coarsening fractions are going to be adjusted in an attempt to reach the maximum number of cells. Be aware though that through proliferation of refinement due to Triangulation::MeshSmoothing, this number is only an indicator. The default value of this argument is to impose no limit on the number of cells.
[in]norm_typeTo determine thresholds, combined errors on subsets of cells are calculated as norms of the criteria on these cells. Different types of norms can be used for this purpose, from which VectorTools::L1_norm and VectorTools::L2_norm are currently supported.

Definition at line 386 of file grid_refinement.cc.

◆ refine_and_coarsen_optimize()

template<int dim, typename Number , int spacedim>
void GridRefinement::refine_and_coarsen_optimize ( Triangulation< dim, spacedim > &  tria,
const Vector< Number > &  criteria,
const unsigned int  order = 2 
)

This function flags cells of a triangulation for refinement with the aim to reach a grid that is optimal with respect to an objective function that tries to balance reducing the error and increasing the numerical cost when the mesh is refined. Specifically, this function makes the assumption that if you refine a cell \(K\) with error indicator \(\eta_K\) provided by the second argument to this function, then the error on the children (for all children together) will only be \(2^{-\text{order}}\eta_K\) where order is the third argument of this function. This makes the assumption that the error is only a local property on a mesh and can be reduced by local refinement – an assumption that is true for the interpolation operator, but not for the usual Galerkin projection, although it is approximately true for elliptic problems where the Greens function decays quickly and the error here is not too much affected by a too coarse mesh somewhere else.

With this, we can define the objective function this function tries to optimize. Let us assume that the mesh currently has \(N_0\) cells. Then, if we refine the \(m\) cells with the largest errors, we expect to get (in \(d\) space dimensions)

\[ N(m) = (N_0-m) + 2^d m = N_0 + (2^d-1)m \]

cells ( \(N_0-m\) are not refined, and each of the \(m\) cells we refine yield \(2^d\) child cells. On the other hand, with refining \(m\) cells, and using the assumptions above, we expect that the error will be

\[ \eta^\text{exp}(m) = \sum_{K, K\; \text{will not be refined}} \eta_K + \sum_{K, K\; \text{will be refined}} 2^{-\text{order}}\eta_K \]

where the first sum extends over \(N_0-m\) cells and the second over the \(m\) cells that will be refined. Note that \(N(m)\) is an increasing function of \(m\) whereas \(\eta^\text{exp}(m)\) is a decreasing function.

This function then tries to find that number \(m\) of cells to mark for refinement for which the objective function

\[ J(m) = N(m)^{\text{order}/d} \eta^\text{exp}(m) \]

is minimal.

The rationale for this function is two-fold. First, compared to the refine_and_coarsen_fixed_fraction() and refine_and_coarsen_fixed_number() functions, this function has the property that if all refinement indicators are the same (i.e., we have achieved a mesh where the error per cell is equilibrated), then the entire mesh is refined. This is based on the observation that a mesh with equilibrated error indicators is the optimal mesh (i.e., has the least overall error) among all meshes with the same number of cells. (For proofs of this, see R. Becker, M. Braack, R. Rannacher: "Numerical simulation of laminar flames at low Mach number with adaptive finite elements", Combustion Theory and Modelling, Vol. 3, Nr. 3, p. 503-534 1999; and W. Bangerth, R. Rannacher: "Adaptive Finite Element Methods for Differential Equations", Birkhauser, 2003.)

Second, the function uses the observation that ideally, the error behaves like \(e \approx c N^{-\alpha}\) with some constant \(\alpha\) that depends on the dimension and the finite element degree. It should - given optimal mesh refinement - not depend so much on the regularity of the solution, as it is based on the idea, that all singularities can be resolved by refinement. Mesh refinement is then based on the idea that we want to make \(c=e N^\alpha\) small. This corresponds to the functional \(J(m)\) above.

Note
This function was originally implemented by Thomas Richter. It follows a strategy described in [167]. See in particular Section 4.3, pp. 42-43.

Definition at line 447 of file grid_refinement.cc.

◆ refine()

template<int dim, typename Number , int spacedim>
void GridRefinement::refine ( Triangulation< dim, spacedim > &  tria,
const Vector< Number > &  criteria,
const double  threshold,
const unsigned int  max_to_mark = numbers::invalid_unsigned_int 
)

Mark all mesh cells for which the value in criteria exceeds threshold for refinement, but only flag up to max_to_mark cells.

The vector criteria contains a nonnegative value for each active cell, ordered in the canonical order of Triangulation::active_cell_iterator.

The cells are only flagged for refinement, they are not actually refined. To do so, you have to call Triangulation::execute_coarsening_and_refinement().

This function does not implement a refinement strategy, it is more a helper function for the actual strategies.

Definition at line 164 of file grid_refinement.cc.

◆ coarsen()

template<int dim, typename Number , int spacedim>
void GridRefinement::coarsen ( Triangulation< dim, spacedim > &  tria,
const Vector< Number > &  criteria,
const double  threshold 
)

Mark all mesh cells for which the value in criteria is less than threshold for coarsening.

The vector criteria contains a nonnegative value for each active cell, ordered in the canonical order of Triangulation::active_cell_iterator.

The cells are only flagged for coarsening, they are not actually coarsened. To do so, you have to call Triangulation::execute_coarsening_and_refinement().

This function does not implement a refinement strategy, it is more a helper function for the actual strategies.

Definition at line 214 of file grid_refinement.cc.