From: Wolfgang Bangerth Date: Thu, 11 Oct 2012 23:58:19 +0000 (+0000) Subject: Properly instantiate a few static variables. X-Git-Tag: v8.0.0~1960 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=fb296a5b26163fbb6583504389da7e4802880c57;p=dealii.git Properly instantiate a few static variables. git-svn-id: https://svn.dealii.org/trunk@27064 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/doc/news/changes.h b/deal.II/doc/news/changes.h index cc9c961a86..31b9d59707 100644 --- a/deal.II/doc/news/changes.h +++ b/deal.II/doc/news/changes.h @@ -87,9 +87,17 @@ DoFHandler, in particular removal of specializations.

Specific improvements

    -
  1. Fixed: handle lucky breakdowns in GMRES/FGMRES. +
  2. Fixed: Several static const member variables of the Accessor +classes were not properly instantiated. This only rarely created +trouble because they are typically only used as template arguments +and the compiler substituted them. However, one would get linker +errors when passing around a reference to them. This is now fixed.
    -(Bärbel Janssen, Timo Heister, 2012/10/09) +(Wolfgang Bangerth, Guido Kanschat, 2012/10/11) + +
  3. Fixed: Handle lucky breakdowns in GMRES/FGMRES. +
    +(Bärbel Janssen, Timo Heister, 2012/10/09)
  4. Fixed: GridTools::find_cells_adjacent_to_vertex got into trouble with anisotropically refined meshes. This is now fixed. diff --git a/deal.II/source/dofs/dof_accessor.cc b/deal.II/source/dofs/dof_accessor.cc index 49e76304e5..39f713a2e5 100644 --- a/deal.II/source/dofs/dof_accessor.cc +++ b/deal.II/source/dofs/dof_accessor.cc @@ -35,6 +35,15 @@ DEAL_II_NAMESPACE_OPEN +/*------------------------- Static variables: DoFAccessor -----------------------*/ + +template +const unsigned int DoFAccessor::dimension; + +template +const unsigned int DoFAccessor::space_dimension; + + /*------------------------- Functions: DoFCellAccessor -----------------------*/ diff --git a/deal.II/source/grid/tria_accessor.cc b/deal.II/source/grid/tria_accessor.cc index cb713aa258..8607915042 100644 --- a/deal.II/source/grid/tria_accessor.cc +++ b/deal.II/source/grid/tria_accessor.cc @@ -971,6 +971,16 @@ namespace +/*------------------------ Static variables: TriaAccessorBase ---------------------------*/ + +template +const unsigned int TriaAccessorBase::dimension; + +template +const unsigned int TriaAccessorBase::space_dimension; + +template +const unsigned int TriaAccessorBase::structure_dimension; /*------------------------ Functions: TriaAccessor ---------------------------*/ diff --git a/tests/deal.II/accessor_01.cc b/tests/deal.II/accessor_01.cc new file mode 100644 index 0000000000..3ab6749f4d --- /dev/null +++ b/tests/deal.II/accessor_01.cc @@ -0,0 +1,84 @@ +//---------------------------------------------------------------------- +// $Id$ +// +// Copyright (C) 2012 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. +// +//---------------------------------------------------------------------- + +// We used to forget to instantiate a few static const member variables. +// Verify that this is now fixed. + + +#include "../tests.h" +#include +#include +#include +#include +#include +#include +#include + + +template +LogStream& +operator << (LogStream& log, const TriaIterator& i) +{ + log << ACCESSOR::dimension << ' ' + << ACCESSOR::structure_dimension << ' ' + << ACCESSOR::space_dimension << ' '; + i.print(log); + return log; +} + + +template +void test_in_dim(const DH& d1, const DH& d2) +{ + typename DH::active_cell_iterator a = d1.begin_active(); + typename DH::cell_iterator l = d1.begin(d1.get_tria().n_levels()-1); + + deallog << "a " << a << std::endl + << "l " << l << std::endl; +} + + +template +void init_tria (Triangulation& tr) +{ + GridGenerator::hyper_cube(tr); + tr.refine_global(4-dim); +} + + +template +void init_dofs (DoFHandler& dof, + const Triangulation& tr, + const FiniteElement& fe) +{ + dof.initialize(tr, fe); +} + + +int main () +{ + initlog(__FILE__, true); + deallog.depth_console(0); + + Triangulation<2> t2; + init_tria (t2); + + FE_Q<2> q21(1); + DoFHandler<2> d21; + init_dofs(d21, t2, q21); + + FE_Q<2> q22(2); + DoFHandler<2> d22; + init_dofs(d22, t2, q22); + + test_in_dim(d21,d22); +} diff --git a/tests/deal.II/accessor_01/cmp/generic b/tests/deal.II/accessor_01/cmp/generic new file mode 100644 index 0000000000..c6f2e2e815 --- /dev/null +++ b/tests/deal.II/accessor_01/cmp/generic @@ -0,0 +1,3 @@ + +DEAL::a 2.0 +DEAL::l 2 2 2 2.0