From: kanschat Date: Sun, 13 Sep 2009 17:37:30 +0000 (+0000) Subject: new class BlockInfo X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=efccc2ccdc9a859d2572568d250bef7bde497c2a;p=dealii-svn.git new class BlockInfo git-svn-id: https://svn.dealii.org/trunk@19449 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 new file mode 100644 index 0000000000..93e1a3c861 --- /dev/null +++ b/deal.II/deal.II/include/dofs/block_info.h @@ -0,0 +1,189 @@ +//--------------------------------------------------------------------------- +// $Id$ +// Version: $Name$ +// +// Copyright (C) 2009 by the deal.II authors +// +// This file is subject to QPL and may not be distributed +// without copyright and license information. Please refer +// to the file deal.II/doc/license.html for the text and +// further information on this license. +// +//--------------------------------------------------------------------------- +#ifndef __deal2__block_info_h +#define __deal2__block_info_h + +#include +#include + +DEAL_II_NAMESPACE_OPEN + +// Forward declarations + +template class DoFHandler; +template class MGDoFHandler; +namespace hp +{ + template class DoFHandler; +} + + +/** + * @brief A small class collecting the different BlockIndices involved in + * global, multilevel and local computations. + * + * In addition to grouping all the block indices into a single + * instance, this class also provides functions to fill them + * automatically. + * + * @todo Extend the functions local() and renumber() to the concept to + * hpDoFHandler. + * + * @author Guido Kanschat, 2009 + */ +class BlockInfo : public Subscriptor +{ + public: + /** + * @brief Fill the object with values + * describing block structure + * of the DoFHandler. + * + * This function will also clear + * the local() indices. + */ + template + void initialize(const DoFHandler&); + + /** + * @brief Fill the object with values + * describing level block + * structure of the MGDoFHandler. + * + * This function will also clear + * the local() indices. + */ + template + void initialize(const MGDoFHandler&); + + /** + * @brief Initialize block structure + * on cells and compute + * renumbering between cell + * dofs and block cell dofs. + */ + template + void initialize_local(const DoFHandler&); + + /** + * Access the BlockIndices + * structure of the global + * system. + */ + const BlockIndices& global() const; + + /** + * Access BlockIndices for the + * local system on a cell. + */ + const BlockIndices& local() const; + + /** + * Access the BlockIndices + * structure of a level in the + * multilevel hierarchy. + */ + const BlockIndices& level(unsigned int level) const; + + /** + * Return the index after local + * renumbering. + * + * The input of this function is + * an index between zero and the + * number of dofs per cell, + * numbered in local block + * ordering, that is first all + * indices of the first system + * block, then all of the second + * block and so forth. The + * function then outpus the index + * in the standard local + * numbering of DoFAccessor. + */ + unsigned int renumber (unsigned int i) const; + private: + /** + * @brief The block structure + * of the global system. + */ + BlockIndices bi_global; + /** + * @brief The multilevel block structure. + */ + std::vector levels; + + /** + * @brief The block structure + * of the cell systems. + */ + BlockIndices bi_local; + + /** + * The base element associated + * with each block. + */ + std::vector base_element; + + /** + * A vector containing the + * renumbering from the + * standard order of degrees of + * freedom on a cell to a + * component wise + * ordering. Filled by + * initialize(). + */ + std::vector local_renumbering; +}; + + + +//----------------------------------------------------------------------// + +inline +const BlockIndices& +BlockInfo::global() const +{ + return bi_global; +} + + +inline +const BlockIndices& +BlockInfo::local() const +{ + return bi_local; +} + + +inline +const BlockIndices& +BlockInfo::level(unsigned int l) const +{ + AssertIndexRange(l, levels.size()); + return levels[l]; +} + + +inline +unsigned int BlockInfo::renumber(unsigned int i) const +{ + AssertIndexRange(i, local_renumbering.size()); + return local_renumbering[i]; +} + + +DEAL_II_NAMESPACE_CLOSE + +#endif diff --git a/deal.II/deal.II/source/dofs/block_info.cc b/deal.II/deal.II/source/dofs/block_info.cc new file mode 100644 index 0000000000..f6b27d9798 --- /dev/null +++ b/deal.II/deal.II/source/dofs/block_info.cc @@ -0,0 +1,77 @@ +//--------------------------------------------------------------------------- +// $Id$ +// Version: $Name$ +// +// Copyright (C) 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009 by the deal.II authors +// +// This file is subject to QPL and may not be distributed +// without copyright and license information. Please refer +// to the file deal.II/doc/license.html for the text and +// further information on this license. +// +//--------------------------------------------------------------------------- + + +#include +#include +#include +#include +#include +#include +#include + +DEAL_II_NAMESPACE_OPEN + + +template +void +BlockInfo::initialize(const DoFHandler& dof) +{ + const FiniteElement& fe = dof.get_fe(); + std::vector sizes(fe.n_blocks()); + DoFTools::count_dofs_per_block(dof, sizes); + bi_global.reinit(sizes); +} + + +template +void +BlockInfo::initialize_local(const DoFHandler& dof) +{ + const FiniteElement& fe = dof.get_fe(); + std::vector sizes(fe.n_blocks()); + + base_element.resize(fe.n_blocks()); + + for (unsigned int i=0;i +void +BlockInfo::initialize(const MGDoFHandler& dof) +{ + initialize((DoFHandler&) dof); + + std::vector > sizes; + MGTools::count_dofs_per_block(dof, sizes); + + levels.resize(sizes.size()); + for (unsigned int i=0;i&); +template void BlockInfo::initialize(const MGDoFHandler&); +template void BlockInfo::initialize_local(const DoFHandler&); + + +DEAL_II_NAMESPACE_CLOSE diff --git a/deal.II/doc/news/changes.h b/deal.II/doc/news/changes.h index 280beefc68..cdf57d8249 100644 --- a/deal.II/doc/news/changes.h +++ b/deal.II/doc/news/changes.h @@ -420,6 +420,19 @@ inconvenience this causes.
    +
  1. New: The class BlockInfo stores and computes all BlockIndices objects related to + DoFHandler and MGDoFHandler based on FESystem. +
    + (GK 2009/09/13) +

    +
  2. + +
  3. Improved: MGDoFHandler now also has a typedef for MGDoFHandler::Container +
    + (GK 2009/09/13) +

    +
  4. +
  5. FETools::compute_block_renumbering() can nor return block sizes instead of start indices.