]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add test case
authorMartin Kronbichler <martin.kronbichler@uni-a.de>
Sat, 27 Aug 2022 16:31:39 +0000 (18:31 +0200)
committerMartin Kronbichler <martin.kronbichler@uni-a.de>
Mon, 29 Aug 2022 20:08:03 +0000 (22:08 +0200)
tests/matrix_free/laplace_operator_04.cc [new file with mode: 0644]
tests/matrix_free/laplace_operator_04.output [new file with mode: 0644]

diff --git a/tests/matrix_free/laplace_operator_04.cc b/tests/matrix_free/laplace_operator_04.cc
new file mode 100644 (file)
index 0000000..01098d2
--- /dev/null
@@ -0,0 +1,106 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 - 2022 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.
+//
+// ---------------------------------------------------------------------
+
+
+
+// check the diagonal of LaplaceOperator for the scalar and vector-valued case
+
+#include <deal.II/base/utilities.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/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+
+#include <deal.II/lac/affine_constraints.h>
+
+#include <deal.II/matrix_free/operators.h>
+
+#include <iostream>
+
+#include "../tests.h"
+
+
+template <int dim, int fe_degree>
+void
+test()
+{
+  using number = double;
+
+  Triangulation<dim> tria;
+  GridGenerator::hyper_cube(tria);
+  tria.refine_global(2);
+
+  FE_Q<dim>       fe_q(fe_degree);
+  FESystem<dim>   fe(fe_q, dim);
+  DoFHandler<dim> dof(tria);
+  dof.distribute_dofs(fe);
+
+  AffineConstraints<double> constraints;
+  constraints.close();
+
+  deallog << "Testing " << dof.get_fe().get_name() << std::endl;
+
+  std::shared_ptr<MatrixFree<dim, number>> mf_data(
+    new MatrixFree<dim, number>());
+  {
+    const QGauss<1>                                  quad(fe_degree + 1);
+    typename MatrixFree<dim, number>::AdditionalData data;
+    data.tasks_parallel_scheme = MatrixFree<dim, number>::AdditionalData::none;
+    data.mapping_update_flags =
+      update_quadrature_points | update_gradients | update_JxW_values;
+    mf_data->reinit(MappingQ1<dim>{}, dof, constraints, quad, data);
+  }
+
+  MatrixFreeOperators::LaplaceOperator<
+    dim,
+    fe_degree,
+    fe_degree + 1,
+    1,
+    LinearAlgebra::distributed::Vector<number>>
+    mf_scalar;
+  MatrixFreeOperators::LaplaceOperator<
+    dim,
+    fe_degree,
+    fe_degree + 1,
+    dim,
+    LinearAlgebra::distributed::Vector<number>>
+    mf_vector;
+  mf_scalar.initialize(mf_data);
+  mf_vector.initialize(mf_data);
+  mf_scalar.compute_diagonal();
+  mf_vector.compute_diagonal();
+  deallog << "Matrix diagonal scalar: " << std::endl;
+  mf_scalar.get_matrix_diagonal()->get_vector().print(
+    deallog.get_file_stream());
+  deallog << "Matrix diagonal vector: " << std::endl;
+  mf_vector.get_matrix_diagonal()->get_vector().print(
+    deallog.get_file_stream());
+}
+
+
+int
+main()
+{
+  initlog();
+
+  deallog << std::setprecision(4);
+
+  test<2, 2>();
+}
diff --git a/tests/matrix_free/laplace_operator_04.output b/tests/matrix_free/laplace_operator_04.output
new file mode 100644 (file)
index 0000000..cef97f8
--- /dev/null
@@ -0,0 +1,12 @@
+
+DEAL::Testing FESystem<2>[FE_Q<2>(2)^2]
+DEAL::Matrix diagonal scalar: 
+Process #0
+Local range: [0, 162), global size: 162
+Vector data:
+6.222e-01 0.000e+00 1.244e+00 0.000e+00 1.244e+00 0.000e+00 2.489e+00 0.000e+00 1.956e+00 0.000e+00 3.911e+00 0.000e+00 1.956e+00 0.000e+00 3.911e+00 0.000e+00 5.689e+00 0.000e+00 1.244e+00 0.000e+00 2.489e+00 0.000e+00 3.911e+00 0.000e+00 1.956e+00 0.000e+00 3.911e+00 0.000e+00 5.689e+00 0.000e+00 1.244e+00 0.000e+00 2.489e+00 0.000e+00 1.956e+00 0.000e+00 3.911e+00 0.000e+00 3.911e+00 0.000e+00 5.689e+00 0.000e+00 2.489e+00 0.000e+00 3.911e+00 0.000e+00 3.911e+00 0.000e+00 5.689e+00 0.000e+00 1.244e+00 0.000e+00 2.489e+00 0.000e+00 3.911e+00 0.000e+00 1.956e+00 0.000e+00 3.911e+00 0.000e+00 5.689e+00 0.000e+00 6.222e-01 0.000e+00 1.244e+00 0.000e+00 1.956e+00 0.000e+00 1.956e+00 0.000e+00 3.911e+00 0.000e+00 5.689e+00 0.000e+00 2.489e+00 0.000e+00 3.911e+00 0.000e+00 3.911e+00 0.000e+00 5.689e+00 0.000e+00 1.244e+00 0.000e+00 1.956e+00 0.000e+00 3.911e+00 0.000e+00 5.689e+00 0.000e+00 1.244e+00 0.000e+00 2.489e+00 0.000e+00 1.956e+00 0.000e+00 3.911e+00 0.000e+00 3.911e+00 0.000e+00 5.689e+00 0.000e+00 2.489e+00 0.000e+00 3.911e+00 0.000e+00 3.911e+00 0.000e+00 5.689e+00 0.000e+00 6.222e-01 0.000e+00 1.244e+00 0.000e+00 1.956e+00 0.000e+00 3.911e+00 0.000e+00 1.956e+00 0.000e+00 5.689e+00 0.000e+00 1.244e+00 0.000e+00 3.911e+00 0.000e+00 1.956e+00 0.000e+00 5.689e+00 0.000e+00 2.489e+00 0.000e+00 3.911e+00 0.000e+00 3.911e+00 0.000e+00 5.689e+00 0.000e+00 1.244e+00 0.000e+00 1.956e+00 0.000e+00 3.911e+00 0.000e+00 5.689e+00 0.000e+00 1.244e+00 0.000e+00 3.911e+00 0.000e+00 1.956e+00 0.000e+00 5.689e+00 0.000e+00 6.222e-01 0.000e+00 1.956e+00 0.000e+00 1.956e+00 0.000e+00 5.689e+00 0.000e+00 
+DEAL::Matrix diagonal vector: 
+Process #0
+Local range: [0, 162), global size: 162
+Vector data:
+6.222e-01 6.222e-01 1.244e+00 1.244e+00 1.244e+00 1.244e+00 2.489e+00 2.489e+00 1.956e+00 1.956e+00 3.911e+00 3.911e+00 1.956e+00 1.956e+00 3.911e+00 3.911e+00 5.689e+00 5.689e+00 1.244e+00 1.244e+00 2.489e+00 2.489e+00 3.911e+00 3.911e+00 1.956e+00 1.956e+00 3.911e+00 3.911e+00 5.689e+00 5.689e+00 1.244e+00 1.244e+00 2.489e+00 2.489e+00 1.956e+00 1.956e+00 3.911e+00 3.911e+00 3.911e+00 3.911e+00 5.689e+00 5.689e+00 2.489e+00 2.489e+00 3.911e+00 3.911e+00 3.911e+00 3.911e+00 5.689e+00 5.689e+00 1.244e+00 1.244e+00 2.489e+00 2.489e+00 3.911e+00 3.911e+00 1.956e+00 1.956e+00 3.911e+00 3.911e+00 5.689e+00 5.689e+00 6.222e-01 6.222e-01 1.244e+00 1.244e+00 1.956e+00 1.956e+00 1.956e+00 1.956e+00 3.911e+00 3.911e+00 5.689e+00 5.689e+00 2.489e+00 2.489e+00 3.911e+00 3.911e+00 3.911e+00 3.911e+00 5.689e+00 5.689e+00 1.244e+00 1.244e+00 1.956e+00 1.956e+00 3.911e+00 3.911e+00 5.689e+00 5.689e+00 1.244e+00 1.244e+00 2.489e+00 2.489e+00 1.956e+00 1.956e+00 3.911e+00 3.911e+00 3.911e+00 3.911e+00 5.689e+00 5.689e+00 2.489e+00 2.489e+00 3.911e+00 3.911e+00 3.911e+00 3.911e+00 5.689e+00 5.689e+00 6.222e-01 6.222e-01 1.244e+00 1.244e+00 1.956e+00 1.956e+00 3.911e+00 3.911e+00 1.956e+00 1.956e+00 5.689e+00 5.689e+00 1.244e+00 1.244e+00 3.911e+00 3.911e+00 1.956e+00 1.956e+00 5.689e+00 5.689e+00 2.489e+00 2.489e+00 3.911e+00 3.911e+00 3.911e+00 3.911e+00 5.689e+00 5.689e+00 1.244e+00 1.244e+00 1.956e+00 1.956e+00 3.911e+00 3.911e+00 5.689e+00 5.689e+00 1.244e+00 1.244e+00 3.911e+00 3.911e+00 1.956e+00 1.956e+00 5.689e+00 5.689e+00 6.222e-01 6.222e-01 1.956e+00 1.956e+00 1.956e+00 1.956e+00 5.689e+00 5.689e+00 

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.