Reference documentation for deal.II version 9.2.0
\(\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\}}\)
Typedefs
hp::Refinement Namespace Reference

Typedefs

template<typename Number >
using ComparisonFunction = std::function< bool(const Number &, const Number &)>
 

Detailed Description

We supply adaptive methods to align computational resources with the complexity of the numerical solution. Error estimates are an appropriate means of determining where adjustments need to be made.

However with hp-adaptivity, we have two ways to realize these adjustments: For irregular solutions, h-adaptive methods which dynamically assign cell sizes tend to reduce the approximation error, while for smooth solutions p-adaptive methods are better suited in which function spaces will be selected dynamically. This namespace collects tools to decide which type of adaptive methods to apply.

Usage

To successfully apply hp-adaptive methods, we recommend the following workflow:

  1. A suitable error estimate is the basis for any kind of adaptive method. Similar to pure grid refinement, we will determine error estimates in the usual way (i.e. KellyErrorEstimator) and mark cells for refinement or coarsening (i.e. GridRefinement).

    Calling Triangulation::execute_coarsening_and_refinement() at this stage will perform pure grid refinement as expected.

  2. Once all refinement and coarsening flags have been distributed on the mesh, we may determine if those qualify for p-adaptive methods. Corresponding functions will set future_fe_indices on top of the refinement and coarsening flags if they fulfil a certain criterion.

    In case of refinement, the superordinate element of the underlying hp::FECollection will be assigned as the future finite element. Correspondingly, the subordinate element will be selected for coarsening.

    Triangulation::execute_coarsening_and_refinement() will now supply both h- and p-adaptive methods independently.

  3. Right now, there may be cells scheduled for both h- and p-adaptation. If we do not want to impose both methods at once, we need to decide which one to pick for each cell individually and unambiguously. Since grid refinement will be imposed by default and we only determine qualification for p-adaptivity on top, we will always decide in favour of p-adaptive methods.

    Calling Triangulation::execute_coarsening_and_refinement() will now perform either h- or p-adaptive methods uniquely on each cell.

  4. Up to this point, each cell knows its destiny in terms of adaptivity. We can now move on to prepare all data structures to be transferred across mesh changes. Previously set refinement and coarsening flags as well as future_fe_indices will be used to update the data accordingly.

As an example, a realisation of pure p-adaptive methods would look like the following:

// step 1: flag cells for refinement or coarsening
Vector<float> estimated_error_per_cell (triangulation.n_active_cells());
hp_dof_handler,
solution,
estimated_error_per_cell);
estimated_error_per_cell,
top_fraction,
bottom_fraction);
// step 2: set future finite element indices on flagged cells
// step 3: decide whether h- or p-adaptive methods will be supplied
// step 4: prepare solutions to be transferred
...
triangulation.execute_coarsening_and_refinement();
Author
Marc Fehling 2019

Typedef Documentation

◆ ComparisonFunction

template<typename Number >
using hp::Refinement::ComparisonFunction = typedef std::function<bool(const Number &, const Number &)>

An alias that defines the characteristics of a function that can be used as a comparison criterion for deciding whether to perform h- or p-adaptation.

Such functions take two numbers as arguments: The first one corresponds to the provided criterion, while the other one conforms to the reference. The result of the comparison will be returned as a boolean.

Definition at line 143 of file refinement.h.

Function Documentation

◆ full_p_adaptivity()

template<int dim, int spacedim>
void hp::Refinement::full_p_adaptivity ( const hp::DoFHandler< dim, spacedim > &  dof_handler)

Each cell flagged for h-refinement will also be flagged for p-refinement. The same applies to coarsening.

Note
Triangulation::prepare_coarsening_and_refinement() may change refine and coarsen flags. Avoid calling it before this particular function.

Setting p adaptivity flags

Definition at line 43 of file refinement.cc.

◆ p_adaptivity_from_flags()

template<int dim, int spacedim>
void hp::Refinement::p_adaptivity_from_flags ( const hp::DoFHandler< dim, spacedim > &  dof_handler,
const std::vector< bool > &  p_flags 
)

Adapt which finite element to use on cells that have been specifically flagged for p-adaptation via the parameter p_flags. Future finite elements will only be assigned if cells have been flagged for refinement and coarsening beforehand.

Each entry of the parameter p_flags needs to correspond to an active cell.

Note
Triangulation::prepare_coarsening_and_refinement() may change refine and coarsen flags. Avoid calling it before this particular function.

Definition at line 55 of file refinement.cc.

◆ p_adaptivity_from_absolute_threshold()

template<int dim, typename Number , int spacedim>
void hp::Refinement::p_adaptivity_from_absolute_threshold ( const hp::DoFHandler< dim, spacedim > &  dof_handler,
const Vector< Number > &  criteria,
const Number  p_refine_threshold,
const Number  p_coarsen_threshold,
const ComparisonFunction< typename identity< Number >::type > &  compare_refine = std::greater_equal<Number>(),
const ComparisonFunction< typename identity< Number >::type > &  compare_coarsen = std::less_equal<Number>() 
)

Adapt which finite element to use on cells whose criteria meet a certain absolute threshold.

For p-refinement and p-coarsening, two separate thresholds need to provided via parameters p_refine_threshold and p_coarsen_threshold.

We consider a cell for p-adaptivity if it is currently flagged for refinement or coarsening and its criterion successfully compares to the corresponding threshold. Let us be more specific on the default case: We consider a cell for p-refinement if it is flagged for refinement and its criterion is larger than or equal to the corresponding threshold. The same applies for p-coarsening, but the cell's criterion must be lower than or equal to the threshold. However, different compare function objects can be supplied via the parameters compare_refine and compare_coarsen to impose different decision strategies.

Each entry of the parameter criteria needs to correspond to an active cell.

Note
Triangulation::prepare_coarsening_and_refinement() may change refine and coarsen flags. Avoid calling it before this particular function.

Definition at line 91 of file refinement.cc.

◆ p_adaptivity_from_relative_threshold()

template<int dim, typename Number , int spacedim>
void hp::Refinement::p_adaptivity_from_relative_threshold ( const hp::DoFHandler< dim, spacedim > &  dof_handler,
const Vector< Number > &  criteria,
const double  p_refine_fraction = 0.5,
const double  p_coarsen_fraction = 0.5,
const ComparisonFunction< typename identity< Number >::type > &  compare_refine = std::greater_equal<Number>(),
const ComparisonFunction< typename identity< Number >::type > &  compare_coarsen = std::less_equal<Number>() 
)

Adapt which finite element to use on cells whose criteria meet a certain threshold relative to the overall range of criterion values.

The threshold will be determined for refined and coarsened cells separately based on the currently set refinement markers. For each class of cells, we determine the maximal and minimal values of all criteria and determine the threshold by linear interpolation between these limits. Parameters p_refine_fraction and p_refine_coarsen are used as interpolation factors, where 0 corresponds to the minimal and 1 to the maximal value. By default, mean values are considered as thresholds.

We consider a cell for p-adaptivity if it is currently flagged for refinement or coarsening and its criterion successfully compares to the corresponding threshold. Let us be more specific on the default case: We consider a cell for p-refinement if it is flagged for refinement and its criterion is larger than or equal to the corresponding threshold. The same applies for p-coarsening, but the cell's criterion must be lower than or equal to the threshold. However, different compare function objects can be supplied via the parameters compare_refine and compare_coarsen to impose different decision strategies.

Each entry of the parameter criteria needs to correspond to an active cell. Parameters p_refine_fraction and p_coarsen_fraction need to be in the interval \([0,1]\).

Note
Triangulation::prepare_coarsening_and_refinement() may change refine and coarsen flags. Avoid calling it before this particular function.

Definition at line 123 of file refinement.cc.

◆ p_adaptivity_fixed_number()

template<int dim, typename Number , int spacedim>
void hp::Refinement::p_adaptivity_fixed_number ( const hp::DoFHandler< dim, spacedim > &  dof_handler,
const Vector< Number > &  criteria,
const double  p_refine_fraction = 0.5,
const double  p_coarsen_fraction = 0.5,
const ComparisonFunction< typename identity< Number >::type > &  compare_refine = std::greater_equal<Number>(),
const ComparisonFunction< typename identity< Number >::type > &  compare_coarsen = std::less_equal<Number>() 
)

Adapt which finite element to use on a given fraction of cells.

Out of all cells flagged for a certain type of adaptation, be it refinement or coarsening, we will determine a fixed number of cells among this subset that will be flagged for the corresponding p-adaptive variant.

For each of both refinement and coarsening subsets, we will determine a threshold based on the provided parameter criteria containing indicators for every active cell. In the default case for refinement, all cells with an indicator larger than or equal to the corresponding threshold will be considered for p-refinement, while for coarsening all cells with an indicator less than or equal to the matching threshold are taken into account. However, different compare function objects can be supplied via the parameters compare_refine and compare_coarsen to impose different decision strategies.

For refinement, the threshold will be associated with the cell that has the p_refine_fraction times Triangulation::n_active_cells() largest indicator, while it is the cell with the p_refine_coarsen times Triangulation::n_active_cells() lowest indicator for coarsening.

Each entry of the parameter criteria needs to correspond to an active cell. Parameters p_refine_fraction and p_coarsen_fraction need to be in the interval \([0,1]\).

Note
Triangulation::prepare_coarsening_and_refinement() may change refine and coarsen flags. Avoid calling it before this particular function.

Definition at line 213 of file refinement.cc.

◆ p_adaptivity_from_regularity()

template<int dim, typename Number , int spacedim>
void hp::Refinement::p_adaptivity_from_regularity ( const hp::DoFHandler< dim, spacedim > &  dof_handler,
const Vector< Number > &  sobolev_indices 
)

Adapt which finite element to use on cells based on the regularity of the (unknown) analytical solution.

With an approximation of the local Sobolev regularity index \(k_K\), we may assess to which finite element space our local solution on cell \(K\) belongs. Since the regularity index is only an estimate, we won't use it to assign the finite element space directly, but rather consider it as an indicator for adaptation. If a cell is flagged for refinement, we will perform p-refinement once it satisfies \(k_K > p_{K,\text{super}}\), where \(p_{K,\text{super}}\) is the polynomial degree of the finite element superordinate to the currently active element on cell \(K\). In case of coarsening, the criterion \(k_K < p_{K,\text{sub}}\) has to be met, with \(p_{K,\text{sub}}\) the degree of the subordinate element.

Each entry of the parameter sobolev_indices needs to correspond to an active cell.

For more theoretical details see [1] .

Note
Triangulation::prepare_coarsening_and_refinement() may change refine and coarsen flags. Avoid calling it before this particular function.

Definition at line 407 of file refinement.cc.

◆ p_adaptivity_from_reference()

template<int dim, typename Number , int spacedim>
void hp::Refinement::p_adaptivity_from_reference ( const hp::DoFHandler< dim, spacedim > &  dof_handler,
const Vector< Number > &  criteria,
const Vector< Number > &  references,
const ComparisonFunction< typename identity< Number >::type > &  compare_refine,
const ComparisonFunction< typename identity< Number >::type > &  compare_coarsen 
)

Adapt which finite element to use on each cell based on how its criterion relates to a reference.

We consider a cell for p-adaptivity if it is currently flagged for refinement or coarsening and its criterion successfully compares to the corresponding reference. Other than functions p_adaptivity_from_absolute_threshold() and p_adaptivity_from_relative_threshold(), compare function objects have to be provided explicitly via the parameters compare_refine and compare_coarsen.

Each entry of the parameters criteria and references needs to correspond to an active cell.

Note
Triangulation::prepare_coarsening_and_refinement() may change refine and coarsen flags. Avoid calling it before this particular function.

Definition at line 458 of file refinement.cc.

◆ predict_error()

template<int dim, typename Number , int spacedim>
void hp::Refinement::predict_error ( const hp::DoFHandler< dim, spacedim > &  dof_handler,
const Vector< Number > &  error_indicators,
Vector< Number > &  predicted_errors,
const double  gamma_p = std::sqrt(0.4),
const double  gamma_h = 2.,
const double  gamma_n = 1. 
)

Predict how the current error_indicators will adapt after refinement and coarsening were to happen on the provided dof_handler, and write its results to predicted_errors. Each entry of error_indicators and predicted_errors corresponds to an active cell on the underlying Triangulation, thus each container has to be of size Triangulation::n_active_cells(). The errors are interpreted to be measured in the energy norm; this assumption enters the rate of convergence that is used in the prediction. The \(l_2\)-norm of the output argument predicted_errors corresponds to the predicted global error after adaptation.

For p-adaptation, the local error is expected to converge exponentially with the polynomial degree of the assigned finite element. Each increase or decrease of the degree will thus change its value by a user-defined control parameter gamma_p.

For h-adaptation, we expect the local error \(\eta_K\) on cell \(K\) to be proportional to \((h_K)^{p_K}\) in the energy norm, where \(h_K\) denotes the cell diameter and \(p_K\) the polynomial degree of the currently assigned finite element on cell \(K\).

During h-coarsening, the finite elements on siblings may be different, and their parent cell will be assigned to their least dominating finite element that belongs to its most general child. Thus, we will always interpolate on an enclosing finite element space. Additionally assuming that the finite elements on the cells to be coarsened are sufficient to represent the solution correctly (e.g. at least quadratic basis functions for a quadratic solution), we are confident to say that the error will not change by sole interpolation on the larger finite element space.

For p-adaptation, the local error is expected to converge exponentially with the polynomial degree of the assigned finite element. Each increase or decrease of the degree will thus change its value by a user-defined control parameter gamma_p. The assumption of exponential convergence is only valid if both h- and p-adaptive methods are combined in a sense that they are both utilitzed throughout a mesh, but do not have to be applied both on a cell simultaneously.

The prediction algorithm is formulated as follows with control parameters gamma_p, gamma_h and gamma_n that may be used to influence prediction for each adaptation type individually. The results for each individual cell are stored in the predicted_errors output argument.

Adaptation type Prediction formula
no adaptation \(\eta_{K,\text{pred}} = \eta_{K} \, \gamma_\text{n}\) \(\gamma_\text{n} \in (0,\infty)\)
p-adaptation \(\eta_{K,\text{pred}} = \eta_{K} \, \gamma_\text{p}^{(p_{K,\text{future}} - p_K)}\) \(\gamma_\text{p} \in (0,1)\)
hp-refinement \(\eta_{K,\text{pred}} = \eta_{K} \, \gamma_\text{h} \, 0.5^{p_{K,\text{future}}} \, \gamma_\text{p}^{(p_{K,\text{future}} - p_{K})}\) \(\gamma_\text{h} \in (0,\infty)\)
hp-coarsening \(\eta_{K,\text{pred}} = \eta_{K} \, (\gamma_\text{h} \, 0.5^{p_{K,\text{future}}})^{-1} \, \gamma_\text{p}^{(p_{K,\text{future}} - p_{K})}\)

On basis of the refinement history, we use the predicted error estimates to decide how cells will be adapted in the next adaptation step. Comparing the predicted error from the previous adaptation step to the error estimates of the current step allows us to justify whether our previous choice of adaptation was justified, and lets us decide how to adapt in the next one.

We thus have to transfer the predicted error from the old to the adapted mesh. When transferring the predicted error to the adapted mesh, make sure to configure your CellDataTransfer object with AdaptationStrategies::Refinement::l2_norm() as a refinement strategy and AdaptationStrategies::Coarsening::l2_norm() as a coarsening strategy. This ensures that the \(l_2\)-norm of the predict errors is preserved on both meshes.

In this context, we assume that the local error on a cell to be h-refined will be divided equally on all of its \(n_{K_c}\) children, whereas local errors on siblings will be summed up on the parent cell in case of h-coarsening. This assumption is often not satisfied in practice: For example, if a cell is at a corner singularity, then the one child cell that ends up closest to the singularity will inherit the majority of the remaining error – but this function can not know where the singularity will be, and consequently assumes equal distribution.

Incorporating the transfer from the old to the adapted mesh, the complete error prediction algorithm reads as follows:

Adaptation type Prediction formula
no adaptation \(\eta_{K,\text{pred}} = \eta_{K} \, \gamma_\text{n}\) \(\gamma_\text{n} \in (0,\infty)\)
p-adaptation \(\eta_{K,\text{pred}} = \eta_{K} \, \gamma_\text{p}^{(p_{K,\text{future}} - p_K)}\) \(\gamma_\text{p} \in (0,1)\)
hp-refinement \(\left( \eta_{K_c,\text{pred}} \right)^2 = n_{K_c}^{-1} \left( \eta_{K_p} \, \gamma_\text{h} \, 0.5^{p_{K_c,\text{future}}} \, \gamma_\text{p}^{(p_{K_c,\text{future}} - p_{K_p})} \right)^2 \quad \forall K_c \text{ children of } K_p\) \(\gamma_\text{h} \in (0,\infty)\)
hp-coarsening \(\left( \eta_{K_p,\text{pred}} \right)^2 = \sum\limits_{K_c} \left( \eta_{K_c} \, (\gamma_\text{h} \, 0.5^{p_{K_p,\text{future}}})^{-1} \, \gamma_\text{p}^{(p_{K_p,\text{future}} - p_{K_c})} \right)^2 \quad \forall K_c \text{ children of } K_p\)

With these predicted error estimates, we are capable of adapting the finite element on cells based on their refinement history or rather the predicted change of their error estimates.

If a cell is flagged for adaptation, we want to perform p-adaptation once the associated error indicators \(\eta_{K}\) on cell \(K\) satisfy \(\eta_{K} < \eta_{K,\text{pred}}\), where the subscript \(\text{pred}\) denotes the predicted error. This corresponds to our assumption of smoothness being correct, else h-adaptation is applied. We achieve this with the function hp::Refinement::p_adaptivity_from_reference() and a function object std::less<Number>() for both comparator parameters.

Also with an alternative strategy, we can determine the fractions of cells to be h- and p-adapted among all cells to be adapted. For this, use hp::Refinement::p_adaptivity_fixed_number() with criteria \((\eta_{K,\text{pred}} - \eta_{K})\).

For the very first adaptation step in either case, the user needs to decide whether h- or p-adaptation is supposed to happen. An h-step will be applied with \(\eta_{K,\text{pred}} = 0\), whereas \(\eta_{K,\text{pred}} = \infty\) ensures a p-step. The latter may be realized with std::numeric_limits::infinity().

The following code snippet demonstrates how to impose hp-adaptivity based on refinement history in an application:

// [initialisation...]
Vector<float> predicted_error_per_cell(triangulation.n_active_cells());
for(unsigned int i = 0; i < triangulation.n_active_cells(); ++i)
predicted_error_per_cell[i] = std::numeric_limits<float>::infinity();
// [during each refinement step...]
// set h-adaptivity flags
Vector<float> estimated_error_per_cell(triangulation.n_active_cells());
GridRefinemet::refine_and_coarsen_fixed_{number|fraction}(...);
// set p-adaptivity flags
hp_dof_handler,
estimated_error_per_cell,
predicted_error_per_cell,
std::less<float>(),
std::less<float>());
hp::Refinement::{choose|force}_p_over_h(hp_dof_handler);
// predict error for the subsequent adaptation
triangulation.prepare_coarsening_and_refinement();
hp_dof_handler,
estimated_error_per_cell,
predicted_error_per_cell);
// perform adaptation
&AdaptationStrategies::Refinement::l2_norm<dim, spacedim, float>,
&AdaptationStrategies::Coarsening::l2_norm<dim, spacedim, float>);
cell_data_transfer.prepare_coarsening_and_refinement();
triangulation.execute_coarsening_and_refinement();
Vector<float> transferred_errors(triangulation.n_active_cells());
cell_data_transfer.unpack(predicted_error_per_cell, transferred_errors);
predicted_error_per_cell = std::move(transferred_errors);

For more theoretical details see [51] , where the default parameters for this function come from as well, i.e. \(\gamma_\text{p}^2 = 0.4\), \(\gamma_\text{h}^2 = 4\), \(\gamma_\text{n}^2 = 1\).

Note
We want to predict the error by how adaptation will actually happen. Thus, this function needs to be called after Triangulation::prepare_coarsening_and_refinement().

Error prediction

Definition at line 494 of file refinement.cc.

◆ force_p_over_h()

template<int dim, int spacedim>
void hp::Refinement::force_p_over_h ( const hp::DoFHandler< dim, spacedim > &  dof_handler)

Choose p-adaptivity over h-adaptivity in any case.

Removes all refine and coarsen flags on cells that have a future_fe_index assigned.

Note
Triangulation::prepare_coarsening_and_refinement() may change refine and coarsen flags. Avoid calling it before this particular function.

Decide between h and p adaptivity

Definition at line 622 of file refinement.cc.

◆ choose_p_over_h()

template<int dim, int spacedim>
void hp::Refinement::choose_p_over_h ( const hp::DoFHandler< dim, spacedim > &  dof_handler)

Choose p-adaptivity over h-adaptivity whenever it is invoked on all related cells.

In case of refinement, information about finite elements will be inherited. Thus we will prefer p-refinement over h-refinement whenever desired, i.e. clear the refine flag and supply a corresponding future_fe_index.

However for coarsening, we follow a different approach. Flagging a cell for h-coarsening does not ultimately mean that it will be coarsened. Only if a cell and all of its siblings are flagged, they will be merged into their parent cell. If we consider p-coarsening on top, we must decide for all siblings together how they will be coarsened. We distinguish between three different cases:

  1. Not all siblings flagged for coarsening: p-coarsening.
    We keep the future_fe_indices and clear the coarsen flags on all siblings.
  2. All siblings flagged for coarsening, but not all for p-adaptation: h-coarsening.
    We keep the coarsen flags and clear all future_fe_indices on all siblings.
  3. All siblings flagged for coarsening and p-adaptation: p-coarsening.
    We keep the future_fe_indices and clear the coarsen flags on all siblings.
Note
The function Triangulation::prepare_coarsening_and_refinement() will clean up all h-coarsening flags if they are not shared among all siblings. In the hp case, we need to bring forward this decision: If the cell will not be coarsened, but qualifies for p-adaptivity, we have to set all flags accordingly. So this function anticipates the decision that Triangulation::prepare_coarsening_and_refinement() would have made later on.
Triangulation::prepare_coarsening_and_refinement() may change refine and coarsen flags. Avoid calling it before this particular function.

Definition at line 636 of file refinement.cc.

CellDataTransfer
Definition: cell_data_transfer.h:108
hp::Refinement::predict_error
void predict_error(const hp::DoFHandler< dim, spacedim > &dof_handler, const Vector< Number > &error_indicators, Vector< Number > &predicted_errors, const double gamma_p=std::sqrt(0.4), const double gamma_h=2., const double gamma_n=1.)
Definition: refinement.cc:494
GridRefinement::refine_and_coarsen_fixed_fraction
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())
Definition: grid_refinement.cc:257
hp::Refinement::full_p_adaptivity
void full_p_adaptivity(const hp::DoFHandler< dim, spacedim > &dof_handler)
Definition: refinement.cc:43
KellyErrorEstimator::estimate
static void estimate(const Mapping< dim, spacedim > &mapping, const DoFHandlerType &dof, const Quadrature< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, typename InputVector::value_type > * > &neumann_bc, const InputVector &solution, Vector< float > &error, const ComponentMask &component_mask=ComponentMask(), const Function< spacedim > *coefficients=nullptr, const unsigned int n_threads=numbers::invalid_unsigned_int, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id, const types::material_id material_id=numbers::invalid_material_id, const Strategy strategy=cell_diameter_over_24)
QGauss
Definition: quadrature_lib.h:40
unsigned int
Particles::Generators::quadrature_points
void quadrature_points(const Triangulation< dim, spacedim > &triangulation, const Quadrature< dim > &quadrature, const std::vector< std::vector< BoundingBox< spacedim >>> &global_bounding_boxes, ParticleHandler< dim, spacedim > &particle_handler, const Mapping< dim, spacedim > &mapping=StaticMappingQ1< dim, spacedim >::mapping)
Definition: generators.cc:442
hp::Refinement
Definition: refinement.h:130
hp::Refinement::p_adaptivity_from_reference
void p_adaptivity_from_reference(const hp::DoFHandler< dim, spacedim > &dof_handler, const Vector< Number > &criteria, const Vector< Number > &references, const ComparisonFunction< typename identity< Number >::type > &compare_refine, const ComparisonFunction< typename identity< Number >::type > &compare_coarsen)
Definition: refinement.cc:458
hp::Refinement::force_p_over_h
void force_p_over_h(const hp::DoFHandler< dim, spacedim > &dof_handler)
Definition: refinement.cc:622
Function
Definition: function.h:151
triangulation
const typename ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
Definition: p4est_wrappers.cc:69
Vector
Definition: mapping_q1_eulerian.h:32