From 59b100d6d0e9911a3862039710c6d2da9dfe45bd Mon Sep 17 00:00:00 2001 From: Wolfgang Bangerth Date: Mon, 11 Jan 2016 15:50:12 -0600 Subject: [PATCH] Add test. --- tests/mpi/data_out_faces_01.cc | 194 ++++++++++++++++++ ...t_faces_01.mpirun=3.with_petsc=true.output | 5 + 2 files changed, 199 insertions(+) create mode 100644 tests/mpi/data_out_faces_01.cc create mode 100644 tests/mpi/data_out_faces_01.mpirun=3.with_petsc=true.output diff --git a/tests/mpi/data_out_faces_01.cc b/tests/mpi/data_out_faces_01.cc new file mode 100644 index 0000000000..f570aac3c8 --- /dev/null +++ b/tests/mpi/data_out_faces_01.cc @@ -0,0 +1,194 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2016 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 at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + +// tests DataOutFaces in parallel + +#include "../tests.h" +#include +#include +#include +#include + +#include +#include +#include +#include + +#include +#include +#include + +#include +#include + +#include +#include +#include +#include +#include + +namespace pdd +{ + using namespace dealii; + + // @sect3{The PDDProblem class template} + + template + class PDDProblem + { + public: + PDDProblem (); + ~PDDProblem (); + void run (); + private: + void setup_system (); + void output_results (); + + MPI_Comm mpi_communicator; + parallel::distributed::Triangulation triangulation; + //Triangulation triangulation; + DoFHandler dof_handler; + FE_Q fe; + + IndexSet locally_owned_dofs; + IndexSet locally_relevant_dofs; + + ConstraintMatrix constraints; + + PETScWrappers::MPI::Vector locally_relevant_solution; + }; + + // The constructor and destructor of the class + template + PDDProblem::PDDProblem () : + mpi_communicator (MPI_COMM_WORLD), + triangulation (mpi_communicator), + dof_handler (triangulation), + fe (1) + {} + + + template + PDDProblem::~PDDProblem () + { + dof_handler.clear (); + } + + // @sect4{PDDProblem::setup_system} + + template + void PDDProblem::setup_system () + { + + // Initialize the mesh + std::vector repetitions; + repetitions.push_back(10); + repetitions.push_back(12); + + GridGenerator::subdivided_hyper_rectangle(triangulation,repetitions,Point (-1.,-1.),Point (1.,1.),true); + // triangulation.refine_global (1); + + // Print out the mesh + dof_handler.distribute_dofs (fe); + + // Initialize the solution + locally_owned_dofs = dof_handler.locally_owned_dofs (); + DoFTools::extract_locally_relevant_dofs (dof_handler,locally_relevant_dofs); + locally_relevant_solution.reinit(locally_owned_dofs,mpi_communicator); + + locally_relevant_solution = 0; + + // Apply a constant function at the boundary + constraints.clear (); + constraints.reinit (locally_relevant_dofs); + DoFTools::make_hanging_node_constraints (dof_handler,constraints); + VectorTools::interpolate_boundary_values (dof_handler, 0, ConstantFunction (1.0), constraints); + constraints.close (); + constraints.distribute(locally_relevant_solution); + } + + // Generate the outputs + template + void PDDProblem::output_results () + { + + // First generate an output for the cells + const unsigned int cycle = 1; + + PETScWrappers::MPI::Vector solution (locally_owned_dofs,locally_relevant_dofs,mpi_communicator); + solution = locally_relevant_solution; + + DataOut data_out; + data_out.attach_dof_handler (dof_handler); + data_out.add_data_vector (solution, "u"); + + // make sure the following works in parallel + data_out.build_patches (); + + + // Then generate the output for the faces. + DataOutFaces data_out_faces; + data_out_faces.attach_dof_handler (dof_handler); + data_out_faces.add_data_vector (solution, "u"); + // make sure the following works in parallel + data_out_faces.build_patches (fe.degree); + } + // @sect4{PDDProblem::run} + + template + void PDDProblem::run () + { + setup_system (); + output_results (); + } +} + + +int main(int argc, char *argv[]) +{ + Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, 1); + MPILogInitAll log; + + try + { + pdd::PDDProblem<2> laplace_problem_2d; + laplace_problem_2d.run (); + + return 0; + } + catch (std::exception &exc) + { + deallog << std::endl << std::endl + << "----------------------------------------------------" + << std::endl; + deallog << "Exception on processing: " << std::endl + << exc.what() << std::endl + << "Aborting!" << std::endl + << "----------------------------------------------------" + << std::endl; + return 1; + } + catch (...) + { + deallog << std::endl << std::endl + << "----------------------------------------------------" + << std::endl; + deallog << "Unknown exception!" << std::endl + << "Aborting!" << std::endl + << "----------------------------------------------------" + << std::endl; + return 1; + }; +} diff --git a/tests/mpi/data_out_faces_01.mpirun=3.with_petsc=true.output b/tests/mpi/data_out_faces_01.mpirun=3.with_petsc=true.output new file mode 100644 index 0000000000..3f2ff2d6cc --- /dev/null +++ b/tests/mpi/data_out_faces_01.mpirun=3.with_petsc=true.output @@ -0,0 +1,5 @@ + + + + + -- 2.39.5