]> https://gitweb.dealii.org/ - dealii.git/commitdiff
DataOutResample for multiple components 14258/head
authorPeter Munch <peterrmuench@gmail.com>
Mon, 12 Sep 2022 17:37:44 +0000 (19:37 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Wed, 14 Sep 2022 07:43:54 +0000 (09:43 +0200)
doc/news/changes/minor/20220912Munch_1 [new file with mode: 0644]
include/deal.II/numerics/vector_tools_evaluate.h
source/numerics/data_out_resample.cc
tests/remote_point_evaluation/data_out_resample_02.cc [new file with mode: 0644]
tests/remote_point_evaluation/data_out_resample_02.with_p4est=true.mpirun=1.output [new file with mode: 0644]
tests/remote_point_evaluation/data_out_resample_02.with_p4est=true.mpirun=4.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20220912Munch_1 b/doc/news/changes/minor/20220912Munch_1
new file mode 100644 (file)
index 0000000..81e4854
--- /dev/null
@@ -0,0 +1,6 @@
+Improved: The functions VectorTools::point_values()/::point_gradients() and 
+the class Utilities::MPI::RemotePointEvaluation now allow to specify the first
+component to be selected. The feature is used in the class DataOutResample
+to output multi-component solutions.
+<br>
+(Peter Munch, Magdalena Schreter, 2022/09/12)
index b8e7132d6bd71c5935917d773f90c64144cdc7cb..0d29e988ba6ba4febe700b689e27b8ca8308fc19 100644 (file)
@@ -111,6 +111,12 @@ namespace VectorTools
    *   cache, dof_handler_2, vector_2);
    * @endcode
    *
+   * The function also works with FiniteElement objects with multiple
+   * components. If one is interested only in a range of components, one can
+   * select these by the parameters @p first_selected_component and
+   * @p n_components. For further details on supported FiniteElement
+   * objects, see the documentation of FEPointEvaluation.
+   *
    * The function can also be used to evaluate cell-data vectors. For this
    * purpose, one passes in a Triangulation instead of a DoFHandler and a
    * vector of size Trinagulation::n_active_cells() or a vector, which
@@ -142,7 +148,8 @@ namespace VectorTools
     const VectorType &                                    vector,
     const std::vector<Point<spacedim>> &                  evaluation_points,
     Utilities::MPI::RemotePointEvaluation<dim, spacedim> &cache,
-    const EvaluationFlags::EvaluationFlags flags = EvaluationFlags::avg);
+    const EvaluationFlags::EvaluationFlags flags = EvaluationFlags::avg,
+    const unsigned int                     first_selected_component = 0);
 
   /**
    * Given a (distributed) solution vector @p vector, evaluate the values at
@@ -171,7 +178,8 @@ namespace VectorTools
     const Utilities::MPI::RemotePointEvaluation<dim, spacedim> &cache,
     const MeshType<dim, spacedim> &                             mesh,
     const VectorType &                                          vector,
-    const EvaluationFlags::EvaluationFlags flags = EvaluationFlags::avg);
+    const EvaluationFlags::EvaluationFlags flags = EvaluationFlags::avg,
+    const unsigned int                     first_selected_component = 0);
 
   /**
    * Given a (distributed) solution vector @p vector, evaluate the gradients at
@@ -200,7 +208,8 @@ namespace VectorTools
     const VectorType &                                    vector,
     const std::vector<Point<spacedim>> &                  evaluation_points,
     Utilities::MPI::RemotePointEvaluation<dim, spacedim> &cache,
-    const EvaluationFlags::EvaluationFlags flags = EvaluationFlags::avg);
+    const EvaluationFlags::EvaluationFlags flags = EvaluationFlags::avg,
+    const unsigned int                     first_selected_component = 0);
 
   /**
    * Given a (distributed) solution vector @p vector, evaluate the gradients at
@@ -228,7 +237,8 @@ namespace VectorTools
     const Utilities::MPI::RemotePointEvaluation<dim, spacedim> &cache,
     const MeshType<dim, spacedim> &                             mesh,
     const VectorType &                                          vector,
-    const EvaluationFlags::EvaluationFlags flags = EvaluationFlags::avg);
+    const EvaluationFlags::EvaluationFlags flags = EvaluationFlags::avg,
+    const unsigned int                     first_selected_component = 0);
 
 
 
@@ -252,11 +262,13 @@ namespace VectorTools
                const VectorType &                  vector,
                const std::vector<Point<spacedim>> &evaluation_points,
                Utilities::MPI::RemotePointEvaluation<dim, spacedim> &cache,
-               const EvaluationFlags::EvaluationFlags                flags)
+               const EvaluationFlags::EvaluationFlags                flags,
+               const unsigned int first_selected_component)
   {
     cache.reinit(evaluation_points, mesh.get_triangulation(), mapping);
 
-    return point_values<n_components>(cache, mesh, vector, flags);
+    return point_values<n_components>(
+      cache, mesh, vector, flags, first_selected_component);
   }
 
 
@@ -277,11 +289,13 @@ namespace VectorTools
                   const VectorType &                  vector,
                   const std::vector<Point<spacedim>> &evaluation_points,
                   Utilities::MPI::RemotePointEvaluation<dim, spacedim> &cache,
-                  const EvaluationFlags::EvaluationFlags                flags)
+                  const EvaluationFlags::EvaluationFlags                flags,
+                  const unsigned int first_selected_component)
   {
     cache.reinit(evaluation_points, mesh.get_triangulation(), mapping);
 
-    return point_gradients<n_components>(cache, mesh, vector, flags);
+    return point_gradients<n_components>(
+      cache, mesh, vector, flags, first_selected_component);
   }
 
 
@@ -396,6 +410,7 @@ namespace VectorTools
       const VectorType &                                          vector,
       const UpdateFlags                                           update_flags,
       const dealii::EvaluationFlags::EvaluationFlags evaluation_flags,
+      const unsigned int                             first_selected_component,
       const std::function<
         value_type(const FEPointEvaluation<n_components,
                                            dim,
@@ -438,7 +453,10 @@ namespace VectorTools
                                              dim,
                                              spacedim,
                                              typename VectorType::value_type>>(
-            cache.get_mapping(), cell->get_fe(), update_flags);
+            cache.get_mapping(),
+            cell->get_fe(),
+            update_flags,
+            first_selected_component);
       auto &evaluator = *evaluators[cell->active_fe_index()];
 
       evaluator.reinit(cell, unit_points);
@@ -532,6 +550,7 @@ namespace VectorTools
       const VectorType &                  vector,
       const UpdateFlags,
       const dealii::EvaluationFlags::EvaluationFlags evaluation_flags,
+      const unsigned int                             first_selected_component,
       const std::function<
         value_type(const FEPointEvaluation<n_components,
                                            dim,
@@ -547,8 +566,9 @@ namespace VectorTools
                                           typename VectorType::value_type>>> &)
     {
       (void)evaluation_flags;
+      (void)first_selected_component;
 
-      Assert(n_components == 1,
+      Assert(n_components == 1 && first_selected_component == 0,
              ExcMessage(
                "A cell-data vector can only have a single component."));
 
@@ -581,7 +601,8 @@ namespace VectorTools
       const MeshType &                                            mesh,
       const VectorType &                                          vector,
       const EvaluationFlags::EvaluationFlags                      flags,
-      const UpdateFlags                                           update_flags,
+      const unsigned int                             first_selected_component,
+      const UpdateFlags                              update_flags,
       const dealii::EvaluationFlags::EvaluationFlags evaluation_flags,
       const std::function<
         value_type(const FEPointEvaluation<n_components,
@@ -627,6 +648,7 @@ namespace VectorTools
               vector,
               update_flags,
               evaluation_flags,
+              first_selected_component,
               process_quadrature_point,
               values,
               solution_values,
@@ -688,7 +710,8 @@ namespace VectorTools
     const Utilities::MPI::RemotePointEvaluation<dim, spacedim> &cache,
     const MeshType<dim, spacedim> &                             mesh,
     const VectorType &                                          vector,
-    const EvaluationFlags::EvaluationFlags                      flags)
+    const EvaluationFlags::EvaluationFlags                      flags,
+    const unsigned int first_selected_component)
   {
     return internal::evaluate_at_points<
       n_components,
@@ -704,6 +727,7 @@ namespace VectorTools
       mesh,
       vector,
       flags,
+      first_selected_component,
       update_values,
       dealii::EvaluationFlags::values,
       [](const auto &evaluator, const auto &q) {
@@ -726,7 +750,8 @@ namespace VectorTools
     const Utilities::MPI::RemotePointEvaluation<dim, spacedim> &cache,
     const MeshType<dim, spacedim> &                             mesh,
     const VectorType &                                          vector,
-    const EvaluationFlags::EvaluationFlags                      flags)
+    const EvaluationFlags::EvaluationFlags                      flags,
+    const unsigned int first_selected_component)
   {
     return internal::evaluate_at_points<
       n_components,
@@ -743,6 +768,7 @@ namespace VectorTools
       mesh,
       vector,
       flags,
+      first_selected_component,
       update_gradients,
       dealii::EvaluationFlags::gradients,
       [](const auto &evaluator, const unsigned &q) {
index 156a7ffbc51a2b8c728b3b52c8528124f4892b09..1c7a35d0985e7b0ff633b7adfe36e18624886331 100644 (file)
@@ -151,32 +151,42 @@ DataOutResample<dim, patch_dim, spacedim>::build_patches(
 
       const auto &dh = *data_ptr->dof_handler;
 
-      // TODO: enable more components
-      AssertDimension(dh.get_fe_collection().n_components(), 1);
-
-      const auto values =
-        VectorTools::point_values<1>(rpe, dh, data_ptr->vector);
-
-      vectors.emplace_back(
-        std::make_shared<LinearAlgebra::distributed::Vector<double>>(
-          partitioner));
-
-      for (unsigned int j = 0; j < values.size(); ++j)
-        vectors.back()->local_element(point_to_local_vector_indices[j]) =
-          values[j];
-
-      vectors.back()->set_ghost_state(true);
-
-      // we can give the vectors arbitrary names ("temp_*") here, since these
-      // are only used internally (by patch_data_out) but not later on during
-      // the actual output to file
-      patch_data_out.add_data_vector(
-        *vectors.back(),
-        std::string("temp_" + std::to_string(counter)),
-        DataOut_DoFData<patch_dim, patch_dim, spacedim, spacedim>::
-          DataVectorType::type_dof_data);
-
-      counter++;
+#ifdef DEBUG
+      for (const auto &fe : dh.get_fe_collection())
+        Assert(
+          fe.n_base_elements() == 1,
+          ExcMessage(
+            "This class currently only supports scalar elements and elements "
+            "with a single base element."));
+#endif
+
+      for (unsigned int comp = 0; comp < dh.get_fe_collection().n_components();
+           ++comp)
+        {
+          const auto values = VectorTools::point_values<1>(
+            rpe, dh, data_ptr->vector, VectorTools::EvaluationFlags::avg, comp);
+
+          vectors.emplace_back(
+            std::make_shared<LinearAlgebra::distributed::Vector<double>>(
+              partitioner));
+
+          for (unsigned int j = 0; j < values.size(); ++j)
+            vectors.back()->local_element(point_to_local_vector_indices[j]) =
+              values[j];
+
+          vectors.back()->set_ghost_state(true);
+
+          // we can give the vectors arbitrary names ("temp_*") here, since
+          // these are only used internally (by patch_data_out) but not later on
+          // during the actual output to file
+          patch_data_out.add_data_vector(
+            *vectors.back(),
+            std::string("temp_" + std::to_string(counter)),
+            DataOut_DoFData<patch_dim, patch_dim, spacedim, spacedim>::
+              DataVectorType::type_dof_data);
+
+          counter++;
+        }
     }
 
   patch_data_out.build_patches(*patch_mapping,
diff --git a/tests/remote_point_evaluation/data_out_resample_02.cc b/tests/remote_point_evaluation/data_out_resample_02.cc
new file mode 100644 (file)
index 0000000..747197f
--- /dev/null
@@ -0,0 +1,124 @@
+// ---------------------------------------------------------------------
+//
+// 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 DataOutResample to create a slice for FESystem with multiple
+// components.
+
+#include <deal.II/base/conditional_ostream.h>
+#include <deal.II/base/mpi_remote_point_evaluation.h>
+
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/dofs/dof_tools.h>
+
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_system.h>
+#include <deal.II/fe/mapping_q1.h>
+
+#include <deal.II/grid/grid_generator.h>
+
+#include <deal.II/numerics/data_out_resample.h>
+#include <deal.II/numerics/vector_tools.h>
+
+#include "../tests.h"
+
+using namespace dealii;
+
+template <int dim>
+class AnalyticalFunction : public Function<dim>
+{
+public:
+  AnalyticalFunction(const unsigned n_components)
+    : Function<dim>(n_components)
+  {}
+
+  virtual double
+  value(const Point<dim> &p, const unsigned int component = 0) const
+  {
+    return p[component] * (component + 1.0) + component;
+  }
+};
+
+
+template <int dim, int spacedim>
+std::shared_ptr<const Utilities::MPI::Partitioner>
+create_partitioner(const DoFHandler<dim, spacedim> &dof_handler)
+{
+  IndexSet locally_relevant_dofs;
+
+  DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant_dofs);
+
+  return std::make_shared<const Utilities::MPI::Partitioner>(
+    dof_handler.locally_owned_dofs(),
+    locally_relevant_dofs,
+    dof_handler.get_communicator());
+}
+
+
+
+int
+main(int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
+  MPILogInitAll                    all;
+
+  const int dim          = 3;
+  const int patch_dim    = 2;
+  const int spacedim     = 3;
+  const int n_components = dim;
+
+  const unsigned int n_refinements_1 = 3;
+  const unsigned int n_refinements_2 = 3;
+  const MPI_Comm     comm            = MPI_COMM_WORLD;
+
+  parallel::distributed::Triangulation<patch_dim, spacedim> tria_slice(comm);
+  GridGenerator::hyper_cube(tria_slice, -1.0, +1.0);
+  tria_slice.refine_global(n_refinements_2);
+
+  MappingQ1<patch_dim, spacedim> mapping_slice;
+
+  parallel::distributed::Triangulation<dim, spacedim> tria_backround(comm);
+  GridGenerator::hyper_cube(tria_backround, -1.0, +1.0);
+  tria_backround.refine_global(n_refinements_1);
+
+  DoFHandler<dim, spacedim> dof_handler(tria_backround);
+  dof_handler.distribute_dofs(
+    FESystem<dim, spacedim>(FE_Q<dim, spacedim>{1}, n_components));
+
+  MappingQ1<dim, spacedim> mapping;
+
+  LinearAlgebra::distributed::Vector<double> vector(
+    create_partitioner(dof_handler));
+
+  VectorTools::interpolate(mapping,
+                           dof_handler,
+                           AnalyticalFunction<dim>(n_components),
+                           vector);
+
+  vector.update_ghost_values();
+
+  DataOutResample<dim, patch_dim, spacedim> data_out(tria_slice, mapping_slice);
+  data_out.add_data_vector(dof_handler, vector, "solution_0");
+  data_out.add_data_vector(dof_handler, vector, "solution_1");
+  data_out.update_mapping(mapping);
+  data_out.build_patches();
+
+#if 1
+  data_out.write_vtk(deallog.get_file_stream());
+#else
+  data_out.write_vtu_with_pvtu_record("./", "data_out_01", 0, comm, 1, 1);
+#endif
+}
diff --git a/tests/remote_point_evaluation/data_out_resample_02.with_p4est=true.mpirun=1.output b/tests/remote_point_evaluation/data_out_resample_02.with_p4est=true.mpirun=1.output
new file mode 100644 (file)
index 0000000..f7cbd89
--- /dev/null
@@ -0,0 +1,351 @@
+
+# vtk DataFile Version 3.0
+#This file was generated 
+ASCII
+DATASET UNSTRUCTURED_GRID
+
+POINTS 256 double
+-1.00000 -1.00000 0.00000
+-0.750000 -1.00000 0.00000
+-1.00000 -0.750000 0.00000
+-0.750000 -0.750000 0.00000
+-0.750000 -1.00000 0.00000
+-0.500000 -1.00000 0.00000
+-0.750000 -0.750000 0.00000
+-0.500000 -0.750000 0.00000
+-1.00000 -0.750000 0.00000
+-0.750000 -0.750000 0.00000
+-1.00000 -0.500000 0.00000
+-0.750000 -0.500000 0.00000
+-0.750000 -0.750000 0.00000
+-0.500000 -0.750000 0.00000
+-0.750000 -0.500000 0.00000
+-0.500000 -0.500000 0.00000
+-0.500000 -1.00000 0.00000
+-0.250000 -1.00000 0.00000
+-0.500000 -0.750000 0.00000
+-0.250000 -0.750000 0.00000
+-0.250000 -1.00000 0.00000
+0.00000 -1.00000 0.00000
+-0.250000 -0.750000 0.00000
+0.00000 -0.750000 0.00000
+-0.500000 -0.750000 0.00000
+-0.250000 -0.750000 0.00000
+-0.500000 -0.500000 0.00000
+-0.250000 -0.500000 0.00000
+-0.250000 -0.750000 0.00000
+0.00000 -0.750000 0.00000
+-0.250000 -0.500000 0.00000
+0.00000 -0.500000 0.00000
+-1.00000 -0.500000 0.00000
+-0.750000 -0.500000 0.00000
+-1.00000 -0.250000 0.00000
+-0.750000 -0.250000 0.00000
+-0.750000 -0.500000 0.00000
+-0.500000 -0.500000 0.00000
+-0.750000 -0.250000 0.00000
+-0.500000 -0.250000 0.00000
+-1.00000 -0.250000 0.00000
+-0.750000 -0.250000 0.00000
+-1.00000 0.00000 0.00000
+-0.750000 0.00000 0.00000
+-0.750000 -0.250000 0.00000
+-0.500000 -0.250000 0.00000
+-0.750000 0.00000 0.00000
+-0.500000 0.00000 0.00000
+-0.500000 -0.500000 0.00000
+-0.250000 -0.500000 0.00000
+-0.500000 -0.250000 0.00000
+-0.250000 -0.250000 0.00000
+-0.250000 -0.500000 0.00000
+0.00000 -0.500000 0.00000
+-0.250000 -0.250000 0.00000
+0.00000 -0.250000 0.00000
+-0.500000 -0.250000 0.00000
+-0.250000 -0.250000 0.00000
+-0.500000 0.00000 0.00000
+-0.250000 0.00000 0.00000
+-0.250000 -0.250000 0.00000
+0.00000 -0.250000 0.00000
+-0.250000 0.00000 0.00000
+0.00000 0.00000 0.00000
+0.00000 -1.00000 0.00000
+0.250000 -1.00000 0.00000
+0.00000 -0.750000 0.00000
+0.250000 -0.750000 0.00000
+0.250000 -1.00000 0.00000
+0.500000 -1.00000 0.00000
+0.250000 -0.750000 0.00000
+0.500000 -0.750000 0.00000
+0.00000 -0.750000 0.00000
+0.250000 -0.750000 0.00000
+0.00000 -0.500000 0.00000
+0.250000 -0.500000 0.00000
+0.250000 -0.750000 0.00000
+0.500000 -0.750000 0.00000
+0.250000 -0.500000 0.00000
+0.500000 -0.500000 0.00000
+0.500000 -1.00000 0.00000
+0.750000 -1.00000 0.00000
+0.500000 -0.750000 0.00000
+0.750000 -0.750000 0.00000
+0.750000 -1.00000 0.00000
+1.00000 -1.00000 0.00000
+0.750000 -0.750000 0.00000
+1.00000 -0.750000 0.00000
+0.500000 -0.750000 0.00000
+0.750000 -0.750000 0.00000
+0.500000 -0.500000 0.00000
+0.750000 -0.500000 0.00000
+0.750000 -0.750000 0.00000
+1.00000 -0.750000 0.00000
+0.750000 -0.500000 0.00000
+1.00000 -0.500000 0.00000
+0.00000 -0.500000 0.00000
+0.250000 -0.500000 0.00000
+0.00000 -0.250000 0.00000
+0.250000 -0.250000 0.00000
+0.250000 -0.500000 0.00000
+0.500000 -0.500000 0.00000
+0.250000 -0.250000 0.00000
+0.500000 -0.250000 0.00000
+0.00000 -0.250000 0.00000
+0.250000 -0.250000 0.00000
+0.00000 0.00000 0.00000
+0.250000 0.00000 0.00000
+0.250000 -0.250000 0.00000
+0.500000 -0.250000 0.00000
+0.250000 0.00000 0.00000
+0.500000 0.00000 0.00000
+0.500000 -0.500000 0.00000
+0.750000 -0.500000 0.00000
+0.500000 -0.250000 0.00000
+0.750000 -0.250000 0.00000
+0.750000 -0.500000 0.00000
+1.00000 -0.500000 0.00000
+0.750000 -0.250000 0.00000
+1.00000 -0.250000 0.00000
+0.500000 -0.250000 0.00000
+0.750000 -0.250000 0.00000
+0.500000 0.00000 0.00000
+0.750000 0.00000 0.00000
+0.750000 -0.250000 0.00000
+1.00000 -0.250000 0.00000
+0.750000 0.00000 0.00000
+1.00000 0.00000 0.00000
+-1.00000 0.00000 0.00000
+-0.750000 0.00000 0.00000
+-1.00000 0.250000 0.00000
+-0.750000 0.250000 0.00000
+-0.750000 0.00000 0.00000
+-0.500000 0.00000 0.00000
+-0.750000 0.250000 0.00000
+-0.500000 0.250000 0.00000
+-1.00000 0.250000 0.00000
+-0.750000 0.250000 0.00000
+-1.00000 0.500000 0.00000
+-0.750000 0.500000 0.00000
+-0.750000 0.250000 0.00000
+-0.500000 0.250000 0.00000
+-0.750000 0.500000 0.00000
+-0.500000 0.500000 0.00000
+-0.500000 0.00000 0.00000
+-0.250000 0.00000 0.00000
+-0.500000 0.250000 0.00000
+-0.250000 0.250000 0.00000
+-0.250000 0.00000 0.00000
+0.00000 0.00000 0.00000
+-0.250000 0.250000 0.00000
+0.00000 0.250000 0.00000
+-0.500000 0.250000 0.00000
+-0.250000 0.250000 0.00000
+-0.500000 0.500000 0.00000
+-0.250000 0.500000 0.00000
+-0.250000 0.250000 0.00000
+0.00000 0.250000 0.00000
+-0.250000 0.500000 0.00000
+0.00000 0.500000 0.00000
+-1.00000 0.500000 0.00000
+-0.750000 0.500000 0.00000
+-1.00000 0.750000 0.00000
+-0.750000 0.750000 0.00000
+-0.750000 0.500000 0.00000
+-0.500000 0.500000 0.00000
+-0.750000 0.750000 0.00000
+-0.500000 0.750000 0.00000
+-1.00000 0.750000 0.00000
+-0.750000 0.750000 0.00000
+-1.00000 1.00000 0.00000
+-0.750000 1.00000 0.00000
+-0.750000 0.750000 0.00000
+-0.500000 0.750000 0.00000
+-0.750000 1.00000 0.00000
+-0.500000 1.00000 0.00000
+-0.500000 0.500000 0.00000
+-0.250000 0.500000 0.00000
+-0.500000 0.750000 0.00000
+-0.250000 0.750000 0.00000
+-0.250000 0.500000 0.00000
+0.00000 0.500000 0.00000
+-0.250000 0.750000 0.00000
+0.00000 0.750000 0.00000
+-0.500000 0.750000 0.00000
+-0.250000 0.750000 0.00000
+-0.500000 1.00000 0.00000
+-0.250000 1.00000 0.00000
+-0.250000 0.750000 0.00000
+0.00000 0.750000 0.00000
+-0.250000 1.00000 0.00000
+0.00000 1.00000 0.00000
+0.00000 0.00000 0.00000
+0.250000 0.00000 0.00000
+0.00000 0.250000 0.00000
+0.250000 0.250000 0.00000
+0.250000 0.00000 0.00000
+0.500000 0.00000 0.00000
+0.250000 0.250000 0.00000
+0.500000 0.250000 0.00000
+0.00000 0.250000 0.00000
+0.250000 0.250000 0.00000
+0.00000 0.500000 0.00000
+0.250000 0.500000 0.00000
+0.250000 0.250000 0.00000
+0.500000 0.250000 0.00000
+0.250000 0.500000 0.00000
+0.500000 0.500000 0.00000
+0.500000 0.00000 0.00000
+0.750000 0.00000 0.00000
+0.500000 0.250000 0.00000
+0.750000 0.250000 0.00000
+0.750000 0.00000 0.00000
+1.00000 0.00000 0.00000
+0.750000 0.250000 0.00000
+1.00000 0.250000 0.00000
+0.500000 0.250000 0.00000
+0.750000 0.250000 0.00000
+0.500000 0.500000 0.00000
+0.750000 0.500000 0.00000
+0.750000 0.250000 0.00000
+1.00000 0.250000 0.00000
+0.750000 0.500000 0.00000
+1.00000 0.500000 0.00000
+0.00000 0.500000 0.00000
+0.250000 0.500000 0.00000
+0.00000 0.750000 0.00000
+0.250000 0.750000 0.00000
+0.250000 0.500000 0.00000
+0.500000 0.500000 0.00000
+0.250000 0.750000 0.00000
+0.500000 0.750000 0.00000
+0.00000 0.750000 0.00000
+0.250000 0.750000 0.00000
+0.00000 1.00000 0.00000
+0.250000 1.00000 0.00000
+0.250000 0.750000 0.00000
+0.500000 0.750000 0.00000
+0.250000 1.00000 0.00000
+0.500000 1.00000 0.00000
+0.500000 0.500000 0.00000
+0.750000 0.500000 0.00000
+0.500000 0.750000 0.00000
+0.750000 0.750000 0.00000
+0.750000 0.500000 0.00000
+1.00000 0.500000 0.00000
+0.750000 0.750000 0.00000
+1.00000 0.750000 0.00000
+0.500000 0.750000 0.00000
+0.750000 0.750000 0.00000
+0.500000 1.00000 0.00000
+0.750000 1.00000 0.00000
+0.750000 0.750000 0.00000
+1.00000 0.750000 0.00000
+0.750000 1.00000 0.00000
+1.00000 1.00000 0.00000
+
+CELLS 64 320
+4      0       1       3       2
+4      4       5       7       6
+4      8       9       11      10
+4      12      13      15      14
+4      16      17      19      18
+4      20      21      23      22
+4      24      25      27      26
+4      28      29      31      30
+4      32      33      35      34
+4      36      37      39      38
+4      40      41      43      42
+4      44      45      47      46
+4      48      49      51      50
+4      52      53      55      54
+4      56      57      59      58
+4      60      61      63      62
+4      64      65      67      66
+4      68      69      71      70
+4      72      73      75      74
+4      76      77      79      78
+4      80      81      83      82
+4      84      85      87      86
+4      88      89      91      90
+4      92      93      95      94
+4      96      97      99      98
+4      100     101     103     102
+4      104     105     107     106
+4      108     109     111     110
+4      112     113     115     114
+4      116     117     119     118
+4      120     121     123     122
+4      124     125     127     126
+4      128     129     131     130
+4      132     133     135     134
+4      136     137     139     138
+4      140     141     143     142
+4      144     145     147     146
+4      148     149     151     150
+4      152     153     155     154
+4      156     157     159     158
+4      160     161     163     162
+4      164     165     167     166
+4      168     169     171     170
+4      172     173     175     174
+4      176     177     179     178
+4      180     181     183     182
+4      184     185     187     186
+4      188     189     191     190
+4      192     193     195     194
+4      196     197     199     198
+4      200     201     203     202
+4      204     205     207     206
+4      208     209     211     210
+4      212     213     215     214
+4      216     217     219     218
+4      220     221     223     222
+4      224     225     227     226
+4      228     229     231     230
+4      232     233     235     234
+4      236     237     239     238
+4      240     241     243     242
+4      244     245     247     246
+4      248     249     251     250
+4      252     253     255     254
+
+CELL_TYPES 64
+ 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9
+POINT_DATA 256
+SCALARS solution_0_0 double 1
+LOOKUP_TABLE default
+-1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 
+SCALARS solution_0_1 double 1
+LOOKUP_TABLE default
+-1.00000 -1.00000 -0.500000 -0.500000 -1.00000 -1.00000 -0.500000 -0.500000 -0.500000 -0.500000 0.00000 0.00000 -0.500000 -0.500000 0.00000 0.00000 -1.00000 -1.00000 -0.500000 -0.500000 -1.00000 -1.00000 -0.500000 -0.500000 -0.500000 -0.500000 0.00000 0.00000 -0.500000 -0.500000 0.00000 0.00000 0.00000 0.00000 0.500000 0.500000 0.00000 0.00000 0.500000 0.500000 0.500000 0.500000 1.00000 1.00000 0.500000 0.500000 1.00000 1.00000 0.00000 0.00000 0.500000 0.500000 0.00000 0.00000 0.500000 0.500000 0.500000 0.500000 1.00000 1.00000 0.500000 0.500000 1.00000 1.00000 -1.00000 -1.00000 -0.500000 -0.500000 -1.00000 -1.00000 -0.500000 -0.500000 -0.500000 -0.500000 0.00000 0.00000 -0.500000 -0.500000 0.00000 0.00000 -1.00000 -1.00000 -0.500000 -0.500000 -1.00000 -1.00000 -0.500000 -0.500000 -0.500000 -0.500000 0.00000 0.00000 -0.500000 -0.500000 0.00000 0.00000 0.00000 0.00000 0.500000 0.500000 0.00000 0.00000 0.500000 0.500000 0.500000 0.500000 1.00000 1.00000 0.500000 0.500000 1.00000 1.00000 0.00000 0.00000 0.500000 0.500000 0.00000 0.00000 0.500000 0.500000 0.500000 0.500000 1.00000 1.00000 0.500000 0.500000 1.00000 1.00000 1.00000 1.00000 1.50000 1.50000 1.00000 1.00000 1.50000 1.50000 1.50000 1.50000 2.00000 2.00000 1.50000 1.50000 2.00000 2.00000 1.00000 1.00000 1.50000 1.50000 1.00000 1.00000 1.50000 1.50000 1.50000 1.50000 2.00000 2.00000 1.50000 1.50000 2.00000 2.00000 2.00000 2.00000 2.50000 2.50000 2.00000 2.00000 2.50000 2.50000 2.50000 2.50000 3.00000 3.00000 2.50000 2.50000 3.00000 3.00000 2.00000 2.00000 2.50000 2.50000 2.00000 2.00000 2.50000 2.50000 2.50000 2.50000 3.00000 3.00000 2.50000 2.50000 3.00000 3.00000 1.00000 1.00000 1.50000 1.50000 1.00000 1.00000 1.50000 1.50000 1.50000 1.50000 2.00000 2.00000 1.50000 1.50000 2.00000 2.00000 1.00000 1.00000 1.50000 1.50000 1.00000 1.00000 1.50000 1.50000 1.50000 1.50000 2.00000 2.00000 1.50000 1.50000 2.00000 2.00000 2.00000 2.00000 2.50000 2.50000 2.00000 2.00000 2.50000 2.50000 2.50000 2.50000 3.00000 3.00000 2.50000 2.50000 3.00000 3.00000 2.00000 2.00000 2.50000 2.50000 2.00000 2.00000 2.50000 2.50000 2.50000 2.50000 3.00000 3.00000 2.50000 2.50000 3.00000 3.00000 
+SCALARS solution_0_2 double 1
+LOOKUP_TABLE default
+2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 
+SCALARS solution_1_0 double 1
+LOOKUP_TABLE default
+-1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 
+SCALARS solution_1_1 double 1
+LOOKUP_TABLE default
+-1.00000 -1.00000 -0.500000 -0.500000 -1.00000 -1.00000 -0.500000 -0.500000 -0.500000 -0.500000 0.00000 0.00000 -0.500000 -0.500000 0.00000 0.00000 -1.00000 -1.00000 -0.500000 -0.500000 -1.00000 -1.00000 -0.500000 -0.500000 -0.500000 -0.500000 0.00000 0.00000 -0.500000 -0.500000 0.00000 0.00000 0.00000 0.00000 0.500000 0.500000 0.00000 0.00000 0.500000 0.500000 0.500000 0.500000 1.00000 1.00000 0.500000 0.500000 1.00000 1.00000 0.00000 0.00000 0.500000 0.500000 0.00000 0.00000 0.500000 0.500000 0.500000 0.500000 1.00000 1.00000 0.500000 0.500000 1.00000 1.00000 -1.00000 -1.00000 -0.500000 -0.500000 -1.00000 -1.00000 -0.500000 -0.500000 -0.500000 -0.500000 0.00000 0.00000 -0.500000 -0.500000 0.00000 0.00000 -1.00000 -1.00000 -0.500000 -0.500000 -1.00000 -1.00000 -0.500000 -0.500000 -0.500000 -0.500000 0.00000 0.00000 -0.500000 -0.500000 0.00000 0.00000 0.00000 0.00000 0.500000 0.500000 0.00000 0.00000 0.500000 0.500000 0.500000 0.500000 1.00000 1.00000 0.500000 0.500000 1.00000 1.00000 0.00000 0.00000 0.500000 0.500000 0.00000 0.00000 0.500000 0.500000 0.500000 0.500000 1.00000 1.00000 0.500000 0.500000 1.00000 1.00000 1.00000 1.00000 1.50000 1.50000 1.00000 1.00000 1.50000 1.50000 1.50000 1.50000 2.00000 2.00000 1.50000 1.50000 2.00000 2.00000 1.00000 1.00000 1.50000 1.50000 1.00000 1.00000 1.50000 1.50000 1.50000 1.50000 2.00000 2.00000 1.50000 1.50000 2.00000 2.00000 2.00000 2.00000 2.50000 2.50000 2.00000 2.00000 2.50000 2.50000 2.50000 2.50000 3.00000 3.00000 2.50000 2.50000 3.00000 3.00000 2.00000 2.00000 2.50000 2.50000 2.00000 2.00000 2.50000 2.50000 2.50000 2.50000 3.00000 3.00000 2.50000 2.50000 3.00000 3.00000 1.00000 1.00000 1.50000 1.50000 1.00000 1.00000 1.50000 1.50000 1.50000 1.50000 2.00000 2.00000 1.50000 1.50000 2.00000 2.00000 1.00000 1.00000 1.50000 1.50000 1.00000 1.00000 1.50000 1.50000 1.50000 1.50000 2.00000 2.00000 1.50000 1.50000 2.00000 2.00000 2.00000 2.00000 2.50000 2.50000 2.00000 2.00000 2.50000 2.50000 2.50000 2.50000 3.00000 3.00000 2.50000 2.50000 3.00000 3.00000 2.00000 2.00000 2.50000 2.50000 2.00000 2.00000 2.50000 2.50000 2.50000 2.50000 3.00000 3.00000 2.50000 2.50000 3.00000 3.00000 
+SCALARS solution_1_2 double 1
+LOOKUP_TABLE default
+2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 
diff --git a/tests/remote_point_evaluation/data_out_resample_02.with_p4est=true.mpirun=4.output b/tests/remote_point_evaluation/data_out_resample_02.with_p4est=true.mpirun=4.output
new file mode 100644 (file)
index 0000000..aae28fc
--- /dev/null
@@ -0,0 +1,447 @@
+
+# vtk DataFile Version 3.0
+#This file was generated 
+ASCII
+DATASET UNSTRUCTURED_GRID
+
+POINTS 64 double
+-1.00000 -1.00000 0.00000
+-0.750000 -1.00000 0.00000
+-1.00000 -0.750000 0.00000
+-0.750000 -0.750000 0.00000
+-0.750000 -1.00000 0.00000
+-0.500000 -1.00000 0.00000
+-0.750000 -0.750000 0.00000
+-0.500000 -0.750000 0.00000
+-1.00000 -0.750000 0.00000
+-0.750000 -0.750000 0.00000
+-1.00000 -0.500000 0.00000
+-0.750000 -0.500000 0.00000
+-0.750000 -0.750000 0.00000
+-0.500000 -0.750000 0.00000
+-0.750000 -0.500000 0.00000
+-0.500000 -0.500000 0.00000
+-0.500000 -1.00000 0.00000
+-0.250000 -1.00000 0.00000
+-0.500000 -0.750000 0.00000
+-0.250000 -0.750000 0.00000
+-0.250000 -1.00000 0.00000
+0.00000 -1.00000 0.00000
+-0.250000 -0.750000 0.00000
+0.00000 -0.750000 0.00000
+-0.500000 -0.750000 0.00000
+-0.250000 -0.750000 0.00000
+-0.500000 -0.500000 0.00000
+-0.250000 -0.500000 0.00000
+-0.250000 -0.750000 0.00000
+0.00000 -0.750000 0.00000
+-0.250000 -0.500000 0.00000
+0.00000 -0.500000 0.00000
+-1.00000 -0.500000 0.00000
+-0.750000 -0.500000 0.00000
+-1.00000 -0.250000 0.00000
+-0.750000 -0.250000 0.00000
+-0.750000 -0.500000 0.00000
+-0.500000 -0.500000 0.00000
+-0.750000 -0.250000 0.00000
+-0.500000 -0.250000 0.00000
+-1.00000 -0.250000 0.00000
+-0.750000 -0.250000 0.00000
+-1.00000 0.00000 0.00000
+-0.750000 0.00000 0.00000
+-0.750000 -0.250000 0.00000
+-0.500000 -0.250000 0.00000
+-0.750000 0.00000 0.00000
+-0.500000 0.00000 0.00000
+-0.500000 -0.500000 0.00000
+-0.250000 -0.500000 0.00000
+-0.500000 -0.250000 0.00000
+-0.250000 -0.250000 0.00000
+-0.250000 -0.500000 0.00000
+0.00000 -0.500000 0.00000
+-0.250000 -0.250000 0.00000
+0.00000 -0.250000 0.00000
+-0.500000 -0.250000 0.00000
+-0.250000 -0.250000 0.00000
+-0.500000 0.00000 0.00000
+-0.250000 0.00000 0.00000
+-0.250000 -0.250000 0.00000
+0.00000 -0.250000 0.00000
+-0.250000 0.00000 0.00000
+0.00000 0.00000 0.00000
+
+CELLS 16 80
+4      0       1       3       2
+4      4       5       7       6
+4      8       9       11      10
+4      12      13      15      14
+4      16      17      19      18
+4      20      21      23      22
+4      24      25      27      26
+4      28      29      31      30
+4      32      33      35      34
+4      36      37      39      38
+4      40      41      43      42
+4      44      45      47      46
+4      48      49      51      50
+4      52      53      55      54
+4      56      57      59      58
+4      60      61      63      62
+
+CELL_TYPES 16
+ 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9
+POINT_DATA 64
+SCALARS solution_0_0 double 1
+LOOKUP_TABLE default
+-1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 
+SCALARS solution_0_1 double 1
+LOOKUP_TABLE default
+-1.00000 -1.00000 -0.500000 -0.500000 -1.00000 -1.00000 -0.500000 -0.500000 -0.500000 -0.500000 0.00000 0.00000 -0.500000 -0.500000 0.00000 0.00000 -1.00000 -1.00000 -0.500000 -0.500000 -1.00000 -1.00000 -0.500000 -0.500000 -0.500000 -0.500000 0.00000 0.00000 -0.500000 -0.500000 0.00000 0.00000 0.00000 0.00000 0.500000 0.500000 0.00000 0.00000 0.500000 0.500000 0.500000 0.500000 1.00000 1.00000 0.500000 0.500000 1.00000 1.00000 0.00000 0.00000 0.500000 0.500000 0.00000 0.00000 0.500000 0.500000 0.500000 0.500000 1.00000 1.00000 0.500000 0.500000 1.00000 1.00000 
+SCALARS solution_0_2 double 1
+LOOKUP_TABLE default
+2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 
+SCALARS solution_1_0 double 1
+LOOKUP_TABLE default
+-1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 
+SCALARS solution_1_1 double 1
+LOOKUP_TABLE default
+-1.00000 -1.00000 -0.500000 -0.500000 -1.00000 -1.00000 -0.500000 -0.500000 -0.500000 -0.500000 0.00000 0.00000 -0.500000 -0.500000 0.00000 0.00000 -1.00000 -1.00000 -0.500000 -0.500000 -1.00000 -1.00000 -0.500000 -0.500000 -0.500000 -0.500000 0.00000 0.00000 -0.500000 -0.500000 0.00000 0.00000 0.00000 0.00000 0.500000 0.500000 0.00000 0.00000 0.500000 0.500000 0.500000 0.500000 1.00000 1.00000 0.500000 0.500000 1.00000 1.00000 0.00000 0.00000 0.500000 0.500000 0.00000 0.00000 0.500000 0.500000 0.500000 0.500000 1.00000 1.00000 0.500000 0.500000 1.00000 1.00000 
+SCALARS solution_1_2 double 1
+LOOKUP_TABLE default
+2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 
+
+# vtk DataFile Version 3.0
+#This file was generated 
+ASCII
+DATASET UNSTRUCTURED_GRID
+
+POINTS 64 double
+0.00000 -1.00000 0.00000
+0.250000 -1.00000 0.00000
+0.00000 -0.750000 0.00000
+0.250000 -0.750000 0.00000
+0.250000 -1.00000 0.00000
+0.500000 -1.00000 0.00000
+0.250000 -0.750000 0.00000
+0.500000 -0.750000 0.00000
+0.00000 -0.750000 0.00000
+0.250000 -0.750000 0.00000
+0.00000 -0.500000 0.00000
+0.250000 -0.500000 0.00000
+0.250000 -0.750000 0.00000
+0.500000 -0.750000 0.00000
+0.250000 -0.500000 0.00000
+0.500000 -0.500000 0.00000
+0.500000 -1.00000 0.00000
+0.750000 -1.00000 0.00000
+0.500000 -0.750000 0.00000
+0.750000 -0.750000 0.00000
+0.750000 -1.00000 0.00000
+1.00000 -1.00000 0.00000
+0.750000 -0.750000 0.00000
+1.00000 -0.750000 0.00000
+0.500000 -0.750000 0.00000
+0.750000 -0.750000 0.00000
+0.500000 -0.500000 0.00000
+0.750000 -0.500000 0.00000
+0.750000 -0.750000 0.00000
+1.00000 -0.750000 0.00000
+0.750000 -0.500000 0.00000
+1.00000 -0.500000 0.00000
+0.00000 -0.500000 0.00000
+0.250000 -0.500000 0.00000
+0.00000 -0.250000 0.00000
+0.250000 -0.250000 0.00000
+0.250000 -0.500000 0.00000
+0.500000 -0.500000 0.00000
+0.250000 -0.250000 0.00000
+0.500000 -0.250000 0.00000
+0.00000 -0.250000 0.00000
+0.250000 -0.250000 0.00000
+0.00000 0.00000 0.00000
+0.250000 0.00000 0.00000
+0.250000 -0.250000 0.00000
+0.500000 -0.250000 0.00000
+0.250000 0.00000 0.00000
+0.500000 0.00000 0.00000
+0.500000 -0.500000 0.00000
+0.750000 -0.500000 0.00000
+0.500000 -0.250000 0.00000
+0.750000 -0.250000 0.00000
+0.750000 -0.500000 0.00000
+1.00000 -0.500000 0.00000
+0.750000 -0.250000 0.00000
+1.00000 -0.250000 0.00000
+0.500000 -0.250000 0.00000
+0.750000 -0.250000 0.00000
+0.500000 0.00000 0.00000
+0.750000 0.00000 0.00000
+0.750000 -0.250000 0.00000
+1.00000 -0.250000 0.00000
+0.750000 0.00000 0.00000
+1.00000 0.00000 0.00000
+
+CELLS 16 80
+4      0       1       3       2
+4      4       5       7       6
+4      8       9       11      10
+4      12      13      15      14
+4      16      17      19      18
+4      20      21      23      22
+4      24      25      27      26
+4      28      29      31      30
+4      32      33      35      34
+4      36      37      39      38
+4      40      41      43      42
+4      44      45      47      46
+4      48      49      51      50
+4      52      53      55      54
+4      56      57      59      58
+4      60      61      63      62
+
+CELL_TYPES 16
+ 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9
+POINT_DATA 64
+SCALARS solution_0_0 double 1
+LOOKUP_TABLE default
+0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 
+SCALARS solution_0_1 double 1
+LOOKUP_TABLE default
+-1.00000 -1.00000 -0.500000 -0.500000 -1.00000 -1.00000 -0.500000 -0.500000 -0.500000 -0.500000 0.00000 0.00000 -0.500000 -0.500000 0.00000 0.00000 -1.00000 -1.00000 -0.500000 -0.500000 -1.00000 -1.00000 -0.500000 -0.500000 -0.500000 -0.500000 0.00000 0.00000 -0.500000 -0.500000 0.00000 0.00000 0.00000 0.00000 0.500000 0.500000 0.00000 0.00000 0.500000 0.500000 0.500000 0.500000 1.00000 1.00000 0.500000 0.500000 1.00000 1.00000 0.00000 0.00000 0.500000 0.500000 0.00000 0.00000 0.500000 0.500000 0.500000 0.500000 1.00000 1.00000 0.500000 0.500000 1.00000 1.00000 
+SCALARS solution_0_2 double 1
+LOOKUP_TABLE default
+2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 
+SCALARS solution_1_0 double 1
+LOOKUP_TABLE default
+0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 
+SCALARS solution_1_1 double 1
+LOOKUP_TABLE default
+-1.00000 -1.00000 -0.500000 -0.500000 -1.00000 -1.00000 -0.500000 -0.500000 -0.500000 -0.500000 0.00000 0.00000 -0.500000 -0.500000 0.00000 0.00000 -1.00000 -1.00000 -0.500000 -0.500000 -1.00000 -1.00000 -0.500000 -0.500000 -0.500000 -0.500000 0.00000 0.00000 -0.500000 -0.500000 0.00000 0.00000 0.00000 0.00000 0.500000 0.500000 0.00000 0.00000 0.500000 0.500000 0.500000 0.500000 1.00000 1.00000 0.500000 0.500000 1.00000 1.00000 0.00000 0.00000 0.500000 0.500000 0.00000 0.00000 0.500000 0.500000 0.500000 0.500000 1.00000 1.00000 0.500000 0.500000 1.00000 1.00000 
+SCALARS solution_1_2 double 1
+LOOKUP_TABLE default
+2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 
+
+
+# vtk DataFile Version 3.0
+#This file was generated 
+ASCII
+DATASET UNSTRUCTURED_GRID
+
+POINTS 64 double
+-1.00000 0.00000 0.00000
+-0.750000 0.00000 0.00000
+-1.00000 0.250000 0.00000
+-0.750000 0.250000 0.00000
+-0.750000 0.00000 0.00000
+-0.500000 0.00000 0.00000
+-0.750000 0.250000 0.00000
+-0.500000 0.250000 0.00000
+-1.00000 0.250000 0.00000
+-0.750000 0.250000 0.00000
+-1.00000 0.500000 0.00000
+-0.750000 0.500000 0.00000
+-0.750000 0.250000 0.00000
+-0.500000 0.250000 0.00000
+-0.750000 0.500000 0.00000
+-0.500000 0.500000 0.00000
+-0.500000 0.00000 0.00000
+-0.250000 0.00000 0.00000
+-0.500000 0.250000 0.00000
+-0.250000 0.250000 0.00000
+-0.250000 0.00000 0.00000
+0.00000 0.00000 0.00000
+-0.250000 0.250000 0.00000
+0.00000 0.250000 0.00000
+-0.500000 0.250000 0.00000
+-0.250000 0.250000 0.00000
+-0.500000 0.500000 0.00000
+-0.250000 0.500000 0.00000
+-0.250000 0.250000 0.00000
+0.00000 0.250000 0.00000
+-0.250000 0.500000 0.00000
+0.00000 0.500000 0.00000
+-1.00000 0.500000 0.00000
+-0.750000 0.500000 0.00000
+-1.00000 0.750000 0.00000
+-0.750000 0.750000 0.00000
+-0.750000 0.500000 0.00000
+-0.500000 0.500000 0.00000
+-0.750000 0.750000 0.00000
+-0.500000 0.750000 0.00000
+-1.00000 0.750000 0.00000
+-0.750000 0.750000 0.00000
+-1.00000 1.00000 0.00000
+-0.750000 1.00000 0.00000
+-0.750000 0.750000 0.00000
+-0.500000 0.750000 0.00000
+-0.750000 1.00000 0.00000
+-0.500000 1.00000 0.00000
+-0.500000 0.500000 0.00000
+-0.250000 0.500000 0.00000
+-0.500000 0.750000 0.00000
+-0.250000 0.750000 0.00000
+-0.250000 0.500000 0.00000
+0.00000 0.500000 0.00000
+-0.250000 0.750000 0.00000
+0.00000 0.750000 0.00000
+-0.500000 0.750000 0.00000
+-0.250000 0.750000 0.00000
+-0.500000 1.00000 0.00000
+-0.250000 1.00000 0.00000
+-0.250000 0.750000 0.00000
+0.00000 0.750000 0.00000
+-0.250000 1.00000 0.00000
+0.00000 1.00000 0.00000
+
+CELLS 16 80
+4      0       1       3       2
+4      4       5       7       6
+4      8       9       11      10
+4      12      13      15      14
+4      16      17      19      18
+4      20      21      23      22
+4      24      25      27      26
+4      28      29      31      30
+4      32      33      35      34
+4      36      37      39      38
+4      40      41      43      42
+4      44      45      47      46
+4      48      49      51      50
+4      52      53      55      54
+4      56      57      59      58
+4      60      61      63      62
+
+CELL_TYPES 16
+ 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9
+POINT_DATA 64
+SCALARS solution_0_0 double 1
+LOOKUP_TABLE default
+-1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 
+SCALARS solution_0_1 double 1
+LOOKUP_TABLE default
+1.00000 1.00000 1.50000 1.50000 1.00000 1.00000 1.50000 1.50000 1.50000 1.50000 2.00000 2.00000 1.50000 1.50000 2.00000 2.00000 1.00000 1.00000 1.50000 1.50000 1.00000 1.00000 1.50000 1.50000 1.50000 1.50000 2.00000 2.00000 1.50000 1.50000 2.00000 2.00000 2.00000 2.00000 2.50000 2.50000 2.00000 2.00000 2.50000 2.50000 2.50000 2.50000 3.00000 3.00000 2.50000 2.50000 3.00000 3.00000 2.00000 2.00000 2.50000 2.50000 2.00000 2.00000 2.50000 2.50000 2.50000 2.50000 3.00000 3.00000 2.50000 2.50000 3.00000 3.00000 
+SCALARS solution_0_2 double 1
+LOOKUP_TABLE default
+2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 
+SCALARS solution_1_0 double 1
+LOOKUP_TABLE default
+-1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -1.00000 -0.750000 -1.00000 -0.750000 -0.750000 -0.500000 -0.750000 -0.500000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 -0.500000 -0.250000 -0.500000 -0.250000 -0.250000 0.00000 -0.250000 0.00000 
+SCALARS solution_1_1 double 1
+LOOKUP_TABLE default
+1.00000 1.00000 1.50000 1.50000 1.00000 1.00000 1.50000 1.50000 1.50000 1.50000 2.00000 2.00000 1.50000 1.50000 2.00000 2.00000 1.00000 1.00000 1.50000 1.50000 1.00000 1.00000 1.50000 1.50000 1.50000 1.50000 2.00000 2.00000 1.50000 1.50000 2.00000 2.00000 2.00000 2.00000 2.50000 2.50000 2.00000 2.00000 2.50000 2.50000 2.50000 2.50000 3.00000 3.00000 2.50000 2.50000 3.00000 3.00000 2.00000 2.00000 2.50000 2.50000 2.00000 2.00000 2.50000 2.50000 2.50000 2.50000 3.00000 3.00000 2.50000 2.50000 3.00000 3.00000 
+SCALARS solution_1_2 double 1
+LOOKUP_TABLE default
+2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 
+
+
+# vtk DataFile Version 3.0
+#This file was generated 
+ASCII
+DATASET UNSTRUCTURED_GRID
+
+POINTS 64 double
+0.00000 0.00000 0.00000
+0.250000 0.00000 0.00000
+0.00000 0.250000 0.00000
+0.250000 0.250000 0.00000
+0.250000 0.00000 0.00000
+0.500000 0.00000 0.00000
+0.250000 0.250000 0.00000
+0.500000 0.250000 0.00000
+0.00000 0.250000 0.00000
+0.250000 0.250000 0.00000
+0.00000 0.500000 0.00000
+0.250000 0.500000 0.00000
+0.250000 0.250000 0.00000
+0.500000 0.250000 0.00000
+0.250000 0.500000 0.00000
+0.500000 0.500000 0.00000
+0.500000 0.00000 0.00000
+0.750000 0.00000 0.00000
+0.500000 0.250000 0.00000
+0.750000 0.250000 0.00000
+0.750000 0.00000 0.00000
+1.00000 0.00000 0.00000
+0.750000 0.250000 0.00000
+1.00000 0.250000 0.00000
+0.500000 0.250000 0.00000
+0.750000 0.250000 0.00000
+0.500000 0.500000 0.00000
+0.750000 0.500000 0.00000
+0.750000 0.250000 0.00000
+1.00000 0.250000 0.00000
+0.750000 0.500000 0.00000
+1.00000 0.500000 0.00000
+0.00000 0.500000 0.00000
+0.250000 0.500000 0.00000
+0.00000 0.750000 0.00000
+0.250000 0.750000 0.00000
+0.250000 0.500000 0.00000
+0.500000 0.500000 0.00000
+0.250000 0.750000 0.00000
+0.500000 0.750000 0.00000
+0.00000 0.750000 0.00000
+0.250000 0.750000 0.00000
+0.00000 1.00000 0.00000
+0.250000 1.00000 0.00000
+0.250000 0.750000 0.00000
+0.500000 0.750000 0.00000
+0.250000 1.00000 0.00000
+0.500000 1.00000 0.00000
+0.500000 0.500000 0.00000
+0.750000 0.500000 0.00000
+0.500000 0.750000 0.00000
+0.750000 0.750000 0.00000
+0.750000 0.500000 0.00000
+1.00000 0.500000 0.00000
+0.750000 0.750000 0.00000
+1.00000 0.750000 0.00000
+0.500000 0.750000 0.00000
+0.750000 0.750000 0.00000
+0.500000 1.00000 0.00000
+0.750000 1.00000 0.00000
+0.750000 0.750000 0.00000
+1.00000 0.750000 0.00000
+0.750000 1.00000 0.00000
+1.00000 1.00000 0.00000
+
+CELLS 16 80
+4      0       1       3       2
+4      4       5       7       6
+4      8       9       11      10
+4      12      13      15      14
+4      16      17      19      18
+4      20      21      23      22
+4      24      25      27      26
+4      28      29      31      30
+4      32      33      35      34
+4      36      37      39      38
+4      40      41      43      42
+4      44      45      47      46
+4      48      49      51      50
+4      52      53      55      54
+4      56      57      59      58
+4      60      61      63      62
+
+CELL_TYPES 16
+ 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9 9
+POINT_DATA 64
+SCALARS solution_0_0 double 1
+LOOKUP_TABLE default
+0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 
+SCALARS solution_0_1 double 1
+LOOKUP_TABLE default
+1.00000 1.00000 1.50000 1.50000 1.00000 1.00000 1.50000 1.50000 1.50000 1.50000 2.00000 2.00000 1.50000 1.50000 2.00000 2.00000 1.00000 1.00000 1.50000 1.50000 1.00000 1.00000 1.50000 1.50000 1.50000 1.50000 2.00000 2.00000 1.50000 1.50000 2.00000 2.00000 2.00000 2.00000 2.50000 2.50000 2.00000 2.00000 2.50000 2.50000 2.50000 2.50000 3.00000 3.00000 2.50000 2.50000 3.00000 3.00000 2.00000 2.00000 2.50000 2.50000 2.00000 2.00000 2.50000 2.50000 2.50000 2.50000 3.00000 3.00000 2.50000 2.50000 3.00000 3.00000 
+SCALARS solution_0_2 double 1
+LOOKUP_TABLE default
+2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 
+SCALARS solution_1_0 double 1
+LOOKUP_TABLE default
+0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.00000 0.250000 0.00000 0.250000 0.250000 0.500000 0.250000 0.500000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 0.500000 0.750000 0.500000 0.750000 0.750000 1.00000 0.750000 1.00000 
+SCALARS solution_1_1 double 1
+LOOKUP_TABLE default
+1.00000 1.00000 1.50000 1.50000 1.00000 1.00000 1.50000 1.50000 1.50000 1.50000 2.00000 2.00000 1.50000 1.50000 2.00000 2.00000 1.00000 1.00000 1.50000 1.50000 1.00000 1.00000 1.50000 1.50000 1.50000 1.50000 2.00000 2.00000 1.50000 1.50000 2.00000 2.00000 2.00000 2.00000 2.50000 2.50000 2.00000 2.00000 2.50000 2.50000 2.50000 2.50000 3.00000 3.00000 2.50000 2.50000 3.00000 3.00000 2.00000 2.00000 2.50000 2.50000 2.00000 2.00000 2.50000 2.50000 2.50000 2.50000 3.00000 3.00000 2.50000 2.50000 3.00000 3.00000 
+SCALARS solution_1_2 double 1
+LOOKUP_TABLE default
+2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 2.00000 
+

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.