From 2204ca401663ce0db9b547c0c03e196d02468850 Mon Sep 17 00:00:00 2001 From: Wolfgang Bangerth Date: Fri, 13 Mar 2015 08:17:41 -0500 Subject: [PATCH] Add GridGenerator::create_triangulation_with_removed_cells(). --- doc/news/changes.h | 6 +++ include/deal.II/grid/grid_generator.h | 38 ++++++++++++++++ source/grid/grid_generator.cc | 45 +++++++++++++++++++ source/grid/grid_generator.inst.in | 6 +++ tests/grid/grid_generator_07.cc | 65 +++++++++++++++++++++++++++ tests/grid/grid_generator_07.output | 32 +++++++++++++ 6 files changed, 192 insertions(+) create mode 100644 tests/grid/grid_generator_07.cc create mode 100644 tests/grid/grid_generator_07.output diff --git a/doc/news/changes.h b/doc/news/changes.h index 909ba4aa08..28ad3ff801 100644 --- a/doc/news/changes.h +++ b/doc/news/changes.h @@ -376,6 +376,12 @@ inconvenience this causes.

Specific improvements

    +
  1. New: GridGenerator::create_triangulation_with_removed_cells() creates + a new mesh out of an existing one by dropping individual cells. +
    + (Wolfgang Bangerth, 2015/03/13) +
  2. +
  3. New: Add MueLu preconditioner from Trilinos. This is a new algebraic multigrid package. The input parameters are almost the same as the ones from ML so that the two preconditioners can be easily swapped. diff --git a/include/deal.II/grid/grid_generator.h b/include/deal.II/grid/grid_generator.h index 11d691e452..6b29f7f018 100644 --- a/include/deal.II/grid/grid_generator.h +++ b/include/deal.II/grid/grid_generator.h @@ -705,6 +705,44 @@ namespace GridGenerator const Triangulation &triangulation_2, Triangulation &result); + /** + * This function creates a triangulation that consists of the same cells + * as are present in the first argument, except those cells that are listed + * in the second argument. The purpose of the function is to generate + * geometries subtractively from the geometry described by an + * existing triangulation. A prototypical case is a 2d domain with + * rectangular holes. This can be achieved by first meshing the entire + * domain and then using this function to get rid of the cells that + * are located at the holes. Likewise, you could create the mesh + * that GridGenerator::hyper_L() produces by starting with a + * GridGenerator::hyper_cube(), refining it once, and then calling + * the current function with a single cell in the second argument. + * + * @param[in] input_triangulation The original triangulation that serves as + * the template from which the new one is to be created. + * @param[in] cells_to_remove A list of cells of the triangulation provided + * as first argument that should be removed (i.e., that should not + * show up in the result. + * @param[out] result The resulting triangulation that consists of the same + * cells as are in @p input_triangulation, with the exception of the cells + * listed in @p cells_to_remove. + * + * @pre Because we cannot create triangulations de novo that contain + * adaptively refined cells, the input triangulation needs to have all + * of its cells on the same level. Oftentimes, this will in fact be + * the coarsest level, but it is allowed to pass in a triangulation + * that has been refined globally a number of times. The output + * triangulation will in that case simply be a mesh with only one + * level that consists of the active cells of the input minus the + * ones listed in the second argument. However, the input triangulation + * must not have been adaptively refined. + */ + template + void + create_triangulation_with_removed_cells (const Triangulation &input_triangulation, + const std::set::active_cell_iterator> &cells_to_remove, + Triangulation &result); + /** * Take a 2d Triangulation that is being extruded in z direction by the diff --git a/source/grid/grid_generator.cc b/source/grid/grid_generator.cc index 6ff40765ff..f6048ba870 100644 --- a/source/grid/grid_generator.cc +++ b/source/grid/grid_generator.cc @@ -3438,6 +3438,51 @@ namespace GridGenerator + template + void + create_triangulation_with_removed_cells (const Triangulation &input_triangulation, + const std::set::active_cell_iterator> &cells_to_remove, + Triangulation &result) + { + // simply copy the vertices; we will later strip those + // that turn out to be unused + std::vector > vertices = input_triangulation.get_vertices(); + + // the loop through the cells and copy stuff, excluding + // the ones we are to remove + std::vector > cells; + for (typename Triangulation::active_cell_iterator + cell = input_triangulation.begin_active(); cell != input_triangulation.end(); ++cell) + if (cells_to_remove.find(cell) == cells_to_remove.end()) + { + Assert (static_cast(cell->level()) == input_triangulation.n_levels()-1, + ExcMessage ("Your input triangulation appears to have " + "adaptively refined cells. This is not allowed. You can " + "only call this function on a triangulation in which " + "all cells are on the same refinement level.")); + + CellData this_cell; + for (unsigned int v=0; v::vertices_per_cell; ++v) + this_cell.vertices[v] = cell->vertex_index(v); + this_cell.material_id = cell->material_id(); + cells.push_back (this_cell); + } + + // throw out duplicated vertices from the two meshes, reorder vertices as + // necessary and create the triangulation + SubCellData subcell_data; + std::vector considered_vertices; + GridTools::delete_duplicated_vertices (vertices, cells, + subcell_data, + considered_vertices); + + // then clear the old triangulation and create the new one + result.clear (); + result.create_triangulation (vertices, cells, subcell_data); + } + + + void extrude_triangulation(const Triangulation<2, 2> &input, const unsigned int n_slices, diff --git a/source/grid/grid_generator.inst.in b/source/grid/grid_generator.inst.in index eb1e25278c..f17220b8d1 100644 --- a/source/grid/grid_generator.inst.in +++ b/source/grid/grid_generator.inst.in @@ -42,6 +42,12 @@ namespace GridGenerator (const Triangulation &triangulation_1, const Triangulation &triangulation_2, Triangulation &result); + + template + void + create_triangulation_with_removed_cells (const Triangulation &input_triangulation, + const std::set::active_cell_iterator> &cells_to_remove, + Triangulation &result); #endif \} diff --git a/tests/grid/grid_generator_07.cc b/tests/grid/grid_generator_07.cc new file mode 100644 index 0000000000..e6c0976755 --- /dev/null +++ b/tests/grid/grid_generator_07.cc @@ -0,0 +1,65 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2005 - 2015 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. +// +// --------------------------------------------------------------------- + + + +// Test GridGenerator::create_triangulation_with_removed_cells + +#include "../tests.h" +#include +#include +#include +#include +#include +#include + +#include +#include + + +template +void test(std::ostream &out) +{ + Triangulation triangulation; + Triangulation tr; + GridGenerator::hyper_cube (triangulation); + triangulation.refine_global (3); + + // remove all cells but the first. this is the hardest case to handle as it + // makes a bunch of vertices unused + std::set::active_cell_iterator> + cells_to_remove; + for (typename Triangulation::active_cell_iterator + cell = ++triangulation.begin_active(); + cell != triangulation.end(); ++cell) + cells_to_remove.insert (cell); + + GridGenerator::create_triangulation_with_removed_cells(triangulation, + cells_to_remove, tr); + GridOut go; + go.write_gnuplot(tr, out); +} + + +int main() +{ + std::ofstream logfile("output"); + deallog.attach(logfile); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + test<2>(logfile); + test<3>(logfile); +} diff --git a/tests/grid/grid_generator_07.output b/tests/grid/grid_generator_07.output new file mode 100644 index 0000000000..e4ce372ce6 --- /dev/null +++ b/tests/grid/grid_generator_07.output @@ -0,0 +1,32 @@ + +0.00000 0.00000 0 0 +0.125000 0.00000 0 0 +0.125000 0.125000 0 0 +0.00000 0.125000 0 0 +0.00000 0.00000 0 0 + + +0.00000 0.00000 0.00000 0 0 +0.125000 0.00000 0.00000 0 0 +0.125000 0.00000 0.125000 0 0 +0.00000 0.00000 0.125000 0 0 +0.00000 0.00000 0.00000 0 0 + +0.00000 0.125000 0.00000 0 0 +0.125000 0.125000 0.00000 0 0 +0.125000 0.125000 0.125000 0 0 +0.00000 0.125000 0.125000 0 0 +0.00000 0.125000 0.00000 0 0 + +0.00000 0.00000 0.00000 0 0 +0.00000 0.125000 0.00000 0 0 + +0.125000 0.00000 0.00000 0 0 +0.125000 0.125000 0.00000 0 0 + +0.125000 0.00000 0.125000 0 0 +0.125000 0.125000 0.125000 0 0 + +0.00000 0.00000 0.125000 0 0 +0.00000 0.125000 0.125000 0 0 + -- 2.39.5