]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Properly instantiate a few static variables.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Thu, 11 Oct 2012 23:58:19 +0000 (23:58 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Thu, 11 Oct 2012 23:58:19 +0000 (23:58 +0000)
git-svn-id: https://svn.dealii.org/trunk@27064 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/changes.h
deal.II/source/dofs/dof_accessor.cc
deal.II/source/grid/tria_accessor.cc
tests/deal.II/accessor_01.cc [new file with mode: 0644]
tests/deal.II/accessor_01/cmp/generic [new file with mode: 0644]

index cc9c961a86e0d0995847724eb48bedc8ce881918..31b9d59707bab9212e01fec83d58945bddc42b51 100644 (file)
@@ -87,9 +87,17 @@ DoFHandler, in particular removal of specializations.
 <h3>Specific improvements</h3>
 
 <ol>
-<li> Fixed: handle lucky breakdowns in GMRES/FGMRES.
+<li> 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.
 <br>
-(Bärbel Janssen, Timo Heister, 2012/10/09)
+(Wolfgang Bangerth, Guido Kanschat, 2012/10/11)
+
+<li> Fixed: Handle lucky breakdowns in GMRES/FGMRES.
+<br>
+(B&auml;rbel Janssen, Timo Heister, 2012/10/09)
 
 <li> Fixed: GridTools::find_cells_adjacent_to_vertex got into
 trouble with anisotropically refined meshes. This is now fixed.
index 49e76304e5a37a462bff8485d73e7ed69508a6fd..39f713a2e5ef705bc8a8b281a49552e23c8ebb53 100644 (file)
 
 DEAL_II_NAMESPACE_OPEN
 
+/*------------------------- Static variables: DoFAccessor -----------------------*/
+
+template <int structdim, class DH>
+const unsigned int DoFAccessor<structdim,DH>::dimension;
+
+template <int structdim, class DH>
+const unsigned int DoFAccessor<structdim,DH>::space_dimension;
+
+
 
 /*------------------------- Functions: DoFCellAccessor -----------------------*/
 
index cb713aa258a712775000cdbde6abf9c341176414..8607915042480a91f1e52a299b45b761a7664ae8 100644 (file)
@@ -971,6 +971,16 @@ namespace
 
 
 
+/*------------------------ Static variables: TriaAccessorBase ---------------------------*/
+
+template <int structdim, int dim, int spacedim>
+const unsigned int TriaAccessorBase<structdim, dim, spacedim>::dimension;
+
+template <int structdim, int dim, int spacedim>
+const unsigned int TriaAccessorBase<structdim, dim, spacedim>::space_dimension;
+
+template <int structdim, int dim, int spacedim>
+const unsigned int TriaAccessorBase<structdim, dim, spacedim>::structure_dimension;
 
 
 /*------------------------ Functions: TriaAccessor ---------------------------*/
diff --git a/tests/deal.II/accessor_01.cc b/tests/deal.II/accessor_01.cc
new file mode 100644 (file)
index 0000000..3ab6749
--- /dev/null
@@ -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 <deal.II/dofs/dof_accessor.h>
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/fe/fe_q.h>
+#include <fstream>
+
+
+template <class ACCESSOR>
+LogStream&
+operator << (LogStream& log, const TriaIterator<ACCESSOR>& i)
+{
+  log << ACCESSOR::dimension << ' '
+      << ACCESSOR::structure_dimension << ' '
+      << ACCESSOR::space_dimension << ' ';
+  i.print(log);
+  return log;
+}
+
+
+template <class DH>
+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<int dim>
+void init_tria (Triangulation<dim>& tr)
+{
+  GridGenerator::hyper_cube(tr);
+  tr.refine_global(4-dim);
+}
+
+
+template<int dim>
+void init_dofs (DoFHandler<dim>& dof,
+               const Triangulation<dim>& tr,
+               const FiniteElement<dim>& 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 (file)
index 0000000..c6f2e2e
--- /dev/null
@@ -0,0 +1,3 @@
+
+DEAL::a 2.0
+DEAL::l 2 2 2 2.0

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.