From: Wolfgang Bangerth Date: Thu, 21 Aug 2014 11:24:06 +0000 (-0500) Subject: Move create_union_triangulation() and extract_boundary_mesh(). X-Git-Tag: v8.2.0-rc1~189^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=f99cae29f1775b5ceab07941c498e4f6f471fda1;p=dealii.git Move create_union_triangulation() and extract_boundary_mesh(). These functions were previously in GridTools, but conceptually they create meshes and so should be in GridGenerator. Move them there but leave the old functions in GridTools with a note that they are now deprecated. Also structure the documentation of the functions in GridGenerator better. --- diff --git a/doc/news/changes.h b/doc/news/changes.h index 859da6936a..aca712d743 100644 --- a/doc/news/changes.h +++ b/doc/news/changes.h @@ -277,6 +277,16 @@ inconvenience this causes.

Specific improvements

    +
  1. Changed: The functions GridTools::extract_boundary_mesh() and + GridTools::create_union_triangulation() have been moved to + GridGenerator::extract_boundary_mesh() and + GridGenerator::create_union_triangulation() since, conceptually, they + generate meshes. The old functions have been retained but are now + deprecated. +
    + (Wolfgang Bangerth, 2014/08/19) +
  2. +
  3. New: TriaAccessor::measure() is now also implemented for faces of 3d cells as long as the face is planar.
    diff --git a/include/deal.II/grid/grid_generator.h b/include/deal.II/grid/grid_generator.h index 14ea1a5924..4bc9a67bbd 100644 --- a/include/deal.II/grid/grid_generator.h +++ b/include/deal.II/grid/grid_generator.h @@ -53,6 +53,11 @@ template class SparseMatrix; */ namespace GridGenerator { + /** + * @name Creating meshes for basic geometries + */ + /*@{*/ + /** * Initialize the given triangulation with a hypercube (line in 1D, * square in 2D, etc) consisting of exactly one cell. The hypercube @@ -640,6 +645,16 @@ namespace GridGenerator const double R, const double r); + /** + * @} + */ + /** + * @name Creating meshes from other meshes + */ + /** + * @{ + */ + /** * Given the two triangulations specified as the first two arguments, create * the triangulation that contains the cells of both triangulation and store @@ -668,7 +683,7 @@ namespace GridGenerator * * @note For a related operation on refined meshes when both meshes are * derived from the same coarse mesh, see - * GridTools::create_union_triangulation . + * GridGenerator::create_union_triangulation(). */ template void @@ -676,6 +691,61 @@ namespace GridGenerator const Triangulation &triangulation_2, Triangulation &result); + /** + * Given the two triangulations + * specified as the first two + * arguments, create the + * triangulation that contains + * the finest cells of both + * triangulation and store it in + * the third parameter. Previous + * content of @p result will be + * deleted. + * + * @note This function is intended + * to create an adaptively refined + * triangulation that contains the + * most refined cells from + * two input triangulations that + * were derived from the same + * coarse grid by adaptive refinement. + * This is an operation sometimes + * needed when one solves for two + * variables of a coupled problem + * on separately refined meshes on + * the same domain (for example + * because these variables have + * boundary layers in different places) + * but then needs to compute something + * that involves both variables or + * wants to output the result into a + * single file. In both cases, in + * order not to lose information, + * the two solutions can not be + * interpolated onto the respectively + * other mesh because that may be + * coarser than the ones on which + * the variable was computed. Rather, + * one needs to have a mesh for the + * domain that is at least as fine + * as each of the two initial meshes. + * This function computes such a mesh. + * + * @note If you want to create + * a mesh that is the merger of + * two other coarse meshes, for + * example in order to compose a mesh + * for a complicated geometry from + * meshes for simpler geometries, + * then this is not the function for you. Instead, consider + * GridGenerator::merge_triangulations(). + */ + template + void + create_union_triangulation (const Triangulation &triangulation_1, + const Triangulation &triangulation_2, + Triangulation &result); + /** * Take a 2d Triangulation that is being extruded in z direction @@ -722,6 +792,134 @@ namespace GridGenerator void flatten_triangulation(const Triangulation &in_tria, Triangulation &out_tria); + /** + * @} + */ + + /** + * @name Creating lower-dimensional meshes from parts of higher-dimensional meshes + */ + /*@{*/ + + +#ifdef _MSC_VER + // Microsoft's VC++ has a bug where it doesn't want to recognize that + // an implementation (definition) of the extract_boundary_mesh function + // matches a declaration. This can apparently only be avoided by + // doing some contortion with the return type using the following + // intermediate type. This is only used when using MS VC++ and uses + // the direct way of doing it otherwise + template