]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
In 1d, VectorTools::interpolate_boundary_values still assumed the old scheme where...
authorbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 14 Feb 2013 20:41:17 +0000 (20:41 +0000)
committerbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 14 Feb 2013 20:41:17 +0000 (20:41 +0000)
git-svn-id: https://svn.dealii.org/trunk@28398 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/changes.h
deal.II/include/deal.II/numerics/vector_tools.templates.h
tests/deal.II/interpolate_boundary_values_02.cc [new file with mode: 0644]
tests/deal.II/interpolate_boundary_values_02/cmp/generic [new file with mode: 0644]

index 29420ce85f8180d4b2f13347b90c457a69219452..f863bdf26baba8bb15e09dbdb37e23ce25589f4d 100644 (file)
@@ -186,6 +186,13 @@ DoFHandler, in particular removal of specializations.
 <h3>Specific improvements</h3>
 
 <ol>
+<li> Fixed: VectorTools::interpolate did not work properly in 1d if
+boundary indicators had been set to anything but the default (i.e.,
+zero at the left and one at the right end of the domain). This was
+a hold-over from the past when these were the only possible values.
+This is now fixed.
+<br>
+(Kevin Dugan, Wolfgang Bangerth, 2013/02/14)
 
 <li> Improved: The iterator class of the deal.II SparseMatrix class
 and SparsityPattern have been revised for performance. Iterating over
index 8a81c86a507bc8a44df620deb273a5b7baff6e9a..ed7f40969d12738808ec932f1fe5836a701cbc96 100644 (file)
@@ -1511,21 +1511,16 @@ namespace VectorTools
       if (function_map.size() == 0)
         return;
 
-      for (typename FunctionMap<spacedim>::type::const_iterator i=function_map.begin();
-           i!=function_map.end(); ++i)
-        Assert (i->first < 2,
-                ExcInvalidBoundaryIndicator());
-
       for (typename DH::active_cell_iterator cell = dof.begin_active();
            cell != dof.end(); ++cell)
         for (unsigned int direction=0;
              direction<GeometryInfo<dim>::faces_per_cell; ++direction)
           if (cell->at_boundary(direction)
               &&
-              (function_map.find(direction) != function_map.end()))
+              (function_map.find(cell->face(direction)->boundary_indicator()) != function_map.end()))
             {
               const Function<DH::space_dimension> &boundary_function
-                = *function_map.find(direction)->second;
+                = *function_map.find(cell->face(direction)->boundary_indicator())->second;
 
               // get the FE corresponding to this
               // cell
diff --git a/tests/deal.II/interpolate_boundary_values_02.cc b/tests/deal.II/interpolate_boundary_values_02.cc
new file mode 100644 (file)
index 0000000..fdb148c
--- /dev/null
@@ -0,0 +1,85 @@
+//----------------------------------------------------------------------
+//    $Id$
+//    Version: $Name$
+//
+//    Copyright (C) 2006, 2007, 2013 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.
+//
+//----------------------------------------------------------------------
+
+
+// VectorTools::interpolate_boundary_values still had bugs in 1d after
+// switching to a scheme where we can assign boundary indicators also in 1d
+
+#include "../tests.h"
+#include <deal.II/base/function_lib.h>
+#include <deal.II/base/logstream.h>
+#include <deal.II/lac/vector.h>
+
+#include <deal.II/grid/tria.h>
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_refinement.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/grid/tria_boundary_lib.h>
+#include <deal.II/dofs/dof_accessor.h>
+#include <deal.II/dofs/dof_tools.h>
+#include <deal.II/numerics/vector_tools.h>
+#include <deal.II/fe/fe_q.h>
+
+#include <fstream>
+#include <iomanip>
+#include <vector>
+
+
+
+template <int dim>
+void test ()
+{
+  deallog << "dim = " << dim << std::endl;
+
+  Triangulation<dim>   triangulation;
+  FE_Q<dim>            fe(1);
+  DoFHandler<dim>      dof_handler(triangulation);
+
+  GridGenerator::hyper_cube (triangulation, 0, 1);
+  triangulation.begin_active()->face(0)->set_boundary_indicator(10);
+  triangulation.begin_active()->face(1)->set_boundary_indicator(20);
+  triangulation.refine_global (1);
+
+  dof_handler.distribute_dofs (fe);
+
+  std::map<unsigned int,double> boundary_values;
+  VectorTools::interpolate_boundary_values (dof_handler,
+                                            10,
+                                            Functions::SquareFunction<dim>(),
+                                            boundary_values);
+  VectorTools::interpolate_boundary_values (dof_handler,
+                                            20,
+                                            Functions::SquareFunction<dim>(),
+                                            boundary_values);
+  deallog << boundary_values.size() << std::endl;
+  for (std::map<unsigned int,double>::const_iterator
+        p = boundary_values.begin();
+       p != boundary_values.end(); ++p)
+    deallog << p->first << ' ' << p->second << std::endl;
+}
+
+
+int main ()
+{
+  std::ofstream logfile("interpolate_boundary_values_02/output");
+  deallog.attach(logfile);
+  deallog.depth_console(0);
+  deallog.threshold_double(1.e-10);
+
+  test<1>();
+  test<2>();
+  test<3>();
+}
+
diff --git a/tests/deal.II/interpolate_boundary_values_02/cmp/generic b/tests/deal.II/interpolate_boundary_values_02/cmp/generic
new file mode 100644 (file)
index 0000000..fc6fb37
--- /dev/null
@@ -0,0 +1,33 @@
+
+DEAL::dim = 1
+DEAL::2
+DEAL::0 0
+DEAL::2 1.00000
+DEAL::dim = 2
+DEAL::6
+DEAL::0 0
+DEAL::2 0.250000
+DEAL::4 1.00000
+DEAL::5 1.25000
+DEAL::6 1.00000
+DEAL::8 2.00000
+DEAL::dim = 3
+DEAL::18
+DEAL::0 0
+DEAL::2 0.250000
+DEAL::4 0.250000
+DEAL::6 0.500000
+DEAL::8 1.00000
+DEAL::9 1.25000
+DEAL::10 1.25000
+DEAL::11 1.50000
+DEAL::12 1.00000
+DEAL::14 1.25000
+DEAL::16 2.00000
+DEAL::17 2.25000
+DEAL::18 1.00000
+DEAL::20 1.25000
+DEAL::22 2.00000
+DEAL::23 2.25000
+DEAL::24 2.00000
+DEAL::26 3.00000

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.