std::vector<unsigned int> target_component)
{
const FiniteElement<dim,spacedim>& fe = dof_handler.get_fe();
- const unsigned int n_components = fe.n_components();
- dofs_per_component.resize (n_components);
- std::fill (dofs_per_component.begin(), dofs_per_component.end(), 0U);
+ std::fill (dofs_per_component.begin(), dofs_per_component.end(), 0U);
+
// If the empty vector was given as
// default argument, set up this
// vector as identity.
if (target_component.size()==0)
{
- target_component.resize(n_components);
- for (unsigned int i=0;i<n_components;++i)
+ target_component.resize(fe.n_components());
+ for (unsigned int i=0; i<fe.n_components(); ++i)
target_component[i] = i;
}
-
- Assert(target_component.size()==n_components,
- ExcDimensionMismatch(target_component.size(),n_components));
+ else
+ Assert (target_component.size()==fe.n_components(),
+ ExcDimensionMismatch(target_component.size(),
+ fe.n_components()));
+
+
+ const unsigned int max_component
+ = *std::max_element (target_component.begin(),
+ target_component.end());
+ const unsigned int n_target_components = max_component + 1;
+ const unsigned int n_components = fe.n_components();
+
+ Assert (dofs_per_component.size() == n_target_components,
+ ExcInternalError());
+
// special case for only one
// component. treat this first
// since it does not require any
std::vector<unsigned int> target_block)
{
const FiniteElement<dim,spacedim>& fe = dof_handler.get_fe();
- const unsigned int n_blocks = fe.n_blocks();
- dofs_per_block.resize (n_blocks);
+
std::fill (dofs_per_block.begin(), dofs_per_block.end(), 0U);
// If the empty vector was given as
// vector as identity.
if (target_block.size()==0)
{
- target_block.resize(n_blocks);
- for (unsigned int i=0;i<n_blocks;++i)
+ target_block.resize(fe.n_blocks());
+ for (unsigned int i=0; i<fe.n_blocks(); ++i)
target_block[i] = i;
}
+ else
+ Assert (target_block.size()==fe.n_blocks(),
+ ExcDimensionMismatch(target_block.size(),
+ fe.n_blocks()));
+
+
+
+ const unsigned int max_block
+ = *std::max_element (target_block.begin(),
+ target_block.end());
+ const unsigned int n_target_blocks = max_block + 1;
+ const unsigned int n_blocks = fe.n_blocks();
- Assert(target_block.size()==n_blocks,
- ExcDimensionMismatch(target_block.size(),n_blocks));
+ Assert (dofs_per_block.size() == n_target_blocks,
+ ExcInternalError());
// special case for only one
// block. treat this first
<h3>deal.II</h3>
<ol>
+ <li>
+ <p>
+ Changed: The DoFTools::count_dofs_per_component and
+ DoFTools::count_dofs_per_block erroneously resized the output argument
+ to the number of components or blocks in the finite element, respectively,
+ even if the target component/block list given as an additional argument
+ needed a different number of output elements. This is now fixed.
+ <br>
+ (WB 2009/02/05)
+ </p>
+
<li>
<p>
New: The GridGenerator::half_hyper_shell function now also exists in 3d.
############################################################
# Makefile,v 1.15 2002/06/13 12:51:13 hartmann Exp
-# Copyright (C) 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008 by the deal.II authors
+# Copyright (C) 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009 by the deal.II authors
############################################################
############################################################
joa_* \
neighboring_cells_at_two_faces \
solution_transfer \
- make_boundary_constraints_*
+ make_boundary_constraints_* \
+ count_dofs_per_block*
# tests for the hp branch:
# fe_collection_*
--- /dev/null
+//---------------------------- count_dofs_per_block_01.cc ---------------------------
+// $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.
+//
+//---------------------------- count_dofs_per_block_01.cc ---------------------------
+
+// the DoFTools::count_dofs_per_{component,block} functions resized the output
+// array to fe.n_components or fe.n_blocks even if a grouping argument was
+// given. this would appear wrong and can lead to all sorts of interesting
+// behavior if not caught by an assertion early enough
+
+
+#include "../tests.h"
+#include <base/logstream.h>
+#include <grid/tria.h>
+#include <grid/tria_iterator.h>
+#include <grid/grid_generator.h>
+#include <dofs/dof_accessor.h>
+#include <dofs/dof_handler.h>
+#include <dofs/dof_tools.h>
+#include <fe/fe_q.h>
+#include <fe/fe_dgq.h>
+#include <fe/fe_dgp.h>
+#include <fe/fe_nedelec.h>
+#include <fe/fe_raviart_thomas.h>
+#include <fe/fe_system.h>
+
+#include <fstream>
+#include <iomanip>
+#include <iomanip>
+#include <string>
+
+
+std::string output_file_name = "count_dofs_per_block_01/output";
+
+
+
+void print (const std::vector<unsigned int> &v)
+{
+ deallog << v.size();
+ for (unsigned int i=0; i<v.size(); ++i)
+ deallog << ' ' << v[i];
+ deallog << std::endl;
+}
+
+
+
+template <int dim>
+void
+check ()
+{
+ // create tria and dofhandler
+ // objects. set different boundary
+ // and sub-domain ids
+ Triangulation<dim> tria;
+ GridGenerator::hyper_cube(tria, 0., 1.);
+ tria.refine_global (1);
+ for (int i=0; i<2; ++i)
+ {
+ tria.begin_active()->set_refine_flag();
+ tria.execute_coarsening_and_refinement ();
+ }
+
+ // Taylor-Hood element
+ FESystem<dim> fe (FE_Q<dim>(1), dim, FE_DGQ<dim>(0), 1);
+ DoFHandler<dim> dof_handler (tria);
+ dof_handler.distribute_dofs (fe);
+
+ // no grouping
+ {
+ std::vector<unsigned int> dpc(dim+1);
+ DoFTools::count_dofs_per_component (dof_handler, dpc);
+ print (dpc);
+ }
+
+ {
+ std::vector<unsigned int> dpc(dim+1);
+ DoFTools::count_dofs_per_block (dof_handler, dpc);
+ print (dpc);
+ }
+
+
+ // grouping into less groups than
+ // components
+ {
+ std::vector<unsigned int> group(dim+1, 0);
+ group[dim] = 1;
+ std::vector<unsigned int> dpc(2);
+ DoFTools::count_dofs_per_component (dof_handler, dpc, false, group);
+ Assert (dpc.size() == 2, ExcInternalError());
+ print (dpc);
+ }
+
+ {
+ std::vector<unsigned int> group(dim+1, 0);
+ group[dim] = 1;
+ std::vector<unsigned int> dpc(2);
+ DoFTools::count_dofs_per_block (dof_handler, dpc, group);
+ Assert (dpc.size() == 2, ExcInternalError());
+ print (dpc);
+ }
+
+ // grouping into more groups than
+ // components
+ {
+ std::vector<unsigned int> group(dim+1, 2*dim);
+ group[dim] = 0;
+ std::vector<unsigned int> dpc(2*dim+1);
+ DoFTools::count_dofs_per_component (dof_handler, dpc, false, group);
+ Assert (dpc.size() == 2*dim+1, ExcInternalError());
+ print (dpc);
+ }
+
+ {
+ std::vector<unsigned int> group(dim+1, 2*dim);
+ group[dim] = 0;
+ std::vector<unsigned int> dpc(2*dim+1);
+ DoFTools::count_dofs_per_block (dof_handler, dpc, group);
+ Assert (dpc.size() == 2*dim+1, ExcInternalError());
+ print (dpc);
+ }
+}
+
+
+
+int main()
+{
+ std::ofstream logfile(output_file_name.c_str());
+ logfile << std::setprecision (2);
+ deallog << std::setprecision (2);
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+ deallog.threshold_double(1.e-10);
+
+ check<1> ();
+ check<2> ();
+ check<3> ();
+}
--- /dev/null
+
+DEAL::2 5 4
+DEAL::2 5 4
+DEAL::2 5 4
+DEAL::2 5 4
+DEAL::3 4 0 5
+DEAL::3 4 0 5
+DEAL::3 18 18 10
+DEAL::3 18 18 10
+DEAL::2 36 10
+DEAL::2 36 10
+DEAL::5 10 0 0 0 36
+DEAL::5 10 0 0 0 36
+DEAL::4 60 60 60 22
+DEAL::4 60 60 60 22
+DEAL::2 180 22
+DEAL::2 180 22
+DEAL::7 22 0 0 0 0 0 180
+DEAL::7 22 0 0 0 0 0 180