]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add tests for deal solvers with TpetraWrappers 16625/head
authorJan Philipp Thiele <thiele@wias-berlin.de>
Sat, 10 Feb 2024 15:00:28 +0000 (16:00 +0100)
committerJan Philipp Thiele <thiele@wias-berlin.de>
Sat, 10 Feb 2024 15:00:28 +0000 (16:00 +0100)
12 files changed:
tests/trilinos_tpetra/deal_solver_01.cc [new file with mode: 0644]
tests/trilinos_tpetra/deal_solver_01.output [new file with mode: 0644]
tests/trilinos_tpetra/deal_solver_02.cc [new file with mode: 0644]
tests/trilinos_tpetra/deal_solver_02.output [new file with mode: 0644]
tests/trilinos_tpetra/deal_solver_03.cc [new file with mode: 0644]
tests/trilinos_tpetra/deal_solver_03.output [new file with mode: 0644]
tests/trilinos_tpetra/deal_solver_04.cc [new file with mode: 0644]
tests/trilinos_tpetra/deal_solver_04.output [new file with mode: 0644]
tests/trilinos_tpetra/deal_solver_05.cc [new file with mode: 0644]
tests/trilinos_tpetra/deal_solver_05.output [new file with mode: 0644]
tests/trilinos_tpetra/deal_solver_06.cc [new file with mode: 0644]
tests/trilinos_tpetra/deal_solver_06.output [new file with mode: 0644]

diff --git a/tests/trilinos_tpetra/deal_solver_01.cc b/tests/trilinos_tpetra/deal_solver_01.cc
new file mode 100644 (file)
index 0000000..fa49c8b
--- /dev/null
@@ -0,0 +1,87 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2004 - 2024 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 the CG solver using the Trilinos matrix and vector classes
+
+
+#include <deal.II/base/utilities.h>
+
+#include <deal.II/lac/precondition.h>
+#include <deal.II/lac/solver_bicgstab.h>
+#include <deal.II/lac/solver_cg.h>
+#include <deal.II/lac/solver_control.h>
+#include <deal.II/lac/solver_gmres.h>
+#include <deal.II/lac/solver_qmrs.h>
+#include <deal.II/lac/solver_richardson.h>
+#include <deal.II/lac/trilinos_tpetra_sparse_matrix.h>
+#include <deal.II/lac/vector_memory.h>
+
+#include <iostream>
+#include <typeinfo>
+
+#include "../tests.h"
+
+#include "../testmatrix.h"
+
+
+int
+main(int argc, char **argv)
+{
+  initlog();
+  deallog << std::setprecision(4);
+
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(
+    argc, argv, testing_max_num_threads());
+
+
+  {
+    SolverControl control(200, 1.e-3);
+
+    const unsigned int size = 32;
+    unsigned int       dim  = (size - 1) * (size - 1);
+
+    deallog << "Size " << size << " Unknowns " << dim << std::endl;
+
+    // Make matrix
+    FDMatrix               testproblem(size, size);
+    DynamicSparsityPattern csp(dim, dim);
+    testproblem.five_point_structure(csp);
+    LinearAlgebra::TpetraWrappers::SparseMatrix<double> A;
+    A.reinit(csp);
+    testproblem.five_point(A);
+
+    LinearAlgebra::TpetraWrappers::Vector<double> f;
+    f.reinit(complete_index_set(dim), MPI_COMM_WORLD);
+    LinearAlgebra::TpetraWrappers::Vector<double> u;
+    u.reinit(complete_index_set(dim), MPI_COMM_WORLD);
+    f.trilinos_rcp()->putScalar(1.0);
+    A.compress(VectorOperation::insert);
+    f.compress(VectorOperation::insert);
+    u.compress(VectorOperation::insert);
+
+    GrowingVectorMemory<LinearAlgebra::TpetraWrappers::Vector<double>> mem;
+    SolverCG<LinearAlgebra::TpetraWrappers::Vector<double>> solver(control,
+                                                                   mem);
+    PreconditionIdentity                                    preconditioner;
+
+    deallog << "Solver type: " << typeid(solver).name() << std::endl;
+
+    check_solver_within_range(solver.solve(A, u, f, preconditioner),
+                              control.last_step(),
+                              42,
+                              44);
+  }
+}
diff --git a/tests/trilinos_tpetra/deal_solver_01.output b/tests/trilinos_tpetra/deal_solver_01.output
new file mode 100644 (file)
index 0000000..9ba3a2d
--- /dev/null
@@ -0,0 +1,4 @@
+
+DEAL::Size 32 Unknowns 961
+DEAL::Solver type: N6dealii8SolverCGINS_13LinearAlgebra14TpetraWrappers6VectorIdNS_11MemorySpace4HostEEEEE
+DEAL::Solver stopped within 42 - 44 iterations
diff --git a/tests/trilinos_tpetra/deal_solver_02.cc b/tests/trilinos_tpetra/deal_solver_02.cc
new file mode 100644 (file)
index 0000000..8418a24
--- /dev/null
@@ -0,0 +1,88 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2004 - 2024 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 the BiCGStab solver using the Trilinos matrix and vector classes
+
+
+
+#include <deal.II/base/utilities.h>
+
+#include <deal.II/lac/precondition.h>
+#include <deal.II/lac/solver_bicgstab.h>
+#include <deal.II/lac/solver_cg.h>
+#include <deal.II/lac/solver_control.h>
+#include <deal.II/lac/solver_gmres.h>
+#include <deal.II/lac/solver_qmrs.h>
+#include <deal.II/lac/solver_richardson.h>
+#include <deal.II/lac/trilinos_tpetra_sparse_matrix.h>
+#include <deal.II/lac/vector_memory.h>
+
+#include <iostream>
+#include <typeinfo>
+
+#include "../tests.h"
+
+#include "../testmatrix.h"
+
+
+int
+main(int argc, char **argv)
+{
+  initlog();
+  deallog << std::setprecision(4);
+
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(
+    argc, argv, testing_max_num_threads());
+
+
+  {
+    SolverControl control(200, 1.3e-10);
+
+    const unsigned int size = 32;
+    unsigned int       dim  = (size - 1) * (size - 1);
+
+    deallog << "Size " << size << " Unknowns " << dim << std::endl;
+
+    // Make matrix
+    FDMatrix               testproblem(size, size);
+    DynamicSparsityPattern csp(dim, dim);
+    testproblem.five_point_structure(csp);
+    LinearAlgebra::TpetraWrappers::SparseMatrix<double> A;
+    A.reinit(csp);
+    testproblem.five_point(A);
+
+    LinearAlgebra::TpetraWrappers::Vector<double> f;
+    f.reinit(complete_index_set(dim), MPI_COMM_WORLD);
+    LinearAlgebra::TpetraWrappers::Vector<double> u;
+    u.reinit(complete_index_set(dim), MPI_COMM_WORLD);
+    f.trilinos_rcp()->putScalar(1.0);
+    A.compress(VectorOperation::insert);
+    f.compress(VectorOperation::insert);
+    u.compress(VectorOperation::insert);
+
+    GrowingVectorMemory<LinearAlgebra::TpetraWrappers::Vector<double>> mem;
+    SolverBicgstab<LinearAlgebra::TpetraWrappers::Vector<double>>      solver(
+      control, mem);
+    PreconditionIdentity preconditioner;
+
+    deallog << "Solver type: " << typeid(solver).name() << std::endl;
+
+    check_solver_within_range(solver.solve(A, u, f, preconditioner),
+                              control.last_step(),
+                              49,
+                              51);
+  }
+}
diff --git a/tests/trilinos_tpetra/deal_solver_02.output b/tests/trilinos_tpetra/deal_solver_02.output
new file mode 100644 (file)
index 0000000..ceafaf1
--- /dev/null
@@ -0,0 +1,4 @@
+
+DEAL::Size 32 Unknowns 961
+DEAL::Solver type: N6dealii14SolverBicgstabINS_13LinearAlgebra14TpetraWrappers6VectorIdNS_11MemorySpace4HostEEEEE
+DEAL::Solver stopped within 49 - 51 iterations
diff --git a/tests/trilinos_tpetra/deal_solver_03.cc b/tests/trilinos_tpetra/deal_solver_03.cc
new file mode 100644 (file)
index 0000000..52cc64f
--- /dev/null
@@ -0,0 +1,90 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2004 - 2024 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 the GMRES solver using the Trilinos matrix and vector classes
+
+
+
+#include <deal.II/base/utilities.h>
+
+#include <deal.II/lac/precondition.h>
+#include <deal.II/lac/solver.h>
+#include <deal.II/lac/solver_bicgstab.h>
+#include <deal.II/lac/solver_cg.h>
+#include <deal.II/lac/solver_control.h>
+#include <deal.II/lac/solver_gmres.h>
+#include <deal.II/lac/solver_qmrs.h>
+#include <deal.II/lac/solver_richardson.h>
+#include <deal.II/lac/trilinos_tpetra_sparse_matrix.h>
+#include <deal.II/lac/trilinos_vector.h>
+#include <deal.II/lac/vector_memory.h>
+
+#include <iostream>
+#include <typeinfo>
+
+#include "../tests.h"
+
+#include "../testmatrix.h"
+
+
+int
+main(int argc, char **argv)
+{
+  initlog();
+  deallog << std::setprecision(2);
+
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(
+    argc, argv, testing_max_num_threads());
+
+
+  {
+    SolverControl control(2000, 1.e-3);
+
+    const unsigned int size = 32;
+    unsigned int       dim  = (size - 1) * (size - 1);
+
+    deallog << "Size " << size << " Unknowns " << dim << std::endl;
+
+    // Make matrix
+    FDMatrix               testproblem(size, size);
+    DynamicSparsityPattern csp(dim, dim);
+    testproblem.five_point_structure(csp);
+    LinearAlgebra::TpetraWrappers::SparseMatrix<double> A;
+    A.reinit(csp);
+    testproblem.five_point(A);
+
+    LinearAlgebra::TpetraWrappers::Vector<double> f;
+    f.reinit(complete_index_set(dim), MPI_COMM_WORLD);
+    LinearAlgebra::TpetraWrappers::Vector<double> u;
+    u.reinit(complete_index_set(dim), MPI_COMM_WORLD);
+    f.trilinos_rcp()->putScalar(1.0);
+    A.compress(VectorOperation::insert);
+    f.compress(VectorOperation::insert);
+    u.compress(VectorOperation::insert);
+
+    GrowingVectorMemory<LinearAlgebra::TpetraWrappers::Vector<double>> mem;
+    SolverGMRES<LinearAlgebra::TpetraWrappers::Vector<double>> solver(control,
+                                                                      mem);
+    PreconditionIdentity                                       preconditioner;
+
+    deallog << "Solver type: " << typeid(solver).name() << std::endl;
+
+    check_solver_within_range(solver.solve(A, u, f, preconditioner),
+                              control.last_step(),
+                              74,
+                              76);
+  }
+}
diff --git a/tests/trilinos_tpetra/deal_solver_03.output b/tests/trilinos_tpetra/deal_solver_03.output
new file mode 100644 (file)
index 0000000..9d72b0e
--- /dev/null
@@ -0,0 +1,4 @@
+
+DEAL::Size 32 Unknowns 961
+DEAL::Solver type: N6dealii11SolverGMRESINS_13LinearAlgebra14TpetraWrappers6VectorIdNS_11MemorySpace4HostEEEEE
+DEAL::Solver stopped within 74 - 76 iterations
diff --git a/tests/trilinos_tpetra/deal_solver_04.cc b/tests/trilinos_tpetra/deal_solver_04.cc
new file mode 100644 (file)
index 0000000..92b55a5
--- /dev/null
@@ -0,0 +1,89 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2004 - 2024 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 the MINRES solver using the Trilinos matrix and vector classes
+
+
+#include <deal.II/base/utilities.h>
+
+#include <deal.II/lac/precondition.h>
+#include <deal.II/lac/solver_bicgstab.h>
+#include <deal.II/lac/solver_cg.h>
+#include <deal.II/lac/solver_control.h>
+#include <deal.II/lac/solver_gmres.h>
+#include <deal.II/lac/solver_minres.h>
+#include <deal.II/lac/solver_qmrs.h>
+#include <deal.II/lac/solver_richardson.h>
+#include <deal.II/lac/trilinos_tpetra_sparse_matrix.h>
+#include <deal.II/lac/trilinos_vector.h>
+#include <deal.II/lac/vector_memory.h>
+
+#include <iostream>
+#include <typeinfo>
+
+#include "../tests.h"
+
+#include "../testmatrix.h"
+
+
+int
+main(int argc, char **argv)
+{
+  initlog();
+  deallog << std::setprecision(4);
+
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(
+    argc, argv, testing_max_num_threads());
+
+
+  {
+    SolverControl control(100, 1.e-3);
+
+    const unsigned int size = 32;
+    unsigned int       dim  = (size - 1) * (size - 1);
+
+    deallog << "Size " << size << " Unknowns " << dim << std::endl;
+
+    // Make matrix
+    FDMatrix               testproblem(size, size);
+    DynamicSparsityPattern csp(dim, dim);
+    testproblem.five_point_structure(csp);
+    LinearAlgebra::TpetraWrappers::SparseMatrix<double> A;
+    A.reinit(csp);
+    testproblem.five_point(A);
+
+    LinearAlgebra::TpetraWrappers::Vector<double> f;
+    f.reinit(complete_index_set(dim), MPI_COMM_WORLD);
+    LinearAlgebra::TpetraWrappers::Vector<double> u;
+    u.reinit(complete_index_set(dim), MPI_COMM_WORLD);
+    f.trilinos_rcp()->putScalar(1.0);
+    A.compress(VectorOperation::insert);
+    f.compress(VectorOperation::insert);
+    u.compress(VectorOperation::insert);
+
+    GrowingVectorMemory<LinearAlgebra::TpetraWrappers::Vector<double>> mem;
+    SolverMinRes<LinearAlgebra::TpetraWrappers::Vector<double>> solver(control,
+                                                                       mem);
+    PreconditionIdentity                                        preconditioner;
+
+    deallog << "Solver type: " << typeid(solver).name() << std::endl;
+
+    check_solver_within_range(solver.solve(A, u, f, preconditioner),
+                              control.last_step(),
+                              0,
+                              42);
+  }
+}
diff --git a/tests/trilinos_tpetra/deal_solver_04.output b/tests/trilinos_tpetra/deal_solver_04.output
new file mode 100644 (file)
index 0000000..662bfbe
--- /dev/null
@@ -0,0 +1,4 @@
+
+DEAL::Size 32 Unknowns 961
+DEAL::Solver type: N6dealii12SolverMinResINS_13LinearAlgebra14TpetraWrappers6VectorIdNS_11MemorySpace4HostEEEEE
+DEAL::Solver stopped after 43 iterations
diff --git a/tests/trilinos_tpetra/deal_solver_05.cc b/tests/trilinos_tpetra/deal_solver_05.cc
new file mode 100644 (file)
index 0000000..b94f4ff
--- /dev/null
@@ -0,0 +1,87 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2004 - 2024 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 the QMRS solver using the Trilinos matrix and vector classes
+
+
+#include <deal.II/base/utilities.h>
+
+#include <deal.II/lac/precondition.h>
+#include <deal.II/lac/solver_bicgstab.h>
+#include <deal.II/lac/solver_cg.h>
+#include <deal.II/lac/solver_control.h>
+#include <deal.II/lac/solver_gmres.h>
+#include <deal.II/lac/solver_qmrs.h>
+#include <deal.II/lac/solver_richardson.h>
+#include <deal.II/lac/trilinos_tpetra_sparse_matrix.h>
+#include <deal.II/lac/trilinos_vector.h>
+#include <deal.II/lac/vector_memory.h>
+
+#include <iostream>
+#include <typeinfo>
+
+#include "../tests.h"
+
+#include "../testmatrix.h"
+
+
+int
+main(int argc, char **argv)
+{
+  initlog();
+  deallog << std::setprecision(4);
+
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(
+    argc, argv, testing_max_num_threads());
+
+
+  {
+    SolverControl control(200, 1.e-3);
+
+    const unsigned int size = 32;
+    unsigned int       dim  = (size - 1) * (size - 1);
+
+    deallog << "Size " << size << " Unknowns " << dim << std::endl;
+
+    // Make matrix
+    FDMatrix               testproblem(size, size);
+    DynamicSparsityPattern csp(dim, dim);
+    testproblem.five_point_structure(csp);
+    LinearAlgebra::TpetraWrappers::SparseMatrix<double> A;
+    A.reinit(csp);
+    testproblem.five_point(A);
+
+    LinearAlgebra::TpetraWrappers::Vector<double> f;
+    f.reinit(complete_index_set(dim), MPI_COMM_WORLD);
+    LinearAlgebra::TpetraWrappers::Vector<double> u;
+    u.reinit(complete_index_set(dim), MPI_COMM_WORLD);
+    f.trilinos_rcp()->putScalar(1.0);
+    A.compress(VectorOperation::insert);
+    f.compress(VectorOperation::insert);
+    u.compress(VectorOperation::insert);
+
+    GrowingVectorMemory<LinearAlgebra::TpetraWrappers::Vector<double>> mem;
+    SolverQMRS<LinearAlgebra::TpetraWrappers::Vector<double>> solver(control,
+                                                                     mem);
+    PreconditionIdentity                                      preconditioner;
+    deallog << "Solver type: " << typeid(solver).name() << std::endl;
+
+    check_solver_within_range(solver.solve(A, u, f, preconditioner),
+                              control.last_step(),
+                              45,
+                              47);
+  }
+}
diff --git a/tests/trilinos_tpetra/deal_solver_05.output b/tests/trilinos_tpetra/deal_solver_05.output
new file mode 100644 (file)
index 0000000..d15b2da
--- /dev/null
@@ -0,0 +1,4 @@
+
+DEAL::Size 32 Unknowns 961
+DEAL::Solver type: N6dealii10SolverQMRSINS_13LinearAlgebra14TpetraWrappers6VectorIdNS_11MemorySpace4HostEEEEE
+DEAL::Solver stopped within 45 - 47 iterations
diff --git a/tests/trilinos_tpetra/deal_solver_06.cc b/tests/trilinos_tpetra/deal_solver_06.cc
new file mode 100644 (file)
index 0000000..08c3f24
--- /dev/null
@@ -0,0 +1,89 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2004 - 2024 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 the Richardson solver using the Trilinos matrix and vector classes
+
+
+#include <deal.II/base/utilities.h>
+
+#include <deal.II/lac/precondition.h>
+#include <deal.II/lac/solver_bicgstab.h>
+#include <deal.II/lac/solver_cg.h>
+#include <deal.II/lac/solver_control.h>
+#include <deal.II/lac/solver_gmres.h>
+#include <deal.II/lac/solver_qmrs.h>
+#include <deal.II/lac/solver_richardson.h>
+#include <deal.II/lac/trilinos_tpetra_sparse_matrix.h>
+#include <deal.II/lac/trilinos_vector.h>
+#include <deal.II/lac/vector_memory.h>
+
+#include <iostream>
+#include <typeinfo>
+
+#include "../tests.h"
+
+#include "../testmatrix.h"
+
+
+int
+main(int argc, char **argv)
+{
+  initlog();
+  deallog << std::setprecision(4);
+
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(
+    argc, argv, testing_max_num_threads());
+
+
+  {
+    const unsigned int size = 32;
+    unsigned int       dim  = (size - 1) * (size - 1);
+    SolverControl      control(10000, 1.e-3);
+
+    deallog << "Size " << size << " Unknowns " << dim << std::endl;
+
+    // Make matrix
+    FDMatrix               testproblem(size, size);
+    DynamicSparsityPattern csp(dim, dim);
+    testproblem.five_point_structure(csp);
+    LinearAlgebra::TpetraWrappers::SparseMatrix<double> A;
+    A.reinit(csp);
+    testproblem.five_point(A);
+
+    LinearAlgebra::TpetraWrappers::Vector<double> f;
+    f.reinit(complete_index_set(dim), MPI_COMM_WORLD);
+    LinearAlgebra::TpetraWrappers::Vector<double> u;
+    u.reinit(complete_index_set(dim), MPI_COMM_WORLD);
+    f.trilinos_rcp()->putScalar(1.0);
+    A.compress(VectorOperation::insert);
+    f.compress(VectorOperation::insert);
+    u.compress(VectorOperation::insert);
+
+    GrowingVectorMemory<LinearAlgebra::TpetraWrappers::Vector<double>> mem;
+    SolverRichardson<LinearAlgebra::TpetraWrappers::Vector<double>>::
+      AdditionalData data(
+        /*omega=*/0.1);
+    SolverRichardson<LinearAlgebra::TpetraWrappers::Vector<double>> solver(
+      control, mem, data);
+    PreconditionIdentity preconditioner;
+    deallog << "Solver type: " << typeid(solver).name() << std::endl;
+
+    check_solver_within_range(solver.solve(A, u, f, preconditioner),
+                              control.last_step(),
+                              5269,
+                              5271);
+  }
+}
diff --git a/tests/trilinos_tpetra/deal_solver_06.output b/tests/trilinos_tpetra/deal_solver_06.output
new file mode 100644 (file)
index 0000000..039bff2
--- /dev/null
@@ -0,0 +1,4 @@
+
+DEAL::Size 32 Unknowns 961
+DEAL::Solver type: N6dealii16SolverRichardsonINS_13LinearAlgebra14TpetraWrappers6VectorIdNS_11MemorySpace4HostEEEEE
+DEAL::Solver stopped within 5269 - 5271 iterations

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.