From 8fcd8d79a89924016b923d67f369940560bf9763 Mon Sep 17 00:00:00 2001 From: Wolfgang Bangerth Date: Wed, 26 Nov 2014 23:33:44 -0600 Subject: [PATCH] Fix an issue reported by David Wells. The problem here is that class DoFHandler only uses pointers to the DoF faces and levels objects in the .h file, so we thought that we don't have to #include the respective header files where these classes are declared. But we also have a serialization function that uses them and if you call it, the compiler will complain about undeclared classes being used. --- doc/news/changes.h | 6 +++ include/deal.II/dofs/dof_faces.h | 5 +-- include/deal.II/dofs/dof_handler.h | 5 +-- include/deal.II/dofs/dof_levels.h | 5 +-- tests/bits/serialize_dof_handler.cc | 58 +++++++++++++++++++++++++ tests/bits/serialize_dof_handler.output | 2 + 6 files changed, 70 insertions(+), 11 deletions(-) create mode 100644 tests/bits/serialize_dof_handler.cc create mode 100644 tests/bits/serialize_dof_handler.output diff --git a/doc/news/changes.h b/doc/news/changes.h index f5d1fa511b..3c1ec70973 100644 --- a/doc/news/changes.h +++ b/doc/news/changes.h @@ -333,6 +333,12 @@ inconvenience this causes.

Specific improvements

    +
  1. Fixed: Serializing an object of type DoFHandler did not work without + including additional header files. +
    + (David Wells, Wolfgang Bangerth, 2014/11/26) +
  2. +
  3. New: The class FEEvaluation with its fast tensor evaluation routines can now be initialized from a mapping, a finite element, a quadrature, and update flags on the fly similar to FEValues. This provides an alternative diff --git a/include/deal.II/dofs/dof_faces.h b/include/deal.II/dofs/dof_faces.h index a228b71056..089412c651 100644 --- a/include/deal.II/dofs/dof_faces.h +++ b/include/deal.II/dofs/dof_faces.h @@ -1,6 +1,6 @@ // --------------------------------------------------------------------- // -// Copyright (C) 2006 - 2013 by the deal.II authors +// Copyright (C) 2006 - 2014 by the deal.II authors // // This file is part of the deal.II library. // @@ -24,9 +24,6 @@ DEAL_II_NAMESPACE_OPEN -template class DoFHandler; -template class MGDoFHandler; - namespace internal { /** diff --git a/include/deal.II/dofs/dof_handler.h b/include/deal.II/dofs/dof_handler.h index c38a675aa7..b3816e5f18 100644 --- a/include/deal.II/dofs/dof_handler.h +++ b/include/deal.II/dofs/dof_handler.h @@ -27,6 +27,8 @@ #include #include #include +#include +#include #include #include @@ -41,9 +43,6 @@ namespace internal { namespace DoFHandler { - template class DoFLevel; - template class DoFFaces; - struct Implementation; namespace Policy diff --git a/include/deal.II/dofs/dof_levels.h b/include/deal.II/dofs/dof_levels.h index db32eeb9d5..eaf9ea1fc9 100644 --- a/include/deal.II/dofs/dof_levels.h +++ b/include/deal.II/dofs/dof_levels.h @@ -1,6 +1,6 @@ // --------------------------------------------------------------------- // -// Copyright (C) 1998 - 2013 by the deal.II authors +// Copyright (C) 1998 - 2014 by the deal.II authors // // This file is part of the deal.II library. // @@ -26,9 +26,6 @@ DEAL_II_NAMESPACE_OPEN -template class DoFHandler; -template class MGDoFHandler; - namespace internal { diff --git a/tests/bits/serialize_dof_handler.cc b/tests/bits/serialize_dof_handler.cc new file mode 100644 index 0000000000..33b992fd4d --- /dev/null +++ b/tests/bits/serialize_dof_handler.cc @@ -0,0 +1,58 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2014 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + + + +// Test an issue reported by David Wells: serializing a DoFHandler object did +// not work right out of the box without manually including additional header +// files + +#include "../tests.h" +#include +#include +#include +#include +#include + +#include + +#include +#include +#include +#include + + +int main () +{ + std::ofstream logfile("output"); + deallog.attach(logfile); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + Triangulation<2> triangulation; + GridGenerator::hyper_cube (triangulation); + triangulation.refine_global (4); + DoFHandler<2> dof_handler (triangulation); + FE_Q<2> finite_element (1); + dof_handler.distribute_dofs (finite_element); + + std::ostringstream out_stream; + boost::archive::text_oarchive archive(out_stream); + + archive << dof_handler; + dof_handler.clear (); + + deallog << "OK" << std::endl; +} diff --git a/tests/bits/serialize_dof_handler.output b/tests/bits/serialize_dof_handler.output new file mode 100644 index 0000000000..0fd8fc12f0 --- /dev/null +++ b/tests/bits/serialize_dof_handler.output @@ -0,0 +1,2 @@ + +DEAL::OK -- 2.39.5