From 28be0592992064559cf89dbb9e951910524d78dd Mon Sep 17 00:00:00 2001 From: guido Date: Tue, 22 Feb 2000 16:37:38 +0000 Subject: [PATCH] parameters exchanged, flux_pattern in MGDoFTools git-svn-id: https://svn.dealii.org/trunk@2463 0785d39b-7218-0410-832d-ea1e28bc413d --- .../deal.II/include/multigrid/mg_dof_tools.h | 42 ++++++---- .../deal.II/source/dofs/dof_renumbering.cc | 2 +- .../deal.II/source/multigrid/mg_dof_tools.cc | 84 +++++++++++++++---- deal.II/deal.II/source/numerics/data_out.cc | 12 ++- 4 files changed, 107 insertions(+), 33 deletions(-) diff --git a/deal.II/deal.II/include/multigrid/mg_dof_tools.h b/deal.II/deal.II/include/multigrid/mg_dof_tools.h index ee586c0aa5..8926c91bb9 100644 --- a/deal.II/deal.II/include/multigrid/mg_dof_tools.h +++ b/deal.II/deal.II/include/multigrid/mg_dof_tools.h @@ -8,6 +8,7 @@ #include #include +//TODO: Consider incorporating these functions in DoFTools /** * This is a collection of functions operating on, and manipulating @@ -18,7 +19,7 @@ * for more information. * * All member functions are static, so there is no need to create an - * object of class #DoFTools#. + * object of class #MGDoFTools#. * * @author Wolfgang Bangerth and others, 1999 */ @@ -28,29 +29,36 @@ class MGDoFTools /** * Write the sparsity structure * of the matrix belonging to the - * specified #level# including - * constrained degrees of freedom - * into the matrix structure. The - * sparsity pattern does not - * include entries introduced by - * the elimination of constrained - * nodes. The sparsity pattern - * is not compressed, since if - * you want to call - * #ConstraintMatrix::condense(1)# - * afterwards, new entries have - * to be added. However, if you - * don't want to call - * #ConstraintMatrix::condense(1)#, + * specified #level#. The sparsity pattern + * is not compressed, so before + * creating the actual matrix * you have to compress the * matrix yourself, using * #SparseMatrixStruct::compress()#. + * + * There is no need to consider + * hanging nodes here, since only + * one level is considered. */ template static void make_sparsity_pattern (const MGDoFHandler &dof_handler, - const unsigned int level, - SparsityPattern &sparsity); + SparsityPattern &sparsity, + const unsigned int level); + + /** + * Make a sparsity pattern including fluxes + * of discontinuous Galerkin methods. + * @see make_sparsity_pattern + * $see DoFTools + */ + template + static void + make_flux_sparsity_pattern (const MGDoFHandler &dof_handler, + SparsityPattern &sparsity, + const unsigned int level); + + }; diff --git a/deal.II/deal.II/source/dofs/dof_renumbering.cc b/deal.II/deal.II/source/dofs/dof_renumbering.cc index 4cb578d796..6394aa15e2 100644 --- a/deal.II/deal.II/source/dofs/dof_renumbering.cc +++ b/deal.II/deal.II/source/dofs/dof_renumbering.cc @@ -217,7 +217,7 @@ void DoFRenumbering::Cuthill_McKee (MGDoFHandler &dof_handler, // make the connection graph SparsityPattern sparsity (dof_handler.n_dofs(level), dof_handler.max_couplings_between_dofs()); - MGDoFTools::make_sparsity_pattern (dof_handler, level, sparsity); + MGDoFTools::make_sparsity_pattern (dof_handler, sparsity, level); const unsigned int n_dofs = sparsity.n_rows(); // store the new dof numbers; invalid_dof_index means diff --git a/deal.II/deal.II/source/multigrid/mg_dof_tools.cc b/deal.II/deal.II/source/multigrid/mg_dof_tools.cc index 3e16b5da04..7cc6d0b601 100644 --- a/deal.II/deal.II/source/multigrid/mg_dof_tools.cc +++ b/deal.II/deal.II/source/multigrid/mg_dof_tools.cc @@ -10,19 +10,21 @@ template -void MGDoFTools::make_sparsity_pattern (const MGDoFHandler &mg_dof_handler, - const unsigned int level, - SparsityPattern &sparsity) +void MGDoFTools::make_sparsity_pattern (const MGDoFHandler &dof, + SparsityPattern &sparsity, + const unsigned int level) { - Assert (sparsity.n_rows() == mg_dof_handler.n_dofs(level), - ExcDimensionMismatch (sparsity.n_rows(), mg_dof_handler.n_dofs(level))); - Assert (sparsity.n_cols() == mg_dof_handler.n_dofs(level), - ExcDimensionMismatch (sparsity.n_cols(), mg_dof_handler.n_dofs(level))); + const unsigned int n_dofs = dof.n_dofs(level); - const unsigned int dofs_per_cell = mg_dof_handler.get_fe().dofs_per_cell; + Assert (sparsity.n_rows() == n_dofs, + ExcDimensionMismatch (sparsity.n_rows(), n_dofs)); + Assert (sparsity.n_cols() == n_dofs, + ExcDimensionMismatch (sparsity.n_cols(), n_dofs)); + + const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell; vector dofs_on_this_cell(dofs_per_cell); - MGDoFHandler::cell_iterator cell = mg_dof_handler.begin(level), - endc = mg_dof_handler.end(level); + MGDoFHandler::cell_iterator cell = dof.begin(level), + endc = dof.end(level); for (; cell!=endc; ++cell) { cell->get_mg_dof_indices (dofs_on_this_cell); @@ -31,13 +33,67 @@ void MGDoFTools::make_sparsity_pattern (const MGDoFHandler &mg_dof_handler, for (unsigned int j=0; j +void +MGDoFTools::make_flux_sparsity_pattern (const MGDoFHandler &dof, + SparsityPattern &sparsity, + const unsigned int level) +{ + const unsigned int n_dofs = dof.n_dofs(level); + + Assert (sparsity.n_rows() == n_dofs, + ExcDimensionMismatch (sparsity.n_rows(), n_dofs)); + Assert (sparsity.n_cols() == n_dofs, + ExcDimensionMismatch (sparsity.n_cols(), n_dofs)); + + const unsigned int dofs_per_cell = dof.get_fe().dofs_per_cell; + vector dofs_on_this_cell(dofs_per_cell); + vector dofs_on_other_cell(dofs_per_cell); + DoFHandler::cell_iterator cell = dof.begin(level), + endc = dof.end(level); + for (; cell!=endc; ++cell) + { + cell->get_dof_indices (dofs_on_this_cell); + // make sparsity pattern for this cell + for (unsigned int i=0; i::faces_per_cell; + ++face) + { + if ( (! cell->at_boundary(face)) && (cell->neighbor_level(face) == level) ) + { + DoFHandler::cell_iterator neighbor = cell->neighbor(face); + neighbor->get_dof_indices (dofs_on_other_cell); + // only add one direction + // The other is taken care of + // by neighbor. + for (unsigned int i=0; i &, - const unsigned int, - SparsityPattern &); + SparsityPattern &, + const unsigned int); +template void MGDoFTools::make_flux_sparsity_pattern (const MGDoFHandler &, + SparsityPattern &, + const unsigned int); diff --git a/deal.II/deal.II/source/numerics/data_out.cc b/deal.II/deal.II/source/numerics/data_out.cc index f521a9e792..c8d83e7f11 100644 --- a/deal.II/deal.II/source/numerics/data_out.cc +++ b/deal.II/deal.II/source/numerics/data_out.cc @@ -85,8 +85,18 @@ template void DataOut_DoFData::add_data_vector (const Vector &vec, const string &name) { +// unsigned int n = dofs->get_fe().n_components (); +// if (n > 1) +// { +// vector names (dofs->get_fe().n_components ()); +// for (unsigned int i=0;i(1,name)); -}; +} -- 2.39.5