From: Giovanni Alzetta Date: Tue, 29 Jan 2019 15:29:15 +0000 (+0100) Subject: Generalizing GridGenerator::general_cell to tria X-Git-Tag: v9.1.0-rc1~388^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=bfd22d07ce53528ad4bdc36069d29cf165169458;p=dealii.git Generalizing GridGenerator::general_cell to tria --- diff --git a/doc/news/changes/minor/20190130GiovanniAlzetta b/doc/news/changes/minor/20190130GiovanniAlzetta new file mode 100644 index 0000000000..a9b64415e3 --- /dev/null +++ b/doc/news/changes/minor/20190130GiovanniAlzetta @@ -0,0 +1,5 @@ +Changed: Function GridGenerator::general_cell() can now generate +a cell of dimension `dim` inside a space of dimension `spacedim` +with `dim <= spacedim` +
+(Giovanni Alzetta, 2019/01/30) diff --git a/include/deal.II/grid/grid_generator.h b/include/deal.II/grid/grid_generator.h index 354037d6f2..0a4d52818f 100644 --- a/include/deal.II/grid/grid_generator.h +++ b/include/deal.II/grid/grid_generator.h @@ -462,11 +462,13 @@ namespace GridGenerator const bool colorize = false); /** - * A general quadrilateral in 2d or a general hexahedron in 3d. It is the - * responsibility of the user to provide the vertices in the right order (see - * the documentation of the GeometryInfo class) because the vertices are - * stored in the same order as they are given. It is also important to make - * sure that the volume of the cell is positive. + * A general @p dim -dimensional cell (a segment if dim is 1, a quadrilateral + * if @p dim is 2, or a hexahedron if @p dim is 3) immersed in a + * @p spacedim -dimensional space. It is the responsibility of the user to + * provide the vertices in the right order (see the documentation of the + * GeometryInfo class) because the vertices are stored in the same order as + * they are given. It is also important to make sure that the volume of the + * cell is positive. * * If the argument @p colorize is false, then all boundary indicators are * set to zero for 2d and 3d. If it is true, the boundary is colorized as in @@ -474,13 +476,17 @@ namespace GridGenerator * @ref GlossColorization "the glossary entry on colorization"). In 1d the * indicators are always colorized, see hyper_rectangle(). * + * @param tria The triangulation that will be created + * @param vertices The 2^dim vertices of the cell + * @param colorize If true, set different boundary ids. + * * @author Bruno Turcksin */ - template + template void - general_cell(Triangulation & tria, - const std::vector> &vertices, - const bool colorize = false); + general_cell(Triangulation & tria, + const std::vector> &vertices, + const bool colorize = false); /** * A parallelogram. The first corner point is the origin. The @p dim diff --git a/source/grid/grid_generator.cc b/source/grid/grid_generator.cc index a51cbe29a3..9738d79372 100644 --- a/source/grid/grid_generator.cc +++ b/source/grid/grid_generator.cc @@ -749,11 +749,11 @@ namespace GridGenerator - template + template void - general_cell(Triangulation & tria, - const std::vector> &vertices, - const bool colorize) + general_cell(Triangulation & tria, + const std::vector> &vertices, + const bool colorize) { Assert(vertices.size() == dealii::GeometryInfo::vertices_per_cell, ExcMessage("Wrong number of vertices.")); @@ -761,7 +761,7 @@ namespace GridGenerator // First create a hyper_rectangle and then deform it. hyper_cube(tria, 0, 1, colorize); - typename Triangulation::active_cell_iterator cell = + typename Triangulation::active_cell_iterator cell = tria.begin_active(); for (unsigned int i = 0; i < GeometryInfo::vertices_per_cell; ++i) cell->vertex(i) = vertices[i]; diff --git a/source/grid/grid_generator.inst.in b/source/grid/grid_generator.inst.in index 0084a28512..2253edee98 100644 --- a/source/grid/grid_generator.inst.in +++ b/source/grid/grid_generator.inst.in @@ -106,6 +106,10 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS) &cells_to_remove, Triangulation &result); + template void + general_cell(Triangulation &, + const std::vector> &, + const bool); #endif \} } @@ -153,11 +157,6 @@ for (deal_II_dimension : DIMENSIONS) const Point &, const bool); - template void - general_cell(Triangulation &, - const std::vector> &, - const bool); - template void simplex( Triangulation &, diff --git a/tests/grid/grid_generator_general_cell.cc b/tests/grid/grid_generator_general_cell_01.cc similarity index 100% rename from tests/grid/grid_generator_general_cell.cc rename to tests/grid/grid_generator_general_cell_01.cc diff --git a/tests/grid/grid_generator_general_cell.output b/tests/grid/grid_generator_general_cell_01.output similarity index 100% rename from tests/grid/grid_generator_general_cell.output rename to tests/grid/grid_generator_general_cell_01.output diff --git a/tests/grid/grid_generator_general_cell_02.cc b/tests/grid/grid_generator_general_cell_02.cc new file mode 100644 index 0000000000..cd60e7aec1 --- /dev/null +++ b/tests/grid/grid_generator_general_cell_02.cc @@ -0,0 +1,95 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2016 - 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. +// +// --------------------------------------------------------------------- + +// Test output for GridGenerator::general_cell() for dim != spacedim + +#include +#include +#include + +#include "../tests.h" + +void +dim_2_3(std::ostream &os) +{ + std::vector> vertices(4); + vertices[0](0) = -1.; + vertices[0](1) = -1.; + vertices[0](0) = 0.; + + vertices[1](0) = 1.; + vertices[1](1) = -1.5; + vertices[1](2) = 0.; + + vertices[2](0) = 1.5; + vertices[2](1) = 1.5; + vertices[2](2) = 0.; + + vertices[3](0) = 2.; + vertices[3](1) = 0.5; + vertices[3](2) = 0.; + + Triangulation<2, 3> tria; + GridGenerator::general_cell<2, 3>(tria, vertices); + + GridOut gout; + gout.write_vtk(tria, os); +} + +void +dim_1_3(std::ostream &os) +{ + std::vector> vertices(2); + vertices[0](0) = -1.; + vertices[0](1) = -1.; + vertices[0](2) = -1.; + vertices[1](0) = 1.; + vertices[1](1) = -1.5; + vertices[1](2) = -1.5; + + Triangulation<1, 3> tria; + GridGenerator::general_cell<1, 3>(tria, vertices); + + GridOut gout; + gout.write_vtk(tria, os); +} + +void +dim_1_2(std::ostream &os) +{ + std::vector> vertices(2); + vertices[0](0) = -1.; + vertices[0](1) = -1.; + + vertices[1](0) = 1.; + vertices[1](1) = -1.5; + + Triangulation<1, 2> tria; + GridGenerator::general_cell<1, 2>(tria, vertices); + + GridOut gout; + gout.write_vtk(tria, os); +} + + +int +main() +{ + initlog(true); + std::ostream &logfile = deallog.get_file_stream(); + dim_2_3(logfile); + dim_1_3(logfile); + dim_1_2(logfile); +} diff --git a/tests/grid/grid_generator_general_cell_02.output b/tests/grid/grid_generator_general_cell_02.output new file mode 100644 index 0000000000..8a44d13862 --- /dev/null +++ b/tests/grid/grid_generator_general_cell_02.output @@ -0,0 +1,87 @@ + +# vtk DataFile Version 3.0 +Triangulation generated with deal.II +ASCII +DATASET UNSTRUCTURED_GRID +POINTS 4 double +0.00000 -1.00000 0.00000 +1.00000 -1.50000 0.00000 +1.50000 1.50000 0.00000 +2.00000 0.500000 0.00000 + +CELLS 1 5 +4 0 1 3 2 + +CELL_TYPES 1 +9 + + + +CELL_DATA 1 +SCALARS MaterialID int 1 +LOOKUP_TABLE default +0 + + + +SCALARS ManifoldID int 1 +LOOKUP_TABLE default +-1 + + +# vtk DataFile Version 3.0 +Triangulation generated with deal.II +ASCII +DATASET UNSTRUCTURED_GRID +POINTS 2 double +-1.00000 -1.00000 -1.00000 +1.00000 -1.50000 -1.50000 + +CELLS 1 3 +2 0 1 + +CELL_TYPES 1 +3 + + + +CELL_DATA 1 +SCALARS MaterialID int 1 +LOOKUP_TABLE default +0 + + + +SCALARS ManifoldID int 1 +LOOKUP_TABLE default +-1 + + +# vtk DataFile Version 3.0 +Triangulation generated with deal.II +ASCII +DATASET UNSTRUCTURED_GRID +POINTS 2 double +-1.00000 -1.00000 0 +1.00000 -1.50000 0 + +CELLS 1 3 +2 0 1 + +CELL_TYPES 1 +3 + + + +CELL_DATA 1 +SCALARS MaterialID int 1 +LOOKUP_TABLE default +0 + + + +SCALARS ManifoldID int 1 +LOOKUP_TABLE default +-1 + +