]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Allow to enforce unique mapping in RemotePointEvaluation 12199/head
authorPeter Munch <peterrmuench@gmail.com>
Wed, 12 May 2021 16:54:06 +0000 (18:54 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Tue, 18 May 2021 09:10:24 +0000 (11:10 +0200)
include/deal.II/base/mpi_remote_point_evaluation.h
include/deal.II/base/mpi_tags.h
include/deal.II/grid/grid_tools.h
source/base/mpi_remote_point_evaluation.cc
source/grid/grid_tools.cc
source/grid/grid_tools.inst.in
tests/remote_point_evaluation/mapping_01.cc [new file with mode: 0644]
tests/remote_point_evaluation/mapping_01.mpirun=1.output [new file with mode: 0644]
tests/remote_point_evaluation/mapping_01.mpirun=4.output [new file with mode: 0644]

index bf1624a6a79f7e463dfff09f16c485fd30e9288a..ebd3f3a5b8a513e27c7cc6020e62900570ca869a 100644 (file)
@@ -53,8 +53,11 @@ namespace Utilities
        *   the tolerance in order to be able to identify a cell.
        *   Floating point arithmetic implies that a point will, in general, not
        *   lie exactly on a vertex, edge, or face.
+       * @param enforce_unique_mapping Enforce unique mapping, i.e.,
+       *   (one-to-one) relation of points and cells.
        */
-      RemotePointEvaluation(const double tolerance = 1e-6);
+      RemotePointEvaluation(const double tolerance              = 1e-6,
+                            const bool   enforce_unique_mapping = false);
 
       /**
        * Destructor.
@@ -168,6 +171,12 @@ namespace Utilities
        */
       const double tolerance;
 
+      /**
+       * Enforce unique mapping, i.e., (one-to-one) relation of points and
+       * cells.
+       */
+      const bool enforce_unique_mapping;
+
       /**
        * Storage for the status of the triangulation signal.
        */
index ae7eff13e75b5c4e0718efc22b13fd70e6d4e34b..a59594d59c2e0639df944b78a1ea7452f3d9e82f 100644 (file)
@@ -146,6 +146,9 @@ namespace Utilities
           // global coarsening transfer
           fine_dof_handler_view_reinit,
 
+          // GridTools::internal::distributed_compute_point_locations
+          distributed_compute_point_locations,
+
         };
       } // namespace Tags
     }   // namespace internal
index 3ce447ad3657933484ecfde0389dd60819b869bc..9b578de60d3e87146b1ba3d62a63836abf7de0c2 100644 (file)
@@ -1149,7 +1149,8 @@ namespace GridTools
       const std::vector<Point<spacedim>> &                   points,
       const std::vector<std::vector<BoundingBox<spacedim>>> &global_bboxes,
       const double                                           tolerance,
-      const bool                                             perform_handshake);
+      const bool                                             perform_handshake,
+      const bool enforce_unique_mapping = false);
 
   } // namespace internal
 
index 8f7a01f4fab1260cc4cfa11472775ccc9873a3b1..f7ceca3aefa006d3bb4223b9cbc0a3d93719026e 100644 (file)
@@ -37,8 +37,10 @@ namespace Utilities
   {
     template <int dim, int spacedim>
     RemotePointEvaluation<dim, spacedim>::RemotePointEvaluation(
-      const double tolerance)
+      const double tolerance,
+      const bool   enforce_unique_mapping)
       : tolerance(tolerance)
+      , enforce_unique_mapping(enforce_unique_mapping)
       , ready_flag(false)
     {}
 
@@ -94,7 +96,12 @@ namespace Utilities
 
       const auto data =
         GridTools::internal::distributed_compute_point_locations(
-          cache, points, global_bboxes, tolerance, true);
+          cache,
+          points,
+          global_bboxes,
+          tolerance,
+          true,
+          enforce_unique_mapping);
 
       this->recv_ranks = data.recv_ranks;
       this->recv_ptrs  = data.recv_ptrs;
index 94abc37f78ab354ea8a7b1f4899bbb3b94553f5a..a30232ab136b3f117be169cf6b52eb2cd7aca2c4 100644 (file)
@@ -5761,6 +5761,11 @@ namespace GridTools
             if (cell.first->is_locally_owned())
               locally_owned_active_cells_around_point.push_back(cell);
         }
+
+      std::sort(locally_owned_active_cells_around_point.begin(),
+                locally_owned_active_cells_around_point.end(),
+                [](const auto &a, const auto &b) { return a.first < b.first; });
+
       return locally_owned_active_cells_around_point;
     }
 
@@ -5773,8 +5778,11 @@ namespace GridTools
       const std::vector<Point<spacedim>> &                   points,
       const std::vector<std::vector<BoundingBox<spacedim>>> &global_bboxes,
       const double                                           tolerance,
-      const bool                                             perform_handshake)
+      const bool                                             perform_handshake,
+      const bool enforce_unique_mapping)
     {
+      Assert(!enforce_unique_mapping || perform_handshake, ExcInternalError());
+
       DistributedComputePointLocationsInternal<dim, spacedim> result;
 
       auto &send_components = result.send_components;
@@ -5860,10 +5868,18 @@ namespace GridTools
                     cell_and_reference_position.second,
                     index_and_point.second,
                     numbers::invalid_unsigned_int);
+
+                  if (enforce_unique_mapping)
+                    break; // in the case of unique mapping, we only need a
+                           // single cell (we take the first)
                 }
 
               if (perform_handshake)
-                request_buffer_temp[i] = cells_and_reference_positions.size();
+                request_buffer_temp[i] =
+                  enforce_unique_mapping ?
+                    std::min<unsigned int>(
+                      1, cells_and_reference_positions.size()) :
+                    cells_and_reference_positions.size();
             }
 
           if (perform_handshake)
@@ -5905,6 +5921,143 @@ namespace GridTools
         process, cache.get_triangulation().get_communicator())
         .run();
 
+      // for unique mapping, we need to modify recv_components and
+      // send_components consistently
+      if (enforce_unique_mapping)
+        {
+          std::vector<unsigned int> mask_recv(recv_components.size());
+          std::vector<unsigned int> mask_send(send_components.size());
+
+          // set up new recv_components and monitor which entries (we keep
+          // the entry of the lowest rank) have been eliminated so that we can
+          // communicate it
+          auto recv_components_copy = recv_components;
+          recv_components.clear();
+
+          for (unsigned int i = 0; i < recv_components_copy.size(); ++i)
+            std::get<2>(recv_components_copy[i]) = i;
+
+          std::sort(recv_components_copy.begin(),
+                    recv_components_copy.end(),
+                    [&](const auto &a, const auto &b) {
+                      if (std::get<0>(a) != std::get<0>(b)) // rank
+                        return std::get<0>(a) < std::get<0>(b);
+
+                      return std::get<2>(a) < std::get<2>(b); // enumeration
+                    });
+
+          std::vector<bool> unique(points.size(), false);
+
+          std::vector<unsigned int> recv_ranks;
+          std::vector<unsigned int> recv_ptrs;
+
+          for (unsigned int i = 0, dummy = numbers::invalid_unsigned_int;
+               i < recv_components_copy.size();
+               ++i)
+            {
+              if (dummy != std::get<0>(recv_components_copy[i]))
+                {
+                  dummy = std::get<0>(recv_components_copy[i]);
+                  recv_ranks.push_back(dummy);
+                  recv_ptrs.push_back(i);
+                }
+
+              if (unique[std::get<1>(recv_components_copy[i])] == false)
+                {
+                  recv_components.emplace_back(recv_components_copy[i]);
+                  mask_recv[i]                                 = 1;
+                  unique[std::get<1>(recv_components_copy[i])] = true;
+                }
+              else
+                {
+                  mask_recv[i] = 0;
+                }
+            }
+          recv_ptrs.push_back(recv_components_copy.size());
+
+          Assert(std::all_of(unique.begin(),
+                             unique.end(),
+                             [](const auto &v) { return v; }),
+                 ExcInternalError());
+
+
+          // prepare send_components so that not needed entries can be
+          // eliminated later on
+          auto send_components_copy = send_components;
+          send_components.clear();
+
+          for (unsigned int i = 0; i < send_components_copy.size(); ++i)
+            std::get<5>(send_components_copy[i]) = i;
+
+          std::sort(send_components_copy.begin(),
+                    send_components_copy.end(),
+                    [&](const auto &a, const auto &b) {
+                      if (std::get<1>(a) != std::get<1>(b)) // rank
+                        return std::get<1>(a) < std::get<1>(b);
+
+                      return std::get<5>(a) < std::get<5>(b); // enumeration
+                    });
+
+          std::vector<unsigned int> send_ranks;
+          std::vector<unsigned int> send_ptrs;
+
+          for (unsigned int i = 0, dummy = numbers::invalid_unsigned_int;
+               i < send_components_copy.size();
+               ++i)
+            {
+              if (dummy != std::get<1>(send_components_copy[i]))
+                {
+                  dummy = std::get<1>(send_components_copy[i]);
+                  send_ranks.push_back(dummy);
+                  send_ptrs.push_back(i);
+                }
+            }
+          send_ptrs.push_back(send_components_copy.size());
+
+          // perform communication
+#ifdef DEAL_II_WITH_MPI
+          std::vector<MPI_Request> req(send_ranks.size() + recv_ranks.size());
+
+          for (unsigned int i = 0; i < send_ranks.size(); ++i)
+            {
+              const auto ierr =
+                MPI_Irecv(mask_send.data() + send_ptrs[i],
+                          send_ptrs[i + 1] - send_ptrs[i],
+                          MPI_UNSIGNED,
+                          send_ranks[i],
+                          Utilities::MPI::internal::Tags::
+                            distributed_compute_point_locations,
+                          cache.get_triangulation().get_communicator(),
+                          &req[i]);
+              AssertThrowMPI(ierr);
+            }
+
+          for (unsigned int i = 0; i < recv_ranks.size(); ++i)
+            {
+              const auto ierr =
+                MPI_Isend(mask_recv.data() + recv_ptrs[i],
+                          recv_ptrs[i + 1] - recv_ptrs[i],
+                          MPI_UNSIGNED,
+                          recv_ranks[i],
+                          Utilities::MPI::internal::Tags::
+                            distributed_compute_point_locations,
+                          cache.get_triangulation().get_communicator(),
+                          &req[i] + send_ranks.size());
+              AssertThrowMPI(ierr);
+            }
+
+          auto ierr = MPI_Waitall(req.size(), req.data(), MPI_STATUSES_IGNORE);
+          AssertThrowMPI(ierr);
+#else
+          mask_send = mask_recv;
+#endif
+
+          // eliminate not needed entries
+          for (unsigned int i = 0; i < send_components_copy.size(); ++i)
+            if (mask_send[i] == 1)
+              send_components.emplace_back(send_components_copy[i]);
+        }
+
       if (true)
         {
           // sort according to rank (and point index and cell) -> make
index 467be82266a76745dd8f50403d25fe899221fd5e..e8fefe6909e0ed783b99489ddd43af73be500090 100644 (file)
@@ -143,6 +143,7 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS)
         const std::vector<Point<deal_II_space_dimension>> &,
         const std::vector<std::vector<BoundingBox<deal_II_space_dimension>>> &,
         const double,
+        const bool,
         const bool);
     \}
 
diff --git a/tests/remote_point_evaluation/mapping_01.cc b/tests/remote_point_evaluation/mapping_01.cc
new file mode 100644 (file)
index 0000000..31e5562
--- /dev/null
@@ -0,0 +1,95 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2021 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.
+//
+// ---------------------------------------------------------------------
+
+// Test Utilities::MPI::RemotePointEvaluation with and without unique mapping.
+
+#include <deal.II/base/mpi.h>
+
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/fe/fe_dgq.h>
+#include <deal.II/fe/mapping_q1.h>
+
+#include <deal.II/grid/grid_generator.h>
+
+#include <deal.II/lac/la_vector.h>
+
+#include <deal.II/numerics/vector_tools_evaluate.h>
+
+#include "../tests.h"
+
+using namespace dealii;
+
+
+template <int dim>
+void
+test(const bool enforce_unique_map)
+{
+  parallel::distributed::Triangulation<dim> tria(MPI_COMM_WORLD);
+  GridGenerator::subdivided_hyper_cube(tria, 2);
+  tria.refine_global(1);
+
+  FE_DGQ<dim>     fe(0);
+  DoFHandler<dim> dof_handler(tria);
+  dof_handler.distribute_dofs(fe);
+
+  Vector<double> vec(dof_handler.n_dofs());
+
+  for (const auto &cell : dof_handler.active_cell_iterators())
+    if (cell->is_locally_owned())
+      vec[cell->global_active_cell_index()] = cell->global_active_cell_index();
+
+  MappingQ1<dim> mapping;
+
+  std::vector<Point<dim>> evaluation_points;
+
+  for (unsigned int j = 0; j <= 4; ++j)
+    for (unsigned int i = 0; i <= 4; ++i)
+      evaluation_points.emplace_back(i * 0.25, j * 0.25);
+
+  Utilities::MPI::RemotePointEvaluation<dim> eval(1e-6, enforce_unique_map);
+
+  const auto result_avg =
+    VectorTools::evaluate_at_points<1>(mapping,
+                                       dof_handler,
+                                       vec,
+                                       evaluation_points,
+                                       eval,
+                                       VectorTools::EvaluationFlags::avg);
+  const auto result_min = VectorTools::evaluate_at_points<1>(
+    eval, dof_handler, vec, VectorTools::EvaluationFlags::min);
+  const auto result_max = VectorTools::evaluate_at_points<1>(
+    eval, dof_handler, vec, VectorTools::EvaluationFlags::max);
+  const auto result_insert = VectorTools::evaluate_at_points<1>(
+    eval, dof_handler, vec, VectorTools::EvaluationFlags::insert);
+
+  for (unsigned int i = 0; i < evaluation_points.size(); ++i)
+    deallog << i << " " << result_avg[i] << " " << result_min[i] << " "
+            << result_max[i] << " " << result_insert[i] << " "
+            << eval.get_point_ptrs()[i + 1] - eval.get_point_ptrs()[i]
+            << std::endl;
+}
+
+
+
+int
+main(int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
+  MPILogInitAll                    all;
+
+  test<2>(false);
+  test<2>(true);
+}
diff --git a/tests/remote_point_evaluation/mapping_01.mpirun=1.output b/tests/remote_point_evaluation/mapping_01.mpirun=1.output
new file mode 100644 (file)
index 0000000..3bc2861
--- /dev/null
@@ -0,0 +1,51 @@
+
+DEAL:0::0 0.00000 0.00000 0.00000 0.00000 1
+DEAL:0::1 0.500000 0.00000 1.00000 0.00000 2
+DEAL:0::2 2.50000 1.00000 4.00000 1.00000 2
+DEAL:0::3 4.50000 4.00000 5.00000 4.00000 2
+DEAL:0::4 5.00000 5.00000 5.00000 5.00000 1
+DEAL:0::5 1.00000 0.00000 2.00000 0.00000 2
+DEAL:0::6 1.50000 0.00000 3.00000 0.00000 4
+DEAL:0::7 3.50000 1.00000 6.00000 1.00000 4
+DEAL:0::8 5.50000 4.00000 7.00000 4.00000 4
+DEAL:0::9 6.00000 5.00000 7.00000 5.00000 2
+DEAL:0::10 5.00000 2.00000 8.00000 2.00000 2
+DEAL:0::11 5.50000 2.00000 9.00000 2.00000 4
+DEAL:0::12 7.50000 3.00000 12.0000 3.00000 4
+DEAL:0::13 9.50000 6.00000 13.0000 6.00000 4
+DEAL:0::14 10.0000 7.00000 13.0000 7.00000 2
+DEAL:0::15 9.00000 8.00000 10.0000 8.00000 2
+DEAL:0::16 9.50000 8.00000 11.0000 8.00000 4
+DEAL:0::17 11.5000 9.00000 14.0000 9.00000 4
+DEAL:0::18 13.5000 12.0000 15.0000 12.0000 4
+DEAL:0::19 14.0000 13.0000 15.0000 13.0000 2
+DEAL:0::20 10.0000 10.0000 10.0000 10.0000 1
+DEAL:0::21 10.5000 10.0000 11.0000 10.0000 2
+DEAL:0::22 12.5000 11.0000 14.0000 11.0000 2
+DEAL:0::23 14.5000 14.0000 15.0000 14.0000 2
+DEAL:0::24 15.0000 15.0000 15.0000 15.0000 1
+DEAL:0::0 0.00000 0.00000 0.00000 0.00000 1
+DEAL:0::1 0.00000 0.00000 0.00000 0.00000 1
+DEAL:0::2 1.00000 1.00000 1.00000 1.00000 1
+DEAL:0::3 4.00000 4.00000 4.00000 4.00000 1
+DEAL:0::4 5.00000 5.00000 5.00000 5.00000 1
+DEAL:0::5 0.00000 0.00000 0.00000 0.00000 1
+DEAL:0::6 0.00000 0.00000 0.00000 0.00000 1
+DEAL:0::7 1.00000 1.00000 1.00000 1.00000 1
+DEAL:0::8 4.00000 4.00000 4.00000 4.00000 1
+DEAL:0::9 5.00000 5.00000 5.00000 5.00000 1
+DEAL:0::10 2.00000 2.00000 2.00000 2.00000 1
+DEAL:0::11 2.00000 2.00000 2.00000 2.00000 1
+DEAL:0::12 3.00000 3.00000 3.00000 3.00000 1
+DEAL:0::13 6.00000 6.00000 6.00000 6.00000 1
+DEAL:0::14 7.00000 7.00000 7.00000 7.00000 1
+DEAL:0::15 8.00000 8.00000 8.00000 8.00000 1
+DEAL:0::16 8.00000 8.00000 8.00000 8.00000 1
+DEAL:0::17 9.00000 9.00000 9.00000 9.00000 1
+DEAL:0::18 12.0000 12.0000 12.0000 12.0000 1
+DEAL:0::19 13.0000 13.0000 13.0000 13.0000 1
+DEAL:0::20 10.0000 10.0000 10.0000 10.0000 1
+DEAL:0::21 10.0000 10.0000 10.0000 10.0000 1
+DEAL:0::22 11.0000 11.0000 11.0000 11.0000 1
+DEAL:0::23 14.0000 14.0000 14.0000 14.0000 1
+DEAL:0::24 15.0000 15.0000 15.0000 15.0000 1
diff --git a/tests/remote_point_evaluation/mapping_01.mpirun=4.output b/tests/remote_point_evaluation/mapping_01.mpirun=4.output
new file mode 100644 (file)
index 0000000..d5a34b5
--- /dev/null
@@ -0,0 +1,207 @@
+
+DEAL:0::0 0.00000 0.00000 0.00000 0.00000 1
+DEAL:0::1 0.500000 0.00000 1.00000 0.00000 2
+DEAL:0::2 2.50000 1.00000 4.00000 1.00000 2
+DEAL:0::3 4.50000 4.00000 5.00000 4.00000 2
+DEAL:0::4 5.00000 5.00000 5.00000 5.00000 1
+DEAL:0::5 1.00000 0.00000 2.00000 0.00000 2
+DEAL:0::6 1.50000 0.00000 3.00000 0.00000 4
+DEAL:0::7 3.50000 1.00000 6.00000 1.00000 4
+DEAL:0::8 5.50000 4.00000 7.00000 4.00000 4
+DEAL:0::9 6.00000 5.00000 7.00000 5.00000 2
+DEAL:0::10 5.00000 2.00000 8.00000 2.00000 2
+DEAL:0::11 5.50000 2.00000 9.00000 2.00000 4
+DEAL:0::12 7.50000 3.00000 12.0000 3.00000 4
+DEAL:0::13 9.50000 6.00000 13.0000 6.00000 4
+DEAL:0::14 10.0000 7.00000 13.0000 7.00000 2
+DEAL:0::15 9.00000 8.00000 10.0000 8.00000 2
+DEAL:0::16 9.50000 8.00000 11.0000 8.00000 4
+DEAL:0::17 11.5000 9.00000 14.0000 9.00000 4
+DEAL:0::18 13.5000 12.0000 15.0000 12.0000 4
+DEAL:0::19 14.0000 13.0000 15.0000 13.0000 2
+DEAL:0::20 10.0000 10.0000 10.0000 10.0000 1
+DEAL:0::21 10.5000 10.0000 11.0000 10.0000 2
+DEAL:0::22 12.5000 11.0000 14.0000 11.0000 2
+DEAL:0::23 14.5000 14.0000 15.0000 14.0000 2
+DEAL:0::24 15.0000 15.0000 15.0000 15.0000 1
+DEAL:0::0 0.00000 0.00000 0.00000 0.00000 1
+DEAL:0::1 0.00000 0.00000 0.00000 0.00000 1
+DEAL:0::2 1.00000 1.00000 1.00000 1.00000 1
+DEAL:0::3 4.00000 4.00000 4.00000 4.00000 1
+DEAL:0::4 5.00000 5.00000 5.00000 5.00000 1
+DEAL:0::5 0.00000 0.00000 0.00000 0.00000 1
+DEAL:0::6 0.00000 0.00000 0.00000 0.00000 1
+DEAL:0::7 1.00000 1.00000 1.00000 1.00000 1
+DEAL:0::8 4.00000 4.00000 4.00000 4.00000 1
+DEAL:0::9 5.00000 5.00000 5.00000 5.00000 1
+DEAL:0::10 2.00000 2.00000 2.00000 2.00000 1
+DEAL:0::11 2.00000 2.00000 2.00000 2.00000 1
+DEAL:0::12 3.00000 3.00000 3.00000 3.00000 1
+DEAL:0::13 6.00000 6.00000 6.00000 6.00000 1
+DEAL:0::14 7.00000 7.00000 7.00000 7.00000 1
+DEAL:0::15 8.00000 8.00000 8.00000 8.00000 1
+DEAL:0::16 8.00000 8.00000 8.00000 8.00000 1
+DEAL:0::17 9.00000 9.00000 9.00000 9.00000 1
+DEAL:0::18 12.0000 12.0000 12.0000 12.0000 1
+DEAL:0::19 13.0000 13.0000 13.0000 13.0000 1
+DEAL:0::20 10.0000 10.0000 10.0000 10.0000 1
+DEAL:0::21 10.0000 10.0000 10.0000 10.0000 1
+DEAL:0::22 11.0000 11.0000 11.0000 11.0000 1
+DEAL:0::23 14.0000 14.0000 14.0000 14.0000 1
+DEAL:0::24 15.0000 15.0000 15.0000 15.0000 1
+
+DEAL:1::0 0.00000 0.00000 0.00000 0.00000 1
+DEAL:1::1 0.500000 0.00000 1.00000 0.00000 2
+DEAL:1::2 2.50000 1.00000 4.00000 1.00000 2
+DEAL:1::3 4.50000 4.00000 5.00000 4.00000 2
+DEAL:1::4 5.00000 5.00000 5.00000 5.00000 1
+DEAL:1::5 1.00000 0.00000 2.00000 0.00000 2
+DEAL:1::6 1.50000 0.00000 3.00000 0.00000 4
+DEAL:1::7 3.50000 1.00000 6.00000 1.00000 4
+DEAL:1::8 5.50000 4.00000 7.00000 4.00000 4
+DEAL:1::9 6.00000 5.00000 7.00000 5.00000 2
+DEAL:1::10 5.00000 2.00000 8.00000 2.00000 2
+DEAL:1::11 5.50000 2.00000 9.00000 2.00000 4
+DEAL:1::12 7.50000 3.00000 12.0000 3.00000 4
+DEAL:1::13 9.50000 6.00000 13.0000 6.00000 4
+DEAL:1::14 10.0000 7.00000 13.0000 7.00000 2
+DEAL:1::15 9.00000 8.00000 10.0000 8.00000 2
+DEAL:1::16 9.50000 8.00000 11.0000 8.00000 4
+DEAL:1::17 11.5000 9.00000 14.0000 9.00000 4
+DEAL:1::18 13.5000 12.0000 15.0000 12.0000 4
+DEAL:1::19 14.0000 13.0000 15.0000 13.0000 2
+DEAL:1::20 10.0000 10.0000 10.0000 10.0000 1
+DEAL:1::21 10.5000 10.0000 11.0000 10.0000 2
+DEAL:1::22 12.5000 11.0000 14.0000 11.0000 2
+DEAL:1::23 14.5000 14.0000 15.0000 14.0000 2
+DEAL:1::24 15.0000 15.0000 15.0000 15.0000 1
+DEAL:1::0 0.00000 0.00000 0.00000 0.00000 1
+DEAL:1::1 0.00000 0.00000 0.00000 0.00000 1
+DEAL:1::2 1.00000 1.00000 1.00000 1.00000 1
+DEAL:1::3 4.00000 4.00000 4.00000 4.00000 1
+DEAL:1::4 5.00000 5.00000 5.00000 5.00000 1
+DEAL:1::5 0.00000 0.00000 0.00000 0.00000 1
+DEAL:1::6 0.00000 0.00000 0.00000 0.00000 1
+DEAL:1::7 1.00000 1.00000 1.00000 1.00000 1
+DEAL:1::8 4.00000 4.00000 4.00000 4.00000 1
+DEAL:1::9 5.00000 5.00000 5.00000 5.00000 1
+DEAL:1::10 2.00000 2.00000 2.00000 2.00000 1
+DEAL:1::11 2.00000 2.00000 2.00000 2.00000 1
+DEAL:1::12 3.00000 3.00000 3.00000 3.00000 1
+DEAL:1::13 6.00000 6.00000 6.00000 6.00000 1
+DEAL:1::14 7.00000 7.00000 7.00000 7.00000 1
+DEAL:1::15 8.00000 8.00000 8.00000 8.00000 1
+DEAL:1::16 8.00000 8.00000 8.00000 8.00000 1
+DEAL:1::17 9.00000 9.00000 9.00000 9.00000 1
+DEAL:1::18 12.0000 12.0000 12.0000 12.0000 1
+DEAL:1::19 13.0000 13.0000 13.0000 13.0000 1
+DEAL:1::20 10.0000 10.0000 10.0000 10.0000 1
+DEAL:1::21 10.0000 10.0000 10.0000 10.0000 1
+DEAL:1::22 11.0000 11.0000 11.0000 11.0000 1
+DEAL:1::23 14.0000 14.0000 14.0000 14.0000 1
+DEAL:1::24 15.0000 15.0000 15.0000 15.0000 1
+
+
+DEAL:2::0 0.00000 0.00000 0.00000 0.00000 1
+DEAL:2::1 0.500000 0.00000 1.00000 0.00000 2
+DEAL:2::2 2.50000 1.00000 4.00000 1.00000 2
+DEAL:2::3 4.50000 4.00000 5.00000 4.00000 2
+DEAL:2::4 5.00000 5.00000 5.00000 5.00000 1
+DEAL:2::5 1.00000 0.00000 2.00000 0.00000 2
+DEAL:2::6 1.50000 0.00000 3.00000 0.00000 4
+DEAL:2::7 3.50000 1.00000 6.00000 1.00000 4
+DEAL:2::8 5.50000 4.00000 7.00000 4.00000 4
+DEAL:2::9 6.00000 5.00000 7.00000 5.00000 2
+DEAL:2::10 5.00000 2.00000 8.00000 2.00000 2
+DEAL:2::11 5.50000 2.00000 9.00000 2.00000 4
+DEAL:2::12 7.50000 3.00000 12.0000 3.00000 4
+DEAL:2::13 9.50000 6.00000 13.0000 6.00000 4
+DEAL:2::14 10.0000 7.00000 13.0000 7.00000 2
+DEAL:2::15 9.00000 8.00000 10.0000 8.00000 2
+DEAL:2::16 9.50000 8.00000 11.0000 8.00000 4
+DEAL:2::17 11.5000 9.00000 14.0000 9.00000 4
+DEAL:2::18 13.5000 12.0000 15.0000 12.0000 4
+DEAL:2::19 14.0000 13.0000 15.0000 13.0000 2
+DEAL:2::20 10.0000 10.0000 10.0000 10.0000 1
+DEAL:2::21 10.5000 10.0000 11.0000 10.0000 2
+DEAL:2::22 12.5000 11.0000 14.0000 11.0000 2
+DEAL:2::23 14.5000 14.0000 15.0000 14.0000 2
+DEAL:2::24 15.0000 15.0000 15.0000 15.0000 1
+DEAL:2::0 0.00000 0.00000 0.00000 0.00000 1
+DEAL:2::1 0.00000 0.00000 0.00000 0.00000 1
+DEAL:2::2 1.00000 1.00000 1.00000 1.00000 1
+DEAL:2::3 4.00000 4.00000 4.00000 4.00000 1
+DEAL:2::4 5.00000 5.00000 5.00000 5.00000 1
+DEAL:2::5 0.00000 0.00000 0.00000 0.00000 1
+DEAL:2::6 0.00000 0.00000 0.00000 0.00000 1
+DEAL:2::7 1.00000 1.00000 1.00000 1.00000 1
+DEAL:2::8 4.00000 4.00000 4.00000 4.00000 1
+DEAL:2::9 5.00000 5.00000 5.00000 5.00000 1
+DEAL:2::10 2.00000 2.00000 2.00000 2.00000 1
+DEAL:2::11 2.00000 2.00000 2.00000 2.00000 1
+DEAL:2::12 3.00000 3.00000 3.00000 3.00000 1
+DEAL:2::13 6.00000 6.00000 6.00000 6.00000 1
+DEAL:2::14 7.00000 7.00000 7.00000 7.00000 1
+DEAL:2::15 8.00000 8.00000 8.00000 8.00000 1
+DEAL:2::16 8.00000 8.00000 8.00000 8.00000 1
+DEAL:2::17 9.00000 9.00000 9.00000 9.00000 1
+DEAL:2::18 12.0000 12.0000 12.0000 12.0000 1
+DEAL:2::19 13.0000 13.0000 13.0000 13.0000 1
+DEAL:2::20 10.0000 10.0000 10.0000 10.0000 1
+DEAL:2::21 10.0000 10.0000 10.0000 10.0000 1
+DEAL:2::22 11.0000 11.0000 11.0000 11.0000 1
+DEAL:2::23 14.0000 14.0000 14.0000 14.0000 1
+DEAL:2::24 15.0000 15.0000 15.0000 15.0000 1
+
+
+DEAL:3::0 0.00000 0.00000 0.00000 0.00000 1
+DEAL:3::1 0.500000 0.00000 1.00000 0.00000 2
+DEAL:3::2 2.50000 1.00000 4.00000 1.00000 2
+DEAL:3::3 4.50000 4.00000 5.00000 4.00000 2
+DEAL:3::4 5.00000 5.00000 5.00000 5.00000 1
+DEAL:3::5 1.00000 0.00000 2.00000 0.00000 2
+DEAL:3::6 1.50000 0.00000 3.00000 0.00000 4
+DEAL:3::7 3.50000 1.00000 6.00000 1.00000 4
+DEAL:3::8 5.50000 4.00000 7.00000 4.00000 4
+DEAL:3::9 6.00000 5.00000 7.00000 5.00000 2
+DEAL:3::10 5.00000 2.00000 8.00000 2.00000 2
+DEAL:3::11 5.50000 2.00000 9.00000 2.00000 4
+DEAL:3::12 7.50000 3.00000 12.0000 3.00000 4
+DEAL:3::13 9.50000 6.00000 13.0000 6.00000 4
+DEAL:3::14 10.0000 7.00000 13.0000 7.00000 2
+DEAL:3::15 9.00000 8.00000 10.0000 8.00000 2
+DEAL:3::16 9.50000 8.00000 11.0000 8.00000 4
+DEAL:3::17 11.5000 9.00000 14.0000 9.00000 4
+DEAL:3::18 13.5000 12.0000 15.0000 12.0000 4
+DEAL:3::19 14.0000 13.0000 15.0000 13.0000 2
+DEAL:3::20 10.0000 10.0000 10.0000 10.0000 1
+DEAL:3::21 10.5000 10.0000 11.0000 10.0000 2
+DEAL:3::22 12.5000 11.0000 14.0000 11.0000 2
+DEAL:3::23 14.5000 14.0000 15.0000 14.0000 2
+DEAL:3::24 15.0000 15.0000 15.0000 15.0000 1
+DEAL:3::0 0.00000 0.00000 0.00000 0.00000 1
+DEAL:3::1 0.00000 0.00000 0.00000 0.00000 1
+DEAL:3::2 1.00000 1.00000 1.00000 1.00000 1
+DEAL:3::3 4.00000 4.00000 4.00000 4.00000 1
+DEAL:3::4 5.00000 5.00000 5.00000 5.00000 1
+DEAL:3::5 0.00000 0.00000 0.00000 0.00000 1
+DEAL:3::6 0.00000 0.00000 0.00000 0.00000 1
+DEAL:3::7 1.00000 1.00000 1.00000 1.00000 1
+DEAL:3::8 4.00000 4.00000 4.00000 4.00000 1
+DEAL:3::9 5.00000 5.00000 5.00000 5.00000 1
+DEAL:3::10 2.00000 2.00000 2.00000 2.00000 1
+DEAL:3::11 2.00000 2.00000 2.00000 2.00000 1
+DEAL:3::12 3.00000 3.00000 3.00000 3.00000 1
+DEAL:3::13 6.00000 6.00000 6.00000 6.00000 1
+DEAL:3::14 7.00000 7.00000 7.00000 7.00000 1
+DEAL:3::15 8.00000 8.00000 8.00000 8.00000 1
+DEAL:3::16 8.00000 8.00000 8.00000 8.00000 1
+DEAL:3::17 9.00000 9.00000 9.00000 9.00000 1
+DEAL:3::18 12.0000 12.0000 12.0000 12.0000 1
+DEAL:3::19 13.0000 13.0000 13.0000 13.0000 1
+DEAL:3::20 10.0000 10.0000 10.0000 10.0000 1
+DEAL:3::21 10.0000 10.0000 10.0000 10.0000 1
+DEAL:3::22 11.0000 11.0000 11.0000 11.0000 1
+DEAL:3::23 14.0000 14.0000 14.0000 14.0000 1
+DEAL:3::24 15.0000 15.0000 15.0000 15.0000 1
+

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.