From: kanschat Date: Thu, 18 Mar 2010 16:16:39 +0000 (+0000) Subject: add BlockInfo to DoFHandler X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=776506f1d98dd995b3bfb6730db6ea98af89476a;p=dealii-svn.git add BlockInfo to DoFHandler git-svn-id: https://svn.dealii.org/trunk@20848 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/deal.II/include/dofs/block_info.h b/deal.II/deal.II/include/dofs/block_info.h index d8bb24f533..26d687c149 100644 --- a/deal.II/deal.II/include/dofs/block_info.h +++ b/deal.II/deal.II/include/dofs/block_info.h @@ -2,7 +2,7 @@ // $Id$ // Version: $Name$ // -// Copyright (C) 2009 by the deal.II authors +// Copyright (C) 2009, 2010 by the deal.II authors // // This file is subject to QPL and may not be distributed // without copyright and license information. Please refer @@ -58,13 +58,17 @@ class BlockInfo : public Subscriptor /** * @brief Fill the object with values * describing level block - * structure of the MGDoFHandler. + * structure of the + * MGDoFHandler. If + * levels_only is false, + * the other initialize() is + * called as well. * * This function will also clear * the local() indices. */ template - void initialize(const MGDoFHandler&); + void initialize(const MGDoFHandler&, bool levels_only = false); /** * @brief Initialize block structure diff --git a/deal.II/deal.II/include/dofs/dof_handler.h b/deal.II/deal.II/include/dofs/dof_handler.h index 21d0199550..9f70cd86f0 100644 --- a/deal.II/deal.II/include/dofs/dof_handler.h +++ b/deal.II/deal.II/include/dofs/dof_handler.h @@ -2,7 +2,7 @@ // $Id$ // Version: $Name$ // -// Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009 by the deal.II authors +// Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009, 2010 by the deal.II authors // // This file is subject to QPL and may not be distributed // without copyright and license information. Please refer @@ -18,6 +18,7 @@ #include #include #include +#include #include #include @@ -250,6 +251,19 @@ class DoFHandler : public Subscriptor virtual void distribute_dofs (const FiniteElement &fe, const unsigned int offset = 0); + /** + * After disatribute_dofs() with + * an FESystem element, the block + * structure of global and level + * vectors is stored in a + * BlockInfo object accessible + * with block_info(). This + * function initializes the local + * block structure on each cell + * in the same object. + */ + void initialize_local_block_info(); + /** * Clear all data of this object and * especially delete the lock this object @@ -918,6 +932,34 @@ class DoFHandler : public Subscriptor unsigned int n_boundary_dofs (const std::set &boundary_indicators) const; + /** + * Access to an object informing + * of the block structure of the + * dof handler. + * + * If an FESystem is used in + * distribute_dofs(), degrees ofd + * freedom naturally split into + * several @ref GlossBlock + * blocks. For each base element + * as many blocks appear as its + * multiplicity. + * + * At the end of + * distribute_dofs(), the number + * of degrees of freedom in each + * block is counted, and stored + * in a BlockInfo object, which + * can be accessed here. In an + * MGDoFHandler, the same is done + * on each level. Additionally, + * the block structure on each + * cell can be generated in this + * object by calling + * initialize_local_block_indices(). + */ + const BlockInfo& block_info() const; + /** * Return a constant reference to * the selected finite element @@ -951,11 +993,16 @@ class DoFHandler : public Subscriptor */ DeclException0 (ExcInvalidTriangulation); /** - * Exception + * No finite element has been + * assigned to this + * DoFHandler. Call the function + * DoFHandler::distribute_dofs() + * before you attemped to do what + * raised this exception. */ DeclException0 (ExcNoFESelected); /** - * Exception + * @todo Replace by ExcInternalError. */ DeclException0 (ExcRenumberingIncomplete); /** @@ -1020,7 +1067,8 @@ class DoFHandler : public Subscriptor * this object as well, though). */ SmartPointer,DoFHandler > selected_fe; - + + BlockInfo block_info_object; private: /** @@ -1137,6 +1185,15 @@ DoFHandler::get_tria () const } +template +inline +const BlockInfo& +DoFHandler::block_info () const +{ + return block_info_object; +} + + #endif // DOXYGEN diff --git a/deal.II/deal.II/source/dofs/block_info.cc b/deal.II/deal.II/source/dofs/block_info.cc index fcbdfef5c3..e3b2991c6f 100644 --- a/deal.II/deal.II/source/dofs/block_info.cc +++ b/deal.II/deal.II/source/dofs/block_info.cc @@ -56,9 +56,10 @@ BlockInfo::initialize_local(const DoFHandler& dof) template void -BlockInfo::initialize(const MGDoFHandler& dof) +BlockInfo::initialize(const MGDoFHandler& dof, bool levels_only) { - initialize((DoFHandler&) dof); + if (!levels_only) + initialize((DoFHandler&) dof); std::vector > sizes (dof.get_tria().n_levels()); for (unsigned int i=0; i& dof) template void BlockInfo::initialize(const DoFHandler&); -template void BlockInfo::initialize(const MGDoFHandler&); +template void BlockInfo::initialize(const MGDoFHandler&, bool); template void BlockInfo::initialize_local(const DoFHandler&); +#if deal_II_dimension==1 || deal_II_dimension==2 +template void BlockInfo::initialize(const DoFHandler&); +template void BlockInfo::initialize(const MGDoFHandler&, bool); +template void BlockInfo::initialize_local(const DoFHandler&); +#endif + DEAL_II_NAMESPACE_CLOSE diff --git a/deal.II/deal.II/source/dofs/dof_handler.cc b/deal.II/deal.II/source/dofs/dof_handler.cc index 36c82dbba3..e776316bbd 100644 --- a/deal.II/deal.II/source/dofs/dof_handler.cc +++ b/deal.II/deal.II/source/dofs/dof_handler.cc @@ -2,7 +2,7 @@ // $Id$ // Version: $Name$ // -// Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009 by the deal.II authors +// Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009, 2010 by the deal.II authors // // This file is subject to QPL and may not be distributed // without copyright and license information. Please refer @@ -1995,9 +1995,17 @@ void DoFHandler:: // finally restore the user flags const_cast &>(*tria).load_user_flags(user_flags); + + block_info_object.initialize(*this); } +template +void DoFHandler::initialize_local_block_info () +{ + block_info_object.initialize_local(*this); +} + template void DoFHandler::clear () diff --git a/deal.II/deal.II/source/dofs/dof_tools.cc b/deal.II/deal.II/source/dofs/dof_tools.cc index fcc25815ae..3366c7b682 100644 --- a/deal.II/deal.II/source/dofs/dof_tools.cc +++ b/deal.II/deal.II/source/dofs/dof_tools.cc @@ -5751,6 +5751,10 @@ DoFTools::count_dofs_per_component ( std::vector&, bool, std::vector); #if deal_II_dimension < 3 +template void DoFTools::extract_level_dofs +(const unsigned int level, const MGDoFHandler&, + const std::vector&, std::vector&, bool); + template void DoFTools::count_dofs_per_component ( @@ -5842,6 +5846,12 @@ DoFTools::map_dofs_to_support_points const DoFHandler&, std::vector >&); +template +void +DoFTools::count_dofs_per_block ( + const DoFHandler&, + std::vector&, std::vector); + #endif template diff --git a/deal.II/deal.II/source/fe/fe_tools.cc b/deal.II/deal.II/source/fe/fe_tools.cc index d9ecc24fa8..8152b48a24 100644 --- a/deal.II/deal.II/source/fe/fe_tools.cc +++ b/deal.II/deal.II/source/fe/fe_tools.cc @@ -2,7 +2,7 @@ // $Id$ // Version: $Name$ // -// Copyright (C) 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009 by the deal.II authors +// Copyright (C) 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009, 2010 by the deal.II authors // // This file is subject to QPL and may not be distributed // without copyright and license information. Please refer @@ -2115,6 +2115,10 @@ void FETools::interpolate #if deal_II_dimension != 3 template +void FETools::compute_block_renumbering ( + const FiniteElement& element, + std::vector&, std::vector&_indices, bool); +template void FETools::interpolate (const DoFHandler &, const Vector &, const DoFHandler &, Vector &); diff --git a/deal.II/deal.II/source/multigrid/mg_dof_handler.cc b/deal.II/deal.II/source/multigrid/mg_dof_handler.cc index 42570597f0..d2819a8947 100644 --- a/deal.II/deal.II/source/multigrid/mg_dof_handler.cc +++ b/deal.II/deal.II/source/multigrid/mg_dof_handler.cc @@ -1992,6 +1992,8 @@ void MGDoFHandler::distribute_dofs (const FiniteElement &>(*(this->tria)) .load_user_flags(user_flags); + + this->block_info_object.initialize(*this, true); } diff --git a/deal.II/deal.II/source/multigrid/mg_tools.cc b/deal.II/deal.II/source/multigrid/mg_tools.cc index a5673cdb2e..2afe759d0c 100644 --- a/deal.II/deal.II/source/multigrid/mg_tools.cc +++ b/deal.II/deal.II/source/multigrid/mg_tools.cc @@ -1024,7 +1024,7 @@ MGTools::count_dofs_per_block ( std::vector >& dofs_per_block, std::vector target_block) { - const FiniteElement& fe = dof_handler.get_fe(); + const FiniteElement& fe = dof_handler.get_fe(); const unsigned int n_blocks = fe.n_blocks(); const unsigned int n_levels = dof_handler.get_tria().n_levels(); @@ -1639,4 +1639,12 @@ MGTools:: extract_inner_interface_dofs (const MGDoFHandler &mg_dof_handler, std::vector > &interface_dofs); +#if deal_II_dimension < 3 + +template void MGTools::count_dofs_per_block ( + const MGDoFHandler&, + std::vector >&, std::vector); + +#endif + DEAL_II_NAMESPACE_CLOSE diff --git a/deal.II/doc/news/changes.h b/deal.II/doc/news/changes.h index d782f6bcff..77558246ec 100644 --- a/deal.II/doc/news/changes.h +++ b/deal.II/doc/news/changes.h @@ -589,6 +589,13 @@ inconvenience this causes.

deal.II

    + +
  1. Improved: DoFHandler now has a BlockInfo object, automatically updated after DoFHandler::distribute_dofs() and accessible by DoFHandler::block_info(). This object can be used to initialize block vectors and obliterates the necessity to count dofs per block (or dofs per component in legacy code) in application programs. +
    +(GK 2010/03/18) +

    +
  2. +
  3. New: The class MGTransferSelect is prepared for use on adaptively refined meshes.