From: Martin Kronbichler Date: Tue, 26 Nov 2019 07:43:44 +0000 (+0100) Subject: New test case X-Git-Tag: v9.2.0-rc1~837^2~1 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=4c7241a4e664eb48d2c3a8b35eccd68dfdf8d92b;p=dealii.git New test case --- diff --git a/tests/matrix_free/face_setup_01.cc b/tests/matrix_free/face_setup_01.cc new file mode 100644 index 0000000000..f6e6519d64 --- /dev/null +++ b/tests/matrix_free/face_setup_01.cc @@ -0,0 +1,88 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2019 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 initialization of MatrixFree with face information on a +// continuous element with hanging nodes + + +#include + +#include + +#include +#include + +#include + +#include + +#include + +#include "../tests.h" + +int +main(int argc, char **argv) +{ + using namespace dealii; + + Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1); + MPILogInitAll log; + + const int dim = 3; + + parallel::distributed::Triangulation tria(MPI_COMM_WORLD); + GridGenerator::hyper_cube(tria); + tria.refine_global(2); + + if (Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0) + { + tria.begin_active()->set_refine_flag(); + } + tria.execute_coarsening_and_refinement(); + + DoFHandler dof_handler(tria); + dof_handler.distribute_dofs(FE_Q(3)); + + MatrixFree matrix_free; + typename MatrixFree::AdditionalData data; + + data.mapping_update_flags = update_values; + data.mapping_update_flags_inner_faces = update_values; // problem + data.mapping_update_flags_boundary_faces = update_values; // problem + + AffineConstraints constraint; + constraint.clear(); + DoFTools::make_hanging_node_constraints(dof_handler, constraint); + constraint.close(); + + matrix_free.reinit(dof_handler, constraint, QGauss<1>(4), data); + + LinearAlgebra::distributed::Vector src; + matrix_free.initialize_dof_vector(src); + + deallog << "main partitioner: size=" + << matrix_free.get_dof_info().vector_partitioner->size() + << " local_size=" + << matrix_free.get_dof_info().vector_partitioner->local_size() + << " n_ghosts=" + << matrix_free.get_dof_info().vector_partitioner->n_ghost_indices() + << std::endl; + for (auto &p : matrix_free.get_dof_info().vector_partitioner_face_variants) + deallog << "partitioner: size=" << p->size() + << " local_size=" << p->local_size() + << " n_ghosts=" << p->n_ghost_indices() << std::endl; +} diff --git a/tests/matrix_free/face_setup_01.with_mpi=true.with_p4est=true.mpirun=1.output b/tests/matrix_free/face_setup_01.with_mpi=true.with_p4est=true.mpirun=1.output new file mode 100644 index 0000000000..c446a4ed14 --- /dev/null +++ b/tests/matrix_free/face_setup_01.with_mpi=true.with_p4est=true.mpirun=1.output @@ -0,0 +1,5 @@ + +DEAL:0::main partitioner: size=2506 local_size=2506 n_ghosts=0 +DEAL:0::partitioner: size=2506 local_size=2506 n_ghosts=0 +DEAL:0::partitioner: size=2506 local_size=2506 n_ghosts=0 +DEAL:0::partitioner: size=2506 local_size=2506 n_ghosts=0 diff --git a/tests/matrix_free/face_setup_01.with_mpi=true.with_p4est=true.mpirun=2.output b/tests/matrix_free/face_setup_01.with_mpi=true.with_p4est=true.mpirun=2.output new file mode 100644 index 0000000000..96b0fa317f --- /dev/null +++ b/tests/matrix_free/face_setup_01.with_mpi=true.with_p4est=true.mpirun=2.output @@ -0,0 +1,11 @@ + +DEAL:0::main partitioner: size=2506 local_size=1492 n_ghosts=273 +DEAL:0::partitioner: size=2506 local_size=1492 n_ghosts=0 +DEAL:0::partitioner: size=2506 local_size=1492 n_ghosts=273 +DEAL:0::partitioner: size=2506 local_size=1492 n_ghosts=273 + +DEAL:1::main partitioner: size=2506 local_size=1014 n_ghosts=442 +DEAL:1::partitioner: size=2506 local_size=1014 n_ghosts=169 +DEAL:1::partitioner: size=2506 local_size=1014 n_ghosts=442 +DEAL:1::partitioner: size=2506 local_size=1014 n_ghosts=442 +