From: hartmann Date: Mon, 3 Feb 2003 10:19:47 +0000 (+0000) Subject: New compute_* functions, that compute and return the renumbering vectors, but do... X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=fa8999d419a9ed038d95cf5cba4e2f8131a2bfba;p=dealii-svn.git New compute_* functions, that compute and return the renumbering vectors, but do not perform the actual renumbering on the dof_handler. git-svn-id: https://svn.dealii.org/trunk@7002 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/deal.II/include/numerics/dof_renumbering.h b/deal.II/deal.II/include/numerics/dof_renumbering.h index 8847d9502a..bbc8d57e6f 100644 --- a/deal.II/deal.II/include/numerics/dof_renumbering.h +++ b/deal.II/deal.II/include/numerics/dof_renumbering.h @@ -216,6 +216,23 @@ class DoFRenumbering const bool use_constraints = false, const std::vector &starting_indices = std::vector()); + /** + * Computes the renumbering + * vector needed by the + * @p{Cuthill_McKee} + * function. Does not perform the + * renumbering on the + * @p{DoFHandler} dofs but + * returns the renumbering + * vector. + */ + template + static std::vector + compute_Cuthill_McKee (DoFHandler &dof_handler, + const bool reversed_numbering = false, + const bool use_constraints = false, + const std::vector &starting_indices = std::vector()); + /** * Renumber the degrees of freedom * according to the Cuthill-McKee method, @@ -294,9 +311,24 @@ class DoFRenumbering */ template static void - component_wise (DoFHandler &dof_handler, + component_wise (DoFHandler &dof_handler, const std::vector &component_order = std::vector()); + /** + * Computes the renumbering + * vector needed by the + * @p{component_wise} + * function. Does not perform the + * renumbering on the + * @p{DoFHandler} dofs but + * returns the renumbering + * vector. + */ + template + static std::vector + compute_component_wise (DoFHandler &dof_handler, + const std::vector &component_order = std::vector()); + /** * Sort the degrees of freedom by * component. It does the same @@ -332,9 +364,23 @@ class DoFRenumbering */ template static void - cell_wise_dg (DoFHandler &dof_handler, + cell_wise_dg (DoFHandler &dof_handler, const std::vector::cell_iterator> &cell_order); - + + /** + * Computes the renumbering + * vector needed by the + * @p{cell_wise_dg} + * function. Does not perform the + * renumbering on the + * @p{DoFHandler} dofs but + * returns the renumbering + * vector. + */ + template + static std::vector + compute_cell_wise_dg (DoFHandler &dof_handler, + const std::vector::cell_iterator> &cell_order); /** * Cell-wise renumbering on one @@ -381,6 +427,21 @@ class DoFRenumbering downstream_dg (DoFHandler &dof_handler, const Point &direction); + /** + * Computes the renumbering + * vector needed by the + * @p{downstream_dg} + * function. Does not perform the + * renumbering on the + * @p{DoFHandler} dofs but + * returns the renumbering + * vector. + */ + template + static std::vector + compute_downstream_dg (DoFHandler &dof_handler, + const Point &direction); + /** * Cell-wise downstream numbering * with respect to a constant @@ -412,6 +473,21 @@ class DoFRenumbering sort_selected_dofs_back (DoFHandler &dof_handler, const std::vector &selected_dofs); + /** + * Computes the renumbering + * vector needed by the + * @p{sort_selected_dofs_back} + * function. Does not perform the + * renumbering on the + * @p{DoFHandler} dofs but + * returns the renumbering + * vector. + */ + template + static std::vector + compute_sort_selected_dofs_back (DoFHandler &dof_handler, + const std::vector &selected_dofs); + /** * Renumber the degrees of * freedom in a random way. @@ -419,7 +495,20 @@ class DoFRenumbering template static void random (DoFHandler &dof_handler); - + + /** + * Computes the renumbering + * vector needed by the + * @p{random} function. Does not + * perform the renumbering on the + * @p{DoFHandler} dofs but + * returns the renumbering + * vector. + */ + template + static std::vector + compute_random (DoFHandler &dof_handler); + /** * Exception */ diff --git a/deal.II/deal.II/source/dofs/dof_renumbering.cc b/deal.II/deal.II/source/dofs/dof_renumbering.cc index 5f6b159ec7..7239aad623 100644 --- a/deal.II/deal.II/source/dofs/dof_renumbering.cc +++ b/deal.II/deal.II/source/dofs/dof_renumbering.cc @@ -50,11 +50,32 @@ extern "C" long int lrand48 (void); #endif + template -void DoFRenumbering::Cuthill_McKee (DoFHandler &dof_handler, - const bool reversed_numbering, - const bool use_constraints, - const std::vector &starting_indices) +void +DoFRenumbering::Cuthill_McKee (DoFHandler &dof_handler, + const bool reversed_numbering, + const bool use_constraints, + const std::vector &starting_indices) +{ + std::vector renumbering + =compute_Cuthill_McKee(dof_handler, reversed_numbering, + use_constraints, starting_indices); + + // actually perform renumbering; + // this is dimension specific and + // thus needs an own function + dof_handler.renumber_dofs (renumbering); +} + + + +template +std::vector +DoFRenumbering::compute_Cuthill_McKee (DoFHandler &dof_handler, + const bool reversed_numbering, + const bool use_constraints, + const std::vector &starting_indices) { // make the connection graph SparsityPattern sparsity (dof_handler.n_dofs(), @@ -253,10 +274,7 @@ void DoFRenumbering::Cuthill_McKee (DoFHandler &dof_handler i!=new_number.end(); ++i) *i = n_dofs-*i-1; - // actually perform renumbering; - // this is dimension specific and - // thus needs an own function - dof_handler.renumber_dofs (new_number); + return new_number; } @@ -432,14 +450,34 @@ template void DoFRenumbering::component_wise (DoFHandler &dof_handler, const std::vector &component_order_arg) +{ + std::vector renumbering + =compute_component_wise(dof_handler, component_order_arg); + + if (renumbering.size()!=0) + dof_handler.renumber_dofs (renumbering); +} + + + +template +std::vector +DoFRenumbering::compute_component_wise (DoFHandler &dof_handler, + const std::vector &component_order_arg) { const unsigned int dofs_per_cell = dof_handler.get_fe().dofs_per_cell; const FiniteElement &fe = dof_handler.get_fe(); + std::vector new_indices (dof_handler.n_dofs(), + DoFHandler::invalid_dof_index); + // do nothing if the FE has only // one component if (fe.n_components() == 1) - return; + { + new_indices.resize(0); + return new_indices; + } std::vector component_order (component_order_arg); if (component_order.size() == 0) @@ -556,8 +594,6 @@ DoFRenumbering::component_wise (DoFHandler &dof_handler, // components in the order the user // desired to see unsigned int next_free_index = 0; - std::vector new_indices (dof_handler.n_dofs(), - DoFHandler::invalid_dof_index); for (unsigned int c=0; c &dof_handler, Assert (next_free_index == dof_handler.n_dofs(), ExcInternalError()); - dof_handler.renumber_dofs (new_indices); + return new_indices; } @@ -735,6 +771,19 @@ template void DoFRenumbering::sort_selected_dofs_back (DoFHandler &dof_handler, const std::vector &selected_dofs) +{ + std::vector renumbering + =compute_sort_selected_dofs_back(dof_handler, selected_dofs); + + dof_handler.renumber_dofs(renumbering); +} + + + +template +std::vector +DoFRenumbering::compute_sort_selected_dofs_back (DoFHandler &dof_handler, + const std::vector &selected_dofs) { const unsigned int n_dofs = dof_handler.n_dofs(); Assert (selected_dofs.size() == n_dofs, @@ -763,15 +812,30 @@ DoFRenumbering::sort_selected_dofs_back (DoFHandler &dof_handler, Assert (next_unselected == n_selected_dofs, ExcInternalError()); Assert (next_selected == n_dofs, ExcInternalError()); - // now perform the renumbering - dof_handler.renumber_dofs (new_dof_indices); + return new_dof_indices; } + template void -DoFRenumbering::cell_wise_dg (DoFHandler& dof, - const typename std::vector::cell_iterator>& cells) +DoFRenumbering::cell_wise_dg ( + DoFHandler& dof, + const typename std::vector::cell_iterator>& cells) +{ + std::vector renumbering + =compute_cell_wise_dg(dof, cells); + + dof.renumber_dofs(renumbering); +} + + + +template +std::vector +DoFRenumbering::compute_cell_wise_dg ( + DoFHandler& dof, + const typename std::vector::cell_iterator>& cells) { Assert(cells.size() == dof.get_tria().n_active_cells(), ExcDimensionMismatch(cells.size(), @@ -816,16 +880,15 @@ DoFRenumbering::cell_wise_dg (DoFHandler& dof, for (unsigned int i=0;i -void -DoFRenumbering::cell_wise_dg (MGDoFHandler& dof, - const unsigned int level, - const typename std::vector::cell_iterator>& cells) +void DoFRenumbering::cell_wise_dg (MGDoFHandler& dof, + const unsigned int level, + const typename std::vector::cell_iterator>& cells) { Assert(cells.size() == dof.get_tria().n_cells(level), ExcDimensionMismatch(cells.size(), @@ -905,11 +968,24 @@ struct CompCells } }; - + template void DoFRenumbering::downstream_dg (DoFHandler& dof, const Point& direction) +{ + std::vector renumbering + =compute_downstream_dg(dof, direction); + + dof.renumber_dofs(renumbering); +} + + + +template +std::vector +DoFRenumbering::compute_downstream_dg (DoFHandler& dof, + const Point& direction) { std::vector::cell_iterator> ordered_cells(dof.get_tria().n_active_cells()); @@ -921,16 +997,15 @@ DoFRenumbering::downstream_dg (DoFHandler& dof, copy (begin, end, ordered_cells.begin()); sort (ordered_cells.begin(), ordered_cells.end(), comparator); - cell_wise_dg(dof, ordered_cells); + return compute_cell_wise_dg(dof, ordered_cells); } #ifdef ENABLE_MULTIGRID template -void -DoFRenumbering::downstream_dg (MGDoFHandler& dof, - const unsigned int level, - const Point& direction) +void DoFRenumbering::downstream_dg (MGDoFHandler& dof, + const unsigned int level, + const Point& direction) { std::vector::cell_iterator> ordered_cells(dof.get_tria().n_cells(level)); @@ -947,9 +1022,21 @@ DoFRenumbering::downstream_dg (MGDoFHandler& dof, #endif + template void DoFRenumbering::random (DoFHandler &dof_handler) +{ + std::vector renumbering + =compute_random(dof_handler); + + dof_handler.renumber_dofs(renumbering); +} + + +template +std::vector +DoFRenumbering::compute_random (DoFHandler &dof_handler) { const unsigned int n_dofs = dof_handler.n_dofs(); @@ -958,7 +1045,7 @@ DoFRenumbering::random (DoFHandler &dof_handler) new_indices[i] = i; std::random_shuffle (new_indices.begin(), new_indices.end()); - dof_handler.renumber_dofs(new_indices); + return new_indices; } @@ -972,11 +1059,25 @@ void DoFRenumbering::Cuthill_McKee const bool, const std::vector&); +template +std::vector +DoFRenumbering::compute_Cuthill_McKee +(DoFHandler&, + const bool, + const bool, + const std::vector&); + template void DoFRenumbering::component_wise (DoFHandler&, const std::vector&); +template +std::vector +DoFRenumbering::compute_component_wise +(DoFHandler&, + const std::vector&); + template void DoFRenumbering::component_wise (MGDoFHandler&, @@ -985,12 +1086,26 @@ void DoFRenumbering::component_wise template -void DoFRenumbering::cell_wise_dg +void +DoFRenumbering::cell_wise_dg (DoFHandler&, const std::vector::cell_iterator>&); template -void DoFRenumbering::downstream_dg +std::vector +DoFRenumbering::compute_cell_wise_dg +(DoFHandler&, + const std::vector::cell_iterator>&); + +template +void +DoFRenumbering::downstream_dg +(DoFHandler&, + const Point&); + +template +std::vector +DoFRenumbering::compute_downstream_dg (DoFHandler&, const Point&); @@ -999,10 +1114,21 @@ void DoFRenumbering::sort_selected_dofs_back (DoFHandler &, const std::vector &); +template +std::vector +DoFRenumbering::compute_sort_selected_dofs_back +(DoFHandler &, + const std::vector &); + template void DoFRenumbering::random (DoFHandler &); +template +std::vector +DoFRenumbering::compute_random +(DoFHandler &); + #ifdef ENABLE_MULTIGRID template void DoFRenumbering::Cuthill_McKee @@ -1010,6 +1136,7 @@ void DoFRenumbering::Cuthill_McKee const unsigned int, const bool, const std::vector&); + template void DoFRenumbering::cell_wise_dg (MGDoFHandler&, diff --git a/deal.II/doc/news/2002/c-3-4.html b/deal.II/doc/news/2002/c-3-4.html index 0798c328ad..74737544d6 100644 --- a/deal.II/doc/news/2002/c-3-4.html +++ b/deal.II/doc/news/2002/c-3-4.html @@ -635,6 +635,18 @@ contributor's names are abbreviated by WB (Wolfgang Bangerth), GK

deal.II

    +
  1. + New: For each of the renumbering functions in the DoFRenumbering class there is now an + additional compute_* + function. These new functions compute and return the + corresponding renumbering vectors but do not perform the actual + renumbering on the DoFHandler + object. The behaviour of the old functions is not changed. +
    + (RH 2003/02/03) +

    +
  2. New: The ThreadMutex classes now have a member class ScopedLock that implements the