]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added Gmsh::create_triangulation_from_boundary_curve
authorLuca Heltai <luca.heltai@sissa.it>
Thu, 3 May 2018 15:43:38 +0000 (17:43 +0200)
committerLuca Heltai <luca.heltai@sissa.it>
Thu, 3 May 2018 19:42:08 +0000 (21:42 +0200)
include/deal.II/gmsh/utilities.h [new file with mode: 0644]
source/CMakeLists.txt
source/gmsh/CMakeLists.txt [new file with mode: 0644]
source/gmsh/utilities.cc [new file with mode: 0644]
source/gmsh/utilities.inst.in [new file with mode: 0644]
tests/gmsh/create_tria_02.cc [new file with mode: 0644]
tests/gmsh/create_tria_02.output [new file with mode: 0644]

diff --git a/include/deal.II/gmsh/utilities.h b/include/deal.II/gmsh/utilities.h
new file mode 100644 (file)
index 0000000..4e3045a
--- /dev/null
@@ -0,0 +1,102 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2018 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+#ifndef dealii_gmsh_parameters_h
+#define dealii_gmsh_parameters_h
+
+#include <deal.II/base/config.h>
+
+#ifdef DEAL_II_WITH_GMSH
+
+#ifdef DEAL_II_WITH_OPENCASCADE
+#include <TopoDS_Shape.hxx>
+#include <TopoDS_Face.hxx>
+#include <TopoDS_Edge.hxx>
+#include <TopoDS_Vertex.hxx>
+#include <TopoDS_Compound.hxx>
+#include <TopoDS_CompSolid.hxx>
+#endif
+
+#include <deal.II/base/parameter_handler.h>
+#include <deal.II/grid/tria.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+/**
+ * A collection of Gmsh related utilities and classes.
+ *
+ * @author Luca Heltai, Dirk Peschka, 2018
+ */
+namespace Gmsh
+{
+
+  /**
+   * A parameter class used to pass options to the Gmsh executable.
+   *
+   * @author Luca Heltai, 2018
+   */
+  class AdditionalParameters
+  {
+  public:
+    /**
+     * Set all additional parameters to their default values.
+     */
+    AdditionalParameters(const double characteristic_length=1.0,
+                         const std::string &output_base_name="");
+
+    /**
+     * Call prm.add_parameter for each member of the AdditionalParameters class.
+     */
+    void add_parameters(ParameterHandler &prm);
+
+    /**
+     * The characteristic length used for the definition of the Gmsh grid.
+     *
+     * Gmsh will try to make sure that the size of each edge is as close as
+     * possible to this value.
+     */
+    double characteristic_length = 1.0;
+
+    /**
+     * Basename for the output files.
+     *
+     * If this is left empty, then temporary files are used, and removed when
+     * not needed any more.
+     */
+    std::string output_base_name = "";
+  };
+
+#ifdef DEAL_II_WITH_OPENCASCADE
+  /**
+   * Given a smooth closed curve, creates a triangulation from it, using
+   * Gmsh.
+   *
+   * The input curve @p boundary should be closed.
+   *
+   * @authors Luca Heltai, Dirk Peschka, 2018
+   */
+  template <int spacedim>
+  void
+  create_triangulation_from_boundary_curve(const TopoDS_Edge &boundary,
+                                           Triangulation<2,spacedim> &tria,
+                                           const AdditionalParameters &prm=AdditionalParameters());
+#endif
+}
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
+#endif
index 0270734e8147107daf878f1d11cde0ae4b2e5513..2704d12c773e15cb552f55a19fb78f5f277eca7c 100644 (file)
@@ -40,6 +40,7 @@ ADD_SUBDIRECTORY(fe)
 ADD_SUBDIRECTORY(dofs)
 ADD_SUBDIRECTORY(lac)
 ADD_SUBDIRECTORY(base)
+ADD_SUBDIRECTORY(gmsh)
 ADD_SUBDIRECTORY(grid)
 ADD_SUBDIRECTORY(hp)
 ADD_SUBDIRECTORY(multigrid)
diff --git a/source/gmsh/CMakeLists.txt b/source/gmsh/CMakeLists.txt
new file mode 100644 (file)
index 0000000..ccf5502
--- /dev/null
@@ -0,0 +1,31 @@
+## ---------------------------------------------------------------------
+##
+## Copyright (C) 2012 - 2017 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 at
+## the top level of the deal.II distribution.
+##
+## ---------------------------------------------------------------------
+
+INCLUDE_DIRECTORIES(BEFORE ${CMAKE_CURRENT_BINARY_DIR})
+
+SET(_src
+  utilities.cc
+  )
+
+SET(_inst
+  utilities.inst.in
+  )
+
+FILE(GLOB _header
+  ${CMAKE_SOURCE_DIR}/include/deal.II/gmsh/*.h
+  )
+
+DEAL_II_ADD_LIBRARY(obj_gmsh OBJECT ${_src} ${_header} ${_inst})
+EXPAND_INSTANTIATIONS(obj_gmsh "${_inst}")
diff --git a/source/gmsh/utilities.cc b/source/gmsh/utilities.cc
new file mode 100644 (file)
index 0000000..bd790da
--- /dev/null
@@ -0,0 +1,92 @@
+#include <deal.II/gmsh/utilities.h>
+#include <deal.II/opencascade/utilities.h>
+
+#include <deal.II/grid/grid_in.h>
+#include <deal.II/grid/grid_out.h>
+
+#include <cstdio>
+
+#ifdef DEAL_II_WITH_GMSH
+
+DEAL_II_NAMESPACE_OPEN
+
+namespace Gmsh
+{
+  AdditionalParameters::AdditionalParameters(const double characteristic_length,
+                                             const std::string &output_base_name):
+    characteristic_length(characteristic_length),
+    output_base_name(output_base_name)
+  {}
+
+
+
+  void AdditionalParameters::add_parameters(ParameterHandler &prm)
+  {
+    prm.add_parameter("Characteristic length", characteristic_length);
+    prm.add_parameter("Intermediate file name base", output_base_name,
+                      "Keep empty, if you want the program to generate "
+                      "temporary files, and then remove them when they "
+                      "are no longer used.");
+  }
+
+
+
+#ifdef DEAL_II_WITH_OPENCASCADE
+  template <int spacedim>
+  void
+  create_triangulation_from_boundary_curve(const TopoDS_Edge &boundary,
+                                           Triangulation<2,spacedim> &tria,
+                                           const AdditionalParameters &prm)
+  {
+    std::string base_name = prm.output_base_name;
+    if (base_name == "")
+      base_name = std::tmpnam(nullptr);
+
+    dealii::OpenCASCADE::write_IGES(boundary, base_name+".iges");
+
+    ofstream geofile;
+    geofile.open( base_name+".geo");
+    geofile << "Merge \"" << base_name << ".iges\";" << std::endl
+            << "Line Loop (2) = {1};" << std::endl
+            << "Plane Surface (3) = {2};" << std::endl
+            << "Characteristic Length { 1 } = " << prm.characteristic_length << ";" << std::endl
+            << "Mesh.RecombineAll = 1;" << std::endl
+            << "Mesh.SubdivisionAlgorithm=1;" << std::endl;
+
+    geofile.close();
+
+    // std::system("gmsh -2 -algo front3d temp_model.geo 1>temp_out.log 2>temp_warn.log");
+    std::stringstream command;
+    command << DEAL_II_GMSH_EXECUTABLE_PATH << " -2 "
+            << base_name << ".geo 1> "
+            << base_name << ".log 2> "
+            << base_name << "_warn.log";
+
+    auto ret_value = std::system(command.str().c_str());
+    AssertThrow(ret_value == 0,
+                ExcMessage("Gmsh failed to run. Check the "+base_name+".log file."));
+
+    std::ifstream grid_file(base_name+".msh");
+    Assert(grid_file, ExcIO());
+
+    GridIn<2,spacedim> gridin;
+    gridin.attach_triangulation(tria);
+    gridin.read_msh(grid_file);
+
+    if (base_name != prm.output_base_name)
+      {
+        std::remove((base_name + ".geo").c_str());
+        std::remove((base_name + ".log").c_str());
+        std::remove((base_name + "_warn.log").c_str());
+        std::remove((base_name + ".msh").c_str());
+        std::remove((base_name + ".iges").c_str());
+      }
+  }
+#endif
+
+#include "utilities.inst"
+}
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
diff --git a/source/gmsh/utilities.inst.in b/source/gmsh/utilities.inst.in
new file mode 100644 (file)
index 0000000..a26d5f8
--- /dev/null
@@ -0,0 +1,27 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2018 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+for (deal_II_dimension : DIMENSIONS)
+{
+#if deal_II_dimension > 1
+
+  template
+  void
+  create_triangulation_from_boundary_curve(const TopoDS_Edge &boundary,
+                                           Triangulation<2,deal_II_dimension> &tria,
+                                           const AdditionalParameters &prm);
+
+#endif
+}
diff --git a/tests/gmsh/create_tria_02.cc b/tests/gmsh/create_tria_02.cc
new file mode 100644 (file)
index 0000000..d75c98c
--- /dev/null
@@ -0,0 +1,49 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2014 - 2018 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+// Create a hyper ball, refine it, extract an iges of the boundary,
+// and create a new mesh using gmsh from that iges file.
+
+#include "../tests.h"
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_out.h>
+
+#include <deal.II/gmsh/utilities.h>
+#include <deal.II/opencascade/utilities.h>
+#include <deal.II/opencascade/boundary_lib.h>
+
+int main ()
+{
+  initlog();
+
+  Triangulation<2> tria;
+  GridGenerator::hyper_ball(tria);
+
+  tria.refine_global(4);
+
+  auto curves = OpenCASCADE::create_curves_from_triangulation_boundary(tria);
+  tria.clear();
+
+  Gmsh::create_triangulation_from_boundary_curve(curves[0], tria);
+
+  GridOut go;
+  std::ofstream ofile("output.svg");
+  go.write_svg(tria, ofile);
+  go.write_svg(tria, deallog.get_file_stream());
+
+  return 0;
+}
diff --git a/tests/gmsh/create_tria_02.output b/tests/gmsh/create_tria_02.output
new file mode 100644 (file)
index 0000000..679a28c
--- /dev/null
@@ -0,0 +1,205 @@
+
+<svg width="1400" height="1000" xmlns="http://www.w3.org/2000/svg" version="1.1">
+
+
+<!-- internal style sheet -->
+<style type="text/css"><![CDATA[
+ rect.background{fill:white}
+ rect{fill:none; stroke:rgb(25,25,25); stroke-width:2}
+ text{font-family:Helvetica; text-anchor:middle; fill:rgb(25,25,25)}
+ line{stroke:rgb(25,25,25); stroke-width:4}
+ path{fill:none; stroke:rgb(25,25,25); stroke-width:2}
+
+ path.p0{fill:rgb(0,102,255); stroke:rgb(25,25,25); stroke-width:2}
+ path.ps0{fill:rgb(0,77,191); stroke:rgb(20,20,20); stroke-width:2}
+ rect.r0{fill:rgb(0,102,255); stroke:rgb(25,25,25); stroke-width:2}
+]]></style>
+
+ <rect class="background" width="1000" height="1000"/>
+  <!-- cells -->
+  <path class="p0" d="M 210 286 L 130 298 L 105 353 L 221 373 L 210 286"/>
+  <text x="167" y="335" style="font-size:15px">0,0</text>
+  <line x1="130" y1="298" x2="105" y2="353"/>
+  <path class="p0" d="M 292 278 L 210 286 L 221 373 L 336 392 L 292 278"/>
+  <text x="265" y="340" style="font-size:15px">0,1</text>
+  <path class="p0" d="M 248 163 L 202 202 L 210 286 L 292 278 L 248 163"/>
+  <text x="238" y="240" style="font-size:15px">0,2</text>
+  <line x1="248" y1="163" x2="202" y2="202"/>
+  <path class="p0" d="M 202 202 L 163 248 L 130 298 L 210 286 L 202 202"/>
+  <text x="176" y="266" style="font-size:15px">0,3</text>
+  <line x1="163" y1="248" x2="130" y2="298"/>
+  <line x1="202" y1="202" x2="163" y2="248"/>
+  <path class="p0" d="M 471 805 L 563 785 L 597 650 L 473 690 L 471 805"/>
+  <text x="526" y="740" style="font-size:15px">0,4</text>
+  <path class="p0" d="M 380 821 L 471 805 L 473 690 L 349 731 L 380 821"/>
+  <text x="418" y="769" style="font-size:15px">0,5</text>
+  <path class="p0" d="M 410 911 L 470 920 L 471 805 L 380 821 L 410 911"/>
+  <text x="433" y="872" style="font-size:15px">0,6</text>
+  <line x1="410" y1="911" x2="470" y2="920"/>
+  <path class="p0" d="M 470 920 L 530 920 L 563 785 L 471 805 L 470 920"/>
+  <text x="509" y="865" style="font-size:15px">0,7</text>
+  <line x1="470" y1="920" x2="530" y2="920"/>
+  <path class="p0" d="M 538 412 L 466 521 L 597 650 L 656 506 L 538 412"/>
+  <text x="564" y="530" style="font-size:16px">0,8</text>
+  <path class="p0" d="M 609 302 L 538 412 L 656 506 L 715 362 L 609 302"/>
+  <text x="629" y="403" style="font-size:15px">0,9</text>
+  <path class="p0" d="M 502 242 L 419 317 L 538 412 L 609 302 L 502 242"/>
+  <text x="517" y="326" style="font-size:15px">0,10</text>
+  <path class="p0" d="M 419 317 L 336 392 L 466 521 L 538 412 L 419 317"/>
+  <text x="440" y="418" style="font-size:15px">0,11</text>
+  <path class="p0" d="M 378 589 L 473 690 L 597 650 L 466 521 L 378 589"/>
+  <text x="479" y="620" style="font-size:15px">0,12</text>
+  <path class="p0" d="M 284 487 L 378 589 L 466 521 L 336 392 L 284 487"/>
+  <text x="366" y="505" style="font-size:15px">0,13</text>
+  <path class="p0" d="M 232 581 L 290 656 L 378 589 L 284 487 L 232 581"/>
+  <text x="296" y="586" style="font-size:15px">0,14</text>
+  <path class="p0" d="M 290 656 L 349 731 L 473 690 L 378 589 L 290 656"/>
+  <text x="373" y="674" style="font-size:15px">0,15</text>
+  <path class="p0" d="M 628 204 L 609 302 L 715 362 L 709 246 L 628 204"/>
+  <text x="665" y="286" style="font-size:15px">0,16</text>
+  <path class="p0" d="M 647 105 L 628 204 L 709 246 L 702 130 L 647 105"/>
+  <text x="671" y="179" style="font-size:15px">0,17</text>
+  <line x1="647" y1="105" x2="702" y2="130"/>
+  <path class="p0" d="M 590 89 L 546 165 L 628 204 L 647 105 L 590 89"/>
+  <text x="603" y="148" style="font-size:15px">0,18</text>
+  <line x1="590" y1="89" x2="647" y2="105"/>
+  <path class="p0" d="M 546 165 L 502 242 L 609 302 L 628 204 L 546 165"/>
+  <text x="571" y="236" style="font-size:15px">0,19</text>
+  <path class="p0" d="M 227 704 L 163 752 L 202 798 L 275 764 L 227 704"/>
+  <text x="217" y="762" style="font-size:15px">0,20</text>
+  <line x1="163" y1="752" x2="202" y2="798"/>
+  <path class="p0" d="M 290 656 L 227 704 L 275 764 L 349 731 L 290 656"/>
+  <text x="285" y="721" style="font-size:15px">0,21</text>
+  <path class="p0" d="M 232 581 L 181 642 L 227 704 L 290 656 L 232 581"/>
+  <text x="233" y="653" style="font-size:15px">0,22</text>
+  <path class="p0" d="M 181 642 L 130 702 L 163 752 L 227 704 L 181 642"/>
+  <text x="175" y="707" style="font-size:15px">0,23</text>
+  <line x1="130" y1="702" x2="163" y2="752"/>
+  <path class="p0" d="M 298 130 L 248 163 L 292 278 L 359 224 L 298 130"/>
+  <text x="299" y="206" style="font-size:15px">0,24</text>
+  <line x1="298" y1="130" x2="248" y2="163"/>
+  <path class="p0" d="M 359 224 L 292 278 L 336 392 L 419 317 L 359 224"/>
+  <text x="351" y="310" style="font-size:15px">0,25</text>
+  <path class="p0" d="M 427 174 L 359 224 L 419 317 L 502 242 L 427 174"/>
+  <text x="427" y="247" style="font-size:15px">0,26</text>
+  <path class="p0" d="M 353 105 L 298 130 L 359 224 L 427 174 L 353 105"/>
+  <text x="359" y="166" style="font-size:15px">0,27</text>
+  <line x1="353" y1="105" x2="298" y2="130"/>
+  <path class="p0" d="M 156 526 L 186 449 L 89 410 L 80 470 L 156 526"/>
+  <text x="128" y="471" style="font-size:15px">0,28</text>
+  <line x1="80" y1="470" x2="89" y2="410"/>
+  <path class="p0" d="M 232 581 L 284 487 L 186 449 L 156 526 L 232 581"/>
+  <text x="215" y="518" style="font-size:15px">0,29</text>
+  <path class="p0" d="M 284 487 L 336 392 L 221 373 L 186 449 L 284 487"/>
+  <text x="257" y="433" style="font-size:15px">0,30</text>
+  <path class="p0" d="M 186 449 L 221 373 L 105 353 L 89 410 L 186 449"/>
+  <text x="150" y="404" style="font-size:15px">0,31</text>
+  <line x1="89" y1="410" x2="105" y2="353"/>
+  <path class="p0" d="M 858 399 L 818 446 L 920 530 L 920 470 L 858 399"/>
+  <text x="879" y="469" style="font-size:15px">0,32</text>
+  <line x1="920" y1="470" x2="920" y2="530"/>
+  <path class="p0" d="M 895 353 L 858 399 L 920 470 L 911 410 L 895 353"/>
+  <text x="896" y="416" style="font-size:15px">0,33</text>
+  <line x1="895" y1="353" x2="911" y2="410"/>
+  <line x1="911" y1="410" x2="920" y2="470"/>
+  <path class="p0" d="M 870 298 L 792 330 L 858 399 L 895 353 L 870 298"/>
+  <text x="854" y="353" style="font-size:15px">0,34</text>
+  <line x1="870" y1="298" x2="895" y2="353"/>
+  <path class="p0" d="M 792 330 L 715 362 L 818 446 L 858 399 L 792 330"/>
+  <text x="796" y="392" style="font-size:15px">0,35</text>
+  <path class="p0" d="M 774 724 L 870 702 L 895 647 L 746 649 L 774 724"/>
+  <text x="821" y="688" style="font-size:15px">0,36</text>
+  <line x1="870" y1="702" x2="895" y2="647"/>
+  <path class="p0" d="M 675 744 L 774 724 L 746 649 L 597 650 L 675 744"/>
+  <text x="698" y="699" style="font-size:15px">0,37</text>
+  <path class="p0" d="M 752 837 L 798 798 L 774 724 L 675 744 L 752 837"/>
+  <text x="750" y="783" style="font-size:15px">0,38</text>
+  <line x1="752" y1="837" x2="798" y2="798"/>
+  <path class="p0" d="M 798 798 L 837 752 L 870 702 L 774 724 L 798 798"/>
+  <text x="820" y="752" style="font-size:15px">0,39</text>
+  <line x1="837" y1="752" x2="870" y2="702"/>
+  <line x1="798" y1="798" x2="837" y2="752"/>
+  <path class="p0" d="M 181 642 L 129 587 L 105 647 L 130 702 L 181 642"/>
+  <text x="136" y="652" style="font-size:15px">0,40</text>
+  <line x1="130" y1="702" x2="105" y2="647"/>
+  <path class="p0" d="M 232 581 L 156 526 L 129 587 L 181 642 L 232 581"/>
+  <text x="175" y="591" style="font-size:15px">0,41</text>
+  <path class="p0" d="M 156 526 L 80 470 L 80 530 L 129 587 L 156 526"/>
+  <text x="111" y="536" style="font-size:15px">0,42</text>
+  <line x1="80" y1="470" x2="80" y2="530"/>
+  <path class="p0" d="M 129 587 L 80 530 L 89 590 L 105 647 L 129 587"/>
+  <text x="101" y="595" style="font-size:14px">0,43</text>
+  <line x1="80" y1="530" x2="89" y2="590"/>
+  <line x1="105" y1="647" x2="89" y2="590"/>
+  <path class="p0" d="M 752 163 L 702 130 L 709 246 L 774 245 L 752 163"/>
+  <text x="734" y="204" style="font-size:15px">0,44</text>
+  <line x1="752" y1="163" x2="702" y2="130"/>
+  <path class="p0" d="M 774 245 L 709 246 L 715 362 L 792 330 L 774 245"/>
+  <text x="748" y="304" style="font-size:15px">0,45</text>
+  <path class="p0" d="M 837 248 L 774 245 L 792 330 L 870 298 L 837 248"/>
+  <text x="818" y="288" style="font-size:15px">0,46</text>
+  <line x1="837" y1="248" x2="870" y2="298"/>
+  <path class="p0" d="M 798 202 L 752 163 L 774 245 L 837 248 L 798 202"/>
+  <text x="790" y="222" style="font-size:15px">0,47</text>
+  <line x1="798" y1="202" x2="837" y2="248"/>
+  <line x1="798" y1="202" x2="752" y2="163"/>
+  <path class="p0" d="M 590 911 L 633 829 L 563 785 L 530 920 L 590 911"/>
+  <text x="579" y="869" style="font-size:15px">0,48</text>
+  <line x1="590" y1="911" x2="530" y2="920"/>
+  <path class="p0" d="M 647 895 L 702 870 L 633 829 L 590 911 L 647 895"/>
+  <text x="643" y="884" style="font-size:15px">0,49</text>
+  <line x1="647" y1="895" x2="590" y2="911"/>
+  <line x1="647" y1="895" x2="702" y2="870"/>
+  <path class="p0" d="M 702 870 L 752 837 L 675 744 L 633 829 L 702 870"/>
+  <text x="690" y="827" style="font-size:15px">0,50</text>
+  <line x1="702" y1="870" x2="752" y2="837"/>
+  <path class="p0" d="M 633 829 L 675 744 L 597 650 L 563 785 L 633 829"/>
+  <text x="617" y="760" style="font-size:15px">0,51</text>
+  <path class="p0" d="M 410 89 L 353 105 L 427 174 L 478 125 L 410 89"/>
+  <text x="417" y="131" style="font-size:15px">0,52</text>
+  <line x1="410" y1="89" x2="353" y2="105"/>
+  <path class="p0" d="M 478 125 L 427 174 L 502 242 L 546 165 L 478 125"/>
+  <text x="488" y="184" style="font-size:15px">0,53</text>
+  <path class="p0" d="M 530 80 L 478 125 L 546 165 L 590 89 L 530 80"/>
+  <text x="536" y="122" style="font-size:15px">0,54</text>
+  <line x1="530" y1="80" x2="590" y2="89"/>
+  <path class="p0" d="M 470 80 L 410 89 L 478 125 L 530 80 L 470 80"/>
+  <text x="472" y="100" style="font-size:14px">0,55</text>
+  <line x1="470" y1="80" x2="530" y2="80"/>
+  <line x1="470" y1="80" x2="410" y2="89"/>
+  <path class="p0" d="M 410 911 L 380 821 L 313 831 L 353 895 L 410 911"/>
+  <text x="364" y="872" style="font-size:15px">0,56</text>
+  <line x1="410" y1="911" x2="353" y2="895"/>
+  <path class="p0" d="M 380 821 L 349 731 L 275 764 L 313 831 L 380 821"/>
+  <text x="329" y="794" style="font-size:15px">0,57</text>
+  <path class="p0" d="M 313 831 L 275 764 L 202 798 L 248 837 L 313 831"/>
+  <text x="260" y="815" style="font-size:15px">0,58</text>
+  <line x1="248" y1="837" x2="202" y2="798"/>
+  <path class="p0" d="M 353 895 L 313 831 L 248 837 L 298 870 L 353 895"/>
+  <text x="303" y="865" style="font-size:14px">0,59</text>
+  <line x1="353" y1="895" x2="298" y2="870"/>
+  <line x1="298" y1="870" x2="248" y2="837"/>
+  <path class="p0" d="M 818 446 L 784 548 L 911 590 L 920 530 L 818 446"/>
+  <text x="858" y="536" style="font-size:15px">0,60</text>
+  <line x1="920" y1="530" x2="911" y2="590"/>
+  <path class="p0" d="M 715 362 L 656 506 L 784 548 L 818 446 L 715 362"/>
+  <text x="743" y="473" style="font-size:15px">0,61</text>
+  <path class="p0" d="M 656 506 L 597 650 L 746 649 L 784 548 L 656 506"/>
+  <text x="696" y="596" style="font-size:15px">0,62</text>
+  <path class="p0" d="M 784 548 L 746 649 L 895 647 L 911 590 L 784 548"/>
+  <text x="834" y="616" style="font-size:15px">0,63</text>
+  <line x1="911" y1="590" x2="895" y2="647"/>
+
+ <!-- legend -->
+ <rect x="1000" y="80" width="320" height="165"/>
+ <text x="1013" y="107" style="text-anchor:start; font-weight:bold; font-size:18px">cell label</text>
+  <text x="1020" y="134" style="text-anchor:start; font-style:oblique; font-size:18px">level_number,</text>
+  <text x="1020" y="161" style="text-anchor:start; font-style:oblique; font-size:18px">cell_index</text>
+  <text x="1000" y="274" style="text-anchor:start; font-size:18px">azimuth: 0°, polar: 0°</text>
+
+ <!-- colorbar -->
+ <text x="1000" y="356" style="text-anchor:start; font-weight:bold; font-size:18px">level_number</text>
+  <rect class="r0" x="1000" y="370" width="25" height="550"/>
+  <text x="1037.50" y="651.000" style="text-anchor:start; font-size:18px; font-weight:bold">0 max min</text>
+
+</svg>
\ No newline at end of file

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.