]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix FEPointEvaluation for different n_q_points with same MappingInfo for different... 14989/head
authorMaximilian Bergbauer <maximilian.bergbauer@tum.de>
Wed, 29 Mar 2023 08:37:55 +0000 (10:37 +0200)
committerMaximilian Bergbauer <maximilian.bergbauer@tum.de>
Wed, 29 Mar 2023 10:47:37 +0000 (12:47 +0200)
include/deal.II/matrix_free/fe_point_evaluation.h
tests/matrix_free/point_evaluation_19.cc [new file with mode: 0644]
tests/matrix_free/point_evaluation_19.output [new file with mode: 0644]

index 724722c4daf63beac53ec3614734d61ae1acea79..03e7592746a14f69bd2943d135ffa8250584de2f 100644 (file)
@@ -977,9 +977,6 @@ FEPointEvaluation<n_components, dim, spacedim, Number>::evaluate(
   const ArrayView<const Number> &         solution_values,
   const EvaluationFlags::EvaluationFlags &evaluation_flag)
 {
-  if (!is_reinitialized)
-    reinit(numbers::invalid_unsigned_int);
-
   if (n_q_points == 0)
     return;
 
@@ -989,6 +986,14 @@ FEPointEvaluation<n_components, dim, spacedim, Number>::evaluate(
       fast_path)
     {
       // fast path with tensor product evaluation
+      const auto unit_points =
+        mapping_info->get_unit_points(current_cell_index, current_face_number);
+
+      // we need to call reinit() here if we reuse the same MappingInfo object
+      // for several FEPointEvaluation objects to resize the data fields
+      if (!is_reinitialized || n_q_points != unit_points.size())
+        reinit(numbers::invalid_unsigned_int);
+
       if (solution_renumbered.size() != dofs_per_component)
         solution_renumbered.resize(dofs_per_component);
       for (unsigned int comp = 0; comp < n_components; ++comp)
@@ -1006,9 +1011,6 @@ FEPointEvaluation<n_components, dim, spacedim, Number>::evaluate(
       unit_gradients.resize(n_q_points,
                             numbers::signaling_nan<gradient_type>());
 
-      const auto unit_points =
-        mapping_info->get_unit_points(current_cell_index, current_face_number);
-
       const std::size_t n_points = unit_points.size();
       const std::size_t n_lanes  = VectorizedArray<Number>::size();
       for (unsigned int i = 0; i < n_points; i += n_lanes)
@@ -1120,9 +1122,6 @@ FEPointEvaluation<n_components, dim, spacedim, Number>::integrate(
   const ArrayView<Number> &               solution_values,
   const EvaluationFlags::EvaluationFlags &integration_flags)
 {
-  if (!is_reinitialized)
-    reinit(numbers::invalid_unsigned_int);
-
   if (n_q_points == 0) // no evaluation points provided
     {
       std::fill(solution_values.begin(), solution_values.end(), 0.0);
@@ -1135,6 +1134,13 @@ FEPointEvaluation<n_components, dim, spacedim, Number>::integrate(
       fast_path)
     {
       // fast path with tensor product integration
+      const auto unit_points =
+        mapping_info->get_unit_points(current_cell_index, current_face_number);
+
+      // we need to call reinit() here if we reuse the same MappingInfo object
+      // for several FEPointEvaluation objects to resize the data fields
+      if (!is_reinitialized || n_q_points != unit_points.size())
+        reinit(numbers::invalid_unsigned_int);
 
       if (integration_flags & EvaluationFlags::values)
         AssertIndexRange(n_q_points, values.size() + 1);
@@ -1150,9 +1156,6 @@ FEPointEvaluation<n_components, dim, spacedim, Number>::integrate(
           n_components,
           VectorizedArray<Number>>::value_type());
 
-      const auto unit_points =
-        mapping_info->get_unit_points(current_cell_index, current_face_number);
-
       const std::size_t n_points = unit_points.size();
       const std::size_t n_lanes  = VectorizedArray<Number>::size();
       for (unsigned int i = 0; i < n_points; i += n_lanes)
diff --git a/tests/matrix_free/point_evaluation_19.cc b/tests/matrix_free/point_evaluation_19.cc
new file mode 100644 (file)
index 0000000..0958782
--- /dev/null
@@ -0,0 +1,146 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2021 - 2022 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.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+// Check FEPointEvaluation::evaluate()/integrate() with precomputed mapping
+// with different number of unit points provided my NonMatching::MappingInfo.
+
+#include <deal.II/base/function_lib.h>
+
+#include <deal.II/dofs/dof_handler.h>
+
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_system.h>
+#include <deal.II/fe/mapping_cartesian.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include <deal.II/lac/vector.h>
+
+#include <deal.II/matrix_free/fe_point_evaluation.h>
+
+#include <deal.II/numerics/vector_tools.h>
+
+#include <iostream>
+
+#include "../tests.h"
+
+
+
+template <int dim>
+class MyFunction : public Function<dim>
+{
+public:
+  MyFunction()
+    : Function<dim>(dim)
+  {}
+
+  double
+  value(const Point<dim> &p, const unsigned int component) const override
+  {
+    return component + p[component];
+  }
+};
+
+
+
+template <int dim>
+void
+test()
+{
+  using namespace dealii;
+  Triangulation<dim> tria;
+
+  GridGenerator::subdivided_hyper_cube(tria, 2, 0, 1);
+
+  MappingCartesian<dim> mapping;
+  deallog << "Cartesian linear mapping" << std::endl;
+
+  std::vector<Point<dim>> unit_points;
+  for (unsigned int i = 0; i < 13; ++i)
+    {
+      Point<dim> p;
+      for (unsigned int d = 0; d < dim; ++d)
+        p[d] = static_cast<double>(i) / 17. + 0.015625 * d;
+      unit_points.push_back(p);
+    }
+
+  std::vector<Point<dim>> unit_points2;
+  for (unsigned int i = 0; i < 15; ++i)
+    {
+      Point<dim> p;
+      for (unsigned int d = 0; d < dim; ++d)
+        p[d] = static_cast<double>(i) / 17. + 0.015625 * d;
+      unit_points2.push_back(p);
+    }
+
+  FESystem<dim>   fe(FE_Q<dim>(2), dim);
+  DoFHandler<dim> dof_handler(tria);
+  dof_handler.distribute_dofs(fe);
+  Vector<double> vector(dof_handler.n_dofs());
+
+  // mapping info object for precomputed mapping
+  NonMatching::MappingInfo<dim, dim> mapping_info(mapping,
+                                                  update_values |
+                                                    update_gradients);
+  // First evaluator to use precomputed mapping
+  FEPointEvaluation<dim, dim> evaluator1(mapping_info, fe);
+
+  VectorTools::interpolate(mapping, dof_handler, MyFunction<dim>(), vector);
+
+  std::vector<double> solution_values(fe.dofs_per_cell);
+
+  for (const auto &cell : dof_handler.active_cell_iterators())
+    {
+      cell->get_dof_values(vector,
+                           solution_values.begin(),
+                           solution_values.end());
+
+      mapping_info.reinit(cell, unit_points);
+
+      evaluator1.evaluate(solution_values,
+                          EvaluationFlags::values | EvaluationFlags::gradients);
+
+      deallog << "Cell with center " << cell->center(true) << std::endl;
+      for (unsigned int i = 0; i < unit_points.size(); ++i)
+        deallog << mapping.transform_unit_to_real_cell(cell, unit_points[i])
+                << ": " << evaluator1.get_value(i) << std::endl;
+      deallog << std::endl;
+
+      mapping_info.reinit(cell, unit_points2);
+
+      evaluator1.evaluate(solution_values,
+                          EvaluationFlags::values | EvaluationFlags::gradients);
+
+      deallog << "Cell with center " << cell->center(true) << std::endl;
+      for (unsigned int i = 0; i < unit_points2.size(); ++i)
+        deallog << mapping.transform_unit_to_real_cell(cell, unit_points2[i])
+                << ": " << evaluator1.get_value(i) << std::endl;
+      deallog << std::endl;
+    }
+}
+
+
+
+int
+main()
+{
+  initlog();
+  deallog << std::setprecision(10);
+
+  test<2>();
+  test<3>();
+}
diff --git a/tests/matrix_free/point_evaluation_19.output b/tests/matrix_free/point_evaluation_19.output
new file mode 100644 (file)
index 0000000..94a5e3d
--- /dev/null
@@ -0,0 +1,387 @@
+
+DEAL::Cartesian linear mapping
+DEAL::Cell with center 0.2500000000 0.2500000000
+DEAL::0.000000000 0.007812500000: 0.000000000 1.007812500
+DEAL::0.02941176471 0.03722426471: 0.02941176471 1.037224265
+DEAL::0.05882352941 0.06663602941: 0.05882352941 1.066636029
+DEAL::0.08823529412 0.09604779412: 0.08823529412 1.096047794
+DEAL::0.1176470588 0.1254595588: 0.1176470588 1.125459559
+DEAL::0.1470588235 0.1548713235: 0.1470588235 1.154871324
+DEAL::0.1764705882 0.1842830882: 0.1764705882 1.184283088
+DEAL::0.2058823529 0.2136948529: 0.2058823529 1.213694853
+DEAL::0.2352941176 0.2431066176: 0.2352941176 1.243106618
+DEAL::0.2647058824 0.2725183824: 0.2647058824 1.272518382
+DEAL::0.2941176471 0.3019301471: 0.2941176471 1.301930147
+DEAL::0.3235294118 0.3313419118: 0.3235294118 1.331341912
+DEAL::0.3529411765 0.3607536765: 0.3529411765 1.360753676
+DEAL::
+DEAL::Cell with center 0.2500000000 0.2500000000
+DEAL::0.000000000 0.007812500000: 0.000000000 1.007812500
+DEAL::0.02941176471 0.03722426471: 0.02941176471 1.037224265
+DEAL::0.05882352941 0.06663602941: 0.05882352941 1.066636029
+DEAL::0.08823529412 0.09604779412: 0.08823529412 1.096047794
+DEAL::0.1176470588 0.1254595588: 0.1176470588 1.125459559
+DEAL::0.1470588235 0.1548713235: 0.1470588235 1.154871324
+DEAL::0.1764705882 0.1842830882: 0.1764705882 1.184283088
+DEAL::0.2058823529 0.2136948529: 0.2058823529 1.213694853
+DEAL::0.2352941176 0.2431066176: 0.2352941176 1.243106618
+DEAL::0.2647058824 0.2725183824: 0.2647058824 1.272518382
+DEAL::0.2941176471 0.3019301471: 0.2941176471 1.301930147
+DEAL::0.3235294118 0.3313419118: 0.3235294118 1.331341912
+DEAL::0.3529411765 0.3607536765: 0.3529411765 1.360753676
+DEAL::0.3823529412 0.3901654412: 0.3823529412 1.390165441
+DEAL::0.4117647059 0.4195772059: 0.4117647059 1.419577206
+DEAL::
+DEAL::Cell with center 0.7500000000 0.2500000000
+DEAL::0.5000000000 0.007812500000: 0.5000000000 1.007812500
+DEAL::0.5294117647 0.03722426471: 0.5294117647 1.037224265
+DEAL::0.5588235294 0.06663602941: 0.5588235294 1.066636029
+DEAL::0.5882352941 0.09604779412: 0.5882352941 1.096047794
+DEAL::0.6176470588 0.1254595588: 0.6176470588 1.125459559
+DEAL::0.6470588235 0.1548713235: 0.6470588235 1.154871324
+DEAL::0.6764705882 0.1842830882: 0.6764705882 1.184283088
+DEAL::0.7058823529 0.2136948529: 0.7058823529 1.213694853
+DEAL::0.7352941176 0.2431066176: 0.7352941176 1.243106618
+DEAL::0.7647058824 0.2725183824: 0.7647058824 1.272518382
+DEAL::0.7941176471 0.3019301471: 0.7941176471 1.301930147
+DEAL::0.8235294118 0.3313419118: 0.8235294118 1.331341912
+DEAL::0.8529411765 0.3607536765: 0.8529411765 1.360753676
+DEAL::
+DEAL::Cell with center 0.7500000000 0.2500000000
+DEAL::0.5000000000 0.007812500000: 0.5000000000 1.007812500
+DEAL::0.5294117647 0.03722426471: 0.5294117647 1.037224265
+DEAL::0.5588235294 0.06663602941: 0.5588235294 1.066636029
+DEAL::0.5882352941 0.09604779412: 0.5882352941 1.096047794
+DEAL::0.6176470588 0.1254595588: 0.6176470588 1.125459559
+DEAL::0.6470588235 0.1548713235: 0.6470588235 1.154871324
+DEAL::0.6764705882 0.1842830882: 0.6764705882 1.184283088
+DEAL::0.7058823529 0.2136948529: 0.7058823529 1.213694853
+DEAL::0.7352941176 0.2431066176: 0.7352941176 1.243106618
+DEAL::0.7647058824 0.2725183824: 0.7647058824 1.272518382
+DEAL::0.7941176471 0.3019301471: 0.7941176471 1.301930147
+DEAL::0.8235294118 0.3313419118: 0.8235294118 1.331341912
+DEAL::0.8529411765 0.3607536765: 0.8529411765 1.360753676
+DEAL::0.8823529412 0.3901654412: 0.8823529412 1.390165441
+DEAL::0.9117647059 0.4195772059: 0.9117647059 1.419577206
+DEAL::
+DEAL::Cell with center 0.2500000000 0.7500000000
+DEAL::0.000000000 0.5078125000: 0.000000000 1.507812500
+DEAL::0.02941176471 0.5372242647: 0.02941176471 1.537224265
+DEAL::0.05882352941 0.5666360294: 0.05882352941 1.566636029
+DEAL::0.08823529412 0.5960477941: 0.08823529412 1.596047794
+DEAL::0.1176470588 0.6254595588: 0.1176470588 1.625459559
+DEAL::0.1470588235 0.6548713235: 0.1470588235 1.654871324
+DEAL::0.1764705882 0.6842830882: 0.1764705882 1.684283088
+DEAL::0.2058823529 0.7136948529: 0.2058823529 1.713694853
+DEAL::0.2352941176 0.7431066176: 0.2352941176 1.743106618
+DEAL::0.2647058824 0.7725183824: 0.2647058824 1.772518382
+DEAL::0.2941176471 0.8019301471: 0.2941176471 1.801930147
+DEAL::0.3235294118 0.8313419118: 0.3235294118 1.831341912
+DEAL::0.3529411765 0.8607536765: 0.3529411765 1.860753676
+DEAL::
+DEAL::Cell with center 0.2500000000 0.7500000000
+DEAL::0.000000000 0.5078125000: 0.000000000 1.507812500
+DEAL::0.02941176471 0.5372242647: 0.02941176471 1.537224265
+DEAL::0.05882352941 0.5666360294: 0.05882352941 1.566636029
+DEAL::0.08823529412 0.5960477941: 0.08823529412 1.596047794
+DEAL::0.1176470588 0.6254595588: 0.1176470588 1.625459559
+DEAL::0.1470588235 0.6548713235: 0.1470588235 1.654871324
+DEAL::0.1764705882 0.6842830882: 0.1764705882 1.684283088
+DEAL::0.2058823529 0.7136948529: 0.2058823529 1.713694853
+DEAL::0.2352941176 0.7431066176: 0.2352941176 1.743106618
+DEAL::0.2647058824 0.7725183824: 0.2647058824 1.772518382
+DEAL::0.2941176471 0.8019301471: 0.2941176471 1.801930147
+DEAL::0.3235294118 0.8313419118: 0.3235294118 1.831341912
+DEAL::0.3529411765 0.8607536765: 0.3529411765 1.860753676
+DEAL::0.3823529412 0.8901654412: 0.3823529412 1.890165441
+DEAL::0.4117647059 0.9195772059: 0.4117647059 1.919577206
+DEAL::
+DEAL::Cell with center 0.7500000000 0.7500000000
+DEAL::0.5000000000 0.5078125000: 0.5000000000 1.507812500
+DEAL::0.5294117647 0.5372242647: 0.5294117647 1.537224265
+DEAL::0.5588235294 0.5666360294: 0.5588235294 1.566636029
+DEAL::0.5882352941 0.5960477941: 0.5882352941 1.596047794
+DEAL::0.6176470588 0.6254595588: 0.6176470588 1.625459559
+DEAL::0.6470588235 0.6548713235: 0.6470588235 1.654871324
+DEAL::0.6764705882 0.6842830882: 0.6764705882 1.684283088
+DEAL::0.7058823529 0.7136948529: 0.7058823529 1.713694853
+DEAL::0.7352941176 0.7431066176: 0.7352941176 1.743106618
+DEAL::0.7647058824 0.7725183824: 0.7647058824 1.772518382
+DEAL::0.7941176471 0.8019301471: 0.7941176471 1.801930147
+DEAL::0.8235294118 0.8313419118: 0.8235294118 1.831341912
+DEAL::0.8529411765 0.8607536765: 0.8529411765 1.860753676
+DEAL::
+DEAL::Cell with center 0.7500000000 0.7500000000
+DEAL::0.5000000000 0.5078125000: 0.5000000000 1.507812500
+DEAL::0.5294117647 0.5372242647: 0.5294117647 1.537224265
+DEAL::0.5588235294 0.5666360294: 0.5588235294 1.566636029
+DEAL::0.5882352941 0.5960477941: 0.5882352941 1.596047794
+DEAL::0.6176470588 0.6254595588: 0.6176470588 1.625459559
+DEAL::0.6470588235 0.6548713235: 0.6470588235 1.654871324
+DEAL::0.6764705882 0.6842830882: 0.6764705882 1.684283088
+DEAL::0.7058823529 0.7136948529: 0.7058823529 1.713694853
+DEAL::0.7352941176 0.7431066176: 0.7352941176 1.743106618
+DEAL::0.7647058824 0.7725183824: 0.7647058824 1.772518382
+DEAL::0.7941176471 0.8019301471: 0.7941176471 1.801930147
+DEAL::0.8235294118 0.8313419118: 0.8235294118 1.831341912
+DEAL::0.8529411765 0.8607536765: 0.8529411765 1.860753676
+DEAL::0.8823529412 0.8901654412: 0.8823529412 1.890165441
+DEAL::0.9117647059 0.9195772059: 0.9117647059 1.919577206
+DEAL::
+DEAL::Cartesian linear mapping
+DEAL::Cell with center 0.2500000000 0.2500000000 0.2500000000
+DEAL::0.000000000 0.007812500000 0.01562500000: 0.000000000 1.007812500 2.015625000
+DEAL::0.02941176471 0.03722426471 0.04503676471: 0.02941176471 1.037224265 2.045036765
+DEAL::0.05882352941 0.06663602941 0.07444852941: 0.05882352941 1.066636029 2.074448529
+DEAL::0.08823529412 0.09604779412 0.1038602941: 0.08823529412 1.096047794 2.103860294
+DEAL::0.1176470588 0.1254595588 0.1332720588: 0.1176470588 1.125459559 2.133272059
+DEAL::0.1470588235 0.1548713235 0.1626838235: 0.1470588235 1.154871324 2.162683824
+DEAL::0.1764705882 0.1842830882 0.1920955882: 0.1764705882 1.184283088 2.192095588
+DEAL::0.2058823529 0.2136948529 0.2215073529: 0.2058823529 1.213694853 2.221507353
+DEAL::0.2352941176 0.2431066176 0.2509191176: 0.2352941176 1.243106618 2.250919118
+DEAL::0.2647058824 0.2725183824 0.2803308824: 0.2647058824 1.272518382 2.280330882
+DEAL::0.2941176471 0.3019301471 0.3097426471: 0.2941176471 1.301930147 2.309742647
+DEAL::0.3235294118 0.3313419118 0.3391544118: 0.3235294118 1.331341912 2.339154412
+DEAL::0.3529411765 0.3607536765 0.3685661765: 0.3529411765 1.360753676 2.368566176
+DEAL::
+DEAL::Cell with center 0.2500000000 0.2500000000 0.2500000000
+DEAL::0.000000000 0.007812500000 0.01562500000: 0.000000000 1.007812500 2.015625000
+DEAL::0.02941176471 0.03722426471 0.04503676471: 0.02941176471 1.037224265 2.045036765
+DEAL::0.05882352941 0.06663602941 0.07444852941: 0.05882352941 1.066636029 2.074448529
+DEAL::0.08823529412 0.09604779412 0.1038602941: 0.08823529412 1.096047794 2.103860294
+DEAL::0.1176470588 0.1254595588 0.1332720588: 0.1176470588 1.125459559 2.133272059
+DEAL::0.1470588235 0.1548713235 0.1626838235: 0.1470588235 1.154871324 2.162683824
+DEAL::0.1764705882 0.1842830882 0.1920955882: 0.1764705882 1.184283088 2.192095588
+DEAL::0.2058823529 0.2136948529 0.2215073529: 0.2058823529 1.213694853 2.221507353
+DEAL::0.2352941176 0.2431066176 0.2509191176: 0.2352941176 1.243106618 2.250919118
+DEAL::0.2647058824 0.2725183824 0.2803308824: 0.2647058824 1.272518382 2.280330882
+DEAL::0.2941176471 0.3019301471 0.3097426471: 0.2941176471 1.301930147 2.309742647
+DEAL::0.3235294118 0.3313419118 0.3391544118: 0.3235294118 1.331341912 2.339154412
+DEAL::0.3529411765 0.3607536765 0.3685661765: 0.3529411765 1.360753676 2.368566176
+DEAL::0.3823529412 0.3901654412 0.3979779412: 0.3823529412 1.390165441 2.397977941
+DEAL::0.4117647059 0.4195772059 0.4273897059: 0.4117647059 1.419577206 2.427389706
+DEAL::
+DEAL::Cell with center 0.7500000000 0.2500000000 0.2500000000
+DEAL::0.5000000000 0.007812500000 0.01562500000: 0.5000000000 1.007812500 2.015625000
+DEAL::0.5294117647 0.03722426471 0.04503676471: 0.5294117647 1.037224265 2.045036765
+DEAL::0.5588235294 0.06663602941 0.07444852941: 0.5588235294 1.066636029 2.074448529
+DEAL::0.5882352941 0.09604779412 0.1038602941: 0.5882352941 1.096047794 2.103860294
+DEAL::0.6176470588 0.1254595588 0.1332720588: 0.6176470588 1.125459559 2.133272059
+DEAL::0.6470588235 0.1548713235 0.1626838235: 0.6470588235 1.154871324 2.162683824
+DEAL::0.6764705882 0.1842830882 0.1920955882: 0.6764705882 1.184283088 2.192095588
+DEAL::0.7058823529 0.2136948529 0.2215073529: 0.7058823529 1.213694853 2.221507353
+DEAL::0.7352941176 0.2431066176 0.2509191176: 0.7352941176 1.243106618 2.250919118
+DEAL::0.7647058824 0.2725183824 0.2803308824: 0.7647058824 1.272518382 2.280330882
+DEAL::0.7941176471 0.3019301471 0.3097426471: 0.7941176471 1.301930147 2.309742647
+DEAL::0.8235294118 0.3313419118 0.3391544118: 0.8235294118 1.331341912 2.339154412
+DEAL::0.8529411765 0.3607536765 0.3685661765: 0.8529411765 1.360753676 2.368566176
+DEAL::
+DEAL::Cell with center 0.7500000000 0.2500000000 0.2500000000
+DEAL::0.5000000000 0.007812500000 0.01562500000: 0.5000000000 1.007812500 2.015625000
+DEAL::0.5294117647 0.03722426471 0.04503676471: 0.5294117647 1.037224265 2.045036765
+DEAL::0.5588235294 0.06663602941 0.07444852941: 0.5588235294 1.066636029 2.074448529
+DEAL::0.5882352941 0.09604779412 0.1038602941: 0.5882352941 1.096047794 2.103860294
+DEAL::0.6176470588 0.1254595588 0.1332720588: 0.6176470588 1.125459559 2.133272059
+DEAL::0.6470588235 0.1548713235 0.1626838235: 0.6470588235 1.154871324 2.162683824
+DEAL::0.6764705882 0.1842830882 0.1920955882: 0.6764705882 1.184283088 2.192095588
+DEAL::0.7058823529 0.2136948529 0.2215073529: 0.7058823529 1.213694853 2.221507353
+DEAL::0.7352941176 0.2431066176 0.2509191176: 0.7352941176 1.243106618 2.250919118
+DEAL::0.7647058824 0.2725183824 0.2803308824: 0.7647058824 1.272518382 2.280330882
+DEAL::0.7941176471 0.3019301471 0.3097426471: 0.7941176471 1.301930147 2.309742647
+DEAL::0.8235294118 0.3313419118 0.3391544118: 0.8235294118 1.331341912 2.339154412
+DEAL::0.8529411765 0.3607536765 0.3685661765: 0.8529411765 1.360753676 2.368566176
+DEAL::0.8823529412 0.3901654412 0.3979779412: 0.8823529412 1.390165441 2.397977941
+DEAL::0.9117647059 0.4195772059 0.4273897059: 0.9117647059 1.419577206 2.427389706
+DEAL::
+DEAL::Cell with center 0.2500000000 0.7500000000 0.2500000000
+DEAL::0.000000000 0.5078125000 0.01562500000: 0.000000000 1.507812500 2.015625000
+DEAL::0.02941176471 0.5372242647 0.04503676471: 0.02941176471 1.537224265 2.045036765
+DEAL::0.05882352941 0.5666360294 0.07444852941: 0.05882352941 1.566636029 2.074448529
+DEAL::0.08823529412 0.5960477941 0.1038602941: 0.08823529412 1.596047794 2.103860294
+DEAL::0.1176470588 0.6254595588 0.1332720588: 0.1176470588 1.625459559 2.133272059
+DEAL::0.1470588235 0.6548713235 0.1626838235: 0.1470588235 1.654871324 2.162683824
+DEAL::0.1764705882 0.6842830882 0.1920955882: 0.1764705882 1.684283088 2.192095588
+DEAL::0.2058823529 0.7136948529 0.2215073529: 0.2058823529 1.713694853 2.221507353
+DEAL::0.2352941176 0.7431066176 0.2509191176: 0.2352941176 1.743106618 2.250919118
+DEAL::0.2647058824 0.7725183824 0.2803308824: 0.2647058824 1.772518382 2.280330882
+DEAL::0.2941176471 0.8019301471 0.3097426471: 0.2941176471 1.801930147 2.309742647
+DEAL::0.3235294118 0.8313419118 0.3391544118: 0.3235294118 1.831341912 2.339154412
+DEAL::0.3529411765 0.8607536765 0.3685661765: 0.3529411765 1.860753676 2.368566176
+DEAL::
+DEAL::Cell with center 0.2500000000 0.7500000000 0.2500000000
+DEAL::0.000000000 0.5078125000 0.01562500000: 0.000000000 1.507812500 2.015625000
+DEAL::0.02941176471 0.5372242647 0.04503676471: 0.02941176471 1.537224265 2.045036765
+DEAL::0.05882352941 0.5666360294 0.07444852941: 0.05882352941 1.566636029 2.074448529
+DEAL::0.08823529412 0.5960477941 0.1038602941: 0.08823529412 1.596047794 2.103860294
+DEAL::0.1176470588 0.6254595588 0.1332720588: 0.1176470588 1.625459559 2.133272059
+DEAL::0.1470588235 0.6548713235 0.1626838235: 0.1470588235 1.654871324 2.162683824
+DEAL::0.1764705882 0.6842830882 0.1920955882: 0.1764705882 1.684283088 2.192095588
+DEAL::0.2058823529 0.7136948529 0.2215073529: 0.2058823529 1.713694853 2.221507353
+DEAL::0.2352941176 0.7431066176 0.2509191176: 0.2352941176 1.743106618 2.250919118
+DEAL::0.2647058824 0.7725183824 0.2803308824: 0.2647058824 1.772518382 2.280330882
+DEAL::0.2941176471 0.8019301471 0.3097426471: 0.2941176471 1.801930147 2.309742647
+DEAL::0.3235294118 0.8313419118 0.3391544118: 0.3235294118 1.831341912 2.339154412
+DEAL::0.3529411765 0.8607536765 0.3685661765: 0.3529411765 1.860753676 2.368566176
+DEAL::0.3823529412 0.8901654412 0.3979779412: 0.3823529412 1.890165441 2.397977941
+DEAL::0.4117647059 0.9195772059 0.4273897059: 0.4117647059 1.919577206 2.427389706
+DEAL::
+DEAL::Cell with center 0.7500000000 0.7500000000 0.2500000000
+DEAL::0.5000000000 0.5078125000 0.01562500000: 0.5000000000 1.507812500 2.015625000
+DEAL::0.5294117647 0.5372242647 0.04503676471: 0.5294117647 1.537224265 2.045036765
+DEAL::0.5588235294 0.5666360294 0.07444852941: 0.5588235294 1.566636029 2.074448529
+DEAL::0.5882352941 0.5960477941 0.1038602941: 0.5882352941 1.596047794 2.103860294
+DEAL::0.6176470588 0.6254595588 0.1332720588: 0.6176470588 1.625459559 2.133272059
+DEAL::0.6470588235 0.6548713235 0.1626838235: 0.6470588235 1.654871324 2.162683824
+DEAL::0.6764705882 0.6842830882 0.1920955882: 0.6764705882 1.684283088 2.192095588
+DEAL::0.7058823529 0.7136948529 0.2215073529: 0.7058823529 1.713694853 2.221507353
+DEAL::0.7352941176 0.7431066176 0.2509191176: 0.7352941176 1.743106618 2.250919118
+DEAL::0.7647058824 0.7725183824 0.2803308824: 0.7647058824 1.772518382 2.280330882
+DEAL::0.7941176471 0.8019301471 0.3097426471: 0.7941176471 1.801930147 2.309742647
+DEAL::0.8235294118 0.8313419118 0.3391544118: 0.8235294118 1.831341912 2.339154412
+DEAL::0.8529411765 0.8607536765 0.3685661765: 0.8529411765 1.860753676 2.368566176
+DEAL::
+DEAL::Cell with center 0.7500000000 0.7500000000 0.2500000000
+DEAL::0.5000000000 0.5078125000 0.01562500000: 0.5000000000 1.507812500 2.015625000
+DEAL::0.5294117647 0.5372242647 0.04503676471: 0.5294117647 1.537224265 2.045036765
+DEAL::0.5588235294 0.5666360294 0.07444852941: 0.5588235294 1.566636029 2.074448529
+DEAL::0.5882352941 0.5960477941 0.1038602941: 0.5882352941 1.596047794 2.103860294
+DEAL::0.6176470588 0.6254595588 0.1332720588: 0.6176470588 1.625459559 2.133272059
+DEAL::0.6470588235 0.6548713235 0.1626838235: 0.6470588235 1.654871324 2.162683824
+DEAL::0.6764705882 0.6842830882 0.1920955882: 0.6764705882 1.684283088 2.192095588
+DEAL::0.7058823529 0.7136948529 0.2215073529: 0.7058823529 1.713694853 2.221507353
+DEAL::0.7352941176 0.7431066176 0.2509191176: 0.7352941176 1.743106618 2.250919118
+DEAL::0.7647058824 0.7725183824 0.2803308824: 0.7647058824 1.772518382 2.280330882
+DEAL::0.7941176471 0.8019301471 0.3097426471: 0.7941176471 1.801930147 2.309742647
+DEAL::0.8235294118 0.8313419118 0.3391544118: 0.8235294118 1.831341912 2.339154412
+DEAL::0.8529411765 0.8607536765 0.3685661765: 0.8529411765 1.860753676 2.368566176
+DEAL::0.8823529412 0.8901654412 0.3979779412: 0.8823529412 1.890165441 2.397977941
+DEAL::0.9117647059 0.9195772059 0.4273897059: 0.9117647059 1.919577206 2.427389706
+DEAL::
+DEAL::Cell with center 0.2500000000 0.2500000000 0.7500000000
+DEAL::0.000000000 0.007812500000 0.5156250000: 0.000000000 1.007812500 2.515625000
+DEAL::0.02941176471 0.03722426471 0.5450367647: 0.02941176471 1.037224265 2.545036765
+DEAL::0.05882352941 0.06663602941 0.5744485294: 0.05882352941 1.066636029 2.574448529
+DEAL::0.08823529412 0.09604779412 0.6038602941: 0.08823529412 1.096047794 2.603860294
+DEAL::0.1176470588 0.1254595588 0.6332720588: 0.1176470588 1.125459559 2.633272059
+DEAL::0.1470588235 0.1548713235 0.6626838235: 0.1470588235 1.154871324 2.662683824
+DEAL::0.1764705882 0.1842830882 0.6920955882: 0.1764705882 1.184283088 2.692095588
+DEAL::0.2058823529 0.2136948529 0.7215073529: 0.2058823529 1.213694853 2.721507353
+DEAL::0.2352941176 0.2431066176 0.7509191176: 0.2352941176 1.243106618 2.750919118
+DEAL::0.2647058824 0.2725183824 0.7803308824: 0.2647058824 1.272518382 2.780330882
+DEAL::0.2941176471 0.3019301471 0.8097426471: 0.2941176471 1.301930147 2.809742647
+DEAL::0.3235294118 0.3313419118 0.8391544118: 0.3235294118 1.331341912 2.839154412
+DEAL::0.3529411765 0.3607536765 0.8685661765: 0.3529411765 1.360753676 2.868566176
+DEAL::
+DEAL::Cell with center 0.2500000000 0.2500000000 0.7500000000
+DEAL::0.000000000 0.007812500000 0.5156250000: 0.000000000 1.007812500 2.515625000
+DEAL::0.02941176471 0.03722426471 0.5450367647: 0.02941176471 1.037224265 2.545036765
+DEAL::0.05882352941 0.06663602941 0.5744485294: 0.05882352941 1.066636029 2.574448529
+DEAL::0.08823529412 0.09604779412 0.6038602941: 0.08823529412 1.096047794 2.603860294
+DEAL::0.1176470588 0.1254595588 0.6332720588: 0.1176470588 1.125459559 2.633272059
+DEAL::0.1470588235 0.1548713235 0.6626838235: 0.1470588235 1.154871324 2.662683824
+DEAL::0.1764705882 0.1842830882 0.6920955882: 0.1764705882 1.184283088 2.692095588
+DEAL::0.2058823529 0.2136948529 0.7215073529: 0.2058823529 1.213694853 2.721507353
+DEAL::0.2352941176 0.2431066176 0.7509191176: 0.2352941176 1.243106618 2.750919118
+DEAL::0.2647058824 0.2725183824 0.7803308824: 0.2647058824 1.272518382 2.780330882
+DEAL::0.2941176471 0.3019301471 0.8097426471: 0.2941176471 1.301930147 2.809742647
+DEAL::0.3235294118 0.3313419118 0.8391544118: 0.3235294118 1.331341912 2.839154412
+DEAL::0.3529411765 0.3607536765 0.8685661765: 0.3529411765 1.360753676 2.868566176
+DEAL::0.3823529412 0.3901654412 0.8979779412: 0.3823529412 1.390165441 2.897977941
+DEAL::0.4117647059 0.4195772059 0.9273897059: 0.4117647059 1.419577206 2.927389706
+DEAL::
+DEAL::Cell with center 0.7500000000 0.2500000000 0.7500000000
+DEAL::0.5000000000 0.007812500000 0.5156250000: 0.5000000000 1.007812500 2.515625000
+DEAL::0.5294117647 0.03722426471 0.5450367647: 0.5294117647 1.037224265 2.545036765
+DEAL::0.5588235294 0.06663602941 0.5744485294: 0.5588235294 1.066636029 2.574448529
+DEAL::0.5882352941 0.09604779412 0.6038602941: 0.5882352941 1.096047794 2.603860294
+DEAL::0.6176470588 0.1254595588 0.6332720588: 0.6176470588 1.125459559 2.633272059
+DEAL::0.6470588235 0.1548713235 0.6626838235: 0.6470588235 1.154871324 2.662683824
+DEAL::0.6764705882 0.1842830882 0.6920955882: 0.6764705882 1.184283088 2.692095588
+DEAL::0.7058823529 0.2136948529 0.7215073529: 0.7058823529 1.213694853 2.721507353
+DEAL::0.7352941176 0.2431066176 0.7509191176: 0.7352941176 1.243106618 2.750919118
+DEAL::0.7647058824 0.2725183824 0.7803308824: 0.7647058824 1.272518382 2.780330882
+DEAL::0.7941176471 0.3019301471 0.8097426471: 0.7941176471 1.301930147 2.809742647
+DEAL::0.8235294118 0.3313419118 0.8391544118: 0.8235294118 1.331341912 2.839154412
+DEAL::0.8529411765 0.3607536765 0.8685661765: 0.8529411765 1.360753676 2.868566176
+DEAL::
+DEAL::Cell with center 0.7500000000 0.2500000000 0.7500000000
+DEAL::0.5000000000 0.007812500000 0.5156250000: 0.5000000000 1.007812500 2.515625000
+DEAL::0.5294117647 0.03722426471 0.5450367647: 0.5294117647 1.037224265 2.545036765
+DEAL::0.5588235294 0.06663602941 0.5744485294: 0.5588235294 1.066636029 2.574448529
+DEAL::0.5882352941 0.09604779412 0.6038602941: 0.5882352941 1.096047794 2.603860294
+DEAL::0.6176470588 0.1254595588 0.6332720588: 0.6176470588 1.125459559 2.633272059
+DEAL::0.6470588235 0.1548713235 0.6626838235: 0.6470588235 1.154871324 2.662683824
+DEAL::0.6764705882 0.1842830882 0.6920955882: 0.6764705882 1.184283088 2.692095588
+DEAL::0.7058823529 0.2136948529 0.7215073529: 0.7058823529 1.213694853 2.721507353
+DEAL::0.7352941176 0.2431066176 0.7509191176: 0.7352941176 1.243106618 2.750919118
+DEAL::0.7647058824 0.2725183824 0.7803308824: 0.7647058824 1.272518382 2.780330882
+DEAL::0.7941176471 0.3019301471 0.8097426471: 0.7941176471 1.301930147 2.809742647
+DEAL::0.8235294118 0.3313419118 0.8391544118: 0.8235294118 1.331341912 2.839154412
+DEAL::0.8529411765 0.3607536765 0.8685661765: 0.8529411765 1.360753676 2.868566176
+DEAL::0.8823529412 0.3901654412 0.8979779412: 0.8823529412 1.390165441 2.897977941
+DEAL::0.9117647059 0.4195772059 0.9273897059: 0.9117647059 1.419577206 2.927389706
+DEAL::
+DEAL::Cell with center 0.2500000000 0.7500000000 0.7500000000
+DEAL::0.000000000 0.5078125000 0.5156250000: 0.000000000 1.507812500 2.515625000
+DEAL::0.02941176471 0.5372242647 0.5450367647: 0.02941176471 1.537224265 2.545036765
+DEAL::0.05882352941 0.5666360294 0.5744485294: 0.05882352941 1.566636029 2.574448529
+DEAL::0.08823529412 0.5960477941 0.6038602941: 0.08823529412 1.596047794 2.603860294
+DEAL::0.1176470588 0.6254595588 0.6332720588: 0.1176470588 1.625459559 2.633272059
+DEAL::0.1470588235 0.6548713235 0.6626838235: 0.1470588235 1.654871324 2.662683824
+DEAL::0.1764705882 0.6842830882 0.6920955882: 0.1764705882 1.684283088 2.692095588
+DEAL::0.2058823529 0.7136948529 0.7215073529: 0.2058823529 1.713694853 2.721507353
+DEAL::0.2352941176 0.7431066176 0.7509191176: 0.2352941176 1.743106618 2.750919118
+DEAL::0.2647058824 0.7725183824 0.7803308824: 0.2647058824 1.772518382 2.780330882
+DEAL::0.2941176471 0.8019301471 0.8097426471: 0.2941176471 1.801930147 2.809742647
+DEAL::0.3235294118 0.8313419118 0.8391544118: 0.3235294118 1.831341912 2.839154412
+DEAL::0.3529411765 0.8607536765 0.8685661765: 0.3529411765 1.860753676 2.868566176
+DEAL::
+DEAL::Cell with center 0.2500000000 0.7500000000 0.7500000000
+DEAL::0.000000000 0.5078125000 0.5156250000: 0.000000000 1.507812500 2.515625000
+DEAL::0.02941176471 0.5372242647 0.5450367647: 0.02941176471 1.537224265 2.545036765
+DEAL::0.05882352941 0.5666360294 0.5744485294: 0.05882352941 1.566636029 2.574448529
+DEAL::0.08823529412 0.5960477941 0.6038602941: 0.08823529412 1.596047794 2.603860294
+DEAL::0.1176470588 0.6254595588 0.6332720588: 0.1176470588 1.625459559 2.633272059
+DEAL::0.1470588235 0.6548713235 0.6626838235: 0.1470588235 1.654871324 2.662683824
+DEAL::0.1764705882 0.6842830882 0.6920955882: 0.1764705882 1.684283088 2.692095588
+DEAL::0.2058823529 0.7136948529 0.7215073529: 0.2058823529 1.713694853 2.721507353
+DEAL::0.2352941176 0.7431066176 0.7509191176: 0.2352941176 1.743106618 2.750919118
+DEAL::0.2647058824 0.7725183824 0.7803308824: 0.2647058824 1.772518382 2.780330882
+DEAL::0.2941176471 0.8019301471 0.8097426471: 0.2941176471 1.801930147 2.809742647
+DEAL::0.3235294118 0.8313419118 0.8391544118: 0.3235294118 1.831341912 2.839154412
+DEAL::0.3529411765 0.8607536765 0.8685661765: 0.3529411765 1.860753676 2.868566176
+DEAL::0.3823529412 0.8901654412 0.8979779412: 0.3823529412 1.890165441 2.897977941
+DEAL::0.4117647059 0.9195772059 0.9273897059: 0.4117647059 1.919577206 2.927389706
+DEAL::
+DEAL::Cell with center 0.7500000000 0.7500000000 0.7500000000
+DEAL::0.5000000000 0.5078125000 0.5156250000: 0.5000000000 1.507812500 2.515625000
+DEAL::0.5294117647 0.5372242647 0.5450367647: 0.5294117647 1.537224265 2.545036765
+DEAL::0.5588235294 0.5666360294 0.5744485294: 0.5588235294 1.566636029 2.574448529
+DEAL::0.5882352941 0.5960477941 0.6038602941: 0.5882352941 1.596047794 2.603860294
+DEAL::0.6176470588 0.6254595588 0.6332720588: 0.6176470588 1.625459559 2.633272059
+DEAL::0.6470588235 0.6548713235 0.6626838235: 0.6470588235 1.654871324 2.662683824
+DEAL::0.6764705882 0.6842830882 0.6920955882: 0.6764705882 1.684283088 2.692095588
+DEAL::0.7058823529 0.7136948529 0.7215073529: 0.7058823529 1.713694853 2.721507353
+DEAL::0.7352941176 0.7431066176 0.7509191176: 0.7352941176 1.743106618 2.750919118
+DEAL::0.7647058824 0.7725183824 0.7803308824: 0.7647058824 1.772518382 2.780330882
+DEAL::0.7941176471 0.8019301471 0.8097426471: 0.7941176471 1.801930147 2.809742647
+DEAL::0.8235294118 0.8313419118 0.8391544118: 0.8235294118 1.831341912 2.839154412
+DEAL::0.8529411765 0.8607536765 0.8685661765: 0.8529411765 1.860753676 2.868566176
+DEAL::
+DEAL::Cell with center 0.7500000000 0.7500000000 0.7500000000
+DEAL::0.5000000000 0.5078125000 0.5156250000: 0.5000000000 1.507812500 2.515625000
+DEAL::0.5294117647 0.5372242647 0.5450367647: 0.5294117647 1.537224265 2.545036765
+DEAL::0.5588235294 0.5666360294 0.5744485294: 0.5588235294 1.566636029 2.574448529
+DEAL::0.5882352941 0.5960477941 0.6038602941: 0.5882352941 1.596047794 2.603860294
+DEAL::0.6176470588 0.6254595588 0.6332720588: 0.6176470588 1.625459559 2.633272059
+DEAL::0.6470588235 0.6548713235 0.6626838235: 0.6470588235 1.654871324 2.662683824
+DEAL::0.6764705882 0.6842830882 0.6920955882: 0.6764705882 1.684283088 2.692095588
+DEAL::0.7058823529 0.7136948529 0.7215073529: 0.7058823529 1.713694853 2.721507353
+DEAL::0.7352941176 0.7431066176 0.7509191176: 0.7352941176 1.743106618 2.750919118
+DEAL::0.7647058824 0.7725183824 0.7803308824: 0.7647058824 1.772518382 2.780330882
+DEAL::0.7941176471 0.8019301471 0.8097426471: 0.7941176471 1.801930147 2.809742647
+DEAL::0.8235294118 0.8313419118 0.8391544118: 0.8235294118 1.831341912 2.839154412
+DEAL::0.8529411765 0.8607536765 0.8685661765: 0.8529411765 1.860753676 2.868566176
+DEAL::0.8823529412 0.8901654412 0.8979779412: 0.8823529412 1.890165441 2.897977941
+DEAL::0.9117647059 0.9195772059 0.9273897059: 0.9117647059 1.919577206 2.927389706
+DEAL::

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.