]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Simplify the logic used and extract all locally relevant DoFs 5604/head
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Thu, 21 Dec 2017 17:27:59 +0000 (18:27 +0100)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Sat, 23 Dec 2017 16:35:09 +0000 (17:35 +0100)
doc/news/changes/minor/20171209DanielArndt
include/deal.II/dofs/dof_tools.h
source/dofs/dof_tools.cc
tests/dofs/dof_tools_04.cc
tests/dofs/dof_tools_04.output
tests/dofs/dof_tools_04a.cc
tests/dofs/dof_tools_04a.with_p4est=true.mpirun=1.output
tests/dofs/dof_tools_04a.with_p4est=true.mpirun=3.output

index 17850e19353346c54b3e43a25a9c6623d9f6a4de..7fb7fa7815849cdba402be4ab2fe17507d58a2eb 100644 (file)
@@ -1,5 +1,5 @@
 Improved: DoFTools::extract_hanging_node_dofs works for 
 parallel::shared::Triangualtion and parallel::distributed::Triangulation
-and reports DoFs belonging to locally owned cells.
+and reports locally relevant DoFs.
 <br>
-(Daniel Arndt, 2017/12/09)
+(Daniel Arndt, 2017/12/21)
index 520ed7c066f6e4b97d7891ea9dfc8c7d70c1ca55..066f8ff2e85259fdd16068b74d1fd42d8bf24031 100644 (file)
@@ -1249,9 +1249,13 @@ namespace DoFTools
    * The size of @p selected_dofs shall equal <tt>dof_handler.n_dofs()</tt>.
    * Previous contents of this array or overwritten.
    *
-   * In case of a parallel::shared::Triangualtion or a
-   * parallel::distributed::Triangulation only dofs belonging to locally owned
-   * cells are reported.
+   * In case of a parallel::shared::Triangulation or a
+   * parallel::distributed::Triangulation only locally relevant dofs are
+   * considered. Note that the vector returned through the second argument still
+   * has size <tt>dof_handler.n_dofs()</tt>. Consequently, it can be very large
+   * for large parallel computations -- in fact, it may be too large to store on
+   * each processor. In that case, you may want to choose the variant of this
+   * function that returns an IndexSet object.
    */
   template <int dim, int spacedim>
   void
index d498af038a61ce693b27a400d56b8623d68e8906..b175999346bc72bcc7aee7e6f79a6523b7abdc8d 100644 (file)
@@ -890,7 +890,7 @@ namespace DoFTools
         // this function is similar to the make_sparsity_pattern function,
         // see there for more information
         for (auto cell : dof_handler.active_cell_iterators())
-          if (cell->is_locally_owned())
+          if (!cell->is_artificial())
             {
               for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
                 if (cell->face(face)->has_children())
@@ -902,8 +902,12 @@ namespace DoFTools
                       selected_dofs.add_index(line->child(0)->vertex_dof_index(1,dof));
 
                     for (unsigned int child=0; child<2; ++child)
-                      for (unsigned int dof=0; dof!=fe.dofs_per_line; ++dof)
-                        selected_dofs.add_index(line->child(child)->dof_index(dof));
+                      {
+                        if (cell->neighbor_child_on_subface (face, child)->is_artificial())
+                          continue;
+                        for (unsigned int dof=0; dof!=fe.dofs_per_line; ++dof)
+                          selected_dofs.add_index(line->child(child)->dof_index(dof));
+                      }
                   }
             }
 
@@ -918,52 +922,35 @@ namespace DoFTools
         const unsigned int dim = 3;
 
         IndexSet selected_dofs (dof_handler.n_dofs());
+        IndexSet unconstrained_dofs (dof_handler.n_dofs());
 
         const FiniteElement<dim,spacedim> &fe = dof_handler.get_fe();
 
-        // this function is similar to the make_sparsity_pattern function,
-        // see there for more information
         for (auto cell : dof_handler.active_cell_iterators())
-          if (cell->is_locally_owned())
+          if (!cell->is_artificial())
             for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
-              if (cell->face(f)->has_children())
-                {
-                  const typename dealii::DoFHandler<dim,spacedim>::face_iterator
-                  face = cell->face(f);
-
-                  for (unsigned int dof=0; dof!=fe.dofs_per_vertex; ++dof)
-                    selected_dofs.add_index(face->child(0)->vertex_dof_index(2,dof));
-
-                  // dof numbers on the centers of the lines bounding this
-                  // face
-                  for (unsigned int line=0; line<4; ++line)
-                    for (unsigned int dof=0; dof!=fe.dofs_per_vertex; ++dof)
-                      selected_dofs.add_index(face->line(line)->child(0)->vertex_dof_index(1,dof));
-
-                  // next the dofs on the lines interior to the face; the
-                  // order of these lines is laid down in the FiniteElement
-                  // class documentation
-                  for (unsigned int dof=0; dof<fe.dofs_per_line; ++dof)
-                    selected_dofs.add_index(face->child(0)->line(1)->dof_index(dof));
-                  for (unsigned int dof=0; dof<fe.dofs_per_line; ++dof)
-                    selected_dofs.add_index(face->child(1)->line(2)->dof_index(dof));
-                  for (unsigned int dof=0; dof<fe.dofs_per_line; ++dof)
-                    selected_dofs.add_index(face->child(2)->line(3)->dof_index(dof));
-                  for (unsigned int dof=0; dof<fe.dofs_per_line; ++dof)
-                    selected_dofs.add_index(face->child(3)->line(0)->dof_index(dof));
-
-                  // dofs on the bordering lines
-                  for (unsigned int line=0; line<4; ++line)
-                    for (unsigned int child=0; child<2; ++child)
-                      for (unsigned int dof=0; dof!=fe.dofs_per_line; ++dof)
-                        selected_dofs.add_index(face->line(line)->child(child)->dof_index(dof));
+              {
+                const typename dealii::DoFHandler<dim,spacedim>::face_iterator
+                face = cell->face(f);
+                if (cell->face(f)->has_children())
+                  {
+                    for (unsigned int child=0; child<4; ++child)
+                      if (!cell->neighbor_child_on_subface (f, child)->is_artificial())
+                        {
+                          // simply take all DoFs that live on this subface
+                          std::vector<types::global_dof_index> ldi(fe.dofs_per_face);
+                          face->child(child)->get_dof_indices(ldi);
+                          selected_dofs.add_indices(ldi.begin(), ldi.end());
+                        }
 
-                  // finally, for the dofs interior to the four child faces
-                  for (unsigned int child=0; child<4; ++child)
-                    for (unsigned int dof=0; dof!=fe.dofs_per_quad; ++dof)
-                      selected_dofs.add_index(face->child(child)->dof_index(dof));
-                }
-        selected_dofs.compress();
+                    // and subtract (in the end) all the indices which a shared
+                    // between this face and its subfaces
+                    for (unsigned int vertex=0; vertex<4; ++vertex)
+                      for (unsigned int dof=0; dof!=fe.dofs_per_vertex; ++dof)
+                        unconstrained_dofs.add_index(face->vertex_dof_index(vertex,dof));
+                  }
+              }
+        selected_dofs.subtract_set(unconstrained_dofs);
         return selected_dofs;
       }
     }
index 78400ff1df4ea7a0e84c824fb8585d9f586c419d..ee89ead282cc105a002e29c8a56aae4206470827 100644 (file)
@@ -27,8 +27,21 @@ template <int dim>
 void
 check_this (const DoFHandler<dim> &dof_handler)
 {
-  std::vector<bool> hanging_node_dofs (dof_handler.n_dofs());
+  const types::global_dof_index n_dofs = dof_handler.n_dofs();
+
+  std::vector<bool> hanging_node_dofs (n_dofs);
   DoFTools::extract_hanging_node_dofs (dof_handler,
                                        hanging_node_dofs);
+
+  ConstraintMatrix constraints;
+  DoFTools::make_hanging_node_constraints(dof_handler, constraints);
+  constraints.close();
+
+  for (types::global_dof_index dof=0; dof<n_dofs; ++dof)
+    if (hanging_node_dofs[dof])
+      AssertThrow(constraints.is_constrained(dof), ExcInternalError());
+
+  AssertThrow (std::count(hanging_node_dofs.begin(), hanging_node_dofs.end(), true)
+               == constraints.n_constraints(), ExcInternalError());
   output_bool_vector (hanging_node_dofs);
 }
index 9dc1779f19878d37e93c5d2d9adc595721c1d8a6..8ae1abb3955d3cb16a80f90eeb82951aa653b5c2 100644 (file)
@@ -4,19 +4,19 @@ DEAL::00000
 DEAL::Checking Q1 in 2d:
 DEAL::000000000000100001
 DEAL::Checking Q1 in 3d:
-DEAL::000000000000000000000000000000000000110111011000000001011011
+DEAL::000000000000000000000000000000000000111111111000000001111111
 DEAL::Checking Q2 in 1d:
 DEAL::000000000
 DEAL::Checking Q2 in 2d:
 DEAL::000000000000000000000000000000100100100000000000001010010
 DEAL::Checking Q2 in 3d:
-DEAL::000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000110001001100010001010110100011011100000000010111000000101111100011011101100000000000000000000000000000000000001001010010001010110100010110000001011100000010111001101110110
+DEAL::000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111001001110010001010110100011111110000000010111100000101111110011011101100000000000000000000000000000000000001101011010001010110100011111000001011110000010111101101110110
 DEAL::Checking Q3 in 1d:
 DEAL::0000000000000
 DEAL::Checking Q3 in 2d:
 DEAL::00000000000000000000000000000000000000000000000000000000100001100000011000000000000000000000000000010011000000110000
 DEAL::Checking Q3 in 3d:
-DEAL::00000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000110000011000011110000000000111100000000000000001001100111100001111000000000000000011011111100000000000000000000000000111100000000111110000000000000000001111000000001111111110000000000111111110000000011111100001111111100000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000100011001100000011110000000000000000100110011110000111100000000000000001011110000000000000000001111000000001111100000000000000000011110000000011111000000111111110000000011111100001111111100000000
+DEAL::00000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111000011000011111100000000111100000000000000001001100111100001111000000000000000011111111111000000000000000000000000111100000000111111100000000000000001111000000001111111111100000000111111110000000011111100001111111100000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000110011001111000011110000000000000000100110011110000111100000000000000001111111100000000000000001111000000001111111000000000000000011110000000011111110000111111110000000011111100001111111100000000
 DEAL::Checking DGQ0 in 1d:
 DEAL::0000
 DEAL::Checking DGQ0 in 2d:
@@ -56,7 +56,7 @@ DEAL::00000000000000000000000000000000000000000000000000000000000000000000000000
 DEAL::Checking Nedelec0 in 2d:
 DEAL::00000000000000001010000000101
 DEAL::Checking Nedelec0 in 3d:
-DEAL::0000000000000000000000000000000000000000000000000000000000000000000010011001011111000001100011110111000000000000000001010010111100011000110111
+DEAL::0000000000000000000000000000000000000000000000000000000000000000000010011101011111100001110011111111000000000000000001011010111110011100111111
 DEAL::Checking RaviartThomas0 in 2d:
 DEAL::00000000000000001010000000101
 DEAL::Checking RaviartThomas1 in 2d:
@@ -88,7 +88,7 @@ DEAL::00000000000000000000000000000000000000000000000000000000000000000000000000
 DEAL::Checking FE_DGP<2>(3)1 in 2d:
 DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
 DEAL::Checking FE_Q<3>(1)3 in 3d:
-DEAL::000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111111000111111111000111111000000000000000000000000111000111111000111111
+DEAL::000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111111111111111111111111111000000000000000000000000111111111111111111111
 DEAL::Checking FE_DGQ<3>(2)2 in 3d:
 DEAL::000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
 DEAL::Checking FE_DGP<3>(3)1 in 3d:
@@ -118,7 +118,7 @@ DEAL::00000000000000000000000000000000000000000000000000000000000000000000000000
 DEAL::Checking FE_DGP<2>(3)1FE_DGP<2>(3)1FE_Q<2>(2)3 in 2d:
 DEAL::00000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111000000111000000000000000000000000001110000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000001110001110000000000000000000000000011100000000000000000000000
 DEAL::Checking FE_DGQ<3>(1)3FE_DGP<3>(3)1FE_Q<3>(1)3 in 3d:
-DEAL::00000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000001111110000000000000000000000000000000000000000000000011100000000000000000000000000000000000000000000111111000000000000000000000000000000000000000000000001110000000000000000000000000000000000000000000011100000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111000000000000000000000000000000000000000000000001110000000000000000000000000000000000000000000011100000000000000000000000000000000000000000000000111000000000000000000000000000000000000000000001110000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
+DEAL::00000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000001111111110000000000000000000000000000000000000000000011100000000000000000000000000000000000000000000111111111000000000000000000000000000000000000000000001110000000000000000000000000000000000000000000011100000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111111000000000000000000000000000000000000000000001110000000000000000000000000000000000000000000011111100000000000000000000000000000000000000000000111000000000000000000000000000000000000000000001110000000000000000000000000000000000000000000000000000000000000000000000000000000000000000
 DEAL::Checking (FESystem<2>(FE_Q<2>(1),3))3FE_DGQ<2>(0)1FE_Q<2>(1)3 in 2d:
 DEAL::0000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000000111111111111000000000000000000000000000000000000000000000000000011111111111100
 DEAL::Checking FE_DGQ<2>(3)1FESystem<2>(FE_DGQ<2>(0),3)1FESystem<2>(FE_Q<2>(2),1, FE_DGQ<2>(0),1)2 in 2d:
index 481daf79a7998b16ac19fe996720485d68ce245a..aa0421b7c1ac854e56f00d0e268d1dbe52aa09b3 100644 (file)
 // ---------------------------------------------------------------------
 
 
-// The same as dof_tools_04 but for a parallel::distributed::parallel::distributed::Triangulation
-// instead of a parallel::distributed::Triangulation:
+// Similar to dof_tools_04 but for a parallel::distributed::Triangulation
+// instead of a (serial) Triangulation:
 // check
 //   DoFTools::extract_hanging_node_constraints
-// uses a slightly different refinement and less different FiniteElements
+// using a slightly different refinement and less different FiniteElements
 
 #include "../tests.h"
 #include <deal.II/base/logstream.h>
@@ -44,72 +44,27 @@ template <int dim>
 void
 check_this (const DoFHandler<dim> &dof_handler)
 {
-  types::global_dof_index n_dofs = dof_handler.n_dofs();
+  const types::global_dof_index n_dofs = dof_handler.n_dofs();
+
+  IndexSet locally_relevant_dofs;
+  DoFTools::extract_locally_relevant_dofs (dof_handler, locally_relevant_dofs);
 
   std::vector<bool> is_hanging_node_constrained (n_dofs);
   DoFTools::extract_hanging_node_dofs (dof_handler,
                                        is_hanging_node_constrained);
 
-  MappingQGeneric<dim> mapping(1);
-  std::map< types::global_dof_index, Point< dim > > support_point_of_dof;
-  DoFTools::map_dofs_to_support_points (mapping, dof_handler, support_point_of_dof);
-
-  std::vector<Point<dim>> constrained_points;
-
-  for (const auto &pair : support_point_of_dof)
-    if (is_hanging_node_constrained[pair.first])
-      constrained_points.push_back(pair.second);
-
-  const int root = 0;
-  const int mylen = constrained_points.size()*dim;
-  const unsigned int n_processes = Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD);
-  const unsigned int my_id = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
-  std::vector<int> recvcounts (n_processes);
-
-  MPI_Gather(&mylen, 1, MPI_INT, recvcounts.data(), 1, MPI_INT, 0, MPI_COMM_WORLD);
-
-  int totlen = 0;
-  std::vector<int> displs;
-  std::vector<Point<dim>> all_points;
+  ConstraintMatrix constraints (locally_relevant_dofs);
+  DoFTools::make_hanging_node_constraints(dof_handler, constraints);
+  constraints.close();
 
-  if (my_id==0)
-    {
-      displs.resize(n_processes);
-
-      displs[0] = 0;
-      totlen += recvcounts[0];
-
-      for (unsigned int i=1; i<n_processes; i++)
-        {
-          totlen += recvcounts[i];
-          displs[i] = displs[i-1] + recvcounts[i-1];
-        }
+  for (const auto &dof : locally_relevant_dofs)
+    if (is_hanging_node_constrained[dof])
+      AssertThrow(constraints.is_constrained(dof), ExcInternalError());
 
-      all_points.resize(totlen/dim);
-    }
+  AssertThrow (std::count(is_hanging_node_constrained.begin(), is_hanging_node_constrained.end(), true)
+               == constraints.n_constraints(), ExcInternalError());
 
-  MPI_Gatherv(constrained_points.data(), mylen, MPI_DOUBLE,
-              all_points.data(), recvcounts.data(), displs.data(), MPI_DOUBLE,
-              0, MPI_COMM_WORLD);
-
-  std::sort(all_points.begin(), all_points.end(),
-            [](const Point<dim> &a, const Point<dim> &b)
-  {
-    for (unsigned int i=0; i<dim; ++i)
-      {
-        if (a(i) < b(i))
-          return true;
-        if (a(i) > b(i))
-          return false;
-      }
-    return false;
-  });
-  all_points.erase(std::unique(all_points.begin(), all_points.end()), all_points.end());
-
-  deallog << all_points.size() << std::endl;
-  for (const auto &point: all_points)
-    deallog << point << std::endl;
-  deallog << std::endl;
+  deallog << "OK" << std::endl;
 }
 
 
@@ -137,6 +92,7 @@ check (const FiniteElement<dim> &fe,
           cell->set_refine_flag();
       tria.execute_coarsening_and_refinement ();
     }
+
   DoFHandler<dim> dof_handler (tria);
   dof_handler.distribute_dofs (fe);
 
@@ -157,6 +113,7 @@ int main (int argc, char **argv)
   Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
   mpi_initlog();
   CHECK_ALL(Q,1);
+
   CHECK_ALL(Q,2);
 
   return 0;
index 6feb20d8b2efee0924ee888f05259e56e8b2c08b..c3d49df7e0907994c49df1e6767154796a999f70 100644 (file)
@@ -1,285 +1,9 @@
 
 DEAL::Checking Q1 in 2d:
-DEAL::4
-DEAL::0.125000 0.500000
-DEAL::0.375000 0.500000
-DEAL::0.500000 0.125000
-DEAL::0.500000 0.375000
-DEAL::
+DEAL::OK
 DEAL::Checking Q1 in 3d:
-DEAL::40
-DEAL::0.00000 0.500000 0.125000
-DEAL::0.00000 0.500000 0.375000
-DEAL::0.00000 0.500000 0.625000
-DEAL::0.00000 0.500000 0.875000
-DEAL::0.125000 0.500000 0.00000
-DEAL::0.125000 0.500000 0.250000
-DEAL::0.125000 0.500000 0.500000
-DEAL::0.125000 0.500000 0.750000
-DEAL::0.125000 0.500000 1.00000
-DEAL::0.250000 0.500000 0.125000
-DEAL::0.250000 0.500000 0.375000
-DEAL::0.250000 0.500000 0.625000
-DEAL::0.250000 0.500000 0.875000
-DEAL::0.375000 0.500000 0.00000
-DEAL::0.375000 0.500000 0.250000
-DEAL::0.375000 0.500000 0.500000
-DEAL::0.375000 0.500000 0.750000
-DEAL::0.375000 0.500000 1.00000
-DEAL::0.500000 0.00000 0.125000
-DEAL::0.500000 0.00000 0.375000
-DEAL::0.500000 0.00000 0.625000
-DEAL::0.500000 0.00000 0.875000
-DEAL::0.500000 0.125000 0.00000
-DEAL::0.500000 0.125000 0.250000
-DEAL::0.500000 0.125000 0.750000
-DEAL::0.500000 0.125000 1.00000
-DEAL::0.500000 0.125000 0.500000
-DEAL::0.500000 0.250000 0.125000
-DEAL::0.500000 0.250000 0.375000
-DEAL::0.500000 0.250000 0.625000
-DEAL::0.500000 0.250000 0.875000
-DEAL::0.500000 0.375000 0.00000
-DEAL::0.500000 0.375000 0.250000
-DEAL::0.500000 0.375000 0.750000
-DEAL::0.500000 0.375000 1.00000
-DEAL::0.500000 0.375000 0.500000
-DEAL::0.500000 0.500000 0.125000
-DEAL::0.500000 0.500000 0.375000
-DEAL::0.500000 0.500000 0.625000
-DEAL::0.500000 0.500000 0.875000
-DEAL::
+DEAL::OK
 DEAL::Checking Q2 in 2d:
-DEAL::12
-DEAL::0.0625000 0.500000
-DEAL::0.125000 0.500000
-DEAL::0.187500 0.500000
-DEAL::0.312500 0.500000
-DEAL::0.375000 0.500000
-DEAL::0.437500 0.500000
-DEAL::0.500000 0.0625000
-DEAL::0.500000 0.125000
-DEAL::0.500000 0.187500
-DEAL::0.500000 0.312500
-DEAL::0.500000 0.375000
-DEAL::0.500000 0.437500
-DEAL::
+DEAL::OK
 DEAL::Checking Q2 in 3d:
-DEAL::216
-DEAL::0.00000 0.500000 0.0625000
-DEAL::0.00000 0.500000 0.125000
-DEAL::0.00000 0.500000 0.187500
-DEAL::0.00000 0.500000 0.312500
-DEAL::0.00000 0.500000 0.375000
-DEAL::0.00000 0.500000 0.437500
-DEAL::0.00000 0.500000 0.562500
-DEAL::0.00000 0.500000 0.625000
-DEAL::0.00000 0.500000 0.687500
-DEAL::0.00000 0.500000 0.812500
-DEAL::0.00000 0.500000 0.875000
-DEAL::0.00000 0.500000 0.937500
-DEAL::0.0625000 0.500000 0.00000
-DEAL::0.0625000 0.500000 0.0625000
-DEAL::0.0625000 0.500000 0.125000
-DEAL::0.0625000 0.500000 0.187500
-DEAL::0.0625000 0.500000 0.250000
-DEAL::0.0625000 0.500000 0.312500
-DEAL::0.0625000 0.500000 0.375000
-DEAL::0.0625000 0.500000 0.437500
-DEAL::0.0625000 0.500000 0.500000
-DEAL::0.0625000 0.500000 0.562500
-DEAL::0.0625000 0.500000 0.625000
-DEAL::0.0625000 0.500000 0.687500
-DEAL::0.0625000 0.500000 0.750000
-DEAL::0.0625000 0.500000 0.812500
-DEAL::0.0625000 0.500000 0.875000
-DEAL::0.0625000 0.500000 0.937500
-DEAL::0.0625000 0.500000 1.00000
-DEAL::0.125000 0.500000 0.00000
-DEAL::0.125000 0.500000 0.250000
-DEAL::0.125000 0.500000 0.500000
-DEAL::0.125000 0.500000 0.750000
-DEAL::0.125000 0.500000 1.00000
-DEAL::0.187500 0.500000 0.00000
-DEAL::0.187500 0.500000 0.0625000
-DEAL::0.187500 0.500000 0.125000
-DEAL::0.187500 0.500000 0.187500
-DEAL::0.187500 0.500000 0.250000
-DEAL::0.187500 0.500000 0.312500
-DEAL::0.187500 0.500000 0.375000
-DEAL::0.187500 0.500000 0.437500
-DEAL::0.187500 0.500000 0.500000
-DEAL::0.187500 0.500000 0.562500
-DEAL::0.187500 0.500000 0.625000
-DEAL::0.187500 0.500000 0.687500
-DEAL::0.187500 0.500000 0.750000
-DEAL::0.187500 0.500000 0.812500
-DEAL::0.187500 0.500000 0.875000
-DEAL::0.187500 0.500000 0.937500
-DEAL::0.187500 0.500000 1.00000
-DEAL::0.250000 0.500000 0.0625000
-DEAL::0.250000 0.500000 0.125000
-DEAL::0.250000 0.500000 0.187500
-DEAL::0.250000 0.500000 0.312500
-DEAL::0.250000 0.500000 0.375000
-DEAL::0.250000 0.500000 0.437500
-DEAL::0.250000 0.500000 0.562500
-DEAL::0.250000 0.500000 0.625000
-DEAL::0.250000 0.500000 0.687500
-DEAL::0.250000 0.500000 0.812500
-DEAL::0.250000 0.500000 0.875000
-DEAL::0.250000 0.500000 0.937500
-DEAL::0.312500 0.500000 0.00000
-DEAL::0.312500 0.500000 0.0625000
-DEAL::0.312500 0.500000 0.125000
-DEAL::0.312500 0.500000 0.187500
-DEAL::0.312500 0.500000 0.250000
-DEAL::0.312500 0.500000 0.312500
-DEAL::0.312500 0.500000 0.375000
-DEAL::0.312500 0.500000 0.437500
-DEAL::0.312500 0.500000 0.500000
-DEAL::0.312500 0.500000 0.562500
-DEAL::0.312500 0.500000 0.625000
-DEAL::0.312500 0.500000 0.687500
-DEAL::0.312500 0.500000 0.750000
-DEAL::0.312500 0.500000 0.812500
-DEAL::0.312500 0.500000 0.875000
-DEAL::0.312500 0.500000 0.937500
-DEAL::0.312500 0.500000 1.00000
-DEAL::0.375000 0.500000 0.00000
-DEAL::0.375000 0.500000 0.250000
-DEAL::0.375000 0.500000 0.500000
-DEAL::0.375000 0.500000 0.750000
-DEAL::0.375000 0.500000 1.00000
-DEAL::0.437500 0.500000 0.00000
-DEAL::0.437500 0.500000 0.0625000
-DEAL::0.437500 0.500000 0.125000
-DEAL::0.437500 0.500000 0.187500
-DEAL::0.437500 0.500000 0.250000
-DEAL::0.437500 0.500000 0.312500
-DEAL::0.437500 0.500000 0.375000
-DEAL::0.437500 0.500000 0.437500
-DEAL::0.437500 0.500000 0.500000
-DEAL::0.437500 0.500000 0.562500
-DEAL::0.437500 0.500000 0.625000
-DEAL::0.437500 0.500000 0.687500
-DEAL::0.437500 0.500000 0.750000
-DEAL::0.437500 0.500000 0.812500
-DEAL::0.437500 0.500000 0.875000
-DEAL::0.437500 0.500000 0.937500
-DEAL::0.437500 0.500000 1.00000
-DEAL::0.500000 0.00000 0.0625000
-DEAL::0.500000 0.00000 0.125000
-DEAL::0.500000 0.00000 0.187500
-DEAL::0.500000 0.00000 0.312500
-DEAL::0.500000 0.00000 0.375000
-DEAL::0.500000 0.00000 0.437500
-DEAL::0.500000 0.00000 0.562500
-DEAL::0.500000 0.00000 0.625000
-DEAL::0.500000 0.00000 0.687500
-DEAL::0.500000 0.00000 0.812500
-DEAL::0.500000 0.00000 0.875000
-DEAL::0.500000 0.00000 0.937500
-DEAL::0.500000 0.0625000 0.00000
-DEAL::0.500000 0.0625000 0.0625000
-DEAL::0.500000 0.0625000 0.187500
-DEAL::0.500000 0.0625000 0.250000
-DEAL::0.500000 0.0625000 0.312500
-DEAL::0.500000 0.0625000 0.437500
-DEAL::0.500000 0.0625000 0.562500
-DEAL::0.500000 0.0625000 0.687500
-DEAL::0.500000 0.0625000 0.750000
-DEAL::0.500000 0.0625000 0.812500
-DEAL::0.500000 0.0625000 0.937500
-DEAL::0.500000 0.0625000 1.00000
-DEAL::0.500000 0.0625000 0.500000
-DEAL::0.500000 0.125000 0.00000
-DEAL::0.500000 0.125000 0.0625000
-DEAL::0.500000 0.125000 0.187500
-DEAL::0.500000 0.125000 0.250000
-DEAL::0.500000 0.125000 0.312500
-DEAL::0.500000 0.125000 0.437500
-DEAL::0.500000 0.125000 0.562500
-DEAL::0.500000 0.125000 0.687500
-DEAL::0.500000 0.125000 0.750000
-DEAL::0.500000 0.125000 0.812500
-DEAL::0.500000 0.125000 0.937500
-DEAL::0.500000 0.125000 1.00000
-DEAL::0.500000 0.125000 0.500000
-DEAL::0.500000 0.187500 0.00000
-DEAL::0.500000 0.187500 0.0625000
-DEAL::0.500000 0.187500 0.187500
-DEAL::0.500000 0.187500 0.250000
-DEAL::0.500000 0.187500 0.312500
-DEAL::0.500000 0.187500 0.437500
-DEAL::0.500000 0.187500 0.687500
-DEAL::0.500000 0.187500 0.750000
-DEAL::0.500000 0.187500 0.812500
-DEAL::0.500000 0.187500 0.937500
-DEAL::0.500000 0.187500 1.00000
-DEAL::0.500000 0.187500 0.562500
-DEAL::0.500000 0.187500 0.500000
-DEAL::0.500000 0.250000 0.0625000
-DEAL::0.500000 0.250000 0.125000
-DEAL::0.500000 0.250000 0.187500
-DEAL::0.500000 0.250000 0.312500
-DEAL::0.500000 0.250000 0.375000
-DEAL::0.500000 0.250000 0.437500
-DEAL::0.500000 0.250000 0.562500
-DEAL::0.500000 0.250000 0.625000
-DEAL::0.500000 0.250000 0.687500
-DEAL::0.500000 0.250000 0.812500
-DEAL::0.500000 0.250000 0.875000
-DEAL::0.500000 0.250000 0.937500
-DEAL::0.500000 0.312500 0.00000
-DEAL::0.500000 0.312500 0.0625000
-DEAL::0.500000 0.312500 0.187500
-DEAL::0.500000 0.312500 0.250000
-DEAL::0.500000 0.312500 0.312500
-DEAL::0.500000 0.312500 0.437500
-DEAL::0.500000 0.312500 0.687500
-DEAL::0.500000 0.312500 0.750000
-DEAL::0.500000 0.312500 0.812500
-DEAL::0.500000 0.312500 0.937500
-DEAL::0.500000 0.312500 1.00000
-DEAL::0.500000 0.312500 0.562500
-DEAL::0.500000 0.312500 0.500000
-DEAL::0.500000 0.375000 0.00000
-DEAL::0.500000 0.375000 0.0625000
-DEAL::0.500000 0.375000 0.187500
-DEAL::0.500000 0.375000 0.250000
-DEAL::0.500000 0.375000 0.312500
-DEAL::0.500000 0.375000 0.687500
-DEAL::0.500000 0.375000 0.750000
-DEAL::0.500000 0.375000 0.812500
-DEAL::0.500000 0.375000 0.937500
-DEAL::0.500000 0.375000 1.00000
-DEAL::0.500000 0.375000 0.437500
-DEAL::0.500000 0.375000 0.562500
-DEAL::0.500000 0.375000 0.500000
-DEAL::0.500000 0.437500 0.00000
-DEAL::0.500000 0.437500 0.0625000
-DEAL::0.500000 0.437500 0.187500
-DEAL::0.500000 0.437500 0.250000
-DEAL::0.500000 0.437500 0.312500
-DEAL::0.500000 0.437500 0.437500
-DEAL::0.500000 0.437500 0.687500
-DEAL::0.500000 0.437500 0.750000
-DEAL::0.500000 0.437500 0.812500
-DEAL::0.500000 0.437500 0.937500
-DEAL::0.500000 0.437500 1.00000
-DEAL::0.500000 0.437500 0.562500
-DEAL::0.500000 0.437500 0.500000
-DEAL::0.500000 0.500000 0.0625000
-DEAL::0.500000 0.500000 0.125000
-DEAL::0.500000 0.500000 0.187500
-DEAL::0.500000 0.500000 0.312500
-DEAL::0.500000 0.500000 0.375000
-DEAL::0.500000 0.500000 0.437500
-DEAL::0.500000 0.500000 0.562500
-DEAL::0.500000 0.500000 0.625000
-DEAL::0.500000 0.500000 0.687500
-DEAL::0.500000 0.500000 0.812500
-DEAL::0.500000 0.500000 0.875000
-DEAL::0.500000 0.500000 0.937500
-DEAL::
+DEAL::OK
index 6feb20d8b2efee0924ee888f05259e56e8b2c08b..c3d49df7e0907994c49df1e6767154796a999f70 100644 (file)
@@ -1,285 +1,9 @@
 
 DEAL::Checking Q1 in 2d:
-DEAL::4
-DEAL::0.125000 0.500000
-DEAL::0.375000 0.500000
-DEAL::0.500000 0.125000
-DEAL::0.500000 0.375000
-DEAL::
+DEAL::OK
 DEAL::Checking Q1 in 3d:
-DEAL::40
-DEAL::0.00000 0.500000 0.125000
-DEAL::0.00000 0.500000 0.375000
-DEAL::0.00000 0.500000 0.625000
-DEAL::0.00000 0.500000 0.875000
-DEAL::0.125000 0.500000 0.00000
-DEAL::0.125000 0.500000 0.250000
-DEAL::0.125000 0.500000 0.500000
-DEAL::0.125000 0.500000 0.750000
-DEAL::0.125000 0.500000 1.00000
-DEAL::0.250000 0.500000 0.125000
-DEAL::0.250000 0.500000 0.375000
-DEAL::0.250000 0.500000 0.625000
-DEAL::0.250000 0.500000 0.875000
-DEAL::0.375000 0.500000 0.00000
-DEAL::0.375000 0.500000 0.250000
-DEAL::0.375000 0.500000 0.500000
-DEAL::0.375000 0.500000 0.750000
-DEAL::0.375000 0.500000 1.00000
-DEAL::0.500000 0.00000 0.125000
-DEAL::0.500000 0.00000 0.375000
-DEAL::0.500000 0.00000 0.625000
-DEAL::0.500000 0.00000 0.875000
-DEAL::0.500000 0.125000 0.00000
-DEAL::0.500000 0.125000 0.250000
-DEAL::0.500000 0.125000 0.750000
-DEAL::0.500000 0.125000 1.00000
-DEAL::0.500000 0.125000 0.500000
-DEAL::0.500000 0.250000 0.125000
-DEAL::0.500000 0.250000 0.375000
-DEAL::0.500000 0.250000 0.625000
-DEAL::0.500000 0.250000 0.875000
-DEAL::0.500000 0.375000 0.00000
-DEAL::0.500000 0.375000 0.250000
-DEAL::0.500000 0.375000 0.750000
-DEAL::0.500000 0.375000 1.00000
-DEAL::0.500000 0.375000 0.500000
-DEAL::0.500000 0.500000 0.125000
-DEAL::0.500000 0.500000 0.375000
-DEAL::0.500000 0.500000 0.625000
-DEAL::0.500000 0.500000 0.875000
-DEAL::
+DEAL::OK
 DEAL::Checking Q2 in 2d:
-DEAL::12
-DEAL::0.0625000 0.500000
-DEAL::0.125000 0.500000
-DEAL::0.187500 0.500000
-DEAL::0.312500 0.500000
-DEAL::0.375000 0.500000
-DEAL::0.437500 0.500000
-DEAL::0.500000 0.0625000
-DEAL::0.500000 0.125000
-DEAL::0.500000 0.187500
-DEAL::0.500000 0.312500
-DEAL::0.500000 0.375000
-DEAL::0.500000 0.437500
-DEAL::
+DEAL::OK
 DEAL::Checking Q2 in 3d:
-DEAL::216
-DEAL::0.00000 0.500000 0.0625000
-DEAL::0.00000 0.500000 0.125000
-DEAL::0.00000 0.500000 0.187500
-DEAL::0.00000 0.500000 0.312500
-DEAL::0.00000 0.500000 0.375000
-DEAL::0.00000 0.500000 0.437500
-DEAL::0.00000 0.500000 0.562500
-DEAL::0.00000 0.500000 0.625000
-DEAL::0.00000 0.500000 0.687500
-DEAL::0.00000 0.500000 0.812500
-DEAL::0.00000 0.500000 0.875000
-DEAL::0.00000 0.500000 0.937500
-DEAL::0.0625000 0.500000 0.00000
-DEAL::0.0625000 0.500000 0.0625000
-DEAL::0.0625000 0.500000 0.125000
-DEAL::0.0625000 0.500000 0.187500
-DEAL::0.0625000 0.500000 0.250000
-DEAL::0.0625000 0.500000 0.312500
-DEAL::0.0625000 0.500000 0.375000
-DEAL::0.0625000 0.500000 0.437500
-DEAL::0.0625000 0.500000 0.500000
-DEAL::0.0625000 0.500000 0.562500
-DEAL::0.0625000 0.500000 0.625000
-DEAL::0.0625000 0.500000 0.687500
-DEAL::0.0625000 0.500000 0.750000
-DEAL::0.0625000 0.500000 0.812500
-DEAL::0.0625000 0.500000 0.875000
-DEAL::0.0625000 0.500000 0.937500
-DEAL::0.0625000 0.500000 1.00000
-DEAL::0.125000 0.500000 0.00000
-DEAL::0.125000 0.500000 0.250000
-DEAL::0.125000 0.500000 0.500000
-DEAL::0.125000 0.500000 0.750000
-DEAL::0.125000 0.500000 1.00000
-DEAL::0.187500 0.500000 0.00000
-DEAL::0.187500 0.500000 0.0625000
-DEAL::0.187500 0.500000 0.125000
-DEAL::0.187500 0.500000 0.187500
-DEAL::0.187500 0.500000 0.250000
-DEAL::0.187500 0.500000 0.312500
-DEAL::0.187500 0.500000 0.375000
-DEAL::0.187500 0.500000 0.437500
-DEAL::0.187500 0.500000 0.500000
-DEAL::0.187500 0.500000 0.562500
-DEAL::0.187500 0.500000 0.625000
-DEAL::0.187500 0.500000 0.687500
-DEAL::0.187500 0.500000 0.750000
-DEAL::0.187500 0.500000 0.812500
-DEAL::0.187500 0.500000 0.875000
-DEAL::0.187500 0.500000 0.937500
-DEAL::0.187500 0.500000 1.00000
-DEAL::0.250000 0.500000 0.0625000
-DEAL::0.250000 0.500000 0.125000
-DEAL::0.250000 0.500000 0.187500
-DEAL::0.250000 0.500000 0.312500
-DEAL::0.250000 0.500000 0.375000
-DEAL::0.250000 0.500000 0.437500
-DEAL::0.250000 0.500000 0.562500
-DEAL::0.250000 0.500000 0.625000
-DEAL::0.250000 0.500000 0.687500
-DEAL::0.250000 0.500000 0.812500
-DEAL::0.250000 0.500000 0.875000
-DEAL::0.250000 0.500000 0.937500
-DEAL::0.312500 0.500000 0.00000
-DEAL::0.312500 0.500000 0.0625000
-DEAL::0.312500 0.500000 0.125000
-DEAL::0.312500 0.500000 0.187500
-DEAL::0.312500 0.500000 0.250000
-DEAL::0.312500 0.500000 0.312500
-DEAL::0.312500 0.500000 0.375000
-DEAL::0.312500 0.500000 0.437500
-DEAL::0.312500 0.500000 0.500000
-DEAL::0.312500 0.500000 0.562500
-DEAL::0.312500 0.500000 0.625000
-DEAL::0.312500 0.500000 0.687500
-DEAL::0.312500 0.500000 0.750000
-DEAL::0.312500 0.500000 0.812500
-DEAL::0.312500 0.500000 0.875000
-DEAL::0.312500 0.500000 0.937500
-DEAL::0.312500 0.500000 1.00000
-DEAL::0.375000 0.500000 0.00000
-DEAL::0.375000 0.500000 0.250000
-DEAL::0.375000 0.500000 0.500000
-DEAL::0.375000 0.500000 0.750000
-DEAL::0.375000 0.500000 1.00000
-DEAL::0.437500 0.500000 0.00000
-DEAL::0.437500 0.500000 0.0625000
-DEAL::0.437500 0.500000 0.125000
-DEAL::0.437500 0.500000 0.187500
-DEAL::0.437500 0.500000 0.250000
-DEAL::0.437500 0.500000 0.312500
-DEAL::0.437500 0.500000 0.375000
-DEAL::0.437500 0.500000 0.437500
-DEAL::0.437500 0.500000 0.500000
-DEAL::0.437500 0.500000 0.562500
-DEAL::0.437500 0.500000 0.625000
-DEAL::0.437500 0.500000 0.687500
-DEAL::0.437500 0.500000 0.750000
-DEAL::0.437500 0.500000 0.812500
-DEAL::0.437500 0.500000 0.875000
-DEAL::0.437500 0.500000 0.937500
-DEAL::0.437500 0.500000 1.00000
-DEAL::0.500000 0.00000 0.0625000
-DEAL::0.500000 0.00000 0.125000
-DEAL::0.500000 0.00000 0.187500
-DEAL::0.500000 0.00000 0.312500
-DEAL::0.500000 0.00000 0.375000
-DEAL::0.500000 0.00000 0.437500
-DEAL::0.500000 0.00000 0.562500
-DEAL::0.500000 0.00000 0.625000
-DEAL::0.500000 0.00000 0.687500
-DEAL::0.500000 0.00000 0.812500
-DEAL::0.500000 0.00000 0.875000
-DEAL::0.500000 0.00000 0.937500
-DEAL::0.500000 0.0625000 0.00000
-DEAL::0.500000 0.0625000 0.0625000
-DEAL::0.500000 0.0625000 0.187500
-DEAL::0.500000 0.0625000 0.250000
-DEAL::0.500000 0.0625000 0.312500
-DEAL::0.500000 0.0625000 0.437500
-DEAL::0.500000 0.0625000 0.562500
-DEAL::0.500000 0.0625000 0.687500
-DEAL::0.500000 0.0625000 0.750000
-DEAL::0.500000 0.0625000 0.812500
-DEAL::0.500000 0.0625000 0.937500
-DEAL::0.500000 0.0625000 1.00000
-DEAL::0.500000 0.0625000 0.500000
-DEAL::0.500000 0.125000 0.00000
-DEAL::0.500000 0.125000 0.0625000
-DEAL::0.500000 0.125000 0.187500
-DEAL::0.500000 0.125000 0.250000
-DEAL::0.500000 0.125000 0.312500
-DEAL::0.500000 0.125000 0.437500
-DEAL::0.500000 0.125000 0.562500
-DEAL::0.500000 0.125000 0.687500
-DEAL::0.500000 0.125000 0.750000
-DEAL::0.500000 0.125000 0.812500
-DEAL::0.500000 0.125000 0.937500
-DEAL::0.500000 0.125000 1.00000
-DEAL::0.500000 0.125000 0.500000
-DEAL::0.500000 0.187500 0.00000
-DEAL::0.500000 0.187500 0.0625000
-DEAL::0.500000 0.187500 0.187500
-DEAL::0.500000 0.187500 0.250000
-DEAL::0.500000 0.187500 0.312500
-DEAL::0.500000 0.187500 0.437500
-DEAL::0.500000 0.187500 0.687500
-DEAL::0.500000 0.187500 0.750000
-DEAL::0.500000 0.187500 0.812500
-DEAL::0.500000 0.187500 0.937500
-DEAL::0.500000 0.187500 1.00000
-DEAL::0.500000 0.187500 0.562500
-DEAL::0.500000 0.187500 0.500000
-DEAL::0.500000 0.250000 0.0625000
-DEAL::0.500000 0.250000 0.125000
-DEAL::0.500000 0.250000 0.187500
-DEAL::0.500000 0.250000 0.312500
-DEAL::0.500000 0.250000 0.375000
-DEAL::0.500000 0.250000 0.437500
-DEAL::0.500000 0.250000 0.562500
-DEAL::0.500000 0.250000 0.625000
-DEAL::0.500000 0.250000 0.687500
-DEAL::0.500000 0.250000 0.812500
-DEAL::0.500000 0.250000 0.875000
-DEAL::0.500000 0.250000 0.937500
-DEAL::0.500000 0.312500 0.00000
-DEAL::0.500000 0.312500 0.0625000
-DEAL::0.500000 0.312500 0.187500
-DEAL::0.500000 0.312500 0.250000
-DEAL::0.500000 0.312500 0.312500
-DEAL::0.500000 0.312500 0.437500
-DEAL::0.500000 0.312500 0.687500
-DEAL::0.500000 0.312500 0.750000
-DEAL::0.500000 0.312500 0.812500
-DEAL::0.500000 0.312500 0.937500
-DEAL::0.500000 0.312500 1.00000
-DEAL::0.500000 0.312500 0.562500
-DEAL::0.500000 0.312500 0.500000
-DEAL::0.500000 0.375000 0.00000
-DEAL::0.500000 0.375000 0.0625000
-DEAL::0.500000 0.375000 0.187500
-DEAL::0.500000 0.375000 0.250000
-DEAL::0.500000 0.375000 0.312500
-DEAL::0.500000 0.375000 0.687500
-DEAL::0.500000 0.375000 0.750000
-DEAL::0.500000 0.375000 0.812500
-DEAL::0.500000 0.375000 0.937500
-DEAL::0.500000 0.375000 1.00000
-DEAL::0.500000 0.375000 0.437500
-DEAL::0.500000 0.375000 0.562500
-DEAL::0.500000 0.375000 0.500000
-DEAL::0.500000 0.437500 0.00000
-DEAL::0.500000 0.437500 0.0625000
-DEAL::0.500000 0.437500 0.187500
-DEAL::0.500000 0.437500 0.250000
-DEAL::0.500000 0.437500 0.312500
-DEAL::0.500000 0.437500 0.437500
-DEAL::0.500000 0.437500 0.687500
-DEAL::0.500000 0.437500 0.750000
-DEAL::0.500000 0.437500 0.812500
-DEAL::0.500000 0.437500 0.937500
-DEAL::0.500000 0.437500 1.00000
-DEAL::0.500000 0.437500 0.562500
-DEAL::0.500000 0.437500 0.500000
-DEAL::0.500000 0.500000 0.0625000
-DEAL::0.500000 0.500000 0.125000
-DEAL::0.500000 0.500000 0.187500
-DEAL::0.500000 0.500000 0.312500
-DEAL::0.500000 0.500000 0.375000
-DEAL::0.500000 0.500000 0.437500
-DEAL::0.500000 0.500000 0.562500
-DEAL::0.500000 0.500000 0.625000
-DEAL::0.500000 0.500000 0.687500
-DEAL::0.500000 0.500000 0.812500
-DEAL::0.500000 0.500000 0.875000
-DEAL::0.500000 0.500000 0.937500
-DEAL::
+DEAL::OK

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.