From 4640cc53351d46a306ee5bf8b7050a3efae2c750 Mon Sep 17 00:00:00 2001 From: Martin Kronbichler Date: Wed, 28 Nov 2018 17:55:31 +0100 Subject: [PATCH] Add test. --- tests/matrix_free/matrix_vector_faces_24.cc | 109 ++++++++++++++++++ ...h_mpi=true.with_p4est=true.mpirun=2.output | 5 + 2 files changed, 114 insertions(+) create mode 100644 tests/matrix_free/matrix_vector_faces_24.cc create mode 100644 tests/matrix_free/matrix_vector_faces_24.with_mpi=true.with_p4est=true.mpirun=2.output diff --git a/tests/matrix_free/matrix_vector_faces_24.cc b/tests/matrix_free/matrix_vector_faces_24.cc new file mode 100644 index 0000000000..206b135be7 --- /dev/null +++ b/tests/matrix_free/matrix_vector_faces_24.cc @@ -0,0 +1,109 @@ +// --------------------------------------------------------------------- +// +// 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.md at +// the top level directory of deal.II. +// +// --------------------------------------------------------------------- + + + +// check the setup of ghost faces and elements in case the flag +// MatrixFree::AdditionalData::hold_all_faces_to_owned_cells is set to true + +#include +#include + +#include + +#include + +#include +#include + +#include + +#include "../tests.h" +#include "create_mesh.h" + +std::ofstream logfile("output"); + +#include "matrix_vector_faces_common.h" + + + +template +void +test() +{ + if (fe_degree > 1) + return; + + // raise element degree by two to test cubic and quartic shape functions + // rather than linears and quadratics according to + // matrix_vector_faces_common.h + + parallel::distributed::Triangulation tria(MPI_COMM_WORLD); + GridGenerator::hyper_cube(tria); + tria.refine_global(1); + tria.last()->set_refine_flag(); + tria.execute_coarsening_and_refinement(); + + FE_DGQ fe(fe_degree); + DoFHandler dof(tria); + dof.distribute_dofs(fe); + AffineConstraints constraints; + constraints.close(); + + deallog << "Testing " << dof.get_fe().get_name(); + deallog << std::endl; + // std::cout << "Number of cells: " << + // dof.get_triangulation().n_active_cells() << std::endl; std::cout << "Number + // of degrees of freedom: " << dof.n_dofs() << std::endl; std::cout << "Number + // of constraints: " << constraints.n_constraints() << std::endl; + + MappingQ mapping(dof.get_fe().degree + 1); + + LinearAlgebra::distributed::Vector in, out; + + MatrixFree mf_data; + const QGauss<1> quad(fe_degree + 1); + typename MatrixFree::AdditionalData data; + data.tasks_parallel_scheme = MatrixFree::AdditionalData::none; + data.tasks_block_size = 3; + data.mapping_update_flags_inner_faces = + (update_gradients | update_JxW_values); + data.mapping_update_flags_boundary_faces = + (update_gradients | update_JxW_values); + data.hold_all_faces_to_owned_cells = true; + + mf_data.reinit(mapping, dof, constraints, quad, data); + + mf_data.initialize_dof_vector(in); + mf_data.initialize_dof_vector(out); + + // Set random seed for reproducibility + Testing::srand(42); + for (unsigned int i = 0; i < in.local_size(); ++i) + { + const double entry = Testing::rand() / (double)RAND_MAX; + in.local_element(i) = entry; + } + + MatrixFreeTest> + mf(mf_data); + mf.vmult(out, in); + + deallog << "Norm of result: " << out.l2_norm() << std::endl; +} diff --git a/tests/matrix_free/matrix_vector_faces_24.with_mpi=true.with_p4est=true.mpirun=2.output b/tests/matrix_free/matrix_vector_faces_24.with_mpi=true.with_p4est=true.mpirun=2.output new file mode 100644 index 0000000000..71f9902dab --- /dev/null +++ b/tests/matrix_free/matrix_vector_faces_24.with_mpi=true.with_p4est=true.mpirun=2.output @@ -0,0 +1,5 @@ + +DEAL:2d::Testing FE_DGQ<2>(1) +DEAL:2d::Norm of result: 4.18 +DEAL:3d::Testing FE_DGQ<3>(1) +DEAL:3d::Norm of result: 3.39 -- 2.39.5