Reference documentation for deal.II version 9.1.1
|
Namespaces | |
boost | |
internal | |
Classes | |
struct | CompareDownstream |
struct | ComparePointwiseDownstream |
Functions | |
template<typename DoFHandlerType > | |
void | Cuthill_McKee (DoFHandlerType &dof_handler, const bool reversed_numbering=false, const bool use_constraints=false, const std::vector< types::global_dof_index > &starting_indices=std::vector< types::global_dof_index >()) |
template<typename DoFHandlerType > | |
void | compute_Cuthill_McKee (std::vector< types::global_dof_index > &new_dof_indices, const DoFHandlerType &, const bool reversed_numbering=false, const bool use_constraints=false, const std::vector< types::global_dof_index > &starting_indices=std::vector< types::global_dof_index >()) |
template<typename DoFHandlerType > | |
void | Cuthill_McKee (DoFHandlerType &dof_handler, const unsigned int level, const bool reversed_numbering=false, const std::vector< types::global_dof_index > &starting_indices=std::vector< types::global_dof_index >()) |
static ::ExceptionBase & | ExcDoFHandlerNotInitialized () |
static ::ExceptionBase & | ExcInvalidComponentOrder () |
static ::ExceptionBase & | ExcNotDGFEM () |
Component-wise numberings | |
template<typename DoFHandlerType > | |
void | component_wise (DoFHandlerType &dof_handler, const std::vector< unsigned int > &target_component=std::vector< unsigned int >()) |
template<typename DoFHandlerType > | |
void | component_wise (DoFHandlerType &dof_handler, const unsigned int level, const std::vector< unsigned int > &target_component=std::vector< unsigned int >()) |
template<int dim, int spacedim, typename CellIterator > | |
types::global_dof_index | compute_component_wise (std::vector< types::global_dof_index > &new_dof_indices, const CellIterator &start, const typename identity< CellIterator >::type &end, const std::vector< unsigned int > &target_component, const bool is_level_operation) |
Block-wise numberings | |
template<int dim, int spacedim> | |
void | block_wise (DoFHandler< dim, spacedim > &dof_handler) |
template<int dim, int spacedim> | |
void | block_wise (DoFHandler< dim, spacedim > &dof_handler, const unsigned int level) |
template<int dim, int spacedim> | |
void | block_wise (hp::DoFHandler< dim, spacedim > &dof_handler) |
template<int dim, int spacedim, class ITERATOR , class ENDITERATOR > | |
types::global_dof_index | compute_block_wise (std::vector< types::global_dof_index > &new_dof_indices, const ITERATOR &start, const ENDITERATOR &end, bool is_level_operation) |
Various cell-wise numberings | |
template<typename DoFHandlerType > | |
void | hierarchical (DoFHandlerType &dof_handler) |
template<typename DoFHandlerType > | |
void | cell_wise (DoFHandlerType &dof_handler, const std::vector< typename DoFHandlerType::active_cell_iterator > &cell_order) |
template<typename DoFHandlerType > | |
void | compute_cell_wise (std::vector< types::global_dof_index > &renumbering, std::vector< types::global_dof_index > &inverse_renumbering, const DoFHandlerType &dof_handler, const std::vector< typename DoFHandlerType::active_cell_iterator > &cell_order) |
template<typename DoFHandlerType > | |
void | cell_wise (DoFHandlerType &dof_handler, const unsigned int level, const std::vector< typename DoFHandlerType::level_cell_iterator > &cell_order) |
template<typename DoFHandlerType > | |
void | compute_cell_wise (std::vector< types::global_dof_index > &renumbering, std::vector< types::global_dof_index > &inverse_renumbering, const DoFHandlerType &dof_handler, const unsigned int level, const std::vector< typename DoFHandlerType::level_cell_iterator > &cell_order) |
Directional numberings | |
template<typename DoFHandlerType > | |
void | downstream (DoFHandlerType &dof_handler, const Tensor< 1, DoFHandlerType::space_dimension > &direction, const bool dof_wise_renumbering=false) |
template<typename DoFHandlerType > | |
void | downstream (DoFHandlerType &dof_handler, const unsigned int level, const Tensor< 1, DoFHandlerType::space_dimension > &direction, const bool dof_wise_renumbering=false) |
template<typename DoFHandlerType > | |
void | compute_downstream (std::vector< types::global_dof_index > &new_dof_indices, std::vector< types::global_dof_index > &reverse, const DoFHandlerType &dof_handler, const Tensor< 1, DoFHandlerType::space_dimension > &direction, const bool dof_wise_renumbering) |
template<typename DoFHandlerType > | |
void | compute_downstream (std::vector< types::global_dof_index > &new_dof_indices, std::vector< types::global_dof_index > &reverse, const DoFHandlerType &dof_handler, const unsigned int level, const Tensor< 1, DoFHandlerType::space_dimension > &direction, const bool dof_wise_renumbering) |
template<typename DoFHandlerType > | |
void | clockwise_dg (DoFHandlerType &dof_handler, const Point< DoFHandlerType::space_dimension > ¢er, const bool counter=false) |
template<typename DoFHandlerType > | |
void | clockwise_dg (DoFHandlerType &dof_handler, const unsigned int level, const Point< DoFHandlerType::space_dimension > ¢er, const bool counter=false) |
template<typename DoFHandlerType > | |
void | compute_clockwise_dg (std::vector< types::global_dof_index > &new_dof_indices, const DoFHandlerType &dof_handler, const Point< DoFHandlerType::space_dimension > ¢er, const bool counter) |
Selective and random numberings | |
template<typename DoFHandlerType > | |
void | sort_selected_dofs_back (DoFHandlerType &dof_handler, const std::vector< bool > &selected_dofs) |
template<typename DoFHandlerType > | |
void | sort_selected_dofs_back (DoFHandlerType &dof_handler, const std::vector< bool > &selected_dofs, const unsigned int level) |
template<typename DoFHandlerType > | |
void | compute_sort_selected_dofs_back (std::vector< types::global_dof_index > &new_dof_indices, const DoFHandlerType &dof_handler, const std::vector< bool > &selected_dofs) |
template<typename DoFHandlerType > | |
void | compute_sort_selected_dofs_back (std::vector< types::global_dof_index > &new_dof_indices, const DoFHandlerType &dof_handler, const std::vector< bool > &selected_dofs, const unsigned int level) |
template<typename DoFHandlerType > | |
void | random (DoFHandlerType &dof_handler) |
template<typename DoFHandlerType > | |
void | random (DoFHandlerType &dof_handler, const unsigned int level) |
template<typename DoFHandlerType > | |
void | compute_random (std::vector< types::global_dof_index > &new_dof_indices, const DoFHandlerType &dof_handler) |
template<typename DoFHandlerType > | |
void | compute_random (std::vector< types::global_dof_index > &new_dof_indices, const DoFHandlerType &dof_handler, const unsigned int level) |
Numberings based on cell attributes | |
template<typename DoFHandlerType > | |
void | subdomain_wise (DoFHandlerType &dof_handler) |
template<typename DoFHandlerType > | |
void | compute_subdomain_wise (std::vector< types::global_dof_index > &new_dof_indices, const DoFHandlerType &dof_handler) |
Implementation of a number of renumbering algorithms for the degrees of freedom on a triangulation. The functions in this namespace compute new indices for each degree of freedom of a DoFHandler or hp::DoFHandler object, and then call DoFHandler::renumber_dofs() or hp::DoFHandler::renumber_dofs().
Within this class, the Cuthill-McKee algorithm is implemented. It starts at a degree of freedom, searches the other DoFs for those which are coupled with the one we started with and numbers these in a certain way. It then finds the second level of DoFs, namely those that couple with those of the previous level (which were those that coupled with the initial DoF) and numbers these. And so on. For the details of the algorithm, especially the numbering within each level, please see H. R. Schwarz: "Methode der finiten Elemente". The reverse Cuthill-McKee algorithm does the same job, but numbers all elements in the reverse order.
These algorithms have one major drawback: they require a good starting point, i.e. the degree of freedom index that will get a new index of zero. The renumbering functions therefore allow the caller to specify such an initial DoF, e.g. by exploiting knowledge of the actual topology of the domain. It is also possible to give several starting indices, which may be used to simulate a simple upstream numbering (by giving the inflow dofs as starting values) or to make preconditioning faster (by letting the Dirichlet boundary indices be starting points).
If no starting index is given, one is chosen automatically, namely one with the smallest coordination number (the coordination number is the number of other dofs this dof couples with). This dof is usually located on the boundary of the domain. There is, however, large ambiguity in this when using the hierarchical meshes used in this library, since in most cases the computational domain is not approximated by tilting and deforming elements and by plugging together variable numbers of elements at vertices, but rather by hierarchical refinement. There is therefore a large number of dofs with equal coordination numbers. The renumbering algorithms will therefore not give optimal results.
In the book of Schwarz (H.R.Schwarz: Methode der finiten Elemente), it is advised to test many starting points, if possible all with the smallest coordination number and also those with slightly higher numbers. However, this seems only possible for meshes with at most several dozen or a few hundred elements found in small engineering problems of the early 1980s (the second edition was published in 1984), but certainly not with those used in this library, featuring several 10,000 to a few 100,000 elements.
The renumbering algorithms need quite a lot of memory, since they have to store for each dof with which other dofs it couples. This is done using a SparsityPattern object used to store the sparsity pattern of matrices. It is not useful for the user to do anything between distributing the dofs and renumbering, i.e. the calls to DoFHandler::distribute_dofs and DoFHandler::renumber_dofs should follow each other immediately. If you try to create a sparsity pattern or anything else in between, these will be invalid afterwards.
The renumbering may take care of dof-to-dof couplings only induced by eliminating constraints. In addition to the memory consumption mentioned above, this also takes quite some computational time, but it may be switched off upon calling the renumber_dofs
function. This will then give inferior results, since knots in the graph (representing dofs) are not found to be neighbors even if they would be after condensation.
The renumbering algorithms work on a purely algebraic basis, due to the isomorphism between the graph theoretical groundwork underlying the algorithms and binary matrices (matrices of which the entries are binary values) represented by the sparsity patterns. In special, the algorithms do not try to exploit topological knowledge (e.g. corner detection) to find appropriate starting points. This way, however, they work in arbitrary space dimension.
If you want to give starting points, you may give a list of dof indices which will form the first step of the renumbering. The dofs of the list will be consecutively numbered starting with zero, i.e. this list is not renumbered according to the coordination number of the nodes. Indices not in the allowed range are deleted. If no index is allowed, the algorithm will search for its own starting point.
The renumbering schemes mentioned above do not lead to optimal results. However, after all there is no algorithm that accomplishes this within reasonable time. There are situations where the lack of optimality even leads to worse results than with the original, crude, levelwise numbering scheme; one of these examples is a mesh of four cells of which always those cells are refined which are neighbors to the center (you may call this mesh a ‘zoom in’ mesh). In one such example the bandwidth was increased by about 50 per cent.
In most other cases, the bandwidth is reduced significantly. The reduction is the better the less structured the grid is. With one grid where the cells were refined according to a random driven algorithm, the bandwidth was reduced by a factor of six.
Using the constraint information usually leads to reductions in bandwidth of 10 or 20 per cent, but may for some very unstructured grids also lead to an increase. You have to weigh the decrease in your case with the time spent to use the constraint information, which usually is several times longer than the ‘pure’ renumbering algorithm.
In almost all cases, the renumbering scheme finds a corner to start with. Since there is more than one corner in most grids and since even an interior degree of freedom may be a better starting point, giving the starting point by the user may be a viable way if you have a simple scheme to derive a suitable point (e.g. by successively taking the third child of the cell top left of the coarsest level, taking its third vertex and the dof index thereof, if you want the top left corner vertex). If you do not know beforehand what your grid will look like (e.g. when using adaptive algorithms), searching a best starting point may be difficult, however, and in many cases will not justify the effort.
For finite elements composed of several base elements using the FESystem class, or for elements which provide several components themselves, it may be of interest to sort the DoF indices by component. This will then bring out the block matrix structure, since otherwise the degrees of freedom are numbered cell-wise without taking into account that they may belong to different components. For example, one may want to sort degree of freedom for a Stokes discretization so that we first get all velocities and then all the pressures so that the resulting matrix naturally decomposes into a \(2\times 2\) system.
This kind of numbering may be obtained by calling the component_wise() function of this class. Since it does not touch the order of indices within each component, it may be worthwhile to first renumber using the Cuthill- McKee or a similar algorithm and afterwards renumbering component-wise. This will bring out the matrix structure and additionally have a good numbering within each block.
The component_wise() function allows not only to honor enumeration based on vector components, but also allows to group together vector components into "blocks" using a defaulted argument to the various DoFRenumbering::component_wise() functions (see GlossComponent vs GlossBlock for a description of the difference). The blocks designated through this argument may, but do not have to be, equal to the blocks that the finite element reports. For example, a typical Stokes element would be
This element has dim+1
vector components and equally many blocks. However, one may want to consider the velocities as one logical block so that all velocity degrees of freedom are enumerated the same way, independent of whether they are \(x\)- or \(y\)-velocities. This is done, for example, in step-20 and step-22 as well as several other tutorial programs.
On the other hand, if you really want to use block structure reported by the finite element itself (a case that is often the case if you have finite elements that have multiple vector components, e.g. the FE_RaviartThomas or FE_Nedelec elements) then you can use the DoFRenumbering::block_wise instead of the DoFRenumbering::component_wise functions.
Given an ordered vector of cells, the function cell_wise() sorts the degrees of freedom such that degrees on earlier cells of this vector will occur before degrees on later cells.
This rule produces a well-defined ordering for discontinuous Galerkin methods (FE_DGP, FE_DGQ). For continuous methods, we use the additional rule that each degree of freedom is ordered according to the first cell in the ordered vector it belongs to.
Applications of this scheme are downstream() and clock_wise_dg(). The first orders the cells according to a downstream direction and then applies cell_wise().
The random() function renumbers degrees of freedom randomly. This function is probably seldom of use, except to check the dependence of solvers (iterative or direct ones) on the numbering of the degrees of freedom.
As a benchmark of comparison, let us consider what the different sparsity patterns produced by the various algorithms when using the \(Q_2^d\times Q_1\) element combination typically employed in the discretization of Stokes equations, when used on the mesh obtained in step-22 after one adaptive mesh refinement in 3d. The space dimension together with the coupled finite element leads to a rather dense system matrix with, on average around 180 nonzero entries per row. After applying each of the reordering strategies shown below, the degrees of freedom are also sorted using DoFRenumbering::component_wise into velocity and pressure groups; this produces the \(2\times 2\) block structure seen below with the large velocity-velocity block at top left, small pressure-pressure block at bottom right, and coupling blocks at top right and bottom left.
The goal of reordering strategies is to improve the preconditioner. In step-22 we use a SparseILU to preconditioner for the velocity-velocity block at the top left. The quality of the preconditioner can then be measured by the number of CG iterations required to solve a linear system with this block. For some of the reordering strategies below we record this number for adaptive refinement cycle 3, with 93176 degrees of freedom; because we solve several linear systems with the same matrix in the Schur complement, the average number of iterations is reported. The lower the number the better the preconditioner and consequently the better the renumbering of degrees of freedom is suited for this task. We also state the run-time of the program, in part determined by the number of iterations needed, for the first 4 cycles on one of our machines. Note that the reported times correspond to the run time of the entire program, not just the affected solver; if a program runs twice as fast with one particular ordering than with another one, then this means that the actual solver is actually several times faster.
Enumeration as produced by deal.II's DoFHandler::distribute_dofs function and no further reordering apart from the component-wise one. With this renumbering, we needed an average of 92.2 iterations for the testcase outlined above, and a runtime of 7min53s. | Random enumeration as produced by applying DoFRenumbering::random after calling DoFHandler::distribute_dofs. This enumeration produces nonzero entries in matrices pretty much everywhere, appearing here as an entirely unstructured matrix. With this renumbering, we needed an average of 71 iterations for the testcase outlined above, and a runtime of 10min55s. The longer runtime despite less iterations compared to the default ordering may be due to the fact that computing and applying the ILU requires us to jump back and forth all through memory due to the lack of localization of matrix entries around the diagonal; this then leads to many cache misses and consequently bad timings. | Cuthill-McKee enumeration as produced by calling the deal.II implementation of the algorithm provided by DoFRenumbering::Cuthill_McKee after DoFHandler::distribute_dofs. With this renumbering, we needed an average of 57.3 iterations for the testcase outlined above, and a runtime of 6min10s. |
Cuthill- McKee enumeration as produced by calling the BOOST implementation of the algorithm provided by DoFRenumbering::boost::Cuthill_McKee after DoFHandler::distribute_dofs. With this renumbering, we needed an average of 51.7 iterations for the testcase outlined above, and a runtime of 5min52s. | King enumeration as produced by calling the BOOST implementation of the algorithm provided by DoFRenumbering::boost::king_ordering after DoFHandler::distribute_dofs. The sparsity pattern appears denser than with BOOST's Cuthill-McKee algorithm; however, this is only an illusion: the number of nonzero entries is the same, they are simply not as well clustered. With this renumbering, we needed an average of 51.0 iterations for the testcase outlined above, and a runtime of 5min03s. Although the number of iterations is only slightly less than with BOOST's Cuthill-McKee implementation, runtime is significantly less. This, again, may be due to cache effects. As a consequence, this is the algorithm best suited to the testcase, and is in fact used in step-22. | Minimum degree enumeration as produced by calling the BOOST implementation of the algorithm provided by DoFRenumbering::boost::minimum_degree after DoFHandler::distribute_dofs. The minimum degree algorithm does not attempt to minimize the bandwidth of a matrix but to minimize the amount of fill-in a LU decomposition would produce, i.e. the number of places in the matrix that would be occupied by elements of an LU decomposition that are not already occupied by elements of the original matrix. The resulting sparsity pattern obviously has an entirely different structure than the ones produced by algorithms trying to minimize the bandwidth. With this renumbering, we needed an average of 58.9 iterations for the testcase outlined above, and a runtime of 6min11s. |
Downstream enumeration using DoFRenumbering::downstream using a direction that points diagonally through the domain. With this renumbering, we needed an average of 90.5 iterations for the testcase outlined above, and a runtime of 7min05s. |
Most of the algorithms listed above also work on multigrid degree of freedom numberings. Refer to the actual function declarations to get more information on this.
void DoFRenumbering::Cuthill_McKee | ( | DoFHandlerType & | dof_handler, |
const bool | reversed_numbering = false , |
||
const bool | use_constraints = false , |
||
const std::vector< types::global_dof_index > & | starting_indices = std::vector<types::global_dof_index>() |
||
) |
Renumber the degrees of freedom according to the Cuthill-McKee method, possibly using the reverse numbering scheme.
See the general documentation of this class for details on the different methods.
As an example of the results of this algorithm, take a look at the comparison of various algorithms in the documentation of the DoFRenumbering namespace.
dof_handler | The DoFHandler or hp::DoFHandler object to work on. |
reversed_numbering | Whether to use the original Cuthill-McKee algorithm, or to reverse the ordering. |
use_constraints | Whether or not to use hanging node constraints in determining the reordering of degrees of freedom. |
starting_indices | A set of degrees of freedom that form the first level of renumbered degrees of freedom. If the set is empty, then a single starting entry is chosen automatically among those that have the smallest number of others that couple with it. |
If the given DoFHandler uses a distributed triangulation (i.e., if dof_handler.locally_owned() is not the complete index set), the renumbering is performed on each processor's degrees of freedom individually, without any communication between processors. In other words, the resulting renumbering is an attempt at minimizing the bandwidth of each diagonal block of the matrix corresponding to one processor separately, without making an attempt at minimizing the bandwidth of the global matrix. Furthermore, the renumbering reuses exactly the same set of DoF indices that each processor used before. In other words, if the previous numbering of DoFs on one processor used a contiguous range of DoF indices, then so will the DoFs on that processor after the renumbering, and they will occupy the same range. The same is true if the previous numbering of DoFs on a processor consisted of a number of index ranges or single indices: after renumbering, the locally owned DoFs on that processor will use the exact same indices, just in a different order.
In addition, if the DoFHandler is built on a parallel triangulation, then on every processor, the starting indices for renumbering need to be a (possibly empty) subset of the locally active degrees of freedom. In general, these starting indices will be different on each processor (unless of course you pass an empty list as is the default), and each processor will use them as starting indices for the local renumbering on that processor.
The starting indices must be locally active degrees of freedom, but the function will only renumber the locally owned subset of the locally owned DoFs. The function accepts starting indices from the largest set of locally active degrees of freedom because a typical renumbering operation with this function starts with indices that are located on the boundary – in the case of the current function, that would be the boundary between processor subdomains. Since the degrees of freedom that are located on subdomain interfaces may be owned by either one of the two processors that own the adjacent subdomains, it is not always easy to identify starting indices that are locally owned. On the other hand, all degrees of freedom on subdomain interfaces are locally active, and so the function accepts them as starting indices even though it can only renumber them on a given processor if they are also locally owned.
Definition at line 366 of file dof_renumbering.cc.
void DoFRenumbering::compute_Cuthill_McKee | ( | std::vector< types::global_dof_index > & | new_dof_indices, |
const DoFHandlerType & | dof_handler, | ||
const bool | reversed_numbering = false , |
||
const bool | use_constraints = false , |
||
const std::vector< types::global_dof_index > & | starting_indices = std::vector<types::global_dof_index>() |
||
) |
Compute the renumbering vector needed by the Cuthill_McKee() function. This function does not perform the renumbering on the DoFHandler DoFs but only returns the renumbering vector.
See the Cuthill_McKee() function for an explanation of the arguments.
Definition at line 390 of file dof_renumbering.cc.
void DoFRenumbering::Cuthill_McKee | ( | DoFHandlerType & | dof_handler, |
const unsigned int | level, | ||
const bool | reversed_numbering = false , |
||
const std::vector< types::global_dof_index > & | starting_indices = std::vector<types::global_dof_index>() |
||
) |
Renumber the degrees of freedom according to the Cuthill-McKee method, eventually using the reverse numbering scheme, in this case for a multigrid numbering of degrees of freedom.
You can give a triangulation level to which this function is to be applied. Since with a level-wise numbering there are no hanging nodes, no constraints can be used, so the respective parameter of the previous function is omitted.
See the general documentation of this class for details on the different methods.
Definition at line 601 of file dof_renumbering.cc.
void DoFRenumbering::component_wise | ( | DoFHandlerType & | dof_handler, |
const std::vector< unsigned int > & | target_component = std::vector<unsigned int>() |
||
) |
Sort the degrees of freedom by vector component. The numbering within each component is not touched, so a degree of freedom with index \(i\), belonging to some component, and another degree of freedom with index \(j\) belonging to the same component will be assigned new indices \(n(i)\) and \(n(j)\) with \(n(i)<n(j)\) if \(i<j\) and \(n(i)>n(j)\) if \(i>j\).
You can specify that the components are ordered in a different way than suggested by the FESystem object you use. To this end, set up the vector target_component
such that the entry at index i
denotes the number of the target component for dofs with component i
in the FESystem. Naming the same target component more than once is possible and results in a blocking of several components into one. This is discussed in step-22. If you omit this argument, the same order as given by the finite element is used.
If one of the base finite elements from which the global finite element under consideration here, is a non-primitive one, i.e. its shape functions have more than one non-zero component, then it is not possible to associate these degrees of freedom with a single vector component. In this case, they are associated with the first vector component to which they belong.
For finite elements with only one component, or a single non-primitive base element, this function is the identity operation.
Definition at line 630 of file dof_renumbering.cc.
void DoFRenumbering::component_wise | ( | DoFHandlerType & | dof_handler, |
const unsigned int | level, | ||
const std::vector< unsigned int > & | target_component = std::vector<unsigned int>() |
||
) |
Sort the degrees of freedom by component. It does the same thing as the above function, only that it does this for one single level of a multilevel discretization. The non-multigrid part of the DoFHandler is not touched.
Definition at line 666 of file dof_renumbering.cc.
types::global_dof_index DoFRenumbering::compute_component_wise | ( | std::vector< types::global_dof_index > & | new_dof_indices, |
const CellIterator & | start, | ||
const typename identity< CellIterator >::type & | end, | ||
const std::vector< unsigned int > & | target_component, | ||
const bool | is_level_operation | ||
) |
Compute the renumbering vector needed by the component_wise() functions. Does not perform the renumbering on the DoFHandler dofs but returns the renumbering vector.
Definition at line 699 of file dof_renumbering.cc.
void DoFRenumbering::block_wise | ( | DoFHandler< dim, spacedim > & | dof_handler | ) |
Sort the degrees of freedom by vector block. The numbering within each block is not touched, so a degree of freedom with index \(i\), belonging to some block, and another degree of freedom with index \(j\) belonging to the same block will be assigned new indices \(n(i)\) and \(n(j)\) with \(n(i)<n(j)\) if \(i<j\) and \(n(i)>n(j)\) if \(i>j\).
Definition at line 961 of file dof_renumbering.cc.
void DoFRenumbering::block_wise | ( | DoFHandler< dim, spacedim > & | dof_handler, |
const unsigned int | level | ||
) |
Sort the degrees of freedom by vector block. It does the same thing as the above function, only that it does this for one single level of a multilevel discretization. The non-multigrid part of the DoFHandler is not touched.
Definition at line 1032 of file dof_renumbering.cc.
void DoFRenumbering::block_wise | ( | hp::DoFHandler< dim, spacedim > & | dof_handler | ) |
Sort the degrees of freedom by block. It does the same thing as the above function.
This function only succeeds if each of the elements in the hp::FECollection attached to the hp::DoFHandler argument has exactly the same number of blocks (see the glossary for more information). Note that this is not always given: while the hp::FECollection class ensures that all of its elements have the same number of vector components, they need not have the same number of blocks. At the same time, this function here needs to match individual blocks across elements and therefore requires that elements have the same number of blocks and that subsequent blocks in one element have the same meaning as in another element.
Definition at line 1001 of file dof_renumbering.cc.
types::global_dof_index DoFRenumbering::compute_block_wise | ( | std::vector< types::global_dof_index > & | new_dof_indices, |
const ITERATOR & | start, | ||
const ENDITERATOR & | end, | ||
bool | is_level_operation | ||
) |
Compute the renumbering vector needed by the block_wise() functions. Does not perform the renumbering on the DoFHandler dofs but returns the renumbering vector.
Definition at line 1068 of file dof_renumbering.cc.
void DoFRenumbering::hierarchical | ( | DoFHandlerType & | dof_handler | ) |
Renumber the degrees cell by cell by traversing the cells in Z order.
There are two reasons to use this function:
For meshes based on parallel::shared::Triangulation, the situation is more complex. Here, the set of locally owned cells is determined by a partitioning algorithm (selected by passing an object of type parallel::shared::Triangulation::Settings to the constructor of the triangulation), and in general these partitioning algorithms may assign cells to subdomains based on decisions that may have nothing to do with the Z order. (Though it is possible to select these flags in a way so that partitioning uses the Z order.) As a consequence, the cells of one subdomain are not contiguous in Z order, and if one renumbered degrees of freedom based on the Z order of cells, one would generally end up with DoF indices that on each processor do not form a contiguous range. This is often inconvenient (for example, because PETSc cannot store vectors and matrices for which the locally owned set of indices is not contiguous), and consequently this function uses the following algorithm for parallel::shared::Triangulation objects:
Definition at line 1373 of file dof_renumbering.cc.
void DoFRenumbering::cell_wise | ( | DoFHandlerType & | dof_handler, |
const std::vector< typename DoFHandlerType::active_cell_iterator > & | cell_order | ||
) |
Renumber degrees of freedom by cell. The function takes a vector of cell iterators (which needs to list all locally owned active cells of the DoF handler objects) and will give degrees of freedom new indices based on where in the given list of cells the cell is on which the degree of freedom is located. Degrees of freedom that exist at the interface between two or more cells will be numbered when they are encountered first.
Degrees of freedom that are encountered first on the same cell retain their original ordering before the renumbering step.
[in,out] | dof_handler | The DoFHandler whose degrees of freedom are to be renumbered. |
[in] | cell_order | A vector that contains the order of the cells that defines the order in which degrees of freedom should be renumbered. |
cell_order
must have size dof_handler.get_triangulation().n_active_cells()
, whereas in case of parallel triangulation its size should be parallel::Triangulation::n_locally_owned_active_cells(). Every active cell iterator of that triangulation needs to be present in cell_order
exactly once. Definition at line 1591 of file dof_renumbering.cc.
void DoFRenumbering::compute_cell_wise | ( | std::vector< types::global_dof_index > & | renumbering, |
std::vector< types::global_dof_index > & | inverse_renumbering, | ||
const DoFHandlerType & | dof_handler, | ||
const std::vector< typename DoFHandlerType::active_cell_iterator > & | cell_order | ||
) |
Compute a renumbering of degrees of freedom by cell. The function takes a vector of cell iterators (which needs to list all locally owned active cells of the DoF handler objects) and will give degrees of freedom new indices based on where in the given list of cells the cell is on which the degree of freedom is located. Degrees of freedom that exist at the interface between two or more cells will be numbered when they are encountered first.
Degrees of freedom that are encountered first on the same cell retain their original ordering before the renumbering step.
[out] | renumbering | A vector of length dof_handler.n_locally_owned_dofs() that contains for each degree of freedom (in their current numbering) their future DoF index. This vector therefore presents a (very particular) permutation of the current DoF indices. |
[out] | inverse_renumbering | The reverse of the permutation returned in the previous argument. In case of parallel::Triangulation the inverse is within locally owned DoFs. |
[in] | dof_handler | The DoFHandler whose degrees of freedom are to be renumbered. |
[in] | cell_order | A vector that contains the order of the cells that defines the order in which degrees of freedom should be renumbered. |
cell_order
must have size dof_handler.get_triangulation().n_active_cells()
, whereas in case of parallel triangulation its size should be parallel::Triangulation::n_locally_owned_active_cells(). Every active cell iterator of that triangulation needs to be present in cell_order
exactly once. i
between zero and dof_handler.n_locally_owned_dofs()
, the condition renumbering[inverse_renumbering[i]] == dof_handler.locally_owned_dofs().nth_index_in_set(i)
will hold. void DoFRenumbering::cell_wise | ( | DoFHandlerType & | dof_handler, |
const unsigned int | level, | ||
const std::vector< typename DoFHandlerType::level_cell_iterator > & | cell_order | ||
) |
Like the other cell_wise() function, but for one level of a multilevel enumeration of degrees of freedom.
void DoFRenumbering::compute_cell_wise | ( | std::vector< types::global_dof_index > & | renumbering, |
std::vector< types::global_dof_index > & | inverse_renumbering, | ||
const DoFHandlerType & | dof_handler, | ||
const unsigned int | level, | ||
const std::vector< typename DoFHandlerType::level_cell_iterator > & | cell_order | ||
) |
Like the other compute_cell_wise() function, but for one level of a multilevel enumeration of degrees of freedom.
void DoFRenumbering::downstream | ( | DoFHandlerType & | dof_handler, |
const Tensor< 1, DoFHandlerType::space_dimension > & | direction, | ||
const bool | dof_wise_renumbering = false |
||
) |
Downstream numbering with respect to a constant flow direction. If the additional argument dof_wise_renumbering
is set to false
, the numbering is performed cell-wise, otherwise it is performed based on the location of the support points.
The cells are sorted such that the centers of cells numbered higher are further downstream with respect to the constant vector direction
than the centers of cells numbered lower. Even if this yields a downstream numbering with respect to the flux on the edges for fairly general grids, this might not be guaranteed for all meshes.
If the dof_wise_renumbering
argument is set to false
, this function produces a downstream ordering of the mesh cells and calls cell_wise(). Therefore, the output only makes sense for Discontinuous Galerkin Finite Elements (all degrees of freedom have to be associated with the interior of the cell in that case) in that case.
If dof_wise_renumbering
is set to true
, the degrees of freedom are renumbered based on the support point location of the individual degrees of freedom (obviously, the finite element needs to define support points for this to work). The numbering of points with the same position in downstream location (e.g. those parallel to the flow direction, or several dofs within a FESystem) will be unaffected.
Definition at line 1757 of file dof_renumbering.cc.
void DoFRenumbering::downstream | ( | DoFHandlerType & | dof_handler, |
const unsigned int | level, | ||
const Tensor< 1, DoFHandlerType::space_dimension > & | direction, | ||
const bool | dof_wise_renumbering = false |
||
) |
Cell-wise downstream numbering with respect to a constant flow direction on one level of a multigrid hierarchy. See the other function with the same name.
Definition at line 1880 of file dof_renumbering.cc.
void DoFRenumbering::compute_downstream | ( | std::vector< types::global_dof_index > & | new_dof_indices, |
std::vector< types::global_dof_index > & | reverse, | ||
const DoFHandlerType & | dof_handler, | ||
const Tensor< 1, DoFHandlerType::space_dimension > & | direction, | ||
const bool | dof_wise_renumbering | ||
) |
Compute the set of renumbering indices needed by the downstream() function. Does not perform the renumbering on the DoFHandler dofs but returns the renumbering vector.
Definition at line 1773 of file dof_renumbering.cc.
void DoFRenumbering::compute_downstream | ( | std::vector< types::global_dof_index > & | new_dof_indices, |
std::vector< types::global_dof_index > & | reverse, | ||
const DoFHandlerType & | dof_handler, | ||
const unsigned int | level, | ||
const Tensor< 1, DoFHandlerType::space_dimension > & | direction, | ||
const bool | dof_wise_renumbering | ||
) |
Compute the set of renumbering indices needed by the downstream() function. Does not perform the renumbering on the DoFHandler dofs but returns the renumbering vector.
Definition at line 1897 of file dof_renumbering.cc.
void DoFRenumbering::clockwise_dg | ( | DoFHandlerType & | dof_handler, |
const Point< DoFHandlerType::space_dimension > & | center, | ||
const bool | counter = false |
||
) |
Cell-wise clockwise numbering.
This function produces a (counter)clockwise ordering of the mesh cells with respect to the hub center
and calls cell_wise(). Therefore, it only works with Discontinuous Galerkin Finite Elements, i.e. all degrees of freedom have to be associated with the interior of the cell.
Definition at line 2054 of file dof_renumbering.cc.
void DoFRenumbering::clockwise_dg | ( | DoFHandlerType & | dof_handler, |
const unsigned int | level, | ||
const Point< DoFHandlerType::space_dimension > & | center, | ||
const bool | counter = false |
||
) |
Cell-wise clockwise numbering on one level of a multigrid hierarchy. See the other function with the same name.
Definition at line 2096 of file dof_renumbering.cc.
void DoFRenumbering::compute_clockwise_dg | ( | std::vector< types::global_dof_index > & | new_dof_indices, |
const DoFHandlerType & | dof_handler, | ||
const Point< DoFHandlerType::space_dimension > & | center, | ||
const bool | counter | ||
) |
Compute the renumbering vector needed by the clockwise_dg() functions. Does not perform the renumbering on the DoFHandler dofs but returns the renumbering vector.
Definition at line 2068 of file dof_renumbering.cc.
void DoFRenumbering::sort_selected_dofs_back | ( | DoFHandlerType & | dof_handler, |
const std::vector< bool > & | selected_dofs | ||
) |
Sort those degrees of freedom which are tagged with true
in the selected_dofs
array to the back of the DoF numbers. The sorting is stable, i.e. the relative order within the tagged degrees of freedom is preserved, as is the relative order within the untagged ones.
selected_dofs
array must have as many elements as the dof_handler
has degrees of freedom. Definition at line 1476 of file dof_renumbering.cc.
void DoFRenumbering::sort_selected_dofs_back | ( | DoFHandlerType & | dof_handler, |
const std::vector< bool > & | selected_dofs, | ||
const unsigned int | level | ||
) |
Sort those degrees of freedom which are tagged with true
in the selected_dofs
array on the level level
to the back of the DoF numbers. The sorting is stable, i.e. the relative order within the tagged degrees of freedom is preserved, as is the relative order within the untagged ones.
selected_dofs
array must have as many elements as the dof_handler
has degrees of freedom on the given level. Definition at line 1490 of file dof_renumbering.cc.
void DoFRenumbering::compute_sort_selected_dofs_back | ( | std::vector< types::global_dof_index > & | new_dof_indices, |
const DoFHandlerType & | dof_handler, | ||
const std::vector< bool > & | selected_dofs | ||
) |
Compute the renumbering vector needed by the sort_selected_dofs_back() function. Does not perform the renumbering on the DoFHandler dofs but returns the renumbering vector.
selected_dofs
array must have as many elements as the dof_handler
has degrees of freedom. Definition at line 1511 of file dof_renumbering.cc.
void DoFRenumbering::compute_sort_selected_dofs_back | ( | std::vector< types::global_dof_index > & | new_dof_indices, |
const DoFHandlerType & | dof_handler, | ||
const std::vector< bool > & | selected_dofs, | ||
const unsigned int | level | ||
) |
This function computes the renumbering vector on each level needed by the sort_selected_dofs_back() function. Does not perform the renumbering on the DoFHandler dofs but only computes the renumbering and returns the renumbering vector.
selected_dofs
array must have as many elements as the dof_handler
has degrees of freedom on the given level. Definition at line 1549 of file dof_renumbering.cc.
void DoFRenumbering::random | ( | DoFHandlerType & | dof_handler | ) |
Renumber the degrees of freedom in a random way. The result of this function is repeatable in that two runs of the same program will yield the same result. This is achieved by creating a new random number generator with a fixed seed every time this function is entered. In particular, the function therefore does not rely on an external random number generator for which it would matter how often it has been called before this function (or, for that matter, whether other threads running concurrently to this function also draw random numbers).
Definition at line 2123 of file dof_renumbering.cc.
void DoFRenumbering::random | ( | DoFHandlerType & | dof_handler, |
const unsigned int | level | ||
) |
Renumber the degrees of freedom in a random way. It does the same thing as the above function, only that it does this for one single level of a multilevel discretization. The non-multigrid part of the DoFHandler is not touched.
Definition at line 2136 of file dof_renumbering.cc.
void DoFRenumbering::compute_random | ( | std::vector< types::global_dof_index > & | new_dof_indices, |
const DoFHandlerType & | dof_handler | ||
) |
Compute the renumbering vector needed by the random() function. See there for more information on the computed random renumbering.
This function does not perform the renumbering on the DoFHandler dofs but returns the renumbering vector.
Definition at line 2154 of file dof_renumbering.cc.
void DoFRenumbering::compute_random | ( | std::vector< types::global_dof_index > & | new_dof_indices, |
const DoFHandlerType & | dof_handler, | ||
const unsigned int | level | ||
) |
Compute the renumbering vector needed by the random() function. Same as the above function but for a single level of a multilevel discretization.
Definition at line 2184 of file dof_renumbering.cc.
void DoFRenumbering::subdomain_wise | ( | DoFHandlerType & | dof_handler | ) |
Renumber the degrees of freedom such that they are associated with the subdomain id of the cells they are living on, i.e. first all degrees of freedom that belong to cells with subdomain zero, then all with subdomain one, etc. This is useful when doing parallel computations with a standard Triangulation after assigning subdomain ids using a partitioner (see the GridTools::partition_triangulation function for this). Calling this function is unnecessary when using a parallel::shared::Triangulation or parallel::distributed::Triangulation, as the degrees of freedom are already enumerated according to the MPI process id. Therefore, if the underlying triangulation is of this type then an error will be thrown.
Note that degrees of freedom associated with faces, edges, and vertices may be associated with multiple subdomains if they are sitting on partition boundaries. It would therefore be undefined with which subdomain they have to be associated. For this, we use what we get from the DoFTools::get_subdomain_association function.
The algorithm is stable, i.e. if two dofs i,j have i<j
and belong to the same subdomain, then they will be in this order also after reordering.
Definition at line 2215 of file dof_renumbering.cc.
void DoFRenumbering::compute_subdomain_wise | ( | std::vector< types::global_dof_index > & | new_dof_indices, |
const DoFHandlerType & | dof_handler | ||
) |
Compute the renumbering vector needed by the subdomain_wise() function. Does not perform the renumbering on the DoFHandler
dofs but returns the renumbering vector.
Definition at line 2236 of file dof_renumbering.cc.