From 2a36e3d5ee64c3d19df5c8fed09160b631e29e03 Mon Sep 17 00:00:00 2001 From: bangerth Date: Sun, 26 Jun 2011 17:33:12 +0000 Subject: [PATCH] New test that currently fails. git-svn-id: https://svn.dealii.org/trunk@23865 0785d39b-7218-0410-832d-ea1e28bc413d --- tests/mpi/solution_transfer_01.cc | 132 ++++++++++++++++++ .../solution_transfer_01/ncpu_10/cmp/generic | 5 + .../solution_transfer_01/ncpu_4/cmp/generic | 5 + 3 files changed, 142 insertions(+) create mode 100644 tests/mpi/solution_transfer_01.cc create mode 100644 tests/mpi/solution_transfer_01/ncpu_10/cmp/generic create mode 100644 tests/mpi/solution_transfer_01/ncpu_4/cmp/generic diff --git a/tests/mpi/solution_transfer_01.cc b/tests/mpi/solution_transfer_01.cc new file mode 100644 index 0000000000..b2a78f7397 --- /dev/null +++ b/tests/mpi/solution_transfer_01.cc @@ -0,0 +1,132 @@ +//--------------------------------------------------------------------------- +// $Id$ +// Version: $Name$ +// +// Copyright (C) 2009, 2010, 2011 by the deal.II authors +// +// This file is subject to QPL and may not be distributed +// without copyright and license information. Please refer +// to the file deal.II/doc/license.html for the text and +// further information on this license. +// +//--------------------------------------------------------------------------- + + +// SolutionTransfer locks up when a process has no locally owned cells + +#include "../tests.h" +#include "coarse_grid_common.h" +#include +#include +#include +#include +#include + +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include + +#include +#include + +#include +#include + + +template +void test() +{ + parallel::distributed::Triangulation<2> tria(MPI_COMM_WORLD, + typename Triangulation<2>::MeshSmoothing + (Triangulation<2>::smoothing_on_refinement | + Triangulation<2>::smoothing_on_coarsening)); + + GridGenerator::hyper_cube (tria,-1.0,1.0); + + DoFHandler<2> dh(tria); + FE_Q<2> fe(1); + + dh.distribute_dofs(fe); + + IndexSet locally_owned_dofs = dh.locally_owned_dofs (); + IndexSet locally_relevant_dofs; + + DoFTools::extract_locally_relevant_dofs (dh,locally_relevant_dofs); + + PETScWrappers::MPI::Vector solution(MPI_COMM_WORLD,locally_owned_dofs,locally_relevant_dofs); + solution = 0; + solution.update_ghost_values(); + + parallel::distributed::SolutionTransfer<2,PETScWrappers::MPI::Vector> soltrans(dh); + + tria.set_all_refine_flags(); + tria.prepare_coarsening_and_refinement(); + + soltrans.prepare_for_coarsening_and_refinement (solution); + + + tria.execute_coarsening_and_refinement(); + + dh.distribute_dofs (fe); + locally_owned_dofs = dh.locally_owned_dofs (); + DoFTools::extract_locally_relevant_dofs (dh,locally_relevant_dofs); + + PETScWrappers::MPI::Vector tmp(MPI_COMM_WORLD,dh.n_dofs(),dh.n_locally_owned_dofs()); + + std::cout<<"I WAS HERE"<(); + deallog.pop(); + deallog.push("3d"); + test<3>(); + deallog.pop(); + } + else + { + test<2>(); + test<3>(); + } + + PetscFinalize(); +} diff --git a/tests/mpi/solution_transfer_01/ncpu_10/cmp/generic b/tests/mpi/solution_transfer_01/ncpu_10/cmp/generic new file mode 100644 index 0000000000..801195966b --- /dev/null +++ b/tests/mpi/solution_transfer_01/ncpu_10/cmp/generic @@ -0,0 +1,5 @@ + +DEAL:0:2d::0 4 +DEAL:0:2d::2.00000 +DEAL:0:3d::0 8 +DEAL:0:3d::2.82843 diff --git a/tests/mpi/solution_transfer_01/ncpu_4/cmp/generic b/tests/mpi/solution_transfer_01/ncpu_4/cmp/generic new file mode 100644 index 0000000000..801195966b --- /dev/null +++ b/tests/mpi/solution_transfer_01/ncpu_4/cmp/generic @@ -0,0 +1,5 @@ + +DEAL:0:2d::0 4 +DEAL:0:2d::2.00000 +DEAL:0:3d::0 8 +DEAL:0:3d::2.82843 -- 2.39.5